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)
 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)
 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)
 
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)
 
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)
 
virtual void v_AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 
- Protected Member Functions inherited from Nektar::SolverUtils::Advection
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_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 50 of file AdvectionWeakDG.cpp.

51  {
52  }

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:291

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 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 244 of file AdvectionWeakDG.cpp.

249  {
250  MultiRegions::ExpListSharedPtr explist = pFields[0];
251  std::shared_ptr<LocalRegions::ExpansionVector> pexp = explist->GetExp();
252  int nTotElmt = (*pexp).size();
253  int nElmtPnt,nElmtCoef;
254 
255  SNekMatSharedPtr tmpGmtx;
256  DNekMatSharedPtr ElmtMat;
257 
258  Array<OneD,NekSingle> GMat_data;
259  Array<OneD,NekDouble> Elmt_data;
260  Array<OneD,NekSingle> Elmt_dataSingle;
261 
262  Array<OneD, DNekMatSharedPtr> mtxPerVar(nTotElmt);
263  Array<OneD, int > elmtpnts(nTotElmt);
264  Array<OneD, int > elmtcoef(nTotElmt);
265  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
266  {
267  nElmtCoef = (*pexp)[nelmt]->GetNcoeffs();
268  nElmtPnt = (*pexp)[nelmt]->GetTotPoints();
269  elmtpnts[nelmt] = nElmtPnt;
270  elmtcoef[nelmt] = nElmtCoef;
271  mtxPerVar[nelmt] =MemoryManager<DNekMat>
272  ::AllocateSharedPtr(nElmtCoef, nElmtPnt);
273  }
274 
275  Array<OneD, DNekMatSharedPtr> mtxPerVarCoeff(nTotElmt);
276  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
277  {
278  nElmtCoef = elmtcoef[nelmt];
279  mtxPerVarCoeff[nelmt] =MemoryManager<DNekMat>
280  ::AllocateSharedPtr(nElmtCoef, nElmtCoef);
281  }
282 
283  for(int m = 0; m < nConvectiveFields; m++)
284  {
285  for(int n = 0; n < nConvectiveFields; n++)
286  {
287  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
288  {
289  (*mtxPerVarCoeff[nelmt]) = 0.0;
290  (*mtxPerVar[nelmt]) = 0.0;
291  }
292 
293  explist->GetMatIpwrtDeriveBase(ElmtJacArray[m][n],
294  mtxPerVar);
295  // Memory can be saved by reusing mtxPerVar
296  explist->AddRightIPTBaseMatrix(mtxPerVar,mtxPerVarCoeff);
297 
298  for(int nelmt = 0; nelmt < nTotElmt; nelmt++)
299  {
300  nElmtCoef = elmtcoef[nelmt];
301  nElmtPnt = elmtpnts[nelmt];
302  int ntotDofs = nElmtCoef*nElmtCoef;
303 
304  if(Elmt_dataSingle.size()<ntotDofs)
305  {
306  Elmt_dataSingle = Array<OneD, NekSingle> (ntotDofs);
307  }
308 
309  tmpGmtx = gmtxarray[m][n]->GetBlock(nelmt,nelmt);
310  ElmtMat = mtxPerVarCoeff[nelmt];
311 
312  GMat_data = tmpGmtx->GetPtr();
313  Elmt_data = ElmtMat->GetPtr();
314 
315  for(int i=0;i<ntotDofs;i++)
316  {
317  Elmt_dataSingle[i] = NekSingle( Elmt_data[i] );
318  }
319 
320  Vmath::Vadd(ntotDofs,GMat_data,1,Elmt_dataSingle,1,GMat_data,1);
321  }
322  }
323  }
324  }
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: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

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 
)
protectedvirtual

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 79 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 
102  nConvectiveFields, fields, advVel, inarray, tmp, time,
103  pFwd, pBwd);
104 
105  // why was this broken in many loops over convective fields?
106  // this is terrible for locality
107  LibUtilities::Timer timer;
108  timer.Start();
109  for (int i = 0; i < nConvectiveFields; ++i)
110  {
111  fields[i]->BwdTrans(tmp[i], outarray[i]);
112  }
113  timer.Stop();
114  timer.AccumulateRegion("AdvWeakDG:_BwdTrans",1);
115  timer1.Stop();
116  timer1.AccumulateRegion("AdvWeakDG:All");
117  }
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)

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 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 119 of file AdvectionWeakDG.cpp.

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

Definition at line 197 of file AdvectionWeakDG.cpp.

206  {
207  boost::ignore_unused(advVel, time);
208  int nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
209 
211  "Riemann solver must be provided for AdvectionWeakDG.");
212 
213  // Store forwards/backwards space along trace space
214  Array<OneD, Array<OneD, NekDouble>> Fwd(nConvectiveFields);
215  Array<OneD, Array<OneD, NekDouble>> Bwd(nConvectiveFields);
216 
217  if (pFwd == NullNekDoubleArrayOfArray ||
219  {
220  for (int i = 0; i < nConvectiveFields; ++i)
221  {
222  Fwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
223  Bwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
224  fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
225  }
226  }
227  else
228  {
229  for (int i = 0; i < nConvectiveFields; ++i)
230  {
231  Fwd[i] = pFwd[i];
232  Bwd[i] = pBwd[i];
233  }
234  }
235 
236  LibUtilities::Timer timer;
237  timer.Start();
238  m_riemann->Solve(m_spaceDim, Fwd, Bwd, TraceFlux);
239  timer.Stop();
240  timer.AccumulateRegion("AdvWeakDG:_Riemann",1);
241 
242  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:223
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 > &  VolumeFlux,
const NekDouble time 
)
inlineprotectedvirtual

Definition at line 88 of file AdvectionWeakDG.h.

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

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 
)
protectedvirtual

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 61 of file AdvectionWeakDG.cpp.

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

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

Member Data Documentation

◆ type

std::string Nektar::SolverUtils::AdvectionWeakDG::type
static
Initial value:
RegisterCreatorFunction("WeakDG", AdvectionWeakDG::create)
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.