Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::SolverUtils::AdvectionWeakDG Class Reference

#include <AdvectionWeakDG.h>

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

Static Public Member Functions

static AdvectionSharedPtr create (std::string advType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 AdvectionWeakDG ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 Initialise AdvectionWeakDG objects and store them before starting the time-stepping. More...
 
virtual void v_Advect (const int nConvective, 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) override
 Compute the advection term at each time-step using the Discontinuous Galerkin approach (DG). More...
 
virtual void v_AdvectCoeffs (const int nConvective, 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) override
 
virtual 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 > &VolumeFlux, const NekDouble &time) override
 Advects Volume Flux. More...
 
virtual void v_AdvectTraceFlux (const int nConvective, 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 >> &TraceFlux, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray) override
 Advects Trace Flux. More...
 
virtual void v_AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::Advection
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...
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Advection
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 Attributes inherited from Nektar::SolverUtils::Advection
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

Definition at line 46 of file AdvectionWeakDG.h.

Constructor & Destructor Documentation

◆ AdvectionWeakDG()

Nektar::SolverUtils::AdvectionWeakDG::AdvectionWeakDG ( )
protected

Definition at line 51 of file AdvectionWeakDG.cpp.

52 {
53 }

Referenced by create().

Member Function Documentation

◆ create()

static AdvectionSharedPtr Nektar::SolverUtils::AdvectionWeakDG::create ( std::string  advType)
inlinestatic

Definition at line 49 of file AdvectionWeakDG.h.

50  {
51  boost::ignore_unused(advType);
52  return AdvectionSharedPtr(new AdvectionWeakDG());
53  }
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:278

References AdvectionWeakDG().

◆ v_AddVolumJacToMat()

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

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 233 of file AdvectionWeakDG.cpp.

238 {
239  MultiRegions::ExpListSharedPtr explist = pFields[0];
240  std::shared_ptr<LocalRegions::ExpansionVector> pexp = explist->GetExp();
241  int nTotElmt = (*pexp).size();
242  int nElmtPnt, nElmtCoef;
243 
244  SNekMatSharedPtr tmpGmtx;
245  DNekMatSharedPtr ElmtMat;
246 
247  Array<OneD, NekSingle> GMat_data;
248  Array<OneD, NekDouble> Elmt_data;
249  Array<OneD, NekSingle> Elmt_dataSingle;
250 
251  Array<OneD, DNekMatSharedPtr> mtxPerVar(nTotElmt);
252  Array<OneD, int> elmtpnts(nTotElmt);
253  Array<OneD, int> elmtcoef(nTotElmt);
254  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
255  {
256  nElmtCoef = (*pexp)[nelmt]->GetNcoeffs();
257  nElmtPnt = (*pexp)[nelmt]->GetTotPoints();
258  elmtpnts[nelmt] = nElmtPnt;
259  elmtcoef[nelmt] = nElmtCoef;
260  mtxPerVar[nelmt] =
261  MemoryManager<DNekMat>::AllocateSharedPtr(nElmtCoef, nElmtPnt);
262  }
263 
264  Array<OneD, DNekMatSharedPtr> mtxPerVarCoeff(nTotElmt);
265  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
266  {
267  nElmtCoef = elmtcoef[nelmt];
268  mtxPerVarCoeff[nelmt] =
269  MemoryManager<DNekMat>::AllocateSharedPtr(nElmtCoef, nElmtCoef);
270  }
271 
272  for (int m = 0; m < nConvectiveFields; m++)
273  {
274  for (int n = 0; n < nConvectiveFields; n++)
275  {
276  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
277  {
278  (*mtxPerVarCoeff[nelmt]) = 0.0;
279  (*mtxPerVar[nelmt]) = 0.0;
280  }
281 
282  explist->GetMatIpwrtDeriveBase(ElmtJacArray[m][n], mtxPerVar);
283  // Memory can be saved by reusing mtxPerVar
284  explist->AddRightIPTBaseMatrix(mtxPerVar, mtxPerVarCoeff);
285 
286  for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
287  {
288  nElmtCoef = elmtcoef[nelmt];
289  nElmtPnt = elmtpnts[nelmt];
290  int ntotDofs = nElmtCoef * nElmtCoef;
291 
292  if (Elmt_dataSingle.size() < ntotDofs)
293  {
294  Elmt_dataSingle = Array<OneD, NekSingle>(ntotDofs);
295  }
296 
297  tmpGmtx = gmtxarray[m][n]->GetBlock(nelmt, nelmt);
298  ElmtMat = mtxPerVarCoeff[nelmt];
299 
300  GMat_data = tmpGmtx->GetPtr();
301  Elmt_data = ElmtMat->GetPtr();
302 
303  for (int i = 0; i < ntotDofs; i++)
304  {
305  Elmt_dataSingle[i] = NekSingle(Elmt_data[i]);
306  }
307 
308  Vmath::Vadd(ntotDofs, GMat_data, 1, Elmt_dataSingle, 1,
309  GMat_data, 1);
310  }
311  }
312  }
313 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< SNekMat > SNekMatSharedPtr
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

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Vmath::Vadd().

◆ v_Advect()

void Nektar::SolverUtils::AdvectionWeakDG::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 
)
overrideprotectedvirtual

Compute the advection term at each time-step using the Discontinuous Galerkin approach (DG).

Parameters
nConvectiveFieldsNumber of fields.
fieldsPointer to fields.
advVelAdvection velocities.
inarraySolution at the previous time-step.
outarrayAdvection term to be passed at the time integration class.

Implements Nektar::SolverUtils::Advection.

Definition at line 80 of file AdvectionWeakDG.cpp.

88 {
89  LibUtilities::Timer timer1;
90  timer1.Start();
91 
92  boost::ignore_unused(advVel, time);
93  size_t nCoeffs = fields[0]->GetNcoeffs();
94 
95  Array<OneD, Array<OneD, NekDouble>> tmp{size_t(nConvectiveFields)};
96  for (int i = 0; i < nConvectiveFields; ++i)
97  {
98  tmp[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
99  }
100  LibUtilities::Timer timer;
101  timer.Start();
102  AdvectionWeakDG::v_AdvectCoeffs(nConvectiveFields, fields, advVel, inarray,
103  tmp, time, pFwd, pBwd);
104  timer.Stop();
105  timer.AccumulateRegion("AdvWeakDG:v_AdvectCoeffs", 2);
106 
107  // why was this broken in many loops over convective fields?
108  // this is terrible for locality
109 
110  timer.Start();
111  for (int i = 0; i < nConvectiveFields; ++i)
112  {
113  fields[i]->BwdTrans(tmp[i], outarray[i]);
114  }
115  timer.Stop();
116  timer.AccumulateRegion("AdvWeakDG:_BwdTrans", 10);
117  timer1.Stop();
118  timer1.AccumulateRegion("AdvWeakDG:All", 10);
119 }
virtual void v_AdvectCoeffs(const int nConvective, 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) override

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and v_AdvectCoeffs().

◆ v_AdvectCoeffs()

void Nektar::SolverUtils::AdvectionWeakDG::v_AdvectCoeffs ( const int  nConvective,
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 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 121 of file AdvectionWeakDG.cpp.

129 {
130  size_t nPointsTot = fields[0]->GetTotPoints();
131  size_t nCoeffs = fields[0]->GetNcoeffs();
132  size_t nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
133 
134  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxvector(
135  nConvectiveFields);
136  // Allocate storage for flux vector F(u).
137  for (int i = 0; i < nConvectiveFields; ++i)
138  {
139  fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>{size_t(m_spaceDim)};
140  for (int j = 0; j < m_spaceDim; ++j)
141  {
142  fluxvector[i][j] = Array<OneD, NekDouble>(nPointsTot, 0.0);
143  }
144  }
145 
146  LibUtilities::Timer timer;
147  timer.Start();
148  v_AdvectVolumeFlux(nConvectiveFields, fields, advVel, inarray, fluxvector,
149  time);
150  timer.Stop();
151  timer.AccumulateRegion("AdvWeakDG:_fluxVector", 3);
152 
153  timer.Start();
154  // Get the advection part (without numerical flux)
155  for (int i = 0; i < nConvectiveFields; ++i)
156  {
157  Vmath::Fill(outarray[i].size(), 0.0, outarray[i], 1);
158  fields[i]->IProductWRTDerivBase(fluxvector[i], outarray[i]);
159  }
160  timer.Stop();
161  timer.AccumulateRegion("AdvWeakDG:_IProductWRTDerivBase", 10);
162 
163  Array<OneD, Array<OneD, NekDouble>> numflux{size_t(nConvectiveFields)};
164 
165  for (int i = 0; i < nConvectiveFields; ++i)
166  {
167  numflux[i] = Array<OneD, NekDouble>{nTracePointsTot, 0.0};
168  }
169 
170  v_AdvectTraceFlux(nConvectiveFields, fields, advVel, inarray, numflux, time,
171  pFwd, pBwd);
172 
173  // Evaulate <\phi, \hat{F}\cdot n> - OutField[i]
174  for (int i = 0; i < nConvectiveFields; ++i)
175  {
176  Vmath::Neg(nCoeffs, outarray[i], 1);
177 
178  timer.Start();
179  fields[i]->AddTraceIntegral(numflux[i], outarray[i]);
180  timer.Stop();
181  timer.AccumulateRegion("AdvWeakDG:_AddTraceIntegral", 10);
182 
183  timer.Start();
184  fields[i]->MultiplyByElmtInvMass(outarray[i], outarray[i]);
185  timer.Stop();
186  timer.AccumulateRegion("AdvWeakDG:_MultiplyByElmtInvMass", 10);
187  }
188 }
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
Definition: Advection.h:215
virtual void v_AdvectTraceFlux(const int nConvective, 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 >> &TraceFlux, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray) override
Advects Trace Flux.
virtual 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 > &VolumeFlux, const NekDouble &time) override
Advects Volume Flux.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

References Nektar::LibUtilities::Timer::AccumulateRegion(), Vmath::Fill(), Nektar::SolverUtils::Advection::m_spaceDim, Vmath::Neg(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), v_AdvectTraceFlux(), and v_AdvectVolumeFlux().

Referenced by v_Advect().

◆ v_AdvectTraceFlux()

void Nektar::SolverUtils::AdvectionWeakDG::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 
)
overrideprotectedvirtual

Advects Trace Flux.

calculate the advection flux in the cell the trace integration

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 190 of file AdvectionWeakDG.cpp.

198 {
199  boost::ignore_unused(advVel, time);
200  int nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
201 
202  ASSERTL1(m_riemann, "Riemann solver must be provided for AdvectionWeakDG.");
203 
204  // Store forwards/backwards space along trace space
205  Array<OneD, Array<OneD, NekDouble>> Fwd(nConvectiveFields);
206  Array<OneD, Array<OneD, NekDouble>> Bwd(nConvectiveFields);
207 
209  {
210  for (int i = 0; i < nConvectiveFields; ++i)
211  {
212  Fwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
213  Bwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
214  fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
215  }
216  }
217  else
218  {
219  for (int i = 0; i < nConvectiveFields; ++i)
220  {
221  Fwd[i] = pFwd[i];
222  Bwd[i] = pBwd[i];
223  }
224  }
225 
226  LibUtilities::Timer timer;
227  timer.Start();
228  m_riemann->Solve(m_spaceDim, Fwd, Bwd, TraceFlux);
229  timer.Stop();
230  timer.AccumulateRegion("AdvWeakDG:_Riemann", 10);
231 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:213
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray

References Nektar::LibUtilities::Timer::AccumulateRegion(), ASSERTL1, Nektar::SolverUtils::Advection::m_riemann, Nektar::SolverUtils::Advection::m_spaceDim, Nektar::NullNekDoubleArrayOfArray, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by v_AdvectCoeffs().

◆ v_AdvectVolumeFlux()

virtual void Nektar::SolverUtils::AdvectionWeakDG::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 
)
inlineoverrideprotectedvirtual

Advects Volume Flux.

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 86 of file AdvectionWeakDG.h.

92  {
93  boost::ignore_unused(nConvectiveFields, fields, advVel, inarray,
94  VolumeFlux, time);
95  m_fluxVector(inarray, VolumeFlux);
96  }
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:211

References Nektar::SolverUtils::Advection::m_fluxVector.

Referenced by v_AdvectCoeffs().

◆ v_InitObject()

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

Initialise AdvectionWeakDG objects and store them before starting the time-stepping.

Parameters
pSessionPointer to session reader.
pFieldsPointer to fields.

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 62 of file AdvectionWeakDG.cpp.

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

References Nektar::SolverUtils::Advection::v_InitObject().

Member Data Documentation

◆ type

std::string Nektar::SolverUtils::AdvectionWeakDG::type
static
Initial value:
=
"WeakDG", AdvectionWeakDG::create, "Weak DG")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static AdvectionSharedPtr create(std::string advType)
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47

Definition at line 55 of file AdvectionWeakDG.h.