Nektar++
Advection.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: Advection.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 advection.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_SOLVERUTILS_ADVECTION
36 #define NEKTAR_SOLVERUTILS_ADVECTION
37 
38 #include <functional>
39 #include <string>
40 
43 #include <MultiRegions/ExpList.h>
46 #include <iomanip>
47 
48 namespace Nektar
49 {
50 namespace SolverUtils
51 {
52 
53 /**
54  * Defines a callback function type which evaluates the flux vector \f$ F(u) \f$
55  * in a conservative advection of the form \f$ \nabla\cdot F(u) \f$.
56  */
57 typedef std::function<void(const Array<OneD, Array<OneD, NekDouble>> &,
58  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &)>
60 
61 /**
62  * @brief An abstract base class encapsulating the concept of advection
63  * of a vector field.
64  *
65  * Subclasses override the Advection::v_InitObject function to
66  * initialise the object and the Advection::v_Advect function to
67  * evaluate the advection of the vector field.
68  */
69 class Advection
70 {
71 public:
73 
74  /// Interface function to initialise the advection object.
78 
79  /// Interface function to advect the vector field.
81  const int nConvectiveFields,
83  const Array<OneD, Array<OneD, NekDouble>> &advVel,
84  const Array<OneD, Array<OneD, NekDouble>> &inarray,
85  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
86  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
88  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
90 
91  /// Interface function to advect the Volume field.
93  const int nConvectiveFields,
95  const Array<OneD, Array<OneD, NekDouble>> &pAdvVel,
96  const Array<OneD, Array<OneD, NekDouble>> &pInarray,
97  TensorOfArray3D<NekDouble> &pVolumeFlux, const NekDouble &pTime)
98  {
99  v_AdvectVolumeFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
100  pVolumeFlux, pTime);
101  }
102 
103  /// Interface function to advect the Trace field.
105  const int nConvectiveFields,
107  const Array<OneD, Array<OneD, NekDouble>> &pAdvVel,
108  const Array<OneD, Array<OneD, NekDouble>> &pInarray,
109  Array<OneD, Array<OneD, NekDouble>> &pTraceFlux, const NekDouble &pTime,
110  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
112  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
114  {
115  v_AdvectTraceFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
116  pTraceFlux, pTime, pFwd, pBwd);
117  }
118 
119  /**
120  * @brief Similar with Advection::Advect(): calculate the advection flux
121  * The difference is in the outarray:
122  * it is the coefficients of basis for AdvectCoeffs()
123  * it is the physics on quadrature points for Advect()
124  */
126  const int nConvectiveFields,
128  const Array<OneD, Array<OneD, NekDouble>> &advVel,
129  const Array<OneD, Array<OneD, NekDouble>> &inarray,
130  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
131  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
133  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
135 
136  /**
137  * @brief Set the flux vector callback function.
138  *
139  * This routine is a utility function to avoid the explicit use of
140  * std::bind. A function and object can be passed to this function
141  * instead.
142  */
143  template <typename FuncPointerT, typename ObjectPointerT>
144  void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
145  {
146  m_fluxVector =
147  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2);
148  }
149 
150  /**
151  * @brief Set a Riemann solver object for this advection object.
152  *
153  * @param riemann The RiemannSolver object.
154  */
156  {
157  m_riemann = riemann;
158  }
159 
160  /**
161  * @brief Set the flux vector callback function.
162  *
163  * @param fluxVector The callback function to override.
164  */
165  inline void SetFluxVector(AdvectionFluxVecCB fluxVector)
166  {
167  m_fluxVector = fluxVector;
168  }
169 
170  /**
171  * @brief Set the base flow used for linearised advection objects.
172  *
173  * @param inarray Vector to use as baseflow
174  */
175  inline void SetBaseFlow(
176  const Array<OneD, Array<OneD, NekDouble>> &inarray,
178  {
179  v_SetBaseFlow(inarray, fields);
180  }
181 
182  template <typename DataType, typename TypeNekBlkMatSharedPtr>
184  const int nConvectiveFields, const int nSpaceDim,
186  const Array<OneD, TypeNekBlkMatSharedPtr> &TracePntJacCons,
188  const Array<OneD, TypeNekBlkMatSharedPtr> &TracePntJacGrad,
189  const Array<OneD, Array<OneD, DataType>> &TracePntJacGradSign);
190 
191  template <typename DataType, typename TypeNekBlkMatSharedPtr>
192  void CalcJacobTraceInteg(
193  const Array<OneD, MultiRegions::ExpListSharedPtr> &pFields, const int m,
194  const int n, const Array<OneD, const TypeNekBlkMatSharedPtr> &PntJac,
195  const Array<OneD, const Array<OneD, DataType>> &PntJacSign,
196  Array<OneD, DNekMatSharedPtr> &TraceJacFwd,
197  Array<OneD, DNekMatSharedPtr> &TraceJacBwd);
198 
201  const int &nConvectiveFields,
202  const TensorOfArray5D<NekDouble> &ElmtJacArray,
204  {
205  v_AddVolumJacToMat(pFields, nConvectiveFields, ElmtJacArray, gmtxarray);
206  }
207 
208 protected:
209  /// Callback function to the flux vector (set when advection is in
210  /// conservative form).
212  /// Riemann solver for DG-type schemes.
214  /// Storage for space dimension. Used for homogeneous extension.
216 
217  /// Initialises the advection object.
218  SOLVER_UTILS_EXPORT virtual void v_InitObject(
221 
222  /// Advects a vector field.
224  const int nConvectiveFields,
226  const Array<OneD, Array<OneD, NekDouble>> &advVel,
227  const Array<OneD, Array<OneD, NekDouble>> &inarray,
228  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
229  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
231  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
233 
234  /// Advects Volume Flux.
236  const int nConvectiveFields,
238  const Array<OneD, Array<OneD, NekDouble>> &advVel,
239  const Array<OneD, Array<OneD, NekDouble>> &inarray,
240  TensorOfArray3D<NekDouble> &pVolumeFlux, const NekDouble &time);
241 
242  /// Advects Trace Flux.
244  const int nConvectiveFields,
246  const Array<OneD, Array<OneD, NekDouble>> &advVel,
247  const Array<OneD, Array<OneD, NekDouble>> &inarray,
248  Array<OneD, Array<OneD, NekDouble>> &pTraceFlux, const NekDouble &time,
249  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
251  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
253 
254  SOLVER_UTILS_EXPORT virtual void v_AdvectCoeffs(
255  const int nConvectiveFields,
257  const Array<OneD, Array<OneD, NekDouble>> &advVel,
258  const Array<OneD, Array<OneD, NekDouble>> &inarray,
259  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
260  const Array<OneD, Array<OneD, NekDouble>> &pFwd =
262  const Array<OneD, Array<OneD, NekDouble>> &pBwd =
264 
265  /// Overrides the base flow used during linearised advection
266  SOLVER_UTILS_EXPORT virtual void v_SetBaseFlow(
267  const Array<OneD, Array<OneD, NekDouble>> &inarray,
269 
272  const int &nConvectiveFields,
273  const TensorOfArray5D<NekDouble> &ElmtJacArray,
275 };
276 
277 /// A shared pointer to an Advection object.
278 typedef std::shared_ptr<Advection> AdvectionSharedPtr;
279 
280 /// Datatype of the NekFactory used to instantiate classes derived
281 /// from the Advection class.
284 
285 /// Gets the factory for initialising advection objects.
287 
288 } // namespace SolverUtils
289 } // namespace Nektar
290 
291 #endif
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:105
An abstract base class encapsulating the concept of advection of a vector field.
Definition: Advection.h:70
void SetBaseFlow(const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Set the base flow used for linearised advection objects.
Definition: Advection.h:175
void CalcJacobTraceInteg(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int m, const int n, const Array< OneD, const TypeNekBlkMatSharedPtr > &PntJac, const Array< OneD, const Array< OneD, DataType >> &PntJacSign, Array< OneD, DNekMatSharedPtr > &TraceJacFwd, Array< OneD, DNekMatSharedPtr > &TraceJacBwd)
Definition: Advection.cpp:243
SOLVER_UTILS_EXPORT void Advect(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Interface function to advect the vector field.
Definition: Advection.cpp:71
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
Definition: Advection.h:215
virtual SOLVER_UTILS_EXPORT void v_AddVolumJacToMat(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
Definition: Advection.cpp:403
virtual SOLVER_UTILS_EXPORT ~Advection()
Definition: Advection.h:72
virtual SOLVER_UTILS_EXPORT void v_AdvectCoeffs(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Definition: Advection.cpp:389
virtual SOLVER_UTILS_EXPORT void v_AdvectTraceFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &pTraceFlux, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Advects Trace Flux.
Definition: Advection.cpp:304
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:211
virtual SOLVER_UTILS_EXPORT void v_Advect(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)=0
Advects a vector field.
SOLVER_UTILS_EXPORT void AdvectTraceFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, Array< OneD, Array< OneD, NekDouble >> &pTraceFlux, const NekDouble &pTime, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Interface function to advect the Trace field.
Definition: Advection.h:104
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:353
SOLVER_UTILS_EXPORT void AdvectCoeffs(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Similar with Advection::Advect(): calculate the advection flux The difference is in the outarray: it ...
Definition: Advection.cpp:331
SOLVER_UTILS_EXPORT void InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Interface function to initialise the advection object.
Definition: Advection.cpp:57
SOLVER_UTILS_EXPORT void AddVolumJacToMat(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
Definition: Advection.h:199
virtual SOLVER_UTILS_EXPORT void v_SetBaseFlow(const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Overrides the base flow used during linearised advection.
Definition: Advection.cpp:380
SOLVER_UTILS_EXPORT void AdvectVolumeFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &pTime)
Interface function to advect the Volume field.
Definition: Advection.h:92
SOLVER_UTILS_EXPORT void AddTraceJacToMat(const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType >> &TracePntJacGradSign)
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:213
void SetFluxVector(AdvectionFluxVecCB fluxVector)
Set the flux vector callback function.
Definition: Advection.h:165
void SetRiemannSolver(RiemannSolverSharedPtr riemann)
Set a Riemann solver object for this advection object.
Definition: Advection.h:155
void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
Set the flux vector callback function.
Definition: Advection.h:144
virtual SOLVER_UTILS_EXPORT void v_AdvectVolumeFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &time)
Advects Volume Flux.
Definition: Advection.cpp:289
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::function< void(const Array< OneD, Array< OneD, NekDouble >> &, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &)> AdvectionFluxVecCB
Definition: Advection.h:59
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
LibUtilities::NekFactory< std::string, Advection, std::string > AdvectionFactory
Datatype of the NekFactory used to instantiate classes derived from the Advection class.
Definition: Advection.h:283
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:278
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble