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 
41 namespace Nektar
42 {
43 namespace SolverUtils
44 {
45 class DiffusionIP : public Diffusion
46 {
47 public:
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 
56  /// Calculate the average of conservative variables on traces
57  template <class T, typename = typename std::enable_if<
58  std::is_floating_point<T>::value ||
60  inline void ConsVarAve(const size_t nConvectiveFields, const T &Bweight,
61  const std::vector<T> &vFwd,
62  const std::vector<T> &vBwd, std::vector<T> &aver)
63  {
64  constexpr int nvelst = 1;
65  int nEngy = nConvectiveFields - 1;
66  int nveled = nEngy;
67 
68  T Fweight = 1.0 - Bweight;
69  for (size_t v = 0; v < nEngy; ++v)
70  {
71  aver[v] = Fweight * vFwd[v] + Bweight * vBwd[v];
72  }
73 
74  T LinternalEngy{};
75  T RinternalEngy{};
76  T AinternalEngy{};
77  for (size_t j = nvelst; j < nveled; ++j)
78  {
79  LinternalEngy += vFwd[j] * vFwd[j];
80  RinternalEngy += vBwd[j] * vBwd[j];
81  AinternalEngy += aver[j] * aver[j];
82  }
83  LinternalEngy *= -0.5 / vFwd[0];
84  RinternalEngy *= -0.5 / vBwd[0];
85  LinternalEngy += vFwd[nEngy];
86  RinternalEngy += vBwd[nEngy];
87 
88  aver[nEngy] = Fweight * LinternalEngy + Bweight * RinternalEngy;
89  aver[nEngy] += AinternalEngy * (0.5 / aver[0]);
90  }
91 
92 protected:
93  DiffusionIP();
94 
98 
102  /// Workspace for v_Diffusion
104  /// Workspace for CallTraceNumFlux
107 
113 
114  /// Get IP penalty factor based on order
115  void GetPenaltyFactor(
117  Array<OneD, NekDouble> &factor);
118  /// Get a constant IP penalty factor
121  Array<OneD, NekDouble> &factor);
122 
123  /// Add symmetric flux integration to field (in coefficient space)
125  const std::size_t nConvectiveFields, const size_t nDim,
126  const size_t nPts, const size_t nTracePts,
128  const Array<OneD, const int> &nonZeroIndex,
129  TensorOfArray3D<NekDouble> &tracflux,
130  Array<OneD, Array<OneD, NekDouble>> &outarray);
131 
132  /// Add symmetric flux integration to field (in physical space)
134  const std::size_t nConvectiveFields, const size_t nDim,
135  const size_t nPts, const size_t nTracePts,
137  const Array<OneD, const int> &nonZeroIndex,
138  TensorOfArray3D<NekDouble> &tracflux,
139  Array<OneD, Array<OneD, NekDouble>> &outarray);
140 
141  /// Calculate symmetric flux on traces
142  void CalcTraceSymFlux(
143  const std::size_t nConvectiveFields, const size_t nDim,
144  const Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
145  Array<OneD, Array<OneD, NekDouble>> &solution_jump,
146  Array<OneD, int> &nonZeroIndexsymm,
147  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &traceSymflux);
148 
149  /// Calculate symmetric flux on traces interface
151  const std::size_t nConvectiveFields,
153  const Array<OneD, Array<OneD, NekDouble>> &inarray,
155  TensorOfArray3D<NekDouble> &VolumeFlux,
156  TensorOfArray3D<NekDouble> &SymmFlux,
157  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
158  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
159  Array<OneD, int> &nonZeroIndex);
160 
161  void v_InitObject(
164 
165  virtual void v_Diffuse(
166  const std::size_t nConvectiveFields,
168  const Array<OneD, Array<OneD, NekDouble>> &inarray,
169  Array<OneD, Array<OneD, NekDouble>> &outarray,
170  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
171  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
172 
173  virtual void v_DiffuseCoeffs(
174  const std::size_t nConvectiveFields,
176  const Array<OneD, Array<OneD, NekDouble>> &inarray,
177  Array<OneD, Array<OneD, NekDouble>> &outarray,
178  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
179  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
180 
181  virtual void v_DiffuseCoeffs(
182  const std::size_t nConvectiveFields,
184  const Array<OneD, Array<OneD, NekDouble>> &inarray,
185  Array<OneD, Array<OneD, NekDouble>> &outarray,
186  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
187  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
189  Array<OneD, int> &nonZeroIndex) override;
190 
191  virtual void v_DiffuseVolumeFlux(
193  const Array<OneD, Array<OneD, NekDouble>> &inarray,
195  TensorOfArray3D<NekDouble> &VolumeFlux,
196  Array<OneD, int> &nonZeroIndex) override;
197 
198  virtual void v_DiffuseTraceFlux(
200  const Array<OneD, Array<OneD, NekDouble>> &inarray,
202  TensorOfArray3D<NekDouble> &VolumeFlux,
203  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
204  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
205  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
206  Array<OneD, int> &nonZeroIndex) override;
207 
209  const int nConvectiveFields,
211  const Array<OneD, Array<OneD, NekDouble>> &inarray,
212  const Array<OneD, const Array<OneD, Array<OneD, NekDouble>>> &qfield,
213  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
214  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
215  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
216  const Array<OneD, const Array<OneD, Array<OneD, NekDouble>>> &qFwd,
217  const Array<OneD, const Array<OneD, Array<OneD, NekDouble>>> &qBwd,
218  const Array<OneD, NekDouble> &MuAVTrace, Array<OneD, int> &nonZeroIndex,
219  const Array<OneD, Array<OneD, NekDouble>> &Aver,
220  const Array<OneD, Array<OneD, NekDouble>> &Jump);
221 
222  virtual void v_AddDiffusionSymmFluxToCoeff(
223  const std::size_t nConvectiveFields,
225  const Array<OneD, Array<OneD, NekDouble>> &inarray,
227  TensorOfArray3D<NekDouble> &VolumeFlux,
228  Array<OneD, Array<OneD, NekDouble>> &outarray,
229  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
230  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
231 
232  virtual void v_AddDiffusionSymmFluxToPhys(
233  const std::size_t nConvectiveFields,
235  const Array<OneD, Array<OneD, NekDouble>> &inarray,
237  TensorOfArray3D<NekDouble> &VolumeFlux,
238  Array<OneD, Array<OneD, NekDouble>> &outarray,
239  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
240  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
241 
242  virtual void v_DiffuseCalcDerivative(
244  const Array<OneD, Array<OneD, NekDouble>> &inarray,
246  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
247  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
248 
249  virtual void v_ConsVarAveJump(
250  const std::size_t nConvectiveFields, const size_t npnts,
251  const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
252  const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
254  Array<OneD, Array<OneD, NekDouble>> &jump) override;
255 
257  {
258  return m_traceNormals;
259  }
260 
261  // /// Calculate numerical flux on traces
262  // void CalcTraceNumFlux(
263  // const std::size_t nConvectiveFields, const size_t nDim,
264  // const size_t nPts, const size_t nTracePts,
265  // const NekDouble PenaltyFactor2,
266  // const Array<OneD, MultiRegions::ExpListSharedPtr> &fields,
267  // const Array<OneD, Array<OneD, NekDouble>> &inarray,
268  // const TensorOfArray3D<NekDouble> &qfield,
269  // const Array<OneD, Array<OneD, NekDouble>> &vFwd,
270  // const Array<OneD, Array<OneD, NekDouble>> &vBwd,
271  // const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &qFwd,
272  // const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &qBwd,
273  // const Array<OneD, NekDouble> &MuVarTrace,
274  // Array<OneD, int> &nonZeroIndexflux,
275  // TensorOfArray3D<NekDouble> &traceflux,
276  // Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
277  // Array<OneD, Array<OneD, NekDouble>> &solution_jump);
278 
279  /// Calculate numerical flux on traces
280  void CalcTraceNumFlux(
281  const NekDouble PenaltyFactor2,
283  const Array<OneD, Array<OneD, NekDouble>> &inarray,
284  const TensorOfArray3D<NekDouble> &qfield,
285  const Array<OneD, Array<OneD, NekDouble>> &vFwd,
286  const Array<OneD, Array<OneD, NekDouble>> &vBwd,
287  Array<OneD, int> &nonZeroIndexflux,
288  TensorOfArray3D<NekDouble> &traceflux,
289  Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
290  Array<OneD, Array<OneD, NekDouble>> &solution_jump);
291 
292  /// Add second derivative term to trace jump (for DDG scheme)
294  const std::size_t nConvectiveFields, const size_t nDim,
295  const size_t nPts, const size_t nTracePts,
296  const NekDouble PenaltyFactor2,
298  const TensorOfArray3D<NekDouble> &qfield,
299  TensorOfArray3D<NekDouble> &numDerivFwd,
300  TensorOfArray3D<NekDouble> &numDerivBwd);
301 
302  void ApplyFluxBndConds(
303  const int nConvectiveFields,
306 };
307 } // namespace SolverUtils
308 } // namespace Nektar
309 
310 #endif
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
Definition: DiffusionIP.h:111
void v_DiffuseTraceFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &qfield, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &qFwd, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &qBwd, const Array< OneD, NekDouble > &MuAVTrace, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &Aver, const Array< OneD, Array< OneD, NekDouble >> &Jump)
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionIP.h:48
Array< OneD, NekDouble > m_tracBwdWeightAver
Definition: DiffusionIP.h:108
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:110
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:60
void GetPenaltyFactorConst(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get a constant IP penalty factor.
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)
virtual void v_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) override
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 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)
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Definition: DiffusionIP.h:101
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
Definition: DiffusionIP.h:105
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 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.
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.
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
Definition: DiffusionIP.h:106
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionIP.h:112
virtual void v_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) override
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:100
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
Array< OneD, NekDouble > m_tracBwdWeightJump
Definition: DiffusionIP.h:109
const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal() override
Definition: DiffusionIP.h:256
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
virtual void v_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) override
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Definition: DiffusionIP.cpp:54
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)
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
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
Workspace for v_Diffusion.
Definition: DiffusionIP.h:103
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionIP.h:99
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.
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.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:550
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble