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)
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetJacobian (Array< OneD, NekDouble > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_DotGradJacobian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_LowerIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_RaiseIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelContravar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateGeomInfo ()
 
- 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:58

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 59 of file MappingXofXZ.h.

63  {
66  pFields);
67  p->InitObject(pFields, pMapping);
68  return p;
69  }
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 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 313 of file MappingXofXZ.cpp.

316 {
317  int physTot = m_fields[0]->GetTotPoints();
318  int nvel = m_nConvectiveFields;
319  Array<OneD, NekDouble> wk(physTot, 0.0);
320 
321  for (int i = 0; i< nvel; i++)
322  {
323  for (int j = 0; j< nvel; j++)
324  {
325  outarray[i*nvel+j] = Array<OneD, NekDouble>(physTot,0.0);
326  }
327  }
328 
329  // Calculate non-zero terms
330 
331  // outarray(0,0) = U1 * fxx/fx + U3 * fxz/fx
332  Vmath::Vdiv(physTot,m_GeometricInfo[2],1,m_GeometricInfo[0],1,wk,1);
333  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+0],1);
334  Vmath::Vdiv(physTot,m_GeometricInfo[3],1,m_GeometricInfo[0],1,wk,1);
335  Vmath::Vvtvp(physTot,wk,1,inarray[2],1,outarray[0*nvel+0],1,
336  outarray[0*nvel+0],1);
337 
338  // outarray(0,2) = U1 * fxz/fx + U3 * fzz/fx
339  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+2],1);
340  Vmath::Vdiv(physTot,m_GeometricInfo[4],1,m_GeometricInfo[0],1,wk,1);
341  Vmath::Vvtvp(physTot,wk,1,inarray[2],1,outarray[0*nvel+2],1,
342  outarray[0*nvel+2],1);
343 
344 }
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:416
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:414
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:408
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:192
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:513
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:257

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 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 346 of file MappingXofXZ.cpp.

349 {
350  int physTot = m_fields[0]->GetTotPoints();
351  int nvel = m_nConvectiveFields;
352  Array<OneD, NekDouble> wk(physTot, 0.0);
353 
354  for (int i = 0; i< nvel; i++)
355  {
356  for (int j = 0; j< nvel; j++)
357  {
358  outarray[i*nvel+j] = Array<OneD, NekDouble>(physTot,0.0);
359  }
360  }
361 
362  // Calculate non-zero terms
363 
364  // outarray(0,0) = U1 * fxx/fx
365  Vmath::Vdiv(physTot,m_GeometricInfo[2],1,m_GeometricInfo[0],1,wk,1);
366  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+0],1);
367 
368  //outarray(0,2) = outarray(2,0) = U1 * fxz/fx
369  Vmath::Vdiv(physTot,m_GeometricInfo[3],1,m_GeometricInfo[0],1,wk,1);
370  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+2],1);
371  Vmath::Vcopy(physTot,outarray[0*nvel+2],1,outarray[2*nvel+0],1);
372 
373  // outarray(2,2) = U1 * fzz/fx
374  Vmath::Vdiv(physTot,m_GeometricInfo[4],1,m_GeometricInfo[0],1,wk,1);
375  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[2*nvel+2],1);
376 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 120 of file MappingXofXZ.cpp.

123 {
124  int physTot = m_fields[0]->GetTotPoints();
125  Array<OneD, NekDouble> wk(physTot, 0.0);
126 
127  // U1 = u1/fx - fz/fx * u3
128  Vmath::Vdiv(physTot, inarray[0], 1,
129  m_GeometricInfo[0], 1, outarray[0], 1);
130  Vmath::Vdiv(physTot, m_GeometricInfo[1], 1,
131  m_GeometricInfo[0], 1, wk, 1);
132  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, wk, 1);
133  Vmath::Vsub(physTot, outarray[0], 1, wk, 1, outarray[0], 1);
134 
135  // U2 = u2
136  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
137 
138  // U3 = u3
139  Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
140 }
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:372

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 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 80 of file MappingXofXZ.cpp.

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

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 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 142 of file MappingXofXZ.cpp.

145 {
146  int physTot = m_fields[0]->GetTotPoints();
147 
148  // U1 = u1*fx
149  Vmath::Vmul(physTot, inarray[0], 1, m_GeometricInfo[0], 1,
150  outarray[0], 1);
151 
152  // U2 = u2
153  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
154 
155  // U3 = u3 + fz*u1
156  Vmath::Vmul(physTot, m_GeometricInfo[1], 1,
157  inarray[0], 1, outarray[2], 1);
158  Vmath::Vadd(physTot, inarray[2], 1, outarray[2], 1, outarray[2], 1);
159 }
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:322

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 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 99 of file MappingXofXZ.cpp.

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

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 
)
protectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 168 of file MappingXofXZ.cpp.

171 {
172  int physTot = m_fields[0]->GetTotPoints();
173 
174  Vmath::Vmul(physTot, m_GeometricInfo[2], 1, inarray[0], 1, outarray, 1);
175  Vmath::Vvtvp(physTot, m_GeometricInfo[3], 1, inarray[2], 1,
176  outarray, 1, outarray,1);
177 }

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)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 214 of file MappingXofXZ.cpp.

