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 <string>
39 #include <functional>
40 
43 #include <MultiRegions/ExpList.h>
46 #include <iomanip>
47 
48 using namespace std;
49 
50 namespace Nektar
51 {
52 namespace SolverUtils
53 {
54 
55 /**
56  * Defines a callback function type which evaluates the flux vector \f$ F(u) \f$
57  * in a conservative advection of the form \f$ \nabla\cdot F(u) \f$.
58  */
59 typedef std::function<void (
60  const Array<OneD, Array<OneD, NekDouble> >&,
61  Array<OneD, Array<OneD, Array<OneD, NekDouble> > >&)>
63 
64 /**
65  * @brief An abstract base class encapsulating the concept of advection
66  * of a vector field.
67  *
68  * Subclasses override the Advection::v_InitObject function to
69  * initialise the object and the Advection::v_Advect function to
70  * evaluate the advection of the vector field.
71  */
72 class Advection
73 {
74 public:
75 
77  {};
78 
79  /// Interface function to initialise the advection object.
80  SOLVER_UTILS_EXPORT void InitObject(
83 
84  /// Interface function to advect the vector field.
85  SOLVER_UTILS_EXPORT void Advect(
86  const int nConvectiveFields,
88  const Array<OneD, Array<OneD, NekDouble> > &advVel,
89  const Array<OneD, Array<OneD, NekDouble> > &inarray,
90  Array<OneD, Array<OneD, NekDouble> > &outarray,
91  const NekDouble &time,
94 
95  /// Interface function to advect the Volume field.
97  const int nConvectiveFields,
99  const Array<OneD, Array<OneD, NekDouble>> &pAdvVel,
100  const Array<OneD, Array<OneD, NekDouble>> &pInarray,
101  TensorOfArray3D<NekDouble> &pVolumeFlux,
102  const NekDouble &pTime)
103  {
104  v_AdvectVolumeFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
105  pVolumeFlux, pTime);
106  }
107 
108  /// Interface function to advect the Trace field.
110  const int nConvectiveFields,
112  const Array<OneD, Array<OneD, NekDouble>> &pAdvVel,
113  const Array<OneD, Array<OneD, NekDouble>> &pInarray,
114  Array<OneD, Array<OneD, NekDouble>> &pTraceFlux,
115  const NekDouble &pTime,
120  {
121  v_AdvectTraceFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
122  pTraceFlux, pTime, pFwd, pBwd);
123  }
124 
125  /**
126  * @brief Similar with Advection::Advect(): calculate the advection flux
127  * The difference is in the outarray:
128  * it is the coefficients of basis for AdvectCoeffs()
129  * it is the physics on quadrature points for Advect()
130  */
131  SOLVER_UTILS_EXPORT void AdvectCoeffs(
132  const int nConvectiveFields,
134  const Array<OneD, Array<OneD, NekDouble> > &advVel,
135  const Array<OneD, Array<OneD, NekDouble> > &inarray,
136  Array<OneD, Array<OneD, NekDouble> > &outarray,
137  const NekDouble &time,
141  &pBwd = NullNekDoubleArrayOfArray);
142 
143  /**
144  * @brief Set the flux vector callback function.
145  *
146  * This routine is a utility function to avoid the explicit use of
147  * std::bind. A function and object can be passed to this function
148  * instead.
149  */
150  template<typename FuncPointerT, typename ObjectPointerT>
151  void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
152  {
153  m_fluxVector = std::bind(
154  func, obj, std::placeholders::_1, std::placeholders::_2);
155  }
156 
157  /**
158  * @brief Set a Riemann solver object for this advection object.
159  *
160  * @param riemann The RiemannSolver object.
161  */
163  {
164  m_riemann = riemann;
165  }
166 
167  /**
168  * @brief Set the flux vector callback function.
169  *
170  * @param fluxVector The callback function to override.
171  */
172  inline void SetFluxVector(AdvectionFluxVecCB fluxVector)
173  {
174  m_fluxVector = fluxVector;
175  }
176 
177  /**
178  * @brief Set the base flow used for linearised advection objects.
179  *
180  * @param inarray Vector to use as baseflow
181  */
182  inline void SetBaseFlow(
183  const Array<OneD, Array<OneD, NekDouble> >& inarray,
185  {
186  v_SetBaseFlow(inarray, fields);
187  }
188 
189  template<typename DataType, typename TypeNekBlkMatSharedPtr>
191  const int nConvectiveFields,
192  const int nSpaceDim,
194  const Array<OneD, TypeNekBlkMatSharedPtr> &TracePntJacCons,
196  const Array<OneD, TypeNekBlkMatSharedPtr> &TracePntJacGrad ,
197  const Array<OneD, Array<OneD, DataType> > &TracePntJacGradSign );
198 
199  template<typename DataType, typename TypeNekBlkMatSharedPtr>
200  void CalcJacobTraceInteg(
202  const int m,
203  const int n,
205  const Array<OneD, const Array<OneD, DataType> > & PntJacSign,
206  Array<OneD, DNekMatSharedPtr> & TraceJacFwd,
207  Array<OneD, DNekMatSharedPtr> & TraceJacBwd);
208 
211  const int &nConvectiveFields,
212  const TensorOfArray5D<NekDouble> &ElmtJacArray,
214  {
215  v_AddVolumJacToMat(pFields,nConvectiveFields,ElmtJacArray,gmtxarray);
216  }
217 
218 protected:
219  /// Callback function to the flux vector (set when advection is in
220  /// conservative form).
222  /// Riemann solver for DG-type schemes.
224  /// Storage for space dimension. Used for homogeneous extension.
226 
227  /// Initialises the advection object.
228  SOLVER_UTILS_EXPORT virtual void v_InitObject(
231 
232  /// Advects a vector field.
234  const int nConvectiveFields,
236  const Array<OneD, Array<OneD, NekDouble> > &advVel,
237  const Array<OneD, Array<OneD, NekDouble> > &inarray,
238  Array<OneD, Array<OneD, NekDouble> > &outarray,
239  const NekDouble &time,
242 
243  /// Advects Volume Flux.
244  SOLVER_UTILS_EXPORT virtual void v_AdvectVolumeFlux(
245  const int nConvectiveFields,
247  const Array<OneD, Array<OneD, NekDouble> > &advVel,
248  const Array<OneD, Array<OneD, NekDouble> > &inarray,
249  TensorOfArray3D<NekDouble> &pVolumeFlux,
250  const NekDouble &time);
251 
252  /// Advects Trace Flux.
253  SOLVER_UTILS_EXPORT virtual void v_AdvectTraceFlux(
254  const int nConvectiveFields,
256  const Array<OneD, Array<OneD, NekDouble> > &advVel,
257  const Array<OneD, Array<OneD, NekDouble> > &inarray,
258  Array<OneD, Array<OneD, NekDouble>> &pTraceFlux,
259  const NekDouble &time,
263  &pBwd = NullNekDoubleArrayOfArray);
264 
265  SOLVER_UTILS_EXPORT virtual void v_AdvectCoeffs(
266  const int nConvectiveFields,
268  const Array<OneD, Array<OneD, NekDouble> > &advVel,
269  const Array<OneD, Array<OneD, NekDouble> > &inarray,
270  Array<OneD, Array<OneD, NekDouble> > &outarray,
271  const NekDouble &time,
275  &pBwd = NullNekDoubleArrayOfArray);
276 
277  /// Overrides the base flow used during linearised advection
278  SOLVER_UTILS_EXPORT virtual void v_SetBaseFlow(
279  const Array<OneD, Array<OneD, NekDouble> > &inarray,
281 
282  SOLVER_UTILS_EXPORT virtual void v_AddVolumJacToMat(
284  const int &nConvectiveFields,
285  const TensorOfArray5D<NekDouble> &ElmtJacArray,
287 
288 };
289 
290 /// A shared pointer to an Advection object.
291 typedef std::shared_ptr<Advection> AdvectionSharedPtr;
292 
293 /// Datatype of the NekFactory used to instantiate classes derived
294 /// from the Advection class.
295 typedef LibUtilities::NekFactory<std::string, Advection,
296 std::string> AdvectionFactory;
297 
298 /// Gets the factory for initialising advection objects.
300 
301 }
302 }
303 
304 #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:73
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
Definition: Advection.h:225
virtual SOLVER_UTILS_EXPORT ~Advection()
Definition: Advection.h:76
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:221
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:109
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 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)
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:209
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:182
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:96
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:223
void SetFluxVector(AdvectionFluxVecCB fluxVector)
Set the flux vector callback function.
Definition: Advection.h:172
void SetRiemannSolver(RiemannSolverSharedPtr riemann)
Set a Riemann solver object for this advection object.
Definition: Advection.h:162
void SetFluxVector(FuncPointerT func, ObjectPointerT obj)
Set the flux vector callback function.
Definition: Advection.h:151
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:62
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:296
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:291
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble