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

#include <MappingXofXZ.h>

Inheritance diagram for Nektar::GlobalMapping::MappingXofXZ:
[legend]

Static Public Member Functions

static GLOBAL_MAPPING_EXPORT MappingSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::GlobalMapping::Mapping
static GLOBAL_MAPPING_EXPORT MappingSharedPtr Load (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Return a pointer to the mapping, creating it on first call. More...
 

Static Public Attributes

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

Protected Member Functions

 MappingXofXZ (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 
virtual GLOBAL_MAPPING_EXPORT void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping) override
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarToCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarToCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_GetJacobian (Array< OneD, NekDouble > &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_DotGradJacobian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_GetMetricTensor (Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_LowerIndex (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_RaiseIndex (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelContravar (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateGeomInfo () override
 
- Protected Member Functions inherited from Nektar::GlobalMapping::Mapping
GLOBAL_MAPPING_EXPORT Mapping (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Constructor. More...
 
GLOBAL_MAPPING_EXPORT void EvaluateFunction (Array< OneD, MultiRegions::ExpListSharedPtr > pFields, LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
 
GLOBAL_MAPPING_EXPORT void EvaluateTimeFunction (LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
 
virtual GLOBAL_MAPPING_EXPORT void v_GetCartesianCoordinates (Array< OneD, NekDouble > &out0, Array< OneD, NekDouble > &out1, Array< OneD, NekDouble > &out2)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetCoordVelocity (Array< OneD, Array< OneD, NekDouble >> &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_Divergence (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_VelocityLaplacian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble alpha)
 
virtual GLOBAL_MAPPING_EXPORT void v_gradgradU (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_CurlCurlField (Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const bool generalized)
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateMapping (const NekDouble time, const Array< OneD, Array< OneD, NekDouble >> &coords=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &coordsVel=NullNekDoubleArrayOfArray)
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateBCs (const NekDouble time)
 

Friends

class MemoryManager< MappingXofXZ >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::GlobalMapping::Mapping
virtual GLOBAL_MAPPING_EXPORT ~Mapping ()
 Destructor. More...
 
GLOBAL_MAPPING_EXPORT void InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
 Initialise the mapping object. More...
 
GLOBAL_MAPPING_EXPORT void ReplaceField (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Replace the Expansion List used by the mapping. More...
 
GLOBAL_MAPPING_EXPORT void Output (LibUtilities::FieldMetaDataMap &fieldMetaDataMap, const std::string &outname)
 Output function called when a chk or fld file is written. More...
 
GLOBAL_MAPPING_EXPORT void ContravarToCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Convert a contravariant vector to the Cartesian system. More...
 
GLOBAL_MAPPING_EXPORT void CovarToCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Convert a covariant vector to the Cartesian system. More...
 
GLOBAL_MAPPING_EXPORT void ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Convert a contravariant vector to the transformed system. More...
 
GLOBAL_MAPPING_EXPORT void CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Convert a covariant vector to the transformed system. More...
 
GLOBAL_MAPPING_EXPORT void GetCartesianCoordinates (Array< OneD, NekDouble > &out0, Array< OneD, NekDouble > &out1, Array< OneD, NekDouble > &out2)
 Get the Cartesian coordinates in the field. More...
 
GLOBAL_MAPPING_EXPORT void GetJacobian (Array< OneD, NekDouble > &outarray)
 Get the Jacobian of the transformation. More...
 
GLOBAL_MAPPING_EXPORT void DotGradJacobian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the dot product with the gradient of the Jacobian. More...
 
GLOBAL_MAPPING_EXPORT void GetMetricTensor (Array< OneD, Array< OneD, NekDouble >> &outarray)
 Get the metric tensor \(g_{ij}\). More...
 
GLOBAL_MAPPING_EXPORT void GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble >> &outarray)
 Get the inverse of metric tensor \(g^{ij}\). More...
 
GLOBAL_MAPPING_EXPORT void LowerIndex (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Lower index of vector: \(v_{i} = g_{ij}*v^{j}\). More...
 
GLOBAL_MAPPING_EXPORT void RaiseIndex (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Raise index of vector: \(v^{i} = g^{ij}*v_{j}\). More...
 
GLOBAL_MAPPING_EXPORT void ApplyChristoffelContravar (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Apply the Christoffel symbols to a contravariant vector. More...
 
GLOBAL_MAPPING_EXPORT void ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Apply the Christoffel symbols to a covariant vector. More...
 
GLOBAL_MAPPING_EXPORT void GetCoordVelocity (Array< OneD, Array< OneD, NekDouble >> &outarray)
 Obtain the velocity of the coordinates. More...
 
GLOBAL_MAPPING_EXPORT void Divergence (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the generalised divergence operator. More...
 
GLOBAL_MAPPING_EXPORT void VelocityLaplacian (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble alpha=0.0)
 Generalised (correction to the) velocity Laplacian operator. More...
 
GLOBAL_MAPPING_EXPORT void gradgradU (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Second order covariant derivatives of a contravariant vector. More...
 
GLOBAL_MAPPING_EXPORT void CurlCurlField (Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const bool generalized)
 CurlCurl calculated on the whole field. More...
 
GLOBAL_MAPPING_EXPORT bool IsTimeDependent ()
 Get flag defining if mapping is time-dependent. More...
 
GLOBAL_MAPPING_EXPORT void SetTimeDependent (const bool value)
 Set flag defining if mapping is time-dependent. More...
 
GLOBAL_MAPPING_EXPORT bool IsFromFunction ()
 Get flag defining if mapping is defined by a function. More...
 
GLOBAL_MAPPING_EXPORT void SetFromFunction (const bool value)
 Set flag defining if mapping is defined by a function. More...
 
GLOBAL_MAPPING_EXPORT bool HasConstantJacobian ()
 Get flag defining if mapping has constant Jacobian. More...
 
GLOBAL_MAPPING_EXPORT bool IsDefined ()
 Get flag determining if the mapping was defined or is trivial. More...
 
GLOBAL_MAPPING_EXPORT void UpdateBCs (const NekDouble time)
 Update the Dirichlet Boundary Conditions when using Mappings. More...
 
GLOBAL_MAPPING_EXPORT void UpdateMapping (const NekDouble time, const Array< OneD, Array< OneD, NekDouble >> &coords=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &coordsVel=NullNekDoubleArrayOfArray)
 Update the Mapping with new coordinates. More...
 
GLOBAL_MAPPING_EXPORT void UpdateGeomInfo ()
 Recompute the metric terms of the Mapping. More...
 
- Protected Attributes inherited from Nektar::GlobalMapping::Mapping
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 
Array< OneD, Array< OneD, NekDouble > > m_coords
 Array with the Cartesian coordinates. More...
 
Array< OneD, Array< OneD, NekDouble > > m_coordsVel
 Array with the velocity of the coordinates. More...
 
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
 Array with metric terms of the mapping. More...
 
int m_nConvectiveFields
 Number of velocity components. More...
 
std::string m_funcName
 Name of the function containing the coordinates. More...
 
std::string m_velFuncName
 Name of the function containing the velocity of the coordinates. More...
 
bool m_constantJacobian
 Flag defining if the Jacobian is constant. More...
 
bool m_timeDependent
 Flag defining if the Mapping is time-dependent. More...
 
bool m_fromFunction
 Flag defining if the Mapping is defined by a function. More...
 
Array< OneD, Array< OneD, NekDouble > > m_wk1
 
Array< OneD, Array< OneD, NekDouble > > m_wk2
 
Array< OneD, Array< OneD, NekDouble > > m_tmp
 
- Static Protected Attributes inherited from Nektar::GlobalMapping::Mapping
static MappingSharedPtr m_mappingPtr = MappingSharedPtr()
 
static bool m_init = false
 
static bool m_isDefined = false
 

Detailed Description

This class implements a mapping defined by a transformation of the type

\[ \bar{x} = \bar{x}(x,z) \]

\[ \bar{y} = y \]

\[ \bar{z} = z \]

where \((\bar{x},\bar{y},\bar{z})\) are the Cartesian (physical) coordinates and \((x,y,z)\) are the transformed (computational) coordinates.

Definition at line 51 of file MappingXofXZ.h.

Constructor & Destructor Documentation

◆ MappingXofXZ()

Nektar::GlobalMapping::MappingXofXZ::MappingXofXZ ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields 
)
protected

Definition at line 57 of file MappingXofXZ.cpp.

60  : Mapping(pSession, pFields)
61 {
62 }
GLOBAL_MAPPING_EXPORT Mapping(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Constructor.
Definition: Mapping.cpp:59

Member Function Documentation

◆ create()

static GLOBAL_MAPPING_EXPORT MappingSharedPtr Nektar::GlobalMapping::MappingXofXZ::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const TiXmlElement *  pMapping 
)
inlinestatic

Creates an instance of this class.

Definition at line 58 of file MappingXofXZ.h.

62  {
65  p->InitObject(pFields, pMapping);
66  return p;
67  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:50

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::GlobalMapping::MappingSharedPtr, and CellMLToNektar.cellml_metadata::p.

◆ v_ApplyChristoffelContravar()

void Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelContravar ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 299 of file MappingXofXZ.cpp.

302 {
303  int physTot = m_fields[0]->GetTotPoints();
304  int nvel = m_nConvectiveFields;
305  Array<OneD, NekDouble> wk(physTot, 0.0);
306 
307  for (int i = 0; i < nvel; i++)
308  {
309  for (int j = 0; j < nvel; j++)
310  {
311  outarray[i * nvel + j] = Array<OneD, NekDouble>(physTot, 0.0);
312  }
313  }
314 
315  // Calculate non-zero terms
316 
317  // outarray(0,0) = U1 * fxx/fx + U3 * fxz/fx
318  Vmath::Vdiv(physTot, m_GeometricInfo[2], 1, m_GeometricInfo[0], 1, wk, 1);
319  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[0 * nvel + 0], 1);
320  Vmath::Vdiv(physTot, m_GeometricInfo[3], 1, m_GeometricInfo[0], 1, wk, 1);
321  Vmath::Vvtvp(physTot, wk, 1, inarray[2], 1, outarray[0 * nvel + 0], 1,
322  outarray[0 * nvel + 0], 1);
323 
324  // outarray(0,2) = U1 * fxz/fx + U3 * fzz/fx
325  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[0 * nvel + 2], 1);
326  Vmath::Vdiv(physTot, m_GeometricInfo[4], 1, m_GeometricInfo[0], 1, wk, 1);
327  Vmath::Vvtvp(physTot, wk, 1, inarray[2], 1, outarray[0 * nvel + 2], 1,
328  outarray[0 * nvel + 2], 1);
329 }
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:414
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:412
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:406
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:209
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:574
void Vdiv(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:284

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_ApplyChristoffelCovar()

void Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelCovar ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 331 of file MappingXofXZ.cpp.

334 {
335  int physTot = m_fields[0]->GetTotPoints();
336  int nvel = m_nConvectiveFields;
337  Array<OneD, NekDouble> wk(physTot, 0.0);
338 
339  for (int i = 0; i < nvel; i++)
340  {
341  for (int j = 0; j < nvel; j++)
342  {
343  outarray[i * nvel + j] = Array<OneD, NekDouble>(physTot, 0.0);
344  }
345  }
346 
347  // Calculate non-zero terms
348 
349  // outarray(0,0) = U1 * fxx/fx
350  Vmath::Vdiv(physTot, m_GeometricInfo[2], 1, m_GeometricInfo[0], 1, wk, 1);
351  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[0 * nvel + 0], 1);
352 
353  // outarray(0,2) = outarray(2,0) = U1 * fxz/fx
354  Vmath::Vdiv(physTot, m_GeometricInfo[3], 1, m_GeometricInfo[0], 1, wk, 1);
355  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[0 * nvel + 2], 1);
356  Vmath::Vcopy(physTot, outarray[0 * nvel + 2], 1, outarray[2 * nvel + 0], 1);
357 
358  // outarray(2,2) = U1 * fzz/fx
359  Vmath::Vdiv(physTot, m_GeometricInfo[4], 1, m_GeometricInfo[0], 1, wk, 1);
360  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[2 * nvel + 2], 1);
361 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, Vmath::Vcopy(), Vmath::Vdiv(), and Vmath::Vmul().

◆ v_ContravarFromCartesian()

void Nektar::GlobalMapping::MappingXofXZ::v_ContravarFromCartesian ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 116 of file MappingXofXZ.cpp.

119 {
120  int physTot = m_fields[0]->GetTotPoints();
121  Array<OneD, NekDouble> wk(physTot, 0.0);
122 
123  // U1 = u1/fx - fz/fx * u3
124  Vmath::Vdiv(physTot, inarray[0], 1, m_GeometricInfo[0], 1, outarray[0], 1);
125  Vmath::Vdiv(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[0], 1, wk, 1);
126  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, wk, 1);
127  Vmath::Vsub(physTot, outarray[0], 1, wk, 1, outarray[0], 1);
128 
129  // U2 = u2
130  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
131 
132  // U3 = u3
133  Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
134 }
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:419

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Vsub().

◆ v_ContravarToCartesian()

void Nektar::GlobalMapping::MappingXofXZ::v_ContravarToCartesian ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 79 of file MappingXofXZ.cpp.

82 {
83  int physTot = m_fields[0]->GetTotPoints();
84 
85  // U1 = fx*u1 + fz*u3
86  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, inarray[0], 1, outarray[0], 1);
87  Vmath::Vvtvp(physTot, m_GeometricInfo[1], 1, inarray[2], 1, outarray[0], 1,
88  outarray[0], 1);
89 
90  // U2 = u2
91  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
92 
93  // U3 = u3
94  Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
95 }

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_CovarFromCartesian()

void Nektar::GlobalMapping::MappingXofXZ::v_CovarFromCartesian ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 136 of file MappingXofXZ.cpp.

139 {
140  int physTot = m_fields[0]->GetTotPoints();
141 
142  // U1 = u1*fx
143  Vmath::Vmul(physTot, inarray[0], 1, m_GeometricInfo[0], 1, outarray[0], 1);
144 
145  // U2 = u2
146  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
147 
148  // U3 = u3 + fz*u1
149  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, inarray[0], 1, outarray[2], 1);
150  Vmath::Vadd(physTot, inarray[2], 1, outarray[2], 1, outarray[2], 1);
151 }
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:359

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vmul().

◆ v_CovarToCartesian()

void Nektar::GlobalMapping::MappingXofXZ::v_CovarToCartesian ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 97 of file MappingXofXZ.cpp.

100 {
101  int physTot = m_fields[0]->GetTotPoints();
102  Array<OneD, NekDouble> wk(physTot, 0.0);
103 
104  // U1 = u1/fx
105  Vmath::Vdiv(physTot, inarray[0], 1, m_GeometricInfo[0], 1, outarray[0], 1);
106 
107  // U2 = u2
108  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
109 
110  // U3 = u3 - fz/fx*u1
111  Vmath::Vdiv(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[0], 1, wk, 1);
112  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, wk, 1);
113  Vmath::Vsub(physTot, inarray[2], 1, wk, 1, outarray[2], 1);
114 }

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Vsub().

◆ v_DotGradJacobian()

void Nektar::GlobalMapping::MappingXofXZ::v_DotGradJacobian ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 159 of file MappingXofXZ.cpp.

162 {
163  int physTot = m_fields[0]->GetTotPoints();
164 
165  Vmath::Vmul(physTot, m_GeometricInfo[2], 1, inarray[0], 1, outarray, 1);
166  Vmath::Vvtvp(physTot, m_GeometricInfo[3], 1, inarray[2], 1, outarray, 1,
167  outarray, 1);
168 }

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_GetInvMetricTensor()

void Nektar::GlobalMapping::MappingXofXZ::v_GetInvMetricTensor ( Array< OneD, Array< OneD, NekDouble >> &  outarray)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 205 of file MappingXofXZ.cpp.

207 {
208  int physTot = m_fields[0]->GetTotPoints();
209  int nvel = m_nConvectiveFields;
210  Array<OneD, NekDouble> wk(physTot, 0.0);
211 
212  for (int i = 0; i < nvel * nvel; i++)
213  {
214  outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
215  }
216  // Fill diagonal with 1.0
217  for (int i = 0; i < nvel; i++)
218  {
219  Vmath::Sadd(physTot, 1.0, outarray[i + nvel * i], 1,
220  outarray[i + nvel * i], 1);
221  }
222 
223  // G^{13} and G^{31} = -fz/fx
224  Vmath::Vdiv(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[0], 1, wk,
225  1); // fz/fx
226  Vmath::Neg(physTot, wk, 1);
227  Vmath::Vcopy(physTot, wk, 1, outarray[0 * nvel + 2], 1);
228  Vmath::Vcopy(physTot, wk, 1, outarray[2 * nvel + 0], 1);
229 
230  // G^{11} = (1+fz^2)/(fx^2)
231  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1, wk,
232  1); // fz^2
233  Vmath::Vadd(physTot, wk, 1, outarray[0 * nvel + 0], 1,
234  outarray[0 * nvel + 0], 1);
235 
236  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1, wk,
237  1); // fx^2
238  Vmath::Vdiv(physTot, outarray[0 * nvel + 0], 1, wk, 1,
239  outarray[0 * nvel + 0], 1);
240 }
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, Vmath::Neg(), Vmath::Sadd(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vdiv(), and Vmath::Vmul().

◆ v_GetJacobian()

void Nektar::GlobalMapping::MappingXofXZ::v_GetJacobian ( Array< OneD, NekDouble > &  outarray)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 153 of file MappingXofXZ.cpp.

154 {
155  int physTot = m_fields[0]->GetTotPoints();
156  Vmath::Vcopy(physTot, m_GeometricInfo[0], 1, outarray, 1);
157 }

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, and Vmath::Vcopy().

◆ v_GetMetricTensor()

void Nektar::GlobalMapping::MappingXofXZ::v_GetMetricTensor ( Array< OneD, Array< OneD, NekDouble >> &  outarray)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 170 of file MappingXofXZ.cpp.

172 {
173  int physTot = m_fields[0]->GetTotPoints();
174  int nvel = m_nConvectiveFields;
175  Array<OneD, NekDouble> wk(physTot, 0.0);
176 
177  for (int i = 0; i < nvel * nvel; i++)
178  {
179  outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
180  }
181  // Fill G^{22} and G^{33} with 1.0
182  for (int i = 1; i < nvel; i++)
183  {
184  Vmath::Sadd(physTot, 1.0, outarray[i + nvel * i], 1,
185  outarray[i + nvel * i], 1);
186  }
187 
188  // G_{13} and G_{31} = fz*fx
189  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[0], 1, wk,
190  1); // fz*fx
191  Vmath::Vcopy(physTot, wk, 1, outarray[0 * nvel + 2], 1);
192  Vmath::Vcopy(physTot, wk, 1, outarray[2 * nvel + 0], 1);
193 
194  // G^{11} = (fx^2)
195  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1,
196  outarray[0 * nvel + 0], 1);
197 
198  // G^{33} = (1+fz^2)
199  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1, wk,
200  1); // fz^2
201  Vmath::Vadd(physTot, wk, 1, outarray[2 * nvel + 2], 1,
202  outarray[2 * nvel + 2], 1);
203 }

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, Vmath::Sadd(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vmul().

◆ v_InitObject()

void Nektar::GlobalMapping::MappingXofXZ::v_InitObject ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const TiXmlElement *  pMapping 
)
overrideprotectedvirtual

This function initialises the Mapping object. It computes the coordinates and velocity coordinates, initialises the workspace variables, and calls UpdateGeomInfo, which will perform the calculations specific for each type of Mapping.

Parameters
pFieldsExpList array used in the mapping
pMappingxml element describing the mapping

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 67 of file MappingXofXZ.cpp.

70 {
71  Mapping::v_InitObject(pFields, pMapping);
72 
73  m_constantJacobian = false;
74 
76  "Mapping X = X(x,z) needs 3 velocity components.");
77 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
virtual GLOBAL_MAPPING_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
Definition: Mapping.cpp:101
bool m_constantJacobian
Flag defining if the Jacobian is constant.
Definition: Mapping.h:423

References ASSERTL0, Nektar::GlobalMapping::Mapping::m_constantJacobian, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, and Nektar::GlobalMapping::Mapping::v_InitObject().

◆ v_LowerIndex()

void Nektar::GlobalMapping::MappingXofXZ::v_LowerIndex ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 242 of file MappingXofXZ.cpp.

245 {
246  int physTot = m_fields[0]->GetTotPoints();
247  Array<OneD, NekDouble> wk(physTot, 0.0);
248 
249  // out[0] = in[0]*fx^2 + in[2] * fz*fx
250  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[0], 1, wk,
251  1); // fz*fx
252  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, outarray[0], 1); // in[2]*fz*fx
253  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[2], 1); // in[0]*fz*fx
254 
255  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1, wk,
256  1); // fx^2
257  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, wk, 1); // in[0]*fx^2
258 
259  Vmath::Vadd(physTot, outarray[0], 1, wk, 1, outarray[0], 1);
260 
261  // out[1] = in[1]
262  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
263 
264  // out[2] = fx*fz*in[0] + (1+fz^2)*in[2]
265  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1, wk,
266  1); // fz^2
267  Vmath::Sadd(physTot, 1.0, wk, 1, wk, 1); // 1+fz^2
268  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, wk, 1); // (1+fz^2)*in[2]
269 
270  Vmath::Vadd(physTot, wk, 1, outarray[2], 1, outarray[2], 1);
271 }

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Sadd(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vmul().

◆ v_RaiseIndex()

void Nektar::GlobalMapping::MappingXofXZ::v_RaiseIndex ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 273 of file MappingXofXZ.cpp.

276 {
277  int physTot = m_fields[0]->GetTotPoints();
278  Array<OneD, NekDouble> wk(physTot, 0.0);
279  Array<OneD, NekDouble> wk_2(physTot, 0.0);
280 
281  // out[2] = in[2] - in[0] * fz/fx
282  Vmath::Vdiv(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[0], 1, wk, 1);
283  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[2], 1);
284  Vmath::Vsub(physTot, inarray[2], 1, outarray[2], 1, outarray[2], 1);
285 
286  // out[0] = in[0]*(1+fz^2)/(fx^2) - in[2] * fz/fx
287  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, outarray[0], 1);
288  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1, wk, 1);
289  Vmath::Sadd(physTot, 1.0, wk, 1, wk, 1);
290  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1, wk_2, 1);
291  Vmath::Vdiv(physTot, wk, 1, wk_2, 1, wk, 1);
292  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, wk, 1);
293  Vmath::Vsub(physTot, wk, 1, outarray[0], 1, outarray[0], 1);
294 
295  // out[1] = in[1]
296  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
297 }

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Sadd(), Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Vsub().

