Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Nektar::SolverUtils::Advection Class Referenceabstract

An abstract base class encapsulating the concept of advection of a vector field. More...

#include <Advection.h>

Inheritance diagram for Nektar::SolverUtils::Advection:
[legend]

Public Member Functions

virtual SOLVER_UTILS_EXPORT ~Advection ()
 
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Interface function to initialise the advection object. More...
 
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. More...
 
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. More...
 
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. More...
 
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 is the coefficients of basis for AdvectCoeffs() it is the physics on quadrature points for Advect() More...
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVector (FuncPointerT func, ObjectPointerT obj)
 Set the flux vector callback function. More...
 
void SetRiemannSolver (RiemannSolverSharedPtr riemann)
 Set a Riemann solver object for this advection object. More...
 
void SetFluxVector (AdvectionFluxVecCB fluxVector)
 Set the flux vector callback function. More...
 
void SetBaseFlow (const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 Set the base flow used for linearised advection objects. More...
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
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)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
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)
 
SOLVER_UTILS_EXPORT void AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
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)
 

Protected Member Functions

virtual SOLVER_UTILS_EXPORT void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialises the advection object. More...
 
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. More...
 
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. More...
 
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. More...
 
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)
 
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. More...
 
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)
 

Protected Attributes

AdvectionFluxVecCB m_fluxVector
 Callback function to the flux vector (set when advection is in conservative form). More...
 
RiemannSolverSharedPtr m_riemann
 Riemann solver for DG-type schemes. More...
 
int m_spaceDim
 Storage for space dimension. Used for homogeneous extension. More...
 

Detailed Description

An abstract base class encapsulating the concept of advection of a vector field.

Subclasses override the Advection::v_InitObject function to initialise the object and the Advection::v_Advect function to evaluate the advection of the vector field.

Definition at line 72 of file Advection.h.

Constructor & Destructor Documentation

◆ ~Advection()

virtual SOLVER_UTILS_EXPORT Nektar::SolverUtils::Advection::~Advection ( )
inlinevirtual

Definition at line 76 of file Advection.h.

77  {};

Member Function Documentation

◆ AddTraceJacToMat() [1/2]

template<typename DataType , typename TypeNekBlkMatSharedPtr >
void Nektar::SolverUtils::Advection::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 
)

Definition at line 89 of file Advection.cpp.

97 {
98  MultiRegions::ExpListSharedPtr tracelist = pFields[0]->GetTrace();
99  std::shared_ptr<LocalRegions::ExpansionVector> traceExp
100  = tracelist->GetExp();
101  int ntotTrac = (*traceExp).size();
102  int nTracePnts = tracelist->GetTotPoints();
103  int nTracPnt,nTracCoef;
104 
105  DNekMatSharedPtr tmp2Add;
106  Array<OneD, DataType> MatData0;
107  Array<OneD, NekDouble> MatData1;
108  Array<OneD, NekDouble> MatData2;
109 
110  Array<OneD, DNekMatSharedPtr> TraceJacFwd(ntotTrac);
111  Array<OneD, DNekMatSharedPtr> TraceJacBwd(ntotTrac);
112 
113  for(int nelmt = 0; nelmt < ntotTrac; nelmt++)
114  {
115  nTracCoef = (*traceExp)[nelmt]->GetNcoeffs();
116  nTracPnt = (*traceExp)[nelmt]->GetTotPoints();
117  TraceJacFwd[nelmt] =MemoryManager<DNekMat>
118  ::AllocateSharedPtr(nTracCoef, nTracPnt,0.0);
119  TraceJacBwd[nelmt] =MemoryManager<DNekMat>
120  ::AllocateSharedPtr(nTracCoef, nTracPnt,0.0);
121  }
122 
123  std::shared_ptr<LocalRegions::ExpansionVector> fieldExp
124  = pFields[0]->GetExp();
125  int nTotElmt = (*fieldExp).size();
126  int nElmtPnt,nElmtCoef;
127 
128  Array<OneD, DNekMatSharedPtr> ElmtJacQuad(nTotElmt);
129  Array<OneD, DNekMatSharedPtr> ElmtJacCoef(nTotElmt);
130 
131  Array<OneD, NekDouble> SymmMatData;
132 
133  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
134  {
135  nElmtCoef = (*fieldExp)[nelmt]->GetNcoeffs();
136  nElmtPnt = (*fieldExp)[nelmt]->GetTotPoints();
137 
138  ElmtJacQuad[nelmt] =MemoryManager<DNekMat>
139  ::AllocateSharedPtr(nElmtCoef, nElmtPnt,0.0);
140  ElmtJacCoef[nelmt] =MemoryManager<DNekMat>
141  ::AllocateSharedPtr(nElmtCoef, nElmtCoef,0.0);
142  }
143 
144  bool TracePntJacGradflag = true;
145 
146  Array<OneD, Array<OneD, DataType> > TraceJacConsSign(2);
147  for(int i = 0; i < 2; i++)
148  {
149  TraceJacConsSign[i] = Array<OneD, DataType>(nTracePnts,1.0);
150  }
151 
152  if(0==TracePntJacGrad.size())
153  {
154  TracePntJacGradflag = false;
155  }
156 
157  for(int m = 0; m < nConvectiveFields; m++)
158  {
159  for(int n = 0; n < nConvectiveFields; n++)
160  {
161  // ElmtJacCons to set 0
162  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
163  {
164  (*ElmtJacCoef[nelmt]) = 0.0;
165  (*ElmtJacQuad[nelmt]) = 0.0;
166  }
167 
168  if(TracePntJacGradflag)
169  {
170  //only BJac is stored here because BJac = - FJac
171  for(int ndir = 0; ndir < nSpaceDim; ndir++)
172  {
173  // ElmtJacGrad to set 0
174  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
175  {
176  (*ElmtJacQuad[nelmt]) = 0.0;
177  }
178  int ngrad = n*nSpaceDim+ndir;
179 
180  CalcJacobTraceInteg(pFields, m,ngrad, TracePntJacGrad,
181  TracePntJacGradSign, TraceJacFwd, TraceJacBwd);
182  pFields[0]->AddTraceJacToElmtJac(TraceJacFwd, TraceJacBwd,
183  ElmtJacQuad);
184  // cost can be lower down by doing 3 directions together
185  pFields[0]->AddRightIPTPhysDerivBase(ndir, ElmtJacQuad,
186  ElmtJacCoef);
187  }
188  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
189  {
190  nElmtCoef = (*fieldExp)[nelmt]->GetNcoeffs();
191  if(SymmMatData.size()<nElmtCoef)
192  {
193  SymmMatData =
194  Array<OneD, NekDouble> (nElmtCoef*nElmtCoef);
195  }
196  MatData1 = ElmtJacCoef[nelmt]->GetPtr();
197  for(int i=0;i<nElmtCoef;i++)
198  {
199  Vmath::Vcopy(nElmtCoef, &MatData1[i*nElmtCoef], 1,
200  &SymmMatData[i], nElmtCoef);
201  }
202  Vmath::Vadd(nElmtCoef*nElmtCoef, SymmMatData, 1,
203  MatData1, 1, MatData1, 1);
204  }
205  }
206 
207  CalcJacobTraceInteg(pFields, m, n, TracePntJacCons,
208  TraceJacConsSign, TraceJacFwd, TraceJacBwd);
209 
210  pFields[0]->AddTraceJacToElmtJac(TraceJacFwd, TraceJacBwd,
211  ElmtJacQuad);
212  //Memory can be save to explore reusing ElmtJacQuad as output
213  pFields[0]->AddRightIPTBaseMatrix(ElmtJacQuad, ElmtJacCoef);
214 
215  // add ElmtJacCons to gmtxarray[m][n]
216  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
217  {
218  tmp2Add = ElmtJacCoef[nelmt];
219  MatData0 = gmtxarray[m][n]->GetBlock(nelmt,nelmt)->GetPtr();
220  MatData1 = tmp2Add->GetPtr();
221  for(int i=0;i<MatData0.size();i++)
222  {
223  MatData0[i] += DataType(MatData1[i]);
224  }
225  }
226  }
227  }
228 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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:250
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:322
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), CalcJacobTraceInteg(), Nektar::Array< OneD, const DataType >::size(), Vmath::Vadd(), and Vmath::Vcopy().

◆ AddTraceJacToMat() [2/2]

template<typename DataType , typename TypeNekBlkMatSharedPtr >
SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Advection::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 
)

◆ AddVolumJacToMat()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Advection::AddVolumJacToMat ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const int &  nConvectiveFields,
const TensorOfArray5D< NekDouble > &  ElmtJacArray,
Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &  gmtxarray 
)
inline

Definition at line 209 of file Advection.h.

214  {
215  v_AddVolumJacToMat(pFields,nConvectiveFields,ElmtJacArray,gmtxarray);
216  }
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:418

◆ Advect()

void Nektar::SolverUtils::Advection::Advect ( 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 > > &  pOutarray,
const NekDouble pTime,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd = NullNekDoubleArrayOfArray,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd = NullNekDoubleArrayOfArray 
)

Interface function to advect the vector field.

Parameters
nConvectiveFieldsNumber of velocity components.
pFieldsExpansion lists for scalar fields.
pAdvVelAdvection velocity.
pInarrayScalar data to advect.
pOutarrayAdvected scalar data.
pTimeSimulation time.

Definition at line 74 of file Advection.cpp.

83 {
84  v_Advect(nConvectiveFields, pFields, pAdvVel, pInarray,
85  pOutarray, pTime, pFwd, pBwd);
86 }
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.

References v_Advect().

◆ AdvectCoeffs()

void Nektar::SolverUtils::Advection::AdvectCoeffs ( 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 > > &  pOutarray,
const NekDouble pTime,
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 is the coefficients of basis for AdvectCoeffs() it is the physics on quadrature points for Advect()

Parameters
nConvectiveFieldsNumber of velocity components.
pFieldsExpansion lists for scalar fields.
pAdvVelAdvection velocity.
pInarrayScalar data to advect.
pOutarrayAdvected scalar data.
pTimeSimulation time.

Definition at line 343 of file Advection.cpp.

352 {
353  v_AdvectCoeffs(nConvectiveFields, pFields, pAdvVel, pInarray,
354  pOutarray, pTime, pFwd, pBwd);
355 }
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:403

References v_AdvectCoeffs().

◆ AdvectTraceFlux()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Advection::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 
)
inline

Interface function to advect the Trace field.

Definition at line 109 of file Advection.h.

120  {
121  v_AdvectTraceFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
122  pTraceFlux, pTime, pFwd, pBwd);
123  }
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:315

◆ AdvectVolumeFlux()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Advection::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 
)
inline

Interface function to advect the Volume field.

Definition at line 96 of file Advection.h.

103  {
104  v_AdvectVolumeFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
105  pVolumeFlux, pTime);
106  }
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:299

◆ CalcJacobTraceInteg()

template<typename DataType , typename TypeNekBlkMatSharedPtr >
void Nektar::SolverUtils::Advection::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 at line 250 of file Advection.cpp.

258 {
259  MultiRegions::ExpListSharedPtr tracelist = pFields[0]->GetTrace();
260  std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
261  tracelist->GetExp();
262  int ntotTrac = (*traceExp).size();
263  int nTracPnt,noffset,pntoffset;
264 
265  Array<OneD, int > tracepnts(ntotTrac);
266  Array<OneD, Array<OneD, NekDouble> > JacFwd(ntotTrac);
267  Array<OneD, Array<OneD, NekDouble> > JacBwd(ntotTrac);
268 
269  for(int nelmt = 0; nelmt < ntotTrac; nelmt++)
270  {
271  nTracPnt = (*traceExp)[nelmt]->GetTotPoints();
272  tracepnts[nelmt] = nTracPnt;
273 
274  JacFwd[nelmt] =Array<OneD, NekDouble>(nTracPnt,0.0);
275  JacBwd[nelmt] =Array<OneD, NekDouble>(nTracPnt,0.0);
276  }
277 
278  // DNekMatSharedPtr FtmpMat,BtmpMat;
279  DataType ftmp,btmp;
280  for(int nelmt = 0; nelmt < ntotTrac; nelmt++)
281  {
282  nTracPnt = tracepnts[nelmt];
283  noffset = tracelist->GetPhys_Offset(nelmt);
284  for(int npnt = 0; npnt < nTracPnt; npnt++)
285  {
286  pntoffset = noffset+npnt;
287  ftmp = (*PntJac[0]->GetBlock(pntoffset,pntoffset))(m,n);
288  JacFwd[nelmt][npnt] = NekDouble(PntJacSign[0][pntoffset]*ftmp);
289 
290  btmp = (*PntJac[1]->GetBlock(pntoffset,pntoffset))(m,n);
291  JacBwd[nelmt][npnt] = NekDouble(PntJacSign[1][pntoffset]*btmp);
292  }
293  }
294  tracelist->GetDiagMatIpwrtBase(JacFwd,TraceJacFwd);
295  tracelist->GetDiagMatIpwrtBase(JacBwd,TraceJacBwd);
296 }
double NekDouble

