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  <
59  std::is_floating_point<T>::value ||
61  >::type
62  >
63  inline void ConsVarAve(
64  const size_t nConvectiveFields,
65  const T& Bweight,
66  const std::vector<T>& vFwd,
67  const std::vector<T>& vBwd,
68  std::vector<T>& aver)
69  {
70  constexpr int nvelst = 1;
71  int nEngy = nConvectiveFields - 1;
72  int nveled = nEngy;
73 
74  T Fweight = 1.0 - Bweight;
75  for (size_t v = 0; v < nEngy; ++v)
76  {
77  aver[v] = Fweight * vFwd[v] + Bweight * vBwd[v];
78  }
79 
80  T LinternalEngy{};
81  T RinternalEngy{};
82  T AinternalEngy{};
83  for (size_t j = nvelst; j < nveled; ++j)
84  {
85  LinternalEngy += vFwd[j] * vFwd[j];
86  RinternalEngy += vBwd[j] * vBwd[j];
87  AinternalEngy += aver[j] * aver[j];
88  }
89  LinternalEngy *= -0.5 / vFwd[0];
90  RinternalEngy *= -0.5 / vBwd[0];
91  LinternalEngy += vFwd[nEngy];
92  RinternalEngy += vBwd[nEngy];
93 
94  aver[nEngy] = Fweight * LinternalEngy + Bweight * RinternalEngy;
95  aver[nEngy] += AinternalEngy * (0.5 / aver[0]);
96 
97  }
98 
99 
100 protected:
101  DiffusionIP();
102 
103  std::string m_shockCaptureType;
107 
112  /// Workspace for v_Diffusion
114  /// Workspace for CallTraceNumFlux
117 
123 
124  /// Get IP penalty factor based on order
125  void GetPenaltyFactor(
127  Array<OneD, NekDouble> &factor);
128  /// Get a constant IP penalty factor
131  Array<OneD, NekDouble> &factor);
132 
133  /// Add symmetric flux integration to field (in coefficient space)
135  const std::size_t nConvectiveFields, const size_t nDim,
136  const size_t nPts, const size_t nTracePts,
138  const Array<OneD, const int> &nonZeroIndex,
139  TensorOfArray3D<NekDouble> &tracflux,
140  Array<OneD, Array<OneD, NekDouble>> &outarray);
141 
142  /// Add symmetric flux integration to field (in physical space)
144  const std::size_t nConvectiveFields, const size_t nDim,
145  const size_t nPts, const size_t nTracePts,
147  const Array<OneD, const int> &nonZeroIndex,
148  TensorOfArray3D<NekDouble> &tracflux,
149  Array<OneD, Array<OneD, NekDouble>> &outarray);
150 
151  /// Calculate symmetric flux on traces
152  void CalcTraceSymFlux(
153  const std::size_t nConvectiveFields, const size_t nDim,
155  const Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
156  Array<OneD, Array<OneD, NekDouble>> &solution_jump,
157  Array<OneD, int> &nonZeroIndexsymm,
158  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &traceSymflux);
159 
160  /// Calculate symmetric flux on traces interface
162  const std::size_t nConvectiveFields,
164  const Array<OneD, Array<OneD, NekDouble>> &inarray,
166  TensorOfArray3D<NekDouble> &VolumeFlux,
167  TensorOfArray3D<NekDouble> &SymmFlux,
168  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
169  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
170  Array<OneD, int> &nonZeroIndex);
171 
172  virtual void v_InitObject(
175 
176  virtual void v_Diffuse(
177  const std::size_t nConvectiveFields,
179  const Array<OneD, Array<OneD, NekDouble>> &inarray,
180  Array<OneD, Array<OneD, NekDouble>> &outarray,
181  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
182  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
183 
184  virtual void v_DiffuseCoeffs(
185  const std::size_t nConvectiveFields,
187  const Array<OneD, Array<OneD, NekDouble>> &inarray,
188  Array<OneD, Array<OneD, NekDouble>> &outarray,
189  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
190  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
191 
192  virtual void v_DiffuseCoeffs(
193  const std::size_t nConvectiveFields,
195  const Array<OneD, Array<OneD, NekDouble>> &inarray,
196  Array<OneD, Array<OneD, NekDouble>> &outarray,
197  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
198  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
200  Array< OneD, int > &nonZeroIndex);
201 
202  virtual void v_DiffuseVolumeFlux(
204  const Array<OneD, Array<OneD, NekDouble>> &inarray,
206  TensorOfArray3D<NekDouble> &VolumeFlux,
207  Array<OneD, int> &nonZeroIndex);
208 
209  virtual void v_DiffuseTraceFlux(
211  const Array<OneD, Array<OneD, NekDouble>> &inarray,
213  TensorOfArray3D<NekDouble> &VolumeFlux,
214  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
215  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
216  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
217  Array<OneD, int> &nonZeroIndex);
218 
220  const int nConvectiveFields,
222  const Array<OneD, Array<OneD, NekDouble>> &inarray,
223  const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &qfield,
224  Array<OneD, Array<OneD, NekDouble> > &TraceFlux,
225  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
226  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
227  const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &qFwd,
228  const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &qBwd,
229  const Array<OneD, NekDouble> &MuAVTrace,
230  Array< OneD, int > &nonZeroIndex,
231  const Array<OneD, Array<OneD, NekDouble>> &Aver,
232  const Array<OneD, Array<OneD, NekDouble>> &Jump);
233 
234  virtual void v_AddDiffusionSymmFluxToCoeff(
235  const std::size_t nConvectiveFields,
237  const Array<OneD, Array<OneD, NekDouble>> &inarray,
239  TensorOfArray3D<NekDouble> &VolumeFlux,
240  Array<OneD, Array<OneD, NekDouble>> &outarray,
241  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
242  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
243 
244  virtual void v_AddDiffusionSymmFluxToPhys(
245  const std::size_t nConvectiveFields,
247  const Array<OneD, Array<OneD, NekDouble>> &inarray,
249  TensorOfArray3D<NekDouble> &VolumeFlux,
250  Array<OneD, Array<OneD, NekDouble>> &outarray,
251  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
252  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
253 
254  virtual void v_DiffuseCalcDerivative(
256  const Array<OneD, Array<OneD, NekDouble>> &inarray,
258  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
259  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
260 
261  virtual void v_ConsVarAveJump(
262  const std::size_t nConvectiveFields, const size_t npnts,
263  const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
264  const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
267 
269  {
270  return m_traceNormals;
271  }
272 
273  // /// Calculate numerical flux on traces
274  // void CalcTraceNumFlux(
275  // const std::size_t nConvectiveFields, const size_t nDim,
276  // const size_t nPts, const size_t nTracePts,
277  // const NekDouble PenaltyFactor2,
278  // const Array<OneD, MultiRegions::ExpListSharedPtr> &fields,
279  // const Array<OneD, Array<OneD, NekDouble>> &inarray,
280  // const TensorOfArray3D<NekDouble> &qfield,
281  // const Array<OneD, Array<OneD, NekDouble>> &vFwd,
282  // const Array<OneD, Array<OneD, NekDouble>> &vBwd,
283  // const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &qFwd,
284  // const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &qBwd,
285  // const Array<OneD, NekDouble> &MuVarTrace,
286  // Array<OneD, int> &nonZeroIndexflux,
287  // TensorOfArray3D<NekDouble> &traceflux,
288  // Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
289  // Array<OneD, Array<OneD, NekDouble>> &solution_jump);
290 
291  /// Calculate numerical flux on traces
292  void CalcTraceNumFlux(
293  const NekDouble PenaltyFactor2,
295  const Array<OneD, Array<OneD, NekDouble>> &inarray,
296  const TensorOfArray3D<NekDouble> &qfield,
297  const Array<OneD, Array<OneD, NekDouble>> &vFwd,
298  const Array<OneD, Array<OneD, NekDouble>> &vBwd,
299  const Array<OneD, NekDouble> &MuVarTrace,
300  Array<OneD, int> &nonZeroIndexflux,
301  TensorOfArray3D<NekDouble> &traceflux,
302  Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
303  Array<OneD, Array<OneD, NekDouble>> &solution_jump);
304 
305  /// Add second derivative term to trace jump (for DDG scheme)
307  const std::size_t nConvectiveFields, const size_t nDim,
308  const size_t nPts, const size_t nTracePts,
309  const NekDouble PenaltyFactor2,
311  const TensorOfArray3D<NekDouble> &qfield,
312  TensorOfArray3D<NekDouble> &numDerivFwd,
313  TensorOfArray3D<NekDouble> &numDerivBwd);
314 
315  void ApplyFluxBndConds(
316  const int nConvectiveFields,
319 };
320 } // namespace SolverUtils
321 } // namespace Nektar
322 
323 #endif
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
Definition: DiffusionIP.h:121
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionIP.h:48
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)
Array< OneD, NekDouble > m_tracBwdWeightAver
Definition: DiffusionIP.h:118
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:120
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)
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:63
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)
Diffusion Volume Flux.
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)
Array< OneD, NekDouble > m_MuVarTrace
Definition: DiffusionIP.h:108
virtual const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: DiffusionIP.h:268
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: DiffusionIP.cpp:53
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:111
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)
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
Definition: DiffusionIP.h:115
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)
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.
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
Definition: DiffusionIP.h:116
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionIP.h:122
void CalcTraceSymFlux(const std::size_t nConvectiveFields, const size_t nDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, 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.
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:110
Array< OneD, NekDouble > m_tracBwdWeightJump
Definition: DiffusionIP.h:119
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)
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)
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)
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, const Array< OneD, NekDouble > &MuVarTrace, 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 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_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)
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
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:627
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble