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 ()
 
void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 Initialise AdvectionWeakDG objects and store them before starting the time-stepping. More...
 
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 42 of file AdvectionWeakDG.h.

Constructor & Destructor Documentation

◆ AdvectionWeakDG()

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

Definition at line 47 of file AdvectionWeakDG.cpp.

48{
49}

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

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

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

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

77 {
78 m_fluxVector(inarray, VolumeFlux);
79 }
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:168

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

Referenced by AdvectCoeffs().

◆ create()

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

Definition at line 45 of file AdvectionWeakDG.h.

46 {
48 }
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:54

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

85{
86 LibUtilities::Timer timer1;
87 timer1.Start();
88
89 size_t nCoeffs = fields[0]->GetNcoeffs();
90
91 Array<OneD, Array<OneD, NekDouble>> tmp{size_t(nConvectiveFields)};
92 for (int i = 0; i < nConvectiveFields; ++i)
93 {
94 tmp[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
95 }
96 LibUtilities::Timer timer;
97 timer.Start();
98 AdvectCoeffs(nConvectiveFields, fields, advVel, inarray, tmp, time, pFwd,
99 pBwd);
100 timer.Stop();
101 timer.AccumulateRegion("AdvectCoeffs", 2);
102
103 timer.Start();
104 for (int i = 0; i < nConvectiveFields; ++i)
105 {
106 fields[i]->BwdTrans(tmp[i], outarray[i]);
107 }
108 timer.Stop();
109 timer.AccumulateRegion("AdvWeakDG:_BwdTrans", 10);
110 timer1.Stop();
111 timer1.AccumulateRegion("AdvWeakDG:All", 10);
112}
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 58 of file AdvectionWeakDG.cpp.

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

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.
static AdvectionSharedPtr create(std::string advType)
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43

Definition at line 50 of file AdvectionWeakDG.h.