Nektar++
Diffusion.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: Diffusion.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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Abstract base class for diffusion.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_SOLVERUTILS_DIFFUSION
36 #define NEKTAR_SOLVERUTILS_DIFFUSION
37 
38 #include <functional>
39 #include <string>
40 
41 #include <boost/core/ignore_unused.hpp>
42 
46 #include <MultiRegions/ExpList.h>
49 
50 namespace Nektar
51 {
52 namespace SolverUtils
53 {
54 typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
55  const TensorOfArray3D<NekDouble> &,
56  TensorOfArray3D<NekDouble> &)>
58 
59 typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
63 
64 typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
68 
69 /**
70  * Parameter list meaning:
71  * 1st: field conservative variables
72  * 2th: Devrivatives of field conservative varialbes
73  * 3rd: the current time for time-dependent boundary
74  * 4th: Fwd of field conservative variables optional
75  * 5th: Fwd of Devrivatives(2nd) optional
76  *
77  * a null pointer need to be passed for optional parameters
78  */
79 typedef std::function<void(
80  const Array<OneD, const Array<OneD, NekDouble>> &,
82  const Array<OneD, const Array<OneD, NekDouble>> &,
83  const Array<OneD, const Array<OneD, Array<OneD, NekDouble>>> &)>
85 
86 /**
87  * Parameter list meaning:
88  * 1st: nvariables
89  * 2nd: nspaceDimension
90  * 3rd: field conservative variables
91  * 4th: Devrivatives of field conservative varialbes
92  * 5th: nonzero flux index array, optional
93  * 6th: normal vectors optional
94  *
95  * a null pointer need to be passed for optional parameters
96  */
97 typedef std::function<void(
98  const int, const Array<OneD, Array<OneD, NekDouble>> &,
102 
103 /**
104  * Parameter list meaning:
105  * 1st: nvariables
106  * 2nd: nspaceDimension
107  * 3rd: trace conservative variables for Diffusion Flux Jacobian
108  * 4th: trace conservative variables( usually the jump of trace value)
109  * 5th: trace symmetric flux
110  * 6th: nonzero flux index array, optional
111  *
112  * a null pointer need to be passed for optional parameters
113  */
114 typedef std::function<void(
115  const int, const Array<OneD, Array<OneD, NekDouble>> &,
119 
120 /**
121  * Parameter list meaning:
122  * 1rd: trace conservative variables
123  */
124 typedef std::function<void(Array<OneD, Array<OneD, NekDouble>> &)>
126 
128 {
129 public:
130  /// Params for Ducros sensor
134 
136 
140 
142  const std::size_t nConvectiveFields,
144  const Array<OneD, Array<OneD, NekDouble>> &inarray,
145  Array<OneD, Array<OneD, NekDouble>> &outarray,
146  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
148  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
150 
152  const std::size_t nConvectiveFields,
154  const Array<OneD, Array<OneD, NekDouble>> &inarray,
155  Array<OneD, Array<OneD, NekDouble>> &outarray,
156  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
158  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
160 
162  const std::size_t nConvectiveFields,
164  const Array<OneD, Array<OneD, NekDouble>> &inarray,
165  Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
166  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
168  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
170 
172  const std::size_t nConvectiveFields,
174  const Array<OneD, Array<OneD, NekDouble>> &inarray,
175  Array<OneD, Array<OneD, NekDouble>> &outarray,
176  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
177  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
178  TensorOfArray3D<NekDouble> &qfield, Array<OneD, int> &nonZeroIndex)
179  {
180  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
181  pBwd, qfield, nonZeroIndex);
182  }
183 
185  const std::size_t nConvectiveFields,
187  const Array<OneD, Array<OneD, NekDouble>> &inarray,
188  Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
189  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
191  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
193 
194  // Diffusion Calculate the physical derivatives
197  const Array<OneD, Array<OneD, NekDouble>> &inarray,
199  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
201  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
203 
204  /// Diffusion Volume FLux
205 
208  const Array<OneD, Array<OneD, NekDouble>> &inarray,
210  TensorOfArray3D<NekDouble> &VolumeFlux,
211  Array<OneD, int> &nonZeroIndex = NullInt1DArray)
212  {
213  v_DiffuseVolumeFlux(fields, inarray, qfields, VolumeFlux, nonZeroIndex);
214  }
215 
216  /// Diffusion term Trace Flux
219  const Array<OneD, Array<OneD, NekDouble>> &inarray,
221  TensorOfArray3D<NekDouble> &VolumeFlux,
222  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
223  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
225  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
227  Array<OneD, int> &nonZeroIndex = NullInt1DArray)
228  {
229  v_DiffuseTraceFlux(fields, inarray, qfields, VolumeFlux, TraceFlux,
230  pFwd, pBwd, nonZeroIndex);
231  }
232 
234  const int nConvectiveFields,
236  const Array<OneD, Array<OneD, NekDouble>> &inarray,
237  const TensorOfArray3D<NekDouble> &qfield,
238  const TensorOfArray3D<NekDouble> &VolumeFlux,
239  TensorOfArray3D<NekDouble> &SymmFlux,
240  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
241  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
242  Array<OneD, int> &nonZeroIndex,
243  Array<OneD, Array<OneD, NekDouble>> &solution_Aver =
245  Array<OneD, Array<OneD, NekDouble>> &solution_jump =
247  {
248  v_DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
249  VolumeFlux, SymmFlux, pFwd, pBwd, nonZeroIndex,
250  solution_Aver, solution_jump);
251  }
252 
253  /// Add symmetric flux to field in coeff space
255  const std::size_t nConvectiveFields,
257  const Array<OneD, Array<OneD, NekDouble>> &inarray,
259  TensorOfArray3D<NekDouble> &VolumeFlux,
260  Array<OneD, Array<OneD, NekDouble>> &outarray,
261  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
262  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
263  {
264  v_AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields, inarray,
265  qfield, VolumeFlux, outarray, pFwd, pBwd);
266  }
267 
268  /// Add symmetric flux to field in coeff physical space
270  const std::size_t nConvectiveFields,
272  const Array<OneD, Array<OneD, NekDouble>> &inarray,
274  TensorOfArray3D<NekDouble> &VolumeFlux,
275  Array<OneD, Array<OneD, NekDouble>> &outarray,
276  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
277  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
278  {
279  v_AddDiffusionSymmFluxToPhys(nConvectiveFields, fields, inarray, qfield,
280  VolumeFlux, outarray, pFwd, pBwd);
281  }
282 
284  const int nConvectiveFields, const int nDim, const int nPts,
285  const int nTracePts,
287  const Array<OneD, const int> &nonZeroIndex,
290  Array<OneD, Array<OneD, NekDouble>> &outarray);
291 
292  template <typename FuncPointerT, typename ObjectPointerT>
293  void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
294  {
295  m_fluxVector = std::bind(func, obj, std::placeholders::_1,
296  std::placeholders::_2, std::placeholders::_3);
297  }
298 
300  {
301  m_fluxVector = fluxVector;
302  }
303 
304  template <typename FuncPointerT, typename ObjectPointerT>
305  void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
306  {
308  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
309  std::placeholders::_3);
310  }
311 
313  {
314  m_fluxVectorNS = fluxVector;
315  }
316 
317  template <typename FuncPointerT, typename ObjectPointerT>
318  void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
319  {
321  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
322  std::placeholders::_3);
323  }
324 
326  {
327  m_fluxPenaltyNS = flux;
328  }
329 
330  template <typename FuncPointerT, typename ObjectPointerT>
331  void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
332  {
334  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
335  std::placeholders::_3, std::placeholders::_4,
336  std::placeholders::_5, std::placeholders::_6);
337  }
338 
340  {
342  }
343 
344  template <typename FuncPointerT, typename ObjectPointerT>
345  void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
346  {
348  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
349  std::placeholders::_3, std::placeholders::_4,
350  std::placeholders::_5, std::placeholders::_6);
351  }
352 
354  {
356  }
357 
358  template <typename FuncPointerT, typename ObjectPointerT>
359  void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
360  {
361  m_SpecialBndTreat = std::bind(func, obj, std::placeholders::_1);
362  }
363 
364  template <typename FuncPointerT, typename ObjectPointerT>
365  void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
366  {
368  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
369  std::placeholders::_3, std::placeholders::_4,
370  std::placeholders::_5, std::placeholders::_6);
371  }
372 
374  {
375  v_SetHomoDerivs(deriv);
376  }
377 
379  {
380  return v_GetFluxTensor();
381  }
382 
383  /// Get the average and jump value of conservative variables on trace
385  const std::size_t nConvectiveFields, const size_t npnts,
386  const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
387  const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
390  {
391  v_ConsVarAveJump(nConvectiveFields, npnts, vFwd, vBwd, aver, jump);
392  }
393 
394  /// Get trace normal
397  {
398  return v_GetTraceNormal();
399  }
400 
401 protected:
409 
411 
415  {
416  boost::ignore_unused(pSession, pFields);
417  };
418 
419  SOLVER_UTILS_EXPORT virtual void v_Diffuse(
420  const std::size_t nConvectiveFields,
422  const Array<OneD, Array<OneD, NekDouble>> &inarray,
423  Array<OneD, Array<OneD, NekDouble>> &outarray,
424  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
425  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
426 
428  const std::size_t nConvectiveFields,
430  const Array<OneD, Array<OneD, NekDouble>> &inarray,
431  Array<OneD, Array<OneD, NekDouble>> &outarray,
432  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
433  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
434 
436  const std::size_t nConvectiveFields,
438  const Array<OneD, Array<OneD, NekDouble>> &inarray,
439  Array<OneD, Array<OneD, NekDouble>> &outarray,
440  const Array<OneD, Array<OneD, NekDouble>> &vFwd,
441  const Array<OneD, Array<OneD, NekDouble>> &vBwd,
442  TensorOfArray3D<NekDouble> &qfield, Array<OneD, int> &nonZeroIndex);
443 
445  const std::size_t nConvectiveFields, const size_t npnts,
446  const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
447  const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
450 
451  /// Diffusion Flux, calculate the physical derivatives
454  const Array<OneD, Array<OneD, NekDouble>> &inarray,
456  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
457  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
458 
459  /// Diffusion Volume Flux
462  const Array<OneD, Array<OneD, NekDouble>> &inarray,
464  TensorOfArray3D<NekDouble> &VolumeFlux, Array<OneD, int> &nonZeroIndex);
465 
466  /// Diffusion term Trace Flux
469  const Array<OneD, Array<OneD, NekDouble>> &inarray,
471  TensorOfArray3D<NekDouble> &VolumeFlux,
472  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
473  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
474  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
475  Array<OneD, int> &nonZeroIndex);
476 
478  const int nConvectiveFields,
480  const Array<OneD, Array<OneD, NekDouble>> &inarray,
481  const TensorOfArray3D<NekDouble> &qfield,
482  const TensorOfArray3D<NekDouble> &VolumeFlux,
483  TensorOfArray3D<NekDouble> &SymmFlux,
484  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
485  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
486  Array<OneD, int> &nonZeroIndex,
487  Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
488  Array<OneD, Array<OneD, NekDouble>> &solution_jump)
489  {
490  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
491  VolumeFlux, SymmFlux, pFwd, pBwd, nonZeroIndex,
492  solution_Aver, solution_jump);
493  ASSERTL0(false, "Not defined function v_DiffuseTraceSymmFlux.");
494  }
495 
497  const std::size_t nConvectiveFields,
499  const Array<OneD, Array<OneD, NekDouble>> &inarray,
501  TensorOfArray3D<NekDouble> &VolumeFlux,
502  Array<OneD, Array<OneD, NekDouble>> &outarray,
503  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
504  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
505  {
506  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
507  VolumeFlux, outarray, pFwd, pBwd);
508  }
510  const std::size_t nConvectiveFields,
512  const Array<OneD, Array<OneD, NekDouble>> &inarray,
514  TensorOfArray3D<NekDouble> &VolumeFlux,
515  Array<OneD, Array<OneD, NekDouble>> &outarray,
516  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
517  const Array<OneD, Array<OneD, NekDouble>> &pBwd)
518  {
519  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
520  VolumeFlux, outarray, pFwd, pBwd);
521  }
522 
524  {
525  boost::ignore_unused(deriv);
526  }
527 
529  {
530  static TensorOfArray3D<NekDouble> tmp;
531  return tmp;
532  }
533 
535  &v_GetTraceNormal();
536 
537  /// Compute primary derivatives
538  SOLVER_UTILS_EXPORT virtual void v_GetPrimVar(
540  const Array<OneD, Array<OneD, NekDouble>> &inarray,
541  Array<OneD, Array<OneD, NekDouble>> &primVar);
542 
543  /// Compute divergence and curl squared
546  const TensorOfArray3D<NekDouble> &pVarDer);
547 };
548 
549 /// A shared pointer to an EquationSystem object
550 typedef std::shared_ptr<SolverUtils::Diffusion> DiffusionSharedPtr;
551 
552 /// Datatype of the NekFactory used to instantiate classes derived
553 /// from the Diffusion class.
557 } // namespace SolverUtils
558 } // namespace Nektar
559 
560 #endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:105
SOLVER_UTILS_EXPORT 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)
Add symmetric flux to field in coeff physical space.
Definition: Diffusion.h:269
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:132
virtual SOLVER_UTILS_EXPORT 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)
Definition: Diffusion.cpp:210
void SetFluxVectorNS(DiffusionFluxVecCBNS fluxVector)
Definition: Diffusion.h:312
virtual SOLVER_UTILS_EXPORT void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
Diffusion Flux, calculate the physical derivatives.
Definition: Diffusion.cpp:221
void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:359
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion Volume FLux.
Definition: Diffusion.h:206
SOLVER_UTILS_EXPORT void 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=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Definition: Diffusion.cpp:59
SOLVER_UTILS_EXPORT 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)
Add symmetric flux to field in coeff space.
Definition: Diffusion.h:254
void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:318
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:406
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceSymmFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble >> &solution_Aver, Array< OneD, Array< OneD, NekDouble >> &solution_jump)
Definition: Diffusion.h:477
SOLVER_UTILS_EXPORT void DiffuseTraceSymmFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble >> &solution_Aver=NullNekDoubleArrayOfArray, Array< OneD, Array< OneD, NekDouble >> &solution_jump=NullNekDoubleArrayOfArray)
Definition: Diffusion.h:233
SOLVER_UTILS_EXPORT void 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=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Similar with Diffusion::Diffuse(): calculate diffusion flux The difference is in the outarray: it is ...
Definition: Diffusion.cpp:76
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Definition: Diffusion.cpp:194
virtual SOLVER_UTILS_EXPORT 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.
Definition: Diffusion.cpp:232
virtual SOLVER_UTILS_EXPORT 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)
Definition: Diffusion.cpp:166
virtual TensorOfArray3D< NekDouble > & GetFluxTensor()
Definition: Diffusion.h:378
void SetFluxPenaltyNS(DiffusionFluxPenaltyNS flux)
Definition: Diffusion.h:325
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:131
SOLVER_UTILS_EXPORT void GetDivCurl(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &pVarDer)
Compute divergence and curl squared.
Definition: Diffusion.cpp:281
SOLVER_UTILS_EXPORT void 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=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion term Trace Flux.
Definition: Diffusion.h:217
virtual SOLVER_UTILS_EXPORT ~Diffusion()
Definition: Diffusion.h:135
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:133
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:403
void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:305
void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:345
void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:293
void SetDiffusionFluxCons(DiffusionFluxCons flux)
Definition: Diffusion.h:339
SOLVER_UTILS_EXPORT void InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.cpp:47
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor()
Definition: Diffusion.h:528
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:405
void SetHomoDerivs(Array< OneD, Array< OneD, NekDouble >> &deriv)
Definition: Diffusion.h:373
virtual SOLVER_UTILS_EXPORT void v_GetPrimVar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &primVar)
Compute primary derivatives.
Definition: Diffusion.cpp:260
void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:365
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:407
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal()
Get trace normal.
Definition: Diffusion.h:396
void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:331
virtual void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble >> &deriv)
Definition: Diffusion.h:523
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:408
SOLVER_UTILS_EXPORT void AddSymmFluxIntegralToOffDiag(const int nConvectiveFields, const int nDim, const int nPts, const int nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &Fwdflux, TensorOfArray3D< NekDouble > &Bwdflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: Diffusion.cpp:111
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:404
SOLVER_UTILS_EXPORT void 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, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
Definition: Diffusion.h:171
virtual SOLVER_UTILS_EXPORT 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)
Definition: Diffusion.h:496
virtual SOLVER_UTILS_EXPORT 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)
Diffusion term Trace Flux.
Definition: Diffusion.cpp:242
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:204
void SetDiffusionFluxConsTrace(DiffusionFluxCons flux)
Definition: Diffusion.h:353
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:402
virtual SOLVER_UTILS_EXPORT 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)
Definition: Diffusion.h:509
virtual SOLVER_UTILS_EXPORT 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)
Definition: Diffusion.cpp:153
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.h:412
SOLVER_UTILS_EXPORT void SetFluxVector(DiffusionFluxVecCB fluxVector)
Definition: Diffusion.h:299
SOLVER_UTILS_EXPORT 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)
Get the average and jump value of conservative variables on trace.
Definition: Diffusion.h:384
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::function< void(Array< OneD, Array< OneD, NekDouble >> &)> SpecialBndTreat
Definition: Diffusion.h:125
std::function< void(const int, const Array< OneD, Array< OneD, NekDouble >> &, const TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &, Array< OneD, int > &, const Array< OneD, Array< OneD, NekDouble >> &)> DiffusionFluxCons
Definition: Diffusion.h:101
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
std::function< void(const Array< OneD, Array< OneD, NekDouble >> &, const Array< OneD, Array< OneD, NekDouble >> &, Array< OneD, Array< OneD, NekDouble >> &)> DiffusionFluxPenaltyNS
Definition: Diffusion.h:67
LibUtilities::NekFactory< std::string, Diffusion, std::string > DiffusionFactory
Datatype of the NekFactory used to instantiate classes derived from the Diffusion class.
Definition: Diffusion.h:555
std::function< void(const Array< OneD, const Array< OneD, NekDouble >> &, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &, NekDouble, const Array< OneD, const Array< OneD, NekDouble >> &, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &)> FunctorDerivBndCond
Definition: Diffusion.h:84
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:550
std::function< void(const Array< OneD, Array< OneD, NekDouble >> &, const TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCB
Definition: Diffusion.h:57
std::function< void(const Array< OneD, Array< OneD, NekDouble >> &, TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCBNS
Definition: Diffusion.h:62
std::function< void(const int, const Array< OneD, Array< OneD, NekDouble >> &, const Array< OneD, Array< OneD, NekDouble >> &, TensorOfArray3D< NekDouble > &, Array< OneD, int > &, const Array< OneD, Array< OneD, NekDouble >> &)> DiffusionSymmFluxCons
Definition: Diffusion.h:118
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static Array< OneD, int > NullInt1DArray
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble