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 69 of file Advection.h.

Constructor & Destructor Documentation

◆ ~Advection()

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

Definition at line 72 of file Advection.h.

72 {};

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 85 of file Advection.cpp.

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

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 199 of file Advection.h.

204  {
205  v_AddVolumJacToMat(pFields, nConvectiveFields, ElmtJacArray, gmtxarray);
206  }
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

References v_AddVolumJacToMat().

◆ 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 71 of file Advection.cpp.

79 {
80  v_Advect(nConvectiveFields, pFields, pAdvVel, pInarray, pOutarray, pTime,
81  pFwd, pBwd);
82 }
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 331 of file Advection.cpp.

339 {
340  v_AdvectCoeffs(nConvectiveFields, pFields, pAdvVel, pInarray, pOutarray,
341  pTime, pFwd, pBwd);
342 }
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

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 104 of file Advection.h.

114  {
115  v_AdvectTraceFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
116  pTraceFlux, pTime, pFwd, pBwd);
117  }
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

References v_AdvectTraceFlux().

◆ 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 92 of file Advection.h.

98  {
99  v_AdvectVolumeFlux(nConvectiveFields, pFields, pAdvVel, pInarray,
100  pVolumeFlux, pTime);
101  }
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

References v_AdvectVolumeFlux().

◆ 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 243 of file Advection.cpp.

249 {
250  MultiRegions::ExpListSharedPtr tracelist = pFields[0]->GetTrace();
251  std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
252  tracelist->GetExp();
253  int ntotTrac = (*traceExp).size();
254  int nTracPnt, noffset, pntoffset;
255 
256  Array<OneD, int> tracepnts(ntotTrac);
257  Array<OneD, Array<OneD, NekDouble>> JacFwd(ntotTrac);
258  Array<OneD, Array<OneD, NekDouble>> JacBwd(ntotTrac);
259 
260  for (int nelmt = 0; nelmt < ntotTrac; nelmt++)
261  {
262  nTracPnt = (*traceExp)[nelmt]->GetTotPoints();
263  tracepnts[nelmt] = nTracPnt;
264 
265  JacFwd[nelmt] = Array<OneD, NekDouble>(nTracPnt, 0.0);
266  JacBwd[nelmt] = Array<OneD, NekDouble>(nTracPnt, 0.0);
267  }
268 
269  // DNekMatSharedPtr FtmpMat,BtmpMat;
270  DataType ftmp, btmp;
271  for (int nelmt = 0; nelmt < ntotTrac; nelmt++)
272  {
273  nTracPnt = tracepnts[nelmt];
274  noffset = tracelist->GetPhys_Offset(nelmt);
275  for (int npnt = 0; npnt < nTracPnt; npnt++)
276  {
277  pntoffset = noffset + npnt;
278  ftmp = (*PntJac[0]->GetBlock(pntoffset, pntoffset))(m, n);
279  JacFwd[nelmt][npnt] = NekDouble(PntJacSign[0][pntoffset] * ftmp);
280 
281  btmp = (*PntJac[1]->GetBlock(pntoffset, pntoffset))(m, n);
282  JacBwd[nelmt][npnt] = NekDouble(PntJacSign[1][pntoffset] * btmp);
283  }
284  }
285  tracelist->GetDiagMatIpwrtBase(JacFwd, TraceJacFwd);
286  tracelist->GetDiagMatIpwrtBase(JacBwd, TraceJacBwd);
287 }
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 57 of file Advection.cpp.

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

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 175 of file Advection.h.

178  {
179  v_SetBaseFlow(inarray, fields);
180  }
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

References v_SetBaseFlow().

◆ 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 165 of file Advection.h.

166  {
167  m_fluxVector = fluxVector;
168  }
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:211

References m_fluxVector.

◆ 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 144 of file Advection.h.

145  {
146  m_fluxVector =
147  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2);
148  }

References m_fluxVector.

◆ 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 155 of file Advection.h.

156  {
157  m_riemann = riemann;
158  }
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:213

References m_riemann.

◆ 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 403 of file Advection.cpp.

408 {
409  boost::ignore_unused(pFields, nConvectiveFields, ElmtJacArray, gmtxarray);
410  ASSERTL0(false, "v_AddVolumJacToMat NOT SPECIFIED");
411  return;
412 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215

References ASSERTL0.

Referenced by AddVolumJacToMat().

◆ 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 389 of file Advection.cpp.

397 {
398  boost::ignore_unused(nConvectiveFields, fields, advVel, inarray, outarray,
399  time, pFwd, pBwd);
400  ASSERTL0(false, "v_AdvectCoeffs not defined");
401 }

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

Reimplemented in Nektar::SolverUtils::AdvectionWeakDG.

Definition at line 304 of file Advection.cpp.

312 {
313  boost::ignore_unused(nConvectiveFields, pFields, pAdvVel, pInarray,
314  pTraceFlux, pTime, pFwd, pBwd);
315  ASSERTL0(false, "Not defined for AdvectTraceFlux.");
316 }

References ASSERTL0.

Referenced by AdvectTraceFlux().

◆ 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.

Reimplemented in Nektar::SolverUtils::AdvectionWeakDG.

Definition at line 289 of file Advection.cpp.

295 {
296  boost::ignore_unused(nConvectiveFields, pFields, pAdvVel, pInarray,
297  pVolumeFlux, pTime);
298  ASSERTL0(false, "Not defined for AdvectVolumeFlux.");
299 }

References ASSERTL0.

Referenced by AdvectVolumeFlux().

◆ 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 353 of file Advection.cpp.

356 {
357  m_spaceDim = pFields[0]->GetCoordim(0);
358 
359  if (pSession->DefinesSolverInfo("HOMOGENEOUS"))
360  {
361  std::string HomoStr = pSession->GetSolverInfo("HOMOGENEOUS");
362  if (HomoStr == "HOMOGENEOUS1D" || HomoStr == "Homogeneous1D" ||
363  HomoStr == "1D" || HomoStr == "Homo1D" ||
364  HomoStr == "HOMOGENEOUS2D" || HomoStr == "Homogeneous2D" ||
365  HomoStr == "2D" || HomoStr == "Homo2D")
366  {
367  m_spaceDim++;
368  }
369  else
370  {
372  "Only 1D homogeneous dimension supported.");
373  }
374  }
375 }
#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:215

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 380 of file Advection.cpp.

383 {
384  boost::ignore_unused(inarray, fields);
386  "A baseflow is not appropriate for this advection type.");
387 }

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

Referenced by SetBaseFlow().

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 211 of file Advection.h.

Referenced by SetFluxVector(), 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