◆ v_UpdateGeomInfo()

void Nektar::GlobalMapping::MappingXofXZ::v_UpdateGeomInfo ( )
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 363 of file MappingXofXZ.cpp.

364 {
365  int phystot = m_fields[0]->GetTotPoints();
366  // Allocation of geometry memory
367  m_GeometricInfo = Array<OneD, Array<OneD, NekDouble>>(5);
368  for (int i = 0; i < m_GeometricInfo.size(); i++)
369  {
370  m_GeometricInfo[i] = Array<OneD, NekDouble>(phystot, 0.0);
371  }
372 
373  bool waveSpace = m_fields[0]->GetWaveSpace();
374  m_fields[0]->SetWaveSpace(false);
375 
376  // Calculate derivatives of transformation
377  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0], m_coords[0],
378  m_GeometricInfo[0]); // f_x
379  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2], m_coords[0],
380  m_GeometricInfo[1]); // f_z
381 
383  m_GeometricInfo[2]); // f_xx
385  m_GeometricInfo[3]); // f_xz
387  m_GeometricInfo[4]); // f_zz
388 
389  m_fields[0]->SetWaveSpace(waveSpace);
390 }
Array< OneD, Array< OneD, NekDouble > > m_coords
Array with the Cartesian coordinates.
Definition: Mapping.h:408
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:91

References Nektar::MultiRegions::DirCartesianMap, Nektar::GlobalMapping::Mapping::m_coords, Nektar::GlobalMapping::Mapping::m_fields, and Nektar::GlobalMapping::Mapping::m_GeometricInfo.

Friends And Related Function Documentation

◆ MemoryManager< MappingXofXZ >

friend class MemoryManager< MappingXofXZ >
friend

Definition at line 1 of file MappingXofXZ.h.

Member Data Documentation

◆ className

std::string Nektar::GlobalMapping::MappingXofXZ::className
static
Initial value:
=
"X = X(x,z)")
static GLOBAL_MAPPING_EXPORT MappingSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
Creates an instance of this class.
Definition: MappingXofXZ.h:58
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
MappingFactory & GetMappingFactory()
Declaration of the mapping factory singleton.
Definition: Mapping.cpp:53

Name of the class.

Definition at line 70 of file MappingXofXZ.h.