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::SolverUtils
42{
43class DiffusionIP : public Diffusion
44{
45public:
46 static DiffusionSharedPtr create([[maybe_unused]] std::string diffType)
47 {
48 return DiffusionSharedPtr(new DiffusionIP());
49 }
50
51 static std::string type;
52
53protected:
55
56 void v_InitObject(
59
60 void v_Diffuse(const std::size_t nConvectiveFields,
62 const Array<OneD, Array<OneD, NekDouble>> &inarray,
64 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
65 const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
66
67 void v_DiffuseCoeffs(
68 const std::size_t nConvectiveFields,
70 const Array<OneD, Array<OneD, NekDouble>> &inarray,
72 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
73 const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
74
77 const Array<OneD, Array<OneD, NekDouble>> &inarray,
80 Array<OneD, int> &nonZeroIndex) override;
81
84 const Array<OneD, Array<OneD, NekDouble>> &inarray,
88 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
89 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
90 Array<OneD, int> &nonZeroIndex) override;
91
94 const Array<OneD, Array<OneD, NekDouble>> &inarray,
96 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
97 const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
98
100 {
101 return m_traceNormals;
102 }
103
104private:
108
112 /// Workspace for v_Diffusion
114 /// Workspace for CallTraceNumFlux
117
123
124 /// Calculate the average of conservative variables on traces
125 template <class T, typename = typename std::enable_if<
126 std::is_floating_point<T>::value ||
128 inline void ConsVarAve(const size_t nConvectiveFields, const T &Bweight,
129 const std::vector<T> &vFwd,
130 const std::vector<T> &vBwd, std::vector<T> &aver)
131 {
132 constexpr int nvelst = 1;
133 int nEngy = nConvectiveFields - 1;
134 int nveled = nEngy;
135
136 T Fweight = 1.0 - Bweight;
137 for (size_t v = 0; v < nEngy; ++v)
138 {
139 aver[v] = Fweight * vFwd[v] + Bweight * vBwd[v];
140 }
141
142 T LinternalEngy{};
143 T RinternalEngy{};
144 T AinternalEngy{};
145 for (size_t j = nvelst; j < nveled; ++j)
146 {
147 LinternalEngy += vFwd[j] * vFwd[j];
148 RinternalEngy += vBwd[j] * vBwd[j];
149 AinternalEngy += aver[j] * aver[j];
150 }
151 LinternalEngy *= -0.5 / vFwd[0];
152 RinternalEngy *= -0.5 / vBwd[0];
153 LinternalEngy += vFwd[nEngy];
154 RinternalEngy += vBwd[nEngy];
155
156 aver[nEngy] = Fweight * LinternalEngy + Bweight * RinternalEngy;
157 aver[nEngy] += AinternalEngy * (0.5 / aver[0]);
158 }
159
160 /// Get IP penalty factor based on order
161 void GetPenaltyFactor(
163 Array<OneD, NekDouble> &factor);
164
165 void ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts,
166 const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
167 const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
170
171 /// Add symmetric flux integration to field (in coefficient space)
173 const std::size_t nConvectiveFields, const size_t nDim,
174 const size_t nPts, const size_t nTracePts,
176 const Array<OneD, const int> &nonZeroIndex,
178 Array<OneD, Array<OneD, NekDouble>> &outarray);
179
180 /// Add symmetric flux integration to field (in physical space)
182 const std::size_t nConvectiveFields, const size_t nDim,
183 const size_t nPts, const size_t nTracePts,
185 const Array<OneD, const int> &nonZeroIndex,
187 Array<OneD, Array<OneD, NekDouble>> &outarray);
188
189 /// Calculate symmetric flux on traces
190 void CalcTraceSymFlux(
191 const std::size_t nConvectiveFields, const size_t nDim,
192 const Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
193 Array<OneD, Array<OneD, NekDouble>> &solution_jump,
194 Array<OneD, int> &nonZeroIndexsymm,
195 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &traceSymflux);
196
198 const std::size_t nConvectiveFields,
200 const Array<OneD, Array<OneD, NekDouble>> &inarray,
202 TensorOfArray3D<NekDouble> &VolumeFlux,
204 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
205 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
206
208 const std::size_t nConvectiveFields,
210 const Array<OneD, Array<OneD, NekDouble>> &inarray,
212 TensorOfArray3D<NekDouble> &VolumeFlux,
214 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
215 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
216
217 /// Calculate symmetric flux on traces interface
219 const std::size_t nConvectiveFields,
221 const Array<OneD, Array<OneD, NekDouble>> &inarray,
223 TensorOfArray3D<NekDouble> &VolumeFlux,
225 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
226 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
227 Array<OneD, int> &nonZeroIndex);
228
229 /// Calculate numerical flux on traces
230 void CalcTraceNumFlux(
231 const NekDouble PenaltyFactor2,
233 const Array<OneD, Array<OneD, NekDouble>> &inarray,
234 const TensorOfArray3D<NekDouble> &qfield,
235 const Array<OneD, Array<OneD, NekDouble>> &vFwd,
236 const Array<OneD, Array<OneD, NekDouble>> &vBwd,
237 Array<OneD, int> &nonZeroIndexflux,
239 Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
240 Array<OneD, Array<OneD, NekDouble>> &solution_jump);
241
242 /// Add second derivative term to trace jump (for DDG scheme)
244 const std::size_t nConvectiveFields, const size_t nDim,
245 const size_t nPts, const size_t nTracePts,
246 const NekDouble PenaltyFactor2,
248 const TensorOfArray3D<NekDouble> &qfield,
249 TensorOfArray3D<NekDouble> &numDerivFwd,
250 TensorOfArray3D<NekDouble> &numDerivBwd);
251
253 const int nConvectiveFields,
256};
257} // namespace Nektar::SolverUtils
258
259#endif
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
Definition: DiffusionIP.h:121
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:46
Array< OneD, NekDouble > m_tracBwdWeightAver
Definition: DiffusionIP.h:118
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:120
const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal() override
Definition: DiffusionIP.h:99
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:128
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.
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:111
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:115
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
Definition: DiffusionIP.h:116
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionIP.h:122
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:110
Array< OneD, NekDouble > m_tracBwdWeightJump
Definition: DiffusionIP.h:119
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)
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Definition: DiffusionIP.cpp:52
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:113
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionIP.h:109
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:55
double NekDouble