Referenced by AddTraceJacToMat().

◆ InitObject()

void Nektar::SolverUtils::Advection::InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)

Interface function to initialise the advection object.

Parameters
pSessionSession configuration data.
pFieldsArray of ExpList objects.

Definition at line 58 of file Advection.cpp.

61 {
62  v_InitObject(pSession, pFields);
63 }
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:366

References v_InitObject().

◆ SetBaseFlow()

void Nektar::SolverUtils::Advection::SetBaseFlow ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields 
)
inline

Set the base flow used for linearised advection objects.

Parameters
inarrayVector to use as baseflow

Definition at line 182 of file Advection.h.

185  {
186  v_SetBaseFlow(inarray, fields);
187  }
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:394

◆ SetFluxVector() [1/2]

void Nektar::SolverUtils::Advection::SetFluxVector ( AdvectionFluxVecCB  fluxVector)
inline

Set the flux vector callback function.

Parameters
fluxVectorThe callback function to override.

Definition at line 172 of file Advection.h.

173  {
174  m_fluxVector = fluxVector;
175  }
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:221

◆ SetFluxVector() [2/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Advection::SetFluxVector ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Set the flux vector callback function.

This routine is a utility function to avoid the explicit use of std::bind. A function and object can be passed to this function instead.

Definition at line 151 of file Advection.h.

152  {
153  m_fluxVector = std::bind(
154  func, obj, std::placeholders::_1, std::placeholders::_2);
155  }

◆ SetRiemannSolver()

void Nektar::SolverUtils::Advection::SetRiemannSolver ( RiemannSolverSharedPtr  riemann)
inline

Set a Riemann solver object for this advection object.

Parameters
riemannThe RiemannSolver object.

Definition at line 162 of file Advection.h.

163  {
164  m_riemann = riemann;
165  }
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:223

◆ v_AddVolumJacToMat()

void Nektar::SolverUtils::Advection::v_AddVolumJacToMat ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const int &  nConvectiveFields,
const TensorOfArray5D< NekDouble > &  ElmtJacArray,
Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &  gmtxarray 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::AdvectionWeakDG.

Definition at line 418 of file Advection.cpp.

423  {
424  boost::ignore_unused(pFields,nConvectiveFields,ElmtJacArray,gmtxarray);
425  ASSERTL0(false,"v_AddVolumJacToMat NOT SPECIFIED");
426  return;
427  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216

References ASSERTL0.

◆ v_Advect()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Advection::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 
)
protectedpure virtual

◆ v_AdvectCoeffs()

void Nektar::SolverUtils::Advection::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 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::AdvectionWeakDG.

Definition at line 403 of file Advection.cpp.

412 {
413  boost::ignore_unused(nConvectiveFields, fields, advVel, inarray, outarray,
414  time, pFwd, pBwd);
415  ASSERTL0(false, "v_AdvectCoeffs not defined");
416 }

References ASSERTL0.

Referenced by AdvectCoeffs().

◆ v_AdvectTraceFlux()

void Nektar::SolverUtils::Advection::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 
)
protectedvirtual

Advects Trace Flux.

calculate the advection flux in the cell the trace integration

Definition at line 315 of file Advection.cpp.

324 {
325  boost::ignore_unused(nConvectiveFields, pFields, pAdvVel, pInarray,
326  pTraceFlux, pTime, pFwd, pBwd);
327  ASSERTL0(false, "Not defined for AdvectTraceFlux.");
328 }

References ASSERTL0.

◆ v_AdvectVolumeFlux()

void Nektar::SolverUtils::Advection::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 
)
protectedvirtual

Advects Volume Flux.

Definition at line 299 of file Advection.cpp.

306 {
307  boost::ignore_unused(nConvectiveFields, pFields, pAdvVel, pInarray,
308  pVolumeFlux, pTime);
309  ASSERTL0(false, "Not defined for AdvectVolumeFlux.");
310 }

References ASSERTL0.

◆ v_InitObject()

void Nektar::SolverUtils::Advection::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Initialises the advection object.

This function should be overridden in derived classes to initialise the specific advection data members. However, this base class function should be called as the first statement of the overridden function to ensure the base class is correctly initialised in order.

Parameters
pSessionSession information.
pFieldsExpansion lists for scalar fields.

Reimplemented in Nektar::SkewSymmetricAdvection, Nektar::NoAdvection, Nektar::NavierStokesAdvection, Nektar::LinearisedAdvection, Nektar::AlternateSkewAdvection, Nektar::SolverUtils::AdvectionWeakDG, Nektar::SolverUtils::AdvectionNonConservative, Nektar::SolverUtils::AdvectionFR, and Nektar::SolverUtils::Advection3DHomogeneous1D.

Definition at line 366 of file Advection.cpp.

369 {
370  m_spaceDim = pFields[0]->GetCoordim(0);
371 
372  if (pSession->DefinesSolverInfo("HOMOGENEOUS"))
373  {
374  std::string HomoStr = pSession->GetSolverInfo("HOMOGENEOUS");
375  if (HomoStr == "HOMOGENEOUS1D" || HomoStr == "Homogeneous1D" ||
376  HomoStr == "1D" || HomoStr == "Homo1D" ||
377  HomoStr == "HOMOGENEOUS2D" || HomoStr == "Homogeneous2D" ||
378  HomoStr == "2D" || HomoStr == "Homo2D")
379  {
380  m_spaceDim++;
381  }
382  else
383  {
385  "Only 1D homogeneous dimension supported.");
386  }
387  }
388 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
Definition: Advection.h:225

References Nektar::ErrorUtil::efatal, m_spaceDim, and NEKERROR.

Referenced by InitObject(), Nektar::SolverUtils::AdvectionFR::v_InitObject(), Nektar::SolverUtils::AdvectionNonConservative::v_InitObject(), Nektar::SolverUtils::AdvectionWeakDG::v_InitObject(), Nektar::LinearisedAdvection::v_InitObject(), Nektar::NavierStokesAdvection::v_InitObject(), and Nektar::SkewSymmetricAdvection::v_InitObject().

◆ v_SetBaseFlow()

void Nektar::SolverUtils::Advection::v_SetBaseFlow ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields 
)
protectedvirtual

Overrides the base flow used during linearised advection.

Reimplemented in Nektar::LinearisedAdvection.

Definition at line 394 of file Advection.cpp.

397 {
398  boost::ignore_unused(inarray, fields);
400  "A baseflow is not appropriate for this advection type.");
401 }

References Nektar::ErrorUtil::efatal, and NEKERROR.

Member Data Documentation

◆ m_fluxVector

AdvectionFluxVecCB Nektar::SolverUtils::Advection::m_fluxVector
protected

Callback function to the flux vector (set when advection is in conservative form).

Definition at line 221 of file Advection.h.

Referenced by Nektar::SolverUtils::AdvectionFR::v_Advect(), Nektar::SolverUtils::Advection3DHomogeneous1D::v_Advect(), and Nektar::SolverUtils::AdvectionWeakDG::v_AdvectVolumeFlux().

◆ m_riemann

RiemannSolverSharedPtr Nektar::SolverUtils::Advection::m_riemann
protected

◆ m_spaceDim

int Nektar::SolverUtils::Advection::m_spaceDim
protected