Nektar++
DiffusionIP.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: DiffusionIP.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: IP diffusion class.
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_SOLVERUTILS_DIFFUSIONWEAKDG
37#define NEKTAR_SOLVERUTILS_DIFFUSIONWEAKDG
38
40
41namespace Nektar
42{
43namespace SolverUtils
44{
45class DiffusionIP : public Diffusion
46{
47public:
48 static DiffusionSharedPtr create(std::string diffType)
49 {
50 boost::ignore_unused(diffType);
51 return DiffusionSharedPtr(new DiffusionIP());
52 }
53
54 static std::string type;
55
56protected:
58
59 virtual void v_InitObject(
62
63 virtual void v_Diffuse(
64 const std::size_t nConvectiveFields,
66 const Array<OneD, Array<OneD, NekDouble>> &inarray,
68 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
69 const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
70
71 virtual void v_DiffuseCoeffs(
72 const std::size_t nConvectiveFields,
74 const Array<OneD, Array<OneD, NekDouble>> &inarray,
76 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
77 const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
78
79 virtual void v_DiffuseVolumeFlux(
81 const Array<OneD, Array<OneD, NekDouble>> &inarray,
84 Array<OneD, int> &nonZeroIndex) override;
85
86 virtual void v_DiffuseTraceFlux(
88 const Array<OneD, Array<OneD, NekDouble>> &inarray,
92 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
93 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
94 Array<OneD, int> &nonZeroIndex) override;
95
96 virtual void v_DiffuseCalcDerivative(
98 const Array<OneD, Array<OneD, NekDouble>> &inarray,
100 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
101 const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
102
104 override
105 {
106 return m_traceNormals;
107 }
108
109private:
113
117 /// Workspace for v_Diffusion
119 /// Workspace for CallTraceNumFlux
122
128
129 /// Calculate the average of conservative variables on traces
130 template <class T, typename = typename std::enable_if<
131 std::is_floating_point<T>::value ||
133 inline void ConsVarAve(const size_t nConvectiveFields, const T &Bweight,
134 const std::vector<T> &vFwd,
135 const std::vector<T> &vBwd, std::vector<T> &aver)
136 {
137 constexpr int nvelst = 1;
138 int nEngy = nConvectiveFields - 1;
139 int nveled = nEngy;
140
141 T Fweight = 1.0 - Bweight;
142 for (size_t v = 0; v < nEngy; ++v)
143 {
144 aver[v] = Fweight * vFwd[v] + Bweight * vBwd[v];
145 }
146
147 T LinternalEngy{};
148 T RinternalEngy{};
149 T AinternalEngy{};
150 for (size_t j = nvelst; j < nveled; ++j)
151 {
152 LinternalEngy += vFwd[j] * vFwd[j];
153 RinternalEngy += vBwd[j] * vBwd[j];
154 AinternalEngy += aver[j] * aver[j];
155 }
156 LinternalEngy *= -0.5 / vFwd[0];
157 RinternalEngy *= -0.5 / vBwd[0];
158 LinternalEngy += vFwd[nEngy];
159 RinternalEngy += vBwd[nEngy];
160
161 aver[nEngy] = Fweight * LinternalEngy + Bweight * RinternalEngy;
162 aver[nEngy] += AinternalEngy * (0.5 / aver[0]);
163 }
164
165 /// Get IP penalty factor based on order
166 void GetPenaltyFactor(
168 Array<OneD, NekDouble> &factor);
169
170 void ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts,
171 const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
172 const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
175
176 /// Add symmetric flux integration to field (in coefficient space)
178 const std::size_t nConvectiveFields, const size_t nDim,
179 const size_t nPts, const size_t nTracePts,
181 const Array<OneD, const int> &nonZeroIndex,
183 Array<OneD, Array<OneD, NekDouble>> &outarray);
184
185 /// Add symmetric flux integration to field (in physical space)
187 const std::size_t nConvectiveFields, const size_t nDim,
188 const size_t nPts, const size_t nTracePts,
190 const Array<OneD, const int> &nonZeroIndex,
192 Array<OneD, Array<OneD, NekDouble>> &outarray);
193
194 /// Calculate symmetric flux on traces
195 void CalcTraceSymFlux(
196 const std::size_t nConvectiveFields, const size_t nDim,
197 const Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
198 Array<OneD, Array<OneD, NekDouble>> &solution_jump,
199 Array<OneD, int> &nonZeroIndexsymm,
200 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &traceSymflux);
201
203 const std::size_t nConvectiveFields,
205 const Array<OneD, Array<OneD, NekDouble>> &inarray,
207 TensorOfArray3D<NekDouble> &VolumeFlux,
209 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
210 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
211
213 const std::size_t nConvectiveFields,
215 const Array<OneD, Array<OneD, NekDouble>> &inarray,
217 TensorOfArray3D<NekDouble> &VolumeFlux,
219 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
220 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
221
222 /// Calculate symmetric flux on traces interface
224 const std::size_t nConvectiveFields,
226 const Array<OneD, Array<OneD, NekDouble>> &inarray,
228 TensorOfArray3D<NekDouble> &VolumeFlux,
230 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
231 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
232 Array<OneD, int> &nonZeroIndex);
233
234 /// Calculate numerical flux on traces
235 void CalcTraceNumFlux(
236 const NekDouble PenaltyFactor2,
238 const Array<OneD, Array<OneD, NekDouble>> &inarray,
239 const TensorOfArray3D<NekDouble> &qfield,
240 const Array<OneD, Array<OneD, NekDouble>> &vFwd,
241 const Array<OneD, Array<OneD, NekDouble>> &vBwd,
242 Array<OneD, int> &nonZeroIndexflux,
244 Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
245 Array<OneD, Array<OneD, NekDouble>> &solution_jump);
246
247 /// Add second derivative term to trace jump (for DDG scheme)
249 const std::size_t nConvectiveFields, const size_t nDim,
250 const size_t nPts, const size_t nTracePts,
251 const NekDouble PenaltyFactor2,
253 const TensorOfArray3D<NekDouble> &qfield,
254 TensorOfArray3D<NekDouble> &numDerivFwd,
255 TensorOfArray3D<NekDouble> &numDerivBwd);
256
258 const int nConvectiveFields,
261};
262} // namespace SolverUtils
263} // namespace Nektar
264
265#endif
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
Definition: DiffusionIP.h:126
virtual void v_DiffuseCoeffs(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionIP.h:48
Array< OneD, NekDouble > m_tracBwdWeightAver
Definition: DiffusionIP.h:123
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:125
void ApplyFluxBndConds(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, Array< OneD, NekDouble > > &flux)
aplly Neuman boundary conditions on flux Currently only consider WallAdiabatic
void GetPenaltyFactor(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get IP penalty factor based on order.
void ConsVarAve(const size_t nConvectiveFields, const T &Bweight, const std::vector< T > &vFwd, const std::vector< T > &vBwd, std::vector< T > &aver)
Calculate the average of conservative variables on traces.
Definition: DiffusionIP.h:133
virtual const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal() override
Definition: DiffusionIP.h:103
void AddSecondDerivToTrace(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &numDerivFwd, TensorOfArray3D< NekDouble > &numDerivBwd)
Add second derivative term to trace jump (for DDG scheme)
void CalcTraceNumFlux(const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble > > &vFwd, const Array< OneD, Array< OneD, NekDouble > > &vBwd, Array< OneD, int > &nonZeroIndexflux, TensorOfArray3D< NekDouble > &traceflux, Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump)
Calculate numerical flux on traces.
virtual void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
Diffusion Flux, calculate the physical derivatives.
void AddDiffusionSymmFluxToPhys(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Definition: DiffusionIP.h:116
void AddSymmFluxIntegralToCoeff(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
Add symmetric flux integration to field (in coefficient space)
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
Definition: DiffusionIP.h:120
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
Definition: DiffusionIP.h:121
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionIP.h:127
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:115
Array< OneD, NekDouble > m_tracBwdWeightJump
Definition: DiffusionIP.h:124
virtual void v_Diffuse(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
void ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble > > &vFwd, const Array< OneD, const Array< OneD, NekDouble > > &vBwd, Array< OneD, Array< OneD, NekDouble > > &aver, Array< OneD, Array< OneD, NekDouble > > &jump)
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Definition: DiffusionIP.cpp:54
virtual void v_DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex) override
Diffusion Volume Flux.
void AddDiffusionSymmFluxToCoeff(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
void DiffuseTraceSymmFlux(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex)
Calculate symmetric flux on traces interface.
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
Workspace for v_Diffusion.
Definition: DiffusionIP.h:118
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionIP.h:114
virtual void v_DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex) override
Diffusion term Trace Flux.
void CalcTraceSymFlux(const std::size_t nConvectiveFields, const size_t nDim, const Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump, Array< OneD, int > &nonZeroIndexsymm, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &traceSymflux)
Calculate symmetric flux on traces.
void AddSymmFluxIntegralToPhys(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
Add symmetric flux integration to field (in physical space)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:58
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble