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)
 
GLOBAL_MAPPING_EXPORT void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping) override
 
GLOBAL_MAPPING_EXPORT void v_ContravarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_CovarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_GetJacobian (Array< OneD, NekDouble > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_DotGradJacobian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_GetMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_LowerIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_RaiseIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelContravar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
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_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)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=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_GetJacobian (Array< OneD, NekDouble > &outarray)=0
 
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)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
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)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
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_UpdateGeomInfo ()=0
 
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 49 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 55 of file MappingXofXZ.cpp.

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

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

60 {
63 p->InitObject(pFields, pMapping);
64 return p;
65 }
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:51

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

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

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

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

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

117{
118 int physTot = m_fields[0]->GetTotPoints();
119 Array<OneD, NekDouble> wk(physTot, 0.0);
120
121 // U1 = u1/fx - fz/fx * u3
122 Vmath::Vdiv(physTot, inarray[0], 1, m_GeometricInfo[0], 1, outarray[0], 1);
123 Vmath::Vdiv(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[0], 1, wk, 1);
124 Vmath::Vmul(physTot, wk, 1, inarray[2], 1, wk, 1);
125 Vmath::Vsub(physTot, outarray[0], 1, wk, 1, outarray[0], 1);
126
127 // U2 = u2
128 Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
129
130 // U3 = u3
131 Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
132}
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.hpp:220

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

Name of the class.

Definition at line 68 of file MappingXofXZ.h.