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

#include <AlternateSkewAdvection.h>

Inheritance diagram for Nektar::AlternateSkewAdvection:
[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

 AlternateSkewAdvection ()
 
 ~AlternateSkewAdvection () 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

int m_advectioncalls
 
bool m_SingleMode
 
bool m_HalfMode
 

Friends

class MemoryManager< AlternateSkewAdvection >
 

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 43 of file AlternateSkewAdvection.h.

Constructor & Destructor Documentation

◆ AlternateSkewAdvection()

Nektar::AlternateSkewAdvection::AlternateSkewAdvection ( )
protected

Constructor. Creates ...

Parameters

param

Definition at line 53 of file AlternateSkewAdvection.cpp.

53 : Advection()
54{
55}

◆ ~AlternateSkewAdvection()

Nektar::AlternateSkewAdvection::~AlternateSkewAdvection ( )
overrideprotected

Definition at line 57 of file AlternateSkewAdvection.cpp.

58{
59}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 50 of file AlternateSkewAdvection.h.

51 {
53 }
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::AlternateSkewAdvection::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 69 of file AlternateSkewAdvection.cpp.

78{
79 // use dimension of Velocity vector to dictate dimension of operation
80 int ndim = advVel.size();
81 int nPointsTot = fields[0]->GetNpoints();
82 Array<OneD, Array<OneD, NekDouble>> velocity(ndim);
83 for (int i = 0; i < ndim; ++i)
84 {
85 if (fields[i]->GetWaveSpace() && !m_SingleMode && !m_HalfMode)
86 {
87 velocity[i] = Array<OneD, NekDouble>(nPointsTot, 0.0);
88 fields[i]->HomogeneousBwdTrans(nPointsTot, advVel[i], velocity[i]);
89 }
90 else
91 {
92 velocity[i] = advVel[i];
93 }
94 }
95 for (int n = 0; n < nConvectiveFields; ++n)
96 {
97 // ToDo: here we should add a check that V has right dimension
98 Array<OneD, NekDouble> gradV0, gradV1, gradV2, tmp, Up;
99
100 gradV0 = Array<OneD, NekDouble>(nPointsTot);
101 tmp = Array<OneD, NekDouble>(nPointsTot);
102
103 // Evaluate V\cdot Grad(u)
104 switch (ndim)
105 {
106 case 1:
107 if (m_advectioncalls % 2 == 0)
108 {
109 fields[0]->PhysDeriv(inarray[n], gradV0);
110 Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1,
111 outarray[n], 1);
112 }
113 else
114 {
115 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1,
116 gradV0, 1);
117 fields[0]->PhysDeriv(gradV0, outarray[n]);
118 }
119 Vmath::Smul(nPointsTot, 0.5, outarray[n], 1, outarray[n],
120 1); // must be mult by 0.5????
121 break;
122 case 2:
123 gradV1 = Array<OneD, NekDouble>(nPointsTot);
124 if (m_advectioncalls % 2 == 0)
125 {
126 fields[0]->PhysDeriv(inarray[n], gradV0, gradV1);
127 Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1,
128 outarray[n], 1);
129 Vmath::Vvtvp(nPointsTot, gradV1, 1, velocity[1], 1,
130 outarray[n], 1, outarray[n], 1);
131 }
132 else
133 {
134 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1,
135 gradV0, 1);
136 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[1], 1,
137 gradV1, 1);
138 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
139 gradV0, outarray[n]);
140 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
141 gradV1, tmp);
142 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
143 1);
144 }
145 Vmath::Smul(nPointsTot, 1.0, outarray[n], 1, outarray[n],
146 1); // must be mult by 0.5????
147 break;
148 case 3:
149 gradV1 = Array<OneD, NekDouble>(nPointsTot);
150 gradV2 = Array<OneD, NekDouble>(nPointsTot);
151
152 // outarray[n] = 1/2(u*du/dx + v*du/dy + w*du/dz + duu/dx +
153 // duv/dy + duw/dz)
154
155 if (fields[0]->GetWaveSpace() == true)
156 {
157 if (m_advectioncalls % 2 == 0)
158 {
159 // vector reused to avoid even more memory requirements
160 // names may be misleading
161 fields[0]->PhysDeriv(inarray[n], gradV0, gradV1,
162 gradV2);
163 fields[0]->HomogeneousBwdTrans(nPointsTot, gradV0, tmp);
164 Vmath::Vmul(nPointsTot, tmp, 1, velocity[0], 1,
165 outarray[n], 1); // + u*du/dx
166 fields[0]->HomogeneousBwdTrans(nPointsTot, gradV1, tmp);
167 Vmath::Vvtvp(nPointsTot, tmp, 1, velocity[1], 1,
168 outarray[n], 1, outarray[n],
169 1); // + v*du/dy
170 fields[0]->HomogeneousBwdTrans(nPointsTot, gradV2, tmp);
171 Vmath::Vvtvp(nPointsTot, tmp, 1, velocity[2], 1,
172 outarray[n], 1, outarray[n],
173 1); // + w*du/dz
174 }
175 else
176 {
177 Up = Array<OneD, NekDouble>(nPointsTot);
178 fields[0]->HomogeneousBwdTrans(nPointsTot, inarray[n],
179 Up);
180 Vmath::Vmul(nPointsTot, Up, 1, velocity[0], 1, gradV0,
181 1);
182 Vmath::Vmul(nPointsTot, Up, 1, velocity[1], 1, gradV1,
183 1);
184 Vmath::Vmul(nPointsTot, Up, 1, velocity[2], 1, gradV2,
185 1);
186
187 fields[0]->SetWaveSpace(false);
188 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
189 gradV0, outarray[n]); // duu/dx
190 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
191 gradV1, tmp); // duv/dy
192 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
193 outarray[n], 1);
194 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
195 gradV2, tmp); // duw/dz
196 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
197 outarray[n], 1);
198 fields[0]->SetWaveSpace(true);
199 }
200
201 Vmath::Smul(nPointsTot, 1.0, outarray[n], 1, tmp,
202 1); // must be mult by 0.5????
203 fields[0]->HomogeneousFwdTrans(nPointsTot, tmp,
204 outarray[n]);
205 }
206 else
207 {
208 if (m_advectioncalls % 2 == 0)
209 {
210 fields[0]->PhysDeriv(inarray[n], gradV0, gradV1,
211 gradV2);
212 Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1,
213 outarray[n], 1);
214 Vmath::Vvtvp(nPointsTot, gradV1, 1, velocity[1], 1,
215 outarray[n], 1, outarray[n], 1);
216 Vmath::Vvtvp(nPointsTot, gradV2, 1, velocity[2], 1,
217 outarray[n], 1, outarray[n], 1);
218 }
219 else
220 {
221 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1,
222 gradV0, 1);
223 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[1], 1,
224 gradV1, 1);
225 Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[2], 1,
226 gradV2, 1);
227 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
228 gradV0, outarray[n]);
229 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
230 gradV1, tmp);
231 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
232 outarray[n], 1);
233 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
234 gradV2, tmp);
235 Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
236 outarray[n], 1);
237 }
238 Vmath::Smul(nPointsTot, 1.0, outarray[n], 1, outarray[n],
239 1); // must be mult by 0.5????
240 }
241 break;
242 default:
243 ASSERTL0(false, "dimension unknown");
244 }
245 }
246}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
const std::vector< NekDouble > velocity
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
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 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, Nektar::MultiRegions::DirCartesianMap, m_advectioncalls, m_HalfMode, m_SingleMode, Vmath::Smul(), Vmath::Vadd(), Nektar::MovementTests::velocity, Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_InitObject()

void Nektar::AlternateSkewAdvection::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 61 of file AlternateSkewAdvection.cpp.

64{
65 pSession->MatchSolverInfo("ModeType", "SingleMode", m_SingleMode, false);
66 pSession->MatchSolverInfo("ModeType", "HalfMode", m_HalfMode, false);
67}

References m_HalfMode, and m_SingleMode.

Friends And Related Function Documentation

◆ MemoryManager< AlternateSkewAdvection >

friend class MemoryManager< AlternateSkewAdvection >
friend

Definition at line 1 of file AlternateSkewAdvection.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 55 of file AlternateSkewAdvection.h.

◆ className2

std::string Nektar::AlternateSkewAdvection::className2
static

Definition at line 56 of file AlternateSkewAdvection.h.

◆ m_advectioncalls

int Nektar::AlternateSkewAdvection::m_advectioncalls
private

Definition at line 79 of file AlternateSkewAdvection.h.

Referenced by v_Advect().

◆ m_HalfMode

bool Nektar::AlternateSkewAdvection::m_HalfMode
private

Definition at line 81 of file AlternateSkewAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_SingleMode

bool Nektar::AlternateSkewAdvection::m_SingleMode
private

Definition at line 80 of file AlternateSkewAdvection.h.

Referenced by v_Advect(), and v_InitObject().