216 {
217  int physTot = m_fields[0]->GetTotPoints();
218  int nvel = m_nConvectiveFields;
219  Array<OneD, NekDouble> wk(physTot, 0.0);
220 
221  for (int i=0; i<nvel*nvel; i++)
222  {
223  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
224  }
225  // Fill diagonal with 1.0
226  for (int i=0; i<nvel; i++)
227  {
228  Vmath::Sadd(physTot, 1.0, outarray[i+nvel*i], 1,
229  outarray[i+nvel*i], 1);
230  }
231 
232  // G^{13} and G^{31} = -fz/fx
233  Vmath::Vdiv(physTot,m_GeometricInfo[1],1,
234  m_GeometricInfo[0],1,wk,1); // fz/fx
235  Vmath::Neg(physTot, wk, 1);
236  Vmath::Vcopy(physTot, wk, 1, outarray[0*nvel+2], 1);
237  Vmath::Vcopy(physTot, wk, 1, outarray[2*nvel+0], 1);
238 
239  // G^{11} = (1+fz^2)/(fx^2)
240  Vmath::Vmul(physTot, m_GeometricInfo[1], 1,
241  m_GeometricInfo[1], 1, wk, 1); // fz^2
242  Vmath::Vadd(physTot, wk, 1, outarray[0*nvel+0], 1,
243  outarray[0*nvel+0], 1);
244 
245  Vmath::Vmul(physTot, m_GeometricInfo[0], 1,
246  m_GeometricInfo[0], 1, wk, 1); // fx^2
247  Vmath::Vdiv(physTot, outarray[0*nvel+0], 1, wk,1,
248  outarray[0*nvel+0], 1);
249 }
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha - x.
Definition: Vmath.cpp:341

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)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 161 of file MappingXofXZ.cpp.

163 {
164  int physTot = m_fields[0]->GetTotPoints();
165  Vmath::Vcopy(physTot, m_GeometricInfo[0], 1, outarray, 1);
166 }

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)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 179 of file MappingXofXZ.cpp.

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

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 
)
protectedvirtual

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 68 of file MappingXofXZ.cpp.

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

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 
)
protectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 251 of file MappingXofXZ.cpp.

254 {
255  int physTot = m_fields[0]->GetTotPoints();
256  Array<OneD, NekDouble> wk(physTot, 0.0);
257 
258  // out[0] = in[0]*fx^2 + in[2] * fz*fx
259  Vmath::Vmul(physTot,m_GeometricInfo[1],1,m_GeometricInfo[0],1,
260  wk,1); // fz*fx
261  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, outarray[0], 1); //in[2]*fz*fx
262  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[2], 1); //in[0]*fz*fx
263 
264  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1,
265  wk, 1); //fx^2
266  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, wk, 1); //in[0]*fx^2
267 
268  Vmath::Vadd(physTot, outarray[0], 1, wk, 1, outarray[0], 1);
269 
270  // out[1] = in[1]
271  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
272 
273  // out[2] = fx*fz*in[0] + (1+fz^2)*in[2]
274  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1,
275  wk, 1); // fz^2
276  Vmath::Sadd(physTot, 1.0, wk, 1, wk, 1); // 1+fz^2
277  Vmath::Vmul(physTot, wk, 1, inarray[2],1, wk, 1); // (1+fz^2)*in[2]
278 
279  Vmath::Vadd(physTot, wk, 1, outarray[2],1, outarray[2], 1);
280 }

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 
)
protectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 282 of file MappingXofXZ.cpp.

285 {
286  int physTot = m_fields[0]->GetTotPoints();
287  Array<OneD, NekDouble> wk(physTot, 0.0);
288  Array<OneD, NekDouble> wk_2(physTot, 0.0);
289 
290  // out[2] = in[2] - in[0] * fz/fx
291  Vmath::Vdiv(physTot,m_GeometricInfo[1],1,m_GeometricInfo[0],1,
292  wk,1);
293  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[2], 1);
294  Vmath::Vsub(physTot, inarray[2], 1, outarray[2], 1, outarray[2], 1);
295 
296  // out[0] = in[0]*(1+fz^2)/(fx^2) - in[2] * fz/fx
297  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, outarray[0], 1);
298  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1,
299  wk, 1);
300  Vmath::Sadd(physTot, 1.0, wk, 1, wk, 1);
301  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1,
302  wk_2, 1);
303  Vmath::Vdiv(physTot, wk, 1, wk_2,1, wk, 1);
304  Vmath::Vmul(physTot, wk, 1, inarray[0],1, wk, 1);
305  Vmath::Vsub(physTot, wk, 1, outarray[0], 1, outarray[0], 1);
306 
307  // out[1] = in[1]
308  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
309 
310 
311 }

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 378 of file MappingXofXZ.cpp.

379 {
380  int phystot = m_fields[0]->GetTotPoints();
381  // Allocation of geometry memory
382  m_GeometricInfo = Array<OneD, Array<OneD, NekDouble> >(5);
383  for (int i = 0; i < m_GeometricInfo.size(); i++)
384  {
385  m_GeometricInfo[i] = Array<OneD, NekDouble>(phystot, 0.0);
386  }
387 
388  bool waveSpace = m_fields[0]->GetWaveSpace();
389  m_fields[0]->SetWaveSpace(false);
390 
391  // Calculate derivatives of transformation
392  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
393  m_coords[0], m_GeometricInfo[0]); //f_x
394  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
395  m_coords[0], m_GeometricInfo[1]); //f_z
396 
397  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
398  m_GeometricInfo[0], m_GeometricInfo[2]); //f_xx
399  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
400  m_GeometricInfo[0], m_GeometricInfo[3]); //f_xz
401  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
402  m_GeometricInfo[1], m_GeometricInfo[4]); //f_zz
403 
404  m_fields[0]->SetWaveSpace(waveSpace);
405 }
Array< OneD, Array< OneD, NekDouble > > m_coords
Array with the Cartesian coordinates.
Definition: Mapping.h:410
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90

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:
=
MappingXofXZ::create, "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:59
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
MappingFactory & GetMappingFactory()
Declaration of the mapping factory singleton.
Definition: Mapping.cpp:52

Name of the class.

Definition at line 72 of file MappingXofXZ.h.