Nektar++
Public Member Functions | 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]

Public Member Functions

SOLVER_UTILS_EXPORT void 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)
 
SOLVER_UTILS_EXPORT void 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)
 
SOLVER_UTILS_EXPORT void AdvectVolumeFlux (const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &VolumeFlux)
 
- 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...
 
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)
 
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)
 

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...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Advection
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_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

- 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

◆ AdvectCoeffs()

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

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 AdvectVolumeFlux(inarray, fluxvector);
149 timer.Stop();
150 timer.AccumulateRegion("AdvWeakDG:_fluxVector", 3);
151
152 timer.Start();
153 // Get the advection part (without numerical flux)
154 for (int i = 0; i < nConvectiveFields; ++i)
155 {
156 Vmath::Fill(outarray[i].size(), 0.0, outarray[i], 1);
157 fields[i]->IProductWRTDerivBase(fluxvector[i], outarray[i]);
158 }
159 timer.Stop();
160 timer.AccumulateRegion("AdvWeakDG:_IProductWRTDerivBase", 10);
161
162 Array<OneD, Array<OneD, NekDouble>> numflux{size_t(nConvectiveFields)};
163
164 for (int i = 0; i < nConvectiveFields; ++i)
165 {
166 numflux[i] = Array<OneD, NekDouble>{nTracePointsTot, 0.0};
167 }
168
169 AdvectTraceFlux(nConvectiveFields, fields, advVel, inarray, numflux, time,
170 pFwd, pBwd);
171
172 // Evaulate <\phi, \hat{F}\cdot n> - OutField[i]
173 for (int i = 0; i < nConvectiveFields; ++i)
174 {
175 Vmath::Neg(nCoeffs, outarray[i], 1);
176
177 timer.Start();
178 fields[i]->AddTraceIntegral(numflux[i], outarray[i]);
179 timer.Stop();
180 timer.AccumulateRegion("AdvWeakDG:_AddTraceIntegral", 10);
181
182 timer.Start();
183 fields[i]->MultiplyByElmtInvMass(outarray[i], outarray[i]);
184 timer.Stop();
185 timer.AccumulateRegion("AdvWeakDG:_MultiplyByElmtInvMass", 10);
186 }
187}
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
Definition: Advection.h:174
SOLVER_UTILS_EXPORT void AdvectVolumeFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &VolumeFlux)
SOLVER_UTILS_EXPORT void 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:513
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43

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

Referenced by v_Advect().

◆ AdvectTraceFlux()

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

Definition at line 192 of file AdvectionWeakDG.cpp.

200{
201 boost::ignore_unused(advVel, time);
202 int nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
203
204 ASSERTL1(m_riemann, "Riemann solver must be provided for AdvectionWeakDG.");
205
206 // Store forwards/backwards space along trace space
207 Array<OneD, Array<OneD, NekDouble>> Fwd(nConvectiveFields);
208 Array<OneD, Array<OneD, NekDouble>> Bwd(nConvectiveFields);
209
211 {
212 for (int i = 0; i < nConvectiveFields; ++i)
213 {
214 Fwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
215 Bwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
216 fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
217 }
218 }
219 else
220 {
221 for (int i = 0; i < nConvectiveFields; ++i)
222 {
223 Fwd[i] = pFwd[i];
224 Bwd[i] = pBwd[i];
225 }
226 }
227
228 LibUtilities::Timer timer;
229 timer.Start();
230 m_riemann->Solve(m_spaceDim, Fwd, Bwd, TraceFlux);
231 timer.Stop();
232 timer.AccumulateRegion("AdvWeakDG:_Riemann", 10);
233}
#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:172
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 AdvectCoeffs().

◆ AdvectVolumeFlux()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::AdvectionWeakDG::AdvectVolumeFlux ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
TensorOfArray3D< NekDouble > &  VolumeFlux 
)
inline

Definition at line 79 of file AdvectionWeakDG.h.

82 {
83 m_fluxVector(inarray, VolumeFlux);
84 }
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:170

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

Referenced by AdvectCoeffs().

◆ 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);
53 }
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:56

References AdvectionWeakDG().

◆ 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 AdvectCoeffs(nConvectiveFields, fields, advVel, inarray, tmp, time, pFwd,
103 pBwd);
104 timer.Stop();
105 timer.AccumulateRegion("AdvectCoeffs", 2);
106
107 timer.Start();
108 for (int i = 0; i < nConvectiveFields; ++i)
109 {
110 fields[i]->BwdTrans(tmp[i], outarray[i]);
111 }
112 timer.Stop();
113 timer.AccumulateRegion("AdvWeakDG:_BwdTrans", 10);
114 timer1.Stop();
115 timer1.AccumulateRegion("AdvWeakDG:All", 10);
116}
SOLVER_UTILS_EXPORT void 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(), AdvectCoeffs(), Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

◆ 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:299

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.