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 <string>
39 #include <functional>
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 (
55  const Array<OneD, Array<OneD, NekDouble> > &,
56  const TensorOfArray3D<NekDouble> &,
57  TensorOfArray3D<NekDouble> &)>
59 
60  typedef std::function<void (
65 
66  typedef std::function<void (
71 
72  typedef std::function<void (
76 
77  /**
78  * Parameter list meaning:
79  * 1st: field conservative variables
80  * 2th: Devrivatives of field conservative varialbes
81  * 3rd: the current time for time-dependent boundary
82  * 4th: Fwd of field conservative variables optional
83  * 5th: Fwd of Devrivatives(2nd) optional
84  *
85  * a null pointer need to be passed for optional parameters
86  */
87  typedef std::function<void (
88  const Array<OneD, const Array<OneD, NekDouble> > &,
89  const Array<OneD, const Array<OneD, Array<OneD, NekDouble> > > &,
90  NekDouble ,
91  const Array<OneD, const Array<OneD, NekDouble> > &,
93 
94  /**
95  * Parameter list meaning:
96  * 1st: nvariables
97  * 2nd: nspaceDimension
98  * 3rd: field conservative variables
99  * 4th: Devrivatives of field conservative varialbes
100  * 5th: nonzero flux index array, optional
101  * 6th: normal vectors optional
102  * 7th: aritificial diffusion facotres optional
103  *
104  * a null pointer need to be passed for optional parameters
105  */
106  typedef std::function<void (
107  const int ,
108  const Array<OneD, Array<OneD, NekDouble> > &,
112  const Array<OneD, Array<OneD, NekDouble> > &,
113  const Array<OneD, NekDouble> &)>
115 
116  /**
117  * Parameter list meaning:
118  * 1st: nvariables
119  * 2nd: nspaceDimension
120  * 3rd: trace conservative variables for Diffusion Flux Jacobian
121  * 4th: trace conservative variables( usually the jump of trace value)
122  * 5th: trace symmetric flux
123  * 6th: nonzero flux index array, optional
124  * 7th: normal vectors optional
125  *
126  * a null pointer need to be passed for optional parameters
127  */
128  typedef std::function<void (
129  const int ,
130  const Array<OneD, Array<OneD, NekDouble> > &,
134  const Array<OneD, Array<OneD, NekDouble> > &)>
136 
137  /**
138  * Parameter list meaning:
139  * 1rd: trace conservative variables
140  */
141  typedef std::function<void (
144 
145  /**
146  * Parameter list meaning:
147  * 1st: trace conservative variables
148  * 2rd: dynamic viscosity
149  */
150  typedef std::function<void (
151  const Array<OneD, const Array<OneD, NekDouble> > &,
154 
155  class Diffusion
156  {
157  public:
158  ///Params for Ducros sensor
162 
164  {};
165 
169 
171  const std::size_t nConvectiveFields,
173  const Array<OneD, Array<OneD, NekDouble> > &inarray,
174  Array<OneD, Array<OneD, NekDouble> > &outarray,
177 
179  const std::size_t nConvectiveFields,
181  const Array<OneD, Array<OneD, NekDouble> > &inarray,
182  Array<OneD, Array<OneD, NekDouble> > &outarray,
186  &pBwd = NullNekDoubleArrayOfArray);
187 
189  const std::size_t nConvectiveFields,
191  const Array<OneD, Array<OneD, NekDouble> > &inarray,
192  Array<OneD, Array<OneD, NekDouble> > &outarray,
193  NekDouble time,
197  &pBwd = NullNekDoubleArrayOfArray);
198 
200  const std::size_t nConvectiveFields,
202  const Array<OneD, Array<OneD, NekDouble> > &inarray,
203  Array<OneD, Array<OneD, NekDouble> > &outarray,
204  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
205  const Array<OneD, Array<OneD, NekDouble> > &pBwd,
207  Array< OneD, int > &nonZeroIndex)
208  {
209  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray,
210  pFwd, pBwd,qfield,nonZeroIndex);
211  }
212 
214  const std::size_t nConvectiveFields,
216  const Array<OneD, Array<OneD, NekDouble> > &inarray,
217  Array<OneD, Array<OneD, NekDouble> > &outarray,
218  NekDouble time,
222  &pBwd = NullNekDoubleArrayOfArray);
223 
224  // Diffusion Calculate the physical derivatives
227  const Array<OneD, Array<OneD, NekDouble>> &inarray,
228  TensorOfArray3D<NekDouble> &qfields,
229  const Array<OneD, Array<OneD, NekDouble>> &pFwd
231  const Array<OneD, Array<OneD, NekDouble>> &pBwd
233 
234  /// Diffusion Volume FLux
235 
238  const Array<OneD, Array<OneD, NekDouble>> &inarray,
240  TensorOfArray3D<NekDouble> &VolumeFlux,
242  &nonZeroIndex = NullInt1DArray)
243  {
244  v_DiffuseVolumeFlux(fields, inarray, qfields,
245  VolumeFlux, nonZeroIndex);
246  }
247 
248  /// Diffusion term Trace Flux
251  const Array<OneD, Array<OneD, NekDouble>> &inarray,
253  TensorOfArray3D<NekDouble> &VolumeFlux,
254  Array<OneD, Array<OneD, NekDouble> > &TraceFlux,
260  &nonZeroIndex = NullInt1DArray)
261  {
262  v_DiffuseTraceFlux(fields, inarray, qfields,
263  VolumeFlux, TraceFlux, pFwd, pBwd, nonZeroIndex);
264  }
265 
267  const int nConvectiveFields,
269  const Array<OneD, Array<OneD, NekDouble>> &inarray,
270  const TensorOfArray3D<NekDouble> &qfield,
271  const TensorOfArray3D<NekDouble> &VolumeFlux,
272  TensorOfArray3D<NekDouble> &SymmFlux,
273  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
274  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
275  Array< OneD, int > &nonZeroIndex,
278  {
279  v_DiffuseTraceSymmFlux(nConvectiveFields,fields,inarray,qfield,
280  VolumeFlux,SymmFlux,pFwd,pBwd,
281  nonZeroIndex,solution_Aver,solution_jump);
282  }
283 
284  /// Add symmetric flux to field in coeff space
286  const std::size_t nConvectiveFields,
288  const Array<OneD, Array<OneD, NekDouble> > &inarray,
290  TensorOfArray3D<NekDouble> &VolumeFlux,
291  Array<OneD, Array<OneD, NekDouble> > &outarray,
292  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
293  const Array<OneD, Array<OneD, NekDouble> > &pBwd)
294  {
295  v_AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields,
296  inarray, qfield, VolumeFlux, outarray, pFwd, pBwd);
297  }
298 
299  /// Add symmetric flux to field in coeff physical space
301  const std::size_t nConvectiveFields,
303  const Array<OneD, Array<OneD, NekDouble> > &inarray,
305  TensorOfArray3D<NekDouble> &VolumeFlux,
306  Array<OneD, Array<OneD, NekDouble> > &outarray,
307  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
308  const Array<OneD, Array<OneD, NekDouble> > &pBwd)
309  {
310  v_AddDiffusionSymmFluxToPhys(nConvectiveFields, fields, inarray,
311  qfield, VolumeFlux, outarray, pFwd, pBwd);
312  }
313 
315  const int nConvectiveFields,
316  const int nDim,
317  const int nPts,
318  const int nTracePts,
320  const Array<OneD, const int > &nonZeroIndex,
323  Array<OneD, Array<OneD, NekDouble> > &outarray);
324 
326  TensorOfArray3D<NekDouble> &fluxvector);
327 
328  template<typename FuncPointerT, typename ObjectPointerT>
329  void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
330  {
331  m_fluxVector = std::bind(
332  func, obj, std::placeholders::_1, std::placeholders::_2,
333  std::placeholders::_3);
334  }
335 
337  DiffusionFluxVecCB fluxVector)
338  {
339  m_fluxVector = fluxVector;
340  }
341 
342  template<typename FuncPointerT, typename ObjectPointerT>
343  void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
344  {
345  m_fluxVectorNS = std::bind(
346  func, obj, std::placeholders::_1, std::placeholders::_2,
347  std::placeholders::_3);
348  }
349 
351  {
352  m_fluxVectorNS = fluxVector;
353  }
354 
355  template<typename FuncPointerT, typename ObjectPointerT>
356  void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
357  {
358  m_fluxPenaltyNS = std::bind(
359  func, obj, std::placeholders::_1, std::placeholders::_2,
360  std::placeholders::_3);
361  }
362 
364  {
365  m_fluxPenaltyNS = flux;
366  }
367 
368  template<typename FuncPointerT, typename ObjectPointerT>
369  void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
370  {
371  m_FunctorDiffusionfluxCons = std::bind(
372  func, obj, std::placeholders::_1, std::placeholders::_2,
373  std::placeholders::_3, std::placeholders::_4,
374  std::placeholders::_5, std::placeholders::_6,
375  std::placeholders::_7);
376  }
377 
379  {
381  }
382 
383  template<typename FuncPointerT, typename ObjectPointerT>
384  void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
385  {
387  func, obj, std::placeholders::_1, std::placeholders::_2,
388  std::placeholders::_3, std::placeholders::_4,
389  std::placeholders::_5, std::placeholders::_6,
390  std::placeholders::_7);
391  }
392 
394  {
396  }
397 
398  template<typename FuncPointerT, typename ObjectPointerT>
400  FuncPointerT func, ObjectPointerT obj)
401  {
402  m_ArtificialDiffusionVector = std::bind(
403  func, obj, std::placeholders::_1, std::placeholders::_2);
404  }
405 
406  template<typename FuncPointerT, typename ObjectPointerT>
407  void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
408  {
409  m_SpecialBndTreat = std::bind(
410  func, obj, std::placeholders::_1);
411  }
412 
413  template<typename FuncPointerT, typename ObjectPointerT>
414  void SetCalcViscosity(FuncPointerT func, ObjectPointerT obj)
415  {
416  m_CalcViscosity = std::bind(
417  func, obj, std::placeholders::_1, std::placeholders::_2);
418  }
419 
420  template<typename FuncPointerT, typename ObjectPointerT>
421  void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
422  {
423  m_FunctorSymmetricfluxCons = std::bind(
424  func, obj, std::placeholders::_1, std::placeholders::_2,
425  std::placeholders::_3, std::placeholders::_4,
426  std::placeholders::_5, std::placeholders::_6);
427  }
428 
430  {
431  v_SetHomoDerivs(deriv);
432  }
433 
435  {
436  return v_GetFluxTensor();
437  }
438  /// Get the mu of artifical viscosity(AV)
441  const Array<OneD, Array<OneD, NekDouble> > &inarray,
443  Array<OneD, NekDouble > &MuVarTrace);
444 
445  /// Get the average and jump value of conservative variables on trace
447  const std::size_t nConvectiveFields,
448  const size_t npnts,
449  const Array<OneD, const Array<OneD, NekDouble> > &vFwd,
450  const Array<OneD, const Array<OneD, NekDouble> > &vBwd,
453  {
454  v_ConsVarAveJump(nConvectiveFields, npnts, vFwd, vBwd, aver,
455  jump);
456  }
457 
458  /// Get trace normal
461  {
462  return v_GetTraceNormal();
463  }
464 
465  protected:
475 
477 
481  {
482  boost::ignore_unused(pSession, pFields);
483  };
484 
485  SOLVER_UTILS_EXPORT virtual void v_Diffuse(
486  const std::size_t nConvectiveFields,
488  const Array<OneD, Array<OneD, NekDouble> > &inarray,
489  Array<OneD, Array<OneD, NekDouble> > &outarray,
490  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
491  const Array<OneD, Array<OneD, NekDouble> > &pBwd);
492 
494  const std::size_t nConvectiveFields,
496  const Array<OneD, Array<OneD, NekDouble> > &inarray,
497  Array<OneD, Array<OneD, NekDouble> > &outarray,
498  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
499  const Array<OneD, Array<OneD, NekDouble> > &pBwd);
500 
502  const std::size_t nConvectiveFields,
504  const Array<OneD, Array<OneD, NekDouble> > &inarray,
505  Array<OneD, Array<OneD, NekDouble> > &outarray,
506  const Array<OneD, Array<OneD, NekDouble> > &vFwd,
507  const Array<OneD, Array<OneD, NekDouble> > &vBwd,
509  Array< OneD, int > &nonZeroIndex);
510 
511 
513  const std::size_t nConvectiveFields,
514  const size_t npnts,
515  const Array<OneD, const Array<OneD, NekDouble> > &vFwd,
516  const Array<OneD, const Array<OneD, NekDouble> > &vBwd,
519 
520  /// Diffusion Flux, calculate the physical derivatives
523  const Array<OneD, Array<OneD, NekDouble> > &inarray,
525  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
526  const Array<OneD, Array<OneD, NekDouble> > &pBwd);
527 
528  /// Diffusion Volume Flux
531  const Array<OneD, Array<OneD, NekDouble>> &inarray,
533  TensorOfArray3D<NekDouble> &VolumeFlux,
535  &nonZeroIndex) ;
536 
537  /// Diffusion term Trace Flux
540  const Array<OneD, Array<OneD, NekDouble>> &inarray,
542  TensorOfArray3D<NekDouble> &VolumeFlux,
543  Array<OneD, Array<OneD, NekDouble> > &TraceFlux,
544  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
545  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
547  &nonZeroIndex);
548 
550  const int nConvectiveFields,
552  const Array<OneD, Array<OneD, NekDouble>> &inarray,
553  const TensorOfArray3D<NekDouble> &qfield,
554  const TensorOfArray3D<NekDouble> &VolumeFlux,
555  TensorOfArray3D<NekDouble> &SymmFlux,
556  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
557  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
558  Array< OneD, int > &nonZeroIndex,
559  Array<OneD, Array<OneD, NekDouble> > &solution_Aver,
560  Array<OneD, Array<OneD, NekDouble> > &solution_jump)
561  {
562  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
563  VolumeFlux, SymmFlux, pFwd, pBwd,
564  nonZeroIndex,solution_Aver,solution_jump);
565  ASSERTL0(false, "Not defined function v_DiffuseTraceSymmFlux.");
566  }
567 
569  const std::size_t nConvectiveFields,
571  const Array<OneD, Array<OneD, NekDouble> > &inarray,
573  TensorOfArray3D<NekDouble> &VolumeFlux,
574  Array<OneD, Array<OneD, NekDouble> > &outarray,
575  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
576  const Array<OneD, Array<OneD, NekDouble> > &pBwd)
577  {
578  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
579  VolumeFlux, outarray, pFwd, pBwd);
580 
581  }
583  const std::size_t nConvectiveFields,
585  const Array<OneD, Array<OneD, NekDouble> > &inarray,
587  TensorOfArray3D<NekDouble> &VolumeFlux,
588  Array<OneD, Array<OneD, NekDouble> > &outarray,
589  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
590  const Array<OneD, Array<OneD, NekDouble> > &pBwd)
591  {
592  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
593  VolumeFlux, outarray, pFwd, pBwd);
594 
595  }
596 
597  virtual void v_SetHomoDerivs(
599  {
600  boost::ignore_unused(deriv);
601  }
602 
604  {
605  static TensorOfArray3D<NekDouble> tmp;
606  return tmp;
607  }
608 
609  SOLVER_UTILS_EXPORT virtual
611  &v_GetTraceNormal();
612 
613  /// Compute primary derivatives
614  SOLVER_UTILS_EXPORT virtual void v_GetPrimVar(
616  const Array<OneD, Array<OneD, NekDouble>> &inarray,
617  Array<OneD, Array<OneD, NekDouble>> &primVar);
618 
619  /// Compute divergence and curl squared
622  const TensorOfArray3D<NekDouble> &pVarDer);
623 
624  };
625 
626  /// A shared pointer to an EquationSystem object
627  typedef std::shared_ptr<SolverUtils::Diffusion> DiffusionSharedPtr;
628 
629  /// Datatype of the NekFactory used to instantiate classes derived
630  /// from the Diffusion class.
634 }
635 }
636 
637 #endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:105
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:549
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:60
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:160
void SetFluxVectorNS(DiffusionFluxVecCBNS fluxVector)
Definition: Diffusion.h:350
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:249
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:77
SOLVER_UTILS_EXPORT void GetAVmu(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &muvar, Array< OneD, NekDouble > &MuVarTrace)
Get the mu of artifical viscosity(AV)
Definition: Diffusion.cpp:115
void SetSpecialBndTreat(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:407
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:266
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:236
void SetFluxPenaltyNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:356
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:471
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:195
void SetArtificialDiffusionVector(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:399
void SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:429
DiffusionArtificialDiffusion m_ArtificialDiffusionVector
Definition: Diffusion.h:469
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:224
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:267
virtual TensorOfArray3D< NekDouble > & GetFluxTensor()
Definition: Diffusion.h:434
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:199
void SetFluxPenaltyNS(DiffusionFluxPenaltyNS flux)
Definition: Diffusion.h:363
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:159
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:582
SOLVER_UTILS_EXPORT void GetDivCurl(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &pVarDer)
Compute divergence and curl squared.
Definition: Diffusion.cpp:319
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:300
virtual SOLVER_UTILS_EXPORT ~Diffusion()
Definition: Diffusion.h:163
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:161
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:467
void SetCalcViscosity(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:414
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:285
void SetFluxVectorNS(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:343
void SetDiffusionFluxConsTrace(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:384
void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:329
void SetDiffusionFluxCons(DiffusionFluxCons flux)
Definition: Diffusion.h:378
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:603
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:470
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:298
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:279
void SetDiffusionSymmFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:421
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:472
virtual void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:597
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:255
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:137
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal()
Get trace normal.
Definition: Diffusion.h:460
void SetDiffusionFluxCons(FuncPointerT func, ObjectPointerT obj)
Definition: Diffusion.h:369
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:474
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:242
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:568
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:468
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:182
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:446
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:236
SOLVER_UTILS_EXPORT void FluxVec(TensorOfArray3D< NekDouble > &fluxvector)
void SetDiffusionFluxConsTrace(DiffusionFluxCons flux)
Definition: Diffusion.h:393
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:466
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.h:478
SOLVER_UTILS_EXPORT void SetFluxVector(DiffusionFluxVecCB fluxVector)
Definition: Diffusion.h:336
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::function< void(Array< OneD, Array< OneD, NekDouble > > &)> SpecialBndTreat
Definition: Diffusion.h:143
std::function< void(const Array< OneD, const Array< OneD, NekDouble > > &, Array< OneD, NekDouble > &)> CalcViscosity
Definition: Diffusion.h:153
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, Array< OneD, NekDouble > &)> DiffusionArtificialDiffusion
Definition: Diffusion.h:75
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, const TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCB
Definition: Diffusion.h:58
LibUtilities::NekFactory< std::string, Diffusion, std::string > DiffusionFactory
Datatype of the NekFactory used to instantiate classes derived from the Diffusion class.
Definition: Diffusion.h:632
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 > > &, const Array< OneD, NekDouble > &)> DiffusionFluxCons
Definition: Diffusion.h:114
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:627
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, const Array< OneD, Array< OneD, NekDouble > > &, Array< OneD, Array< OneD, NekDouble > > &)> DiffusionFluxPenaltyNS
Definition: Diffusion.h:70
std::function< void(const Array< OneD, Array< OneD, NekDouble > > &, TensorOfArray3D< NekDouble > &, TensorOfArray3D< NekDouble > &)> DiffusionFluxVecCBNS
Definition: Diffusion.h:64
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:92
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:135
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
static Array< OneD, int > NullInt1DArray
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble