Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
Nektar::SkewSymmetricAdvection Class Reference

#include <SkewSymmetricAdvection.h>

Inheritance diagram for Nektar::SkewSymmetricAdvection:
[legend]

Static Public Member Functions

static SolverUtils::AdvectionSharedPtr create (std::string)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
static std::string className2
 

Protected Member Functions

 SkewSymmetricAdvection ()
 
 ~SkewSymmetricAdvection () override
 
void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 Initialises the advection object. More...
 
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) override
 Advects a vector field. 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...
 

Private Attributes

bool m_homogen_dealiasing
 
bool m_SingleMode
 
bool m_HalfMode
 

Friends

class MemoryManager< SkewSymmetricAdvection >
 

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...
 
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)
 
- 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 44 of file SkewSymmetricAdvection.h.

Constructor & Destructor Documentation

◆ SkewSymmetricAdvection()

Nektar::SkewSymmetricAdvection::SkewSymmetricAdvection ( )
protected

Definition at line 48 of file SkewSymmetricAdvection.cpp.

48 : Advection()
49
50{
51}

◆ ~SkewSymmetricAdvection()

Nektar::SkewSymmetricAdvection::~SkewSymmetricAdvection ( )
overrideprotected

Definition at line 56 of file SkewSymmetricAdvection.cpp.

57{
58}

Member Function Documentation

◆ create()

static SolverUtils::AdvectionSharedPtr Nektar::SkewSymmetricAdvection::create ( std::string  )
inlinestatic

Creates an instance of this class.

Definition at line 51 of file SkewSymmetricAdvection.h.

52 {
54 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ v_Advect()

void Nektar::SkewSymmetricAdvection::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

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 77 of file SkewSymmetricAdvection.cpp.

86{
87 // use dimension of Velocity vector to dictate dimension of operation
88 int ndim = advVel.size();
89 int nqtot = fields[0]->GetTotPoints();
90 ASSERTL1(nConvectiveFields == inarray.size(),
91 "Number of convective fields and Inarray are not compatible");
92
93 Array<OneD, Array<OneD, NekDouble>> velocity(ndim);
94 for (int i = 0; i < ndim; ++i)
95 {
96 if (fields[i]->GetWaveSpace() && !m_SingleMode && !m_HalfMode)
97 {
98 velocity[i] = Array<OneD, NekDouble>(nqtot, 0.0);
99 fields[i]->HomogeneousBwdTrans(nqtot, advVel[i], velocity[i]);
100 }
101 else
102 {
103 velocity[i] = advVel[i];
104 }
105 }
106
107 for (int n = 0; n < nConvectiveFields; ++n)
108 {
109 // ToDo: here we should add a check that V has right dimension
110
111 int nPointsTot = fields[0]->GetNpoints();
112 Array<OneD, NekDouble> gradV0, gradV1, gradV2, tmp, Up;
113
114 gradV0 = Array<OneD, NekDouble>(nPointsTot);
115 tmp = Array<OneD, NekDouble>(nPointsTot);
116
117 // Evaluate V\cdot Grad(u)
118 switch (ndim)
119 {
120 case 1:
121 fields[0]->PhysDeriv(inarray[n], gradV0);
122 Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1, outarray[n],
123 1);
124 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1, gradV0,
125 1);
126 fields[0]->PhysDeriv(gradV0, tmp);
127 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n], 1);
128 Vmath::Smul(nPointsTot, 0.5, outarray[n], 1, outarray[n], 1);
129 break;
130 case 2:
131 gradV1 = Array<OneD, NekDouble>(nPointsTot);
132 fields[0]->PhysDeriv(inarray[n], gradV0, gradV1);
133 Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1, outarray[n],
134 1);
135 Vmath::Vvtvp(nPointsTot, gradV1, 1, velocity[1], 1, outarray[n],
136 1, outarray[n], 1);
137 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1, gradV0,
138 1);
139 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[1], 1, gradV1,
140 1);
141 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0], gradV0,
142 tmp);
143 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n], 1);
144 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1], gradV1,
145 tmp);
146 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n], 1);
147 Vmath::Smul(nPointsTot, 0.5, outarray[n], 1, outarray[n], 1);
148 break;
149 case 3:
150 gradV1 = Array<OneD, NekDouble>(nPointsTot);
151 gradV2 = Array<OneD, NekDouble>(nPointsTot);
152
153 fields[0]->PhysDeriv(inarray[n], gradV0, gradV1, gradV2);
154
155 // outarray[n] = 1/2(u*du/dx + v*du/dy + w*du/dz + duu/dx +
156 // duv/dy + duw/dz)
157
158 if (m_homogen_dealiasing == true &&
159 fields[0]->GetWaveSpace() == false)
160 {
161 fields[0]->DealiasedProd(nPointsTot, velocity[0], gradV0,
162 gradV0);
163 fields[0]->DealiasedProd(nPointsTot, velocity[1], gradV1,
164 gradV1);
165 fields[0]->DealiasedProd(nPointsTot, velocity[2], gradV2,
166 gradV2);
167 Vmath::Vadd(nPointsTot, gradV0, 1, gradV1, 1, outarray[n],
168 1);
169 Vmath::Vadd(nPointsTot, gradV2, 1, outarray[n], 1,
170 outarray[n], 1);
171 fields[0]->DealiasedProd(nPointsTot, inarray[n],
172 velocity[0], gradV0);
173 fields[0]->DealiasedProd(nPointsTot, inarray[n],
174 velocity[1], gradV1);
175 fields[0]->DealiasedProd(nPointsTot, inarray[n],
176 velocity[2], gradV2);
177 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
178 gradV0, tmp);
179 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
180 1);
181 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
182 gradV1, tmp);
183 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
184 1);
185 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
186 gradV2, tmp);
187 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
188 1);
189 Vmath::Smul(nPointsTot, 0.5, outarray[n], 1, outarray[n],
190 1);
191 }
192 else if (fields[0]->GetWaveSpace() == true &&
193 m_homogen_dealiasing == false)
194 {
195 Up = Array<OneD, NekDouble>(nPointsTot);
196 // vector reused to avoid even more memory requirements
197 // names may be misleading
198 fields[0]->HomogeneousBwdTrans(nPointsTot, gradV0, tmp);
199 Vmath::Vmul(nPointsTot, tmp, 1, velocity[0], 1, outarray[n],
200 1); // + u*du/dx
201 fields[0]->HomogeneousBwdTrans(nPointsTot, gradV1, tmp);
202 Vmath::Vvtvp(nPointsTot, tmp, 1, velocity[1], 1,
203 outarray[n], 1, outarray[n], 1); // + v*du/dy
204 fields[0]->HomogeneousBwdTrans(nPointsTot, gradV2, tmp);
205 Vmath::Vvtvp(nPointsTot, tmp, 1, velocity[2], 1,
206 outarray[n], 1, outarray[n], 1); // + w*du/dz
207
208 fields[0]->HomogeneousBwdTrans(nPointsTot, inarray[n], Up);
209 Vmath::Vmul(nPointsTot, Up, 1, velocity[0], 1, gradV0, 1);
210 Vmath::Vmul(nPointsTot, Up, 1, velocity[1], 1, gradV1, 1);
211 Vmath::Vmul(nPointsTot, Up, 1, velocity[2], 1, gradV2, 1);
212
213 fields[0]->SetWaveSpace(false);
214 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
215 gradV0, tmp); // duu/dx
216 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
217 1);
218 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
219 gradV1, tmp); // duv/dy
220 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
221 1);
222 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
223 gradV2, tmp); // duw/dz
224 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
225 1);
226 fields[0]->SetWaveSpace(true);
227
228 Vmath::Smul(nPointsTot, 0.5, outarray[n], 1, tmp, 1);
229 fields[0]->HomogeneousFwdTrans(nPointsTot, tmp,
230 outarray[n]);
231 }
232 else if (fields[0]->GetWaveSpace() == false &&
233 m_homogen_dealiasing == false)
234 {
235 Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1,
236 outarray[n], 1);
237 Vmath::Vvtvp(nPointsTot, gradV1, 1, velocity[1], 1,
238 outarray[n], 1, outarray[n], 1);
239 Vmath::Vvtvp(nPointsTot, gradV2, 1, velocity[2], 1,
240 outarray[n], 1, outarray[n], 1);
241 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1,
242 gradV0, 1);
243 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[1], 1,
244 gradV1, 1);
245 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[2], 1,
246 gradV2, 1);
247 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
248 gradV0, tmp);
249 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
250 1);
251 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
252 gradV1, tmp);
253 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
254 1);
255 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
256 gradV2, tmp);
257 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
258 1);
259 Vmath::Smul(nPointsTot, 0.5, outarray[n], 1, outarray[n],
260 1);
261 }
262 else
263 {
264 ASSERTL0(false,
265 "Dealiasing is not allowed in combination "
266 "with the Skew-Symmetric advection form for "
267 "efficiency reasons.");
268 }
269 break;
270 default:
271 ASSERTL0(false, "dimension unknown");
272 }
273
274 Vmath::Neg(nqtot, outarray[n], 1);
275 }
276}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
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.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::DirCartesianMap, m_HalfMode, m_homogen_dealiasing, m_SingleMode, Vmath::Neg(), Vmath::Smul(), Vmath::Vadd(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_InitObject()

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

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 from Nektar::SolverUtils::Advection.

Definition at line 63 of file SkewSymmetricAdvection.cpp.

66{
67 Advection::v_InitObject(pSession, pFields);
68
69 m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
70 pSession->MatchSolverInfo("ModeType", "SingleMode", m_SingleMode, false);
71 pSession->MatchSolverInfo("ModeType", "HalfMode", m_HalfMode, false);
72}
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:295

References m_HalfMode, m_homogen_dealiasing, m_SingleMode, and Nektar::SolverUtils::Advection::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< SkewSymmetricAdvection >

friend class MemoryManager< SkewSymmetricAdvection >
friend

Definition at line 1 of file SkewSymmetricAdvection.h.

Member Data Documentation

◆ className

string Nektar::SkewSymmetricAdvection::className
static
Initial value:
=
"SkewSymmetric", SkewSymmetricAdvection::create, "Skew Symmetric")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static SolverUtils::AdvectionSharedPtr create(std::string)
Creates an instance of this class.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43

Name of class.

Definition at line 57 of file SkewSymmetricAdvection.h.

◆ className2

std::string Nektar::SkewSymmetricAdvection::className2
static

Definition at line 58 of file SkewSymmetricAdvection.h.

◆ m_HalfMode

bool Nektar::SkewSymmetricAdvection::m_HalfMode
private

Definition at line 83 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_homogen_dealiasing

bool Nektar::SkewSymmetricAdvection::m_homogen_dealiasing
private

Definition at line 81 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_SingleMode

bool Nektar::SkewSymmetricAdvection::m_SingleMode
private

Definition at line 82 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().