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 ()
 
virtual ~AlternateSkewAdvection ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 Initialises the advection object. More...
 
virtual 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 ( )
protectedvirtual

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 71 of file AlternateSkewAdvection.cpp.

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

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, m_advectioncalls, m_HalfMode, m_SingleMode, Vmath::Smul(), Vmath::Vadd(), 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 boost::ignore_unused(fields);
66
67 pSession->MatchSolverInfo("ModeType", "SingleMode", m_SingleMode, false);
68 pSession->MatchSolverInfo("ModeType", "HalfMode", m_HalfMode, false);
69}

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:198
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47

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().