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

#include <MappingXYofXY.h>

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

void CalculateMetricTensor ()
 
void CalculateChristoffel ()
 
 MappingXYofXY (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_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_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)
 

Protected Attributes

Array< OneD, Array< OneD, NekDouble > > m_metricTensor
 
Array< OneD, Array< OneD, NekDouble > > m_Christoffel
 
- 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
 

Friends

class MemoryManager< MappingXYofXY >
 

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...
 
- 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 the transformation

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

\[ \bar{y} = \bar{y}(x,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 MappingXYofXY.h.

Constructor & Destructor Documentation

◆ MappingXYofXY()

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

Definition at line 55 of file MappingXYofXY.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

◆ CalculateChristoffel()

void Nektar::GlobalMapping::MappingXYofXY::CalculateChristoffel ( )
protected

Definition at line 400 of file MappingXYofXY.cpp.

401{
402 int physTot = m_fields[0]->GetTotPoints();
403 int nvel = m_nConvectiveFields;
404
405 Array<OneD, Array<OneD, NekDouble>> G(nvel * nvel);
406 Array<OneD, Array<OneD, NekDouble>> G_inv(nvel * nvel);
407 Array<OneD, Array<OneD, NekDouble>> gradG(2 * 2 * 2);
408 Array<OneD, Array<OneD, NekDouble>> tmp(2 * 2 * 2);
409 m_Christoffel = Array<OneD, Array<OneD, NekDouble>>(6);
410 // Allocate memory
411 for (int i = 0; i < gradG.size(); i++)
412 {
413 gradG[i] = Array<OneD, NekDouble>(physTot, 0.0);
414 tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
415 }
416 for (int i = 0; i < G.size(); i++)
417 {
418 G[i] = Array<OneD, NekDouble>(physTot, 0.0);
419 G_inv[i] = Array<OneD, NekDouble>(physTot, 0.0);
420 }
421
422 // Get the metric tensor and its inverse
424 GetInvMetricTensor(G_inv);
425
426 bool waveSpace = m_fields[0]->GetWaveSpace();
427 m_fields[0]->SetWaveSpace(false);
428 // Calculate gradients of g
429 // consider only 2 dimensions, since the 3rd is trivial
430 for (int i = 0; i < 2; i++)
431 {
432 for (int j = 0; j < 2; j++)
433 {
434 for (int k = 0; k < 2; k++)
435 {
437 G[i * nvel + j],
438 gradG[i * 2 * 2 + j * 2 + k]);
439 }
440 }
441 }
442
443 // Calculate tmp[p,j,k] = 1/2( gradG[pj,k]+ gradG[pk,j]-gradG[jk,p])
444 for (int p = 0; p < 2; p++)
445 {
446 for (int j = 0; j < 2; j++)
447 {
448 for (int k = 0; k < 2; k++)
449 {
450 Vmath::Vadd(physTot, gradG[p * 2 * 2 + j * 2 + k], 1,
451 gradG[p * 2 * 2 + k * 2 + j], 1,
452 tmp[p * 2 * 2 + j * 2 + k], 1);
453 Vmath::Vsub(physTot, tmp[p * 2 * 2 + j * 2 + k], 1,
454 gradG[j * 2 * 2 + k * 2 + p], 1,
455 tmp[p * 2 * 2 + j * 2 + k], 1);
456 Vmath::Smul(physTot, 0.5, tmp[p * 2 * 2 + j * 2 + k], 1,
457 tmp[p * 2 * 2 + j * 2 + k], 1);
458 }
459 }
460 }
461
462 // Calculate Christoffel symbols = g^ip tmp[p,j,k]
463 int n = 0;
464 for (int i = 0; i < 2; i++)
465 {
466 for (int j = 0; j < 2; j++)
467 {
468 for (int k = 0; k <= j; k++)
469 {
470 m_Christoffel[n] = Array<OneD, NekDouble>(physTot, 0.0);
471 for (int p = 0; p < 2; p++)
472 {
473 Vmath::Vvtvp(physTot, G_inv[i * nvel + p], 1,
474 tmp[p * 2 * 2 + j * 2 + k], 1,
475 m_Christoffel[n], 1, m_Christoffel[n], 1);
476 }
477 n = n + 1;
478 }
479 }
480 }
481
482 m_fields[0]->SetWaveSpace(waveSpace);
483}
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:412
GLOBAL_MAPPING_EXPORT void GetInvMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)
Get the inverse of metric tensor .
Definition: Mapping.h:182
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:404
GLOBAL_MAPPING_EXPORT void GetMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)
Get the metric tensor .
Definition: Mapping.h:175
Array< OneD, Array< OneD, NekDouble > > m_Christoffel
Definition: MappingXYofXY.h:78
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
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::MultiRegions::DirCartesianMap, Nektar::GlobalMapping::Mapping::GetInvMetricTensor(), Nektar::GlobalMapping::Mapping::GetMetricTensor(), m_Christoffel, Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, CellMLToNektar.cellml_metadata::p, Vmath::Smul(), Vmath::Vadd(), Vmath::Vsub(), and Vmath::Vvtvp().

Referenced by v_UpdateGeomInfo().

◆ CalculateMetricTensor()

void Nektar::GlobalMapping::MappingXYofXY::CalculateMetricTensor ( )
protected

Definition at line 374 of file MappingXYofXY.cpp.

375{
376 int physTot = m_fields[0]->GetTotPoints();
377 // Allocate memory
378 m_metricTensor = Array<OneD, Array<OneD, NekDouble>>(3);
379 for (int i = 0; i < m_metricTensor.size(); i++)
380 {
381 m_metricTensor[i] = Array<OneD, NekDouble>(physTot, 0.0);
382 }
383 // g_{1,1} = fx^2+gx^2
384 Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1,
385 m_metricTensor[0], 1);
386 Vmath::Vvtvp(physTot, m_GeometricInfo[2], 1, m_GeometricInfo[2], 1,
387 m_metricTensor[0], 1, m_metricTensor[0], 1);
388 // g_{2,2} = fy^2+gy^2
389 Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1,
390 m_metricTensor[1], 1);
391 Vmath::Vvtvp(physTot, m_GeometricInfo[3], 1, m_GeometricInfo[3], 1,
392 m_metricTensor[1], 1, m_metricTensor[1], 1);
393 // g_{1,2} = g_{2,1} = fy*fx+gx*gy
394 Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[1], 1,
395 m_metricTensor[2], 1);
396 Vmath::Vvtvp(physTot, m_GeometricInfo[2], 1, m_GeometricInfo[3], 1,
397 m_metricTensor[2], 1, m_metricTensor[2], 1);
398}
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:410
Array< OneD, Array< OneD, NekDouble > > m_metricTensor
Definition: MappingXYofXY.h:77
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

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

Referenced by v_UpdateGeomInfo().

◆ create()

static GLOBAL_MAPPING_EXPORT MappingSharedPtr Nektar::GlobalMapping::MappingXYofXY::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 MappingXYofXY.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::MappingXYofXY::v_ApplyChristoffelContravar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
overrideprotectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 253 of file MappingXYofXY.cpp.

256{
257 int physTot = m_fields[0]->GetTotPoints();
258 int nvel = m_nConvectiveFields;
259
260 for (int i = 0; i < nvel; i++)
261 {
262 for (int j = 0; j < nvel; j++)
263 {
264 outarray[i * nvel + j] = Array<OneD, NekDouble>(physTot, 0.0);
265 }
266 }
267
268 // Calculate non-zero terms
269
270 // outarray(0,0) = U1 * m_Christoffel[0] + U2 * m_Christoffel[1]
271 Vmath::Vmul(physTot, m_Christoffel[0], 1, inarray[0], 1,
272 outarray[0 * nvel + 0], 1);
273 Vmath::Vvtvp(physTot, m_Christoffel[1], 1, inarray[1], 1,
274 outarray[0 * nvel + 0], 1, outarray[0 * nvel + 0], 1);
275
276 // outarray(0,1) = U1 * m_Christoffel[1] + U2 * m_Christoffel[2]
277 Vmath::Vmul(physTot, m_Christoffel[1], 1, inarray[0], 1,
278 outarray[0 * nvel + 1], 1);
279 Vmath::Vvtvp(physTot, m_Christoffel[2], 1, inarray[1], 1,
280 outarray[0 * nvel + 1], 1, outarray[0 * nvel + 1], 1);
281
282 // outarray(1,0) = U1 * m_Christoffel[3] + U2 * m_Christoffel[4]
283 Vmath::Vmul(physTot, m_Christoffel[3], 1, inarray[0], 1,
284 outarray[1 * nvel + 0], 1);
285 Vmath::Vvtvp(physTot, m_Christoffel[4], 1, inarray[1], 1,
286 outarray[1 * nvel + 0], 1, outarray[1 * nvel + 0], 1);
287
288 // outarray(1,1) = U1 * m_Christoffel[4] + U2 * m_Christoffel[5]
289 Vmath::Vmul(physTot, m_Christoffel[4], 1, inarray[0], 1,
290 outarray[1 * nvel + 1], 1);
291 Vmath::Vvtvp(physTot, m_Christoffel[5], 1, inarray[1], 1,
292 outarray[1 * nvel + 1], 1, outarray[1 * nvel + 1], 1);
293}

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

◆ v_ApplyChristoffelCovar()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 295 of file MappingXYofXY.cpp.

298{
299 int physTot = m_fields[0]->GetTotPoints();
300 int nvel = m_nConvectiveFields;
301
302 for (int i = 0; i < nvel; i++)
303 {
304 for (int j = 0; j < nvel; j++)
305 {
306 outarray[i * nvel + j] = Array<OneD, NekDouble>(physTot, 0.0);
307 }
308 }
309
310 // Calculate non-zero terms
311
312 // outarray(0,0) = U1 * m_Christoffel[0] + U2 * m_Christoffel[3]
313 Vmath::Vmul(physTot, m_Christoffel[0], 1, inarray[0], 1,
314 outarray[0 * nvel + 0], 1);
315 Vmath::Vvtvp(physTot, m_Christoffel[3], 1, inarray[1], 1,
316 outarray[0 * nvel + 0], 1, outarray[0 * nvel + 0], 1);
317
318 // outarray(0,1) = U1 * m_Christoffel[1] + U2 * m_Christoffel[4]
319 Vmath::Vmul(physTot, m_Christoffel[1], 1, inarray[0], 1,
320 outarray[0 * nvel + 1], 1);
321 Vmath::Vvtvp(physTot, m_Christoffel[4], 1, inarray[1], 1,
322 outarray[0 * nvel + 1], 1, outarray[0 * nvel + 1], 1);
323
324 // outarray(1,0) = U1 * m_Christoffel[1] + U2 * m_Christoffel[4]
325 Vmath::Vmul(physTot, m_Christoffel[1], 1, inarray[0], 1,
326 outarray[1 * nvel + 0], 1);
327 Vmath::Vvtvp(physTot, m_Christoffel[4], 1, inarray[1], 1,
328 outarray[1 * nvel + 0], 1, outarray[1 * nvel + 0], 1);
329
330 // outarray(1,1) = U1 * m_Christoffel[2] + U2 * m_Christoffel[5]
331 Vmath::Vmul(physTot, m_Christoffel[2], 1, inarray[0], 1,
332 outarray[1 * nvel + 1], 1);
333 Vmath::Vvtvp(physTot, m_Christoffel[5], 1, inarray[1], 1,
334 outarray[1 * nvel + 1], 1, outarray[1 * nvel + 1], 1);
335}

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

◆ v_ContravarFromCartesian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 126 of file MappingXYofXY.cpp.

129{
130 int physTot = m_fields[0]->GetTotPoints();
131 Array<OneD, NekDouble> wk(physTot, 0.0);
132
133 // U1 = [gy*u1-fy*u2]/(fx*gy-gx*fy)
134 Vmath::Vmul(physTot, inarray[1], 1, m_GeometricInfo[1], 1, outarray[0], 1);
135 Vmath::Vvtvm(physTot, inarray[0], 1, m_GeometricInfo[3], 1, outarray[0], 1,
136 outarray[0], 1);
137 Vmath::Vdiv(physTot, outarray[0], 1, m_GeometricInfo[4], 1, outarray[0], 1);
138
139 // U2 = [fx*u2-gx*u1]/(fx*gy-gx*fy)
140 Vmath::Vmul(physTot, inarray[0], 1, m_GeometricInfo[2], 1, outarray[1], 1);
141 Vmath::Vvtvm(physTot, inarray[1], 1, m_GeometricInfo[0], 1, outarray[1], 1,
142 outarray[1], 1);
143 Vmath::Vdiv(physTot, outarray[1], 1, m_GeometricInfo[4], 1, outarray[1], 1);
144
145 // U3 = u3
146 if (m_nConvectiveFields == 3)
147 {
148 Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
149 }
150}
void Vvtvm(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)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.hpp:381
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
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(), Vmath::Vmul(), and Vmath::Vvtvm().

◆ v_ContravarToCartesian()

void Nektar::GlobalMapping::MappingXYofXY::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 MappingXYofXY.cpp.

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

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

◆ v_CovarFromCartesian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 152 of file MappingXYofXY.cpp.

155{
156 int physTot = m_fields[0]->GetTotPoints();
157
158 // U1 = u1*fx +gx*u2
159 Vmath::Vmul(physTot, m_GeometricInfo[0], 1, inarray[0], 1, outarray[0], 1);
160 Vmath::Vvtvp(physTot, m_GeometricInfo[2], 1, inarray[1], 1, outarray[0], 1,
161 outarray[0], 1);
162
163 // U2 = fy*u1 + gy*u2
164 Vmath::Vmul(physTot, m_GeometricInfo[1], 1, inarray[0], 1, outarray[1], 1);
165 Vmath::Vvtvp(physTot, m_GeometricInfo[3], 1, inarray[1], 1, outarray[1], 1,
166 outarray[1], 1);
167
168 // U3 = u3
169 if (m_nConvectiveFields == 3)
170 {
171 Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
172 }
173}

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

◆ v_CovarToCartesian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 100 of file MappingXYofXY.cpp.

103{
104 int physTot = m_fields[0]->GetTotPoints();
105 Array<OneD, NekDouble> wk(physTot, 0.0);
106
107 // U1 = [gy*u1-gx*u2]/(fx*gy-gx*fy)
108 Vmath::Vmul(physTot, inarray[1], 1, m_GeometricInfo[2], 1, outarray[0], 1);
109 Vmath::Vvtvm(physTot, inarray[0], 1, m_GeometricInfo[3], 1, outarray[0], 1,
110 outarray[0], 1);
111 Vmath::Vdiv(physTot, outarray[0], 1, m_GeometricInfo[4], 1, outarray[0], 1);
112
113 // U2 = [fx*u2 - fy*u1]/(fx*gy-gx*fy)
114 Vmath::Vmul(physTot, inarray[0], 1, m_GeometricInfo[1], 1, outarray[1], 1);
115 Vmath::Vvtvm(physTot, inarray[1], 1, m_GeometricInfo[0], 1, outarray[1], 1,
116 outarray[1], 1);
117 Vmath::Vdiv(physTot, outarray[1], 1, m_GeometricInfo[4], 1, outarray[1], 1);
118
119 // U3 = u3
120 if (m_nConvectiveFields == 3)
121 {
122 Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
123 }
124}

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

◆ v_GetInvMetricTensor()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 210 of file MappingXYofXY.cpp.

212{
213 int physTot = m_fields[0]->GetTotPoints();
214 int nvel = m_nConvectiveFields;
215
216 for (int i = 0; i < nvel * nvel; i++)
217 {
218 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
219 }
220
221 // Get Jacobian
222 Array<OneD, NekDouble> Jac(physTot, 0.0);
223 GetJacobian(Jac);
224
225 // Get Jacobian squared
226 Array<OneD, NekDouble> wk(physTot, 0.0);
227 Vmath::Vmul(physTot, Jac, 1, Jac, 1, wk, 1);
228 // G^{1,1} = m_metricTensor[1]/Jac^2
229 Vmath::Vcopy(physTot, m_metricTensor[1], 1, outarray[0 * nvel + 0], 1);
230 Vmath::Vdiv(physTot, outarray[0 * nvel + 0], 1, wk, 1,
231 outarray[0 * nvel + 0], 1);
232
233 // G^{2,2} = m_metricTensor[0]/Jac^2
234 Vmath::Vcopy(physTot, m_metricTensor[0], 1, outarray[1 * nvel + 1], 1);
235 Vmath::Vdiv(physTot, outarray[1 * nvel + 1], 1, wk, 1,
236 outarray[1 * nvel + 1], 1);
237
238 // G^{1,2} = G^{2,1} = -m_metricTensor[2]/Jac^2
239 Vmath::Vcopy(physTot, m_metricTensor[2], 1, outarray[0 * nvel + 1], 1);
240 Vmath::Neg(physTot, outarray[0 * nvel + 1], 1);
241 Vmath::Vdiv(physTot, outarray[0 * nvel + 1], 1, wk, 1,
242 outarray[0 * nvel + 1], 1);
243 Vmath::Vcopy(physTot, outarray[0 * nvel + 1], 1, outarray[1 * nvel + 0], 1);
244
245 // G^{3,3} = 1
246 if (m_nConvectiveFields == 3)
247 {
248 Vmath::Sadd(physTot, 1.0, outarray[2 * nvel + 2], 1,
249 outarray[2 * nvel + 2], 1);
250 }
251}
GLOBAL_MAPPING_EXPORT void GetJacobian(Array< OneD, NekDouble > &outarray)
Get the Jacobian of the transformation.
Definition: Mapping.h:153
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::GetJacobian(), Nektar::GlobalMapping::Mapping::m_fields, m_metricTensor, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, Vmath::Neg(), Vmath::Sadd(), Vmath::Vcopy(), Vmath::Vdiv(), and Vmath::Vmul().

◆ v_GetJacobian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 175 of file MappingXYofXY.cpp.

176{
177 int physTot = m_fields[0]->GetTotPoints();
178 Vmath::Vabs(physTot, m_GeometricInfo[4], 1, outarray, 1);
179}
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.hpp:352

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

◆ v_GetMetricTensor()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 181 of file MappingXYofXY.cpp.

183{
184 int physTot = m_fields[0]->GetTotPoints();
185 int nvel = m_nConvectiveFields;
186
187 for (int i = 0; i < nvel * nvel; i++)
188 {
189 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
190 }
191
192 // g_{1,1} = m_metricTensor[0]
193 Vmath::Vcopy(physTot, m_metricTensor[0], 1, outarray[0 * nvel + 0], 1);
194
195 // g_{2,2} = m_metricTensor[1]
196 Vmath::Vcopy(physTot, m_metricTensor[1], 1, outarray[1 * nvel + 1], 1);
197
198 // g_{1,2}=g{2,1} = m_metricTensor[2]
199 Vmath::Vcopy(physTot, m_metricTensor[2], 1, outarray[0 * nvel + 1], 1);
200 Vmath::Vcopy(physTot, m_metricTensor[2], 1, outarray[1 * nvel + 0], 1);
201
202 // g_{3,3} = 1
203 if (m_nConvectiveFields == 3)
204 {
205 Vmath::Sadd(physTot, 1.0, outarray[2 * nvel + 2], 1,
206 outarray[2 * nvel + 2], 1);
207 }
208}

References Nektar::GlobalMapping::Mapping::m_fields, m_metricTensor, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, Vmath::Sadd(), and Vmath::Vcopy().

◆ v_InitObject()

void Nektar::GlobalMapping::MappingXYofXY::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 MappingXYofXY.cpp.

68{
69 Mapping::v_InitObject(pFields, pMapping);
70
71 m_constantJacobian = false;
72
74 "Mapping X = X(x,y), Y = Y(x,y) needs 2 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_UpdateGeomInfo()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 337 of file MappingXYofXY.cpp.

338{
339 int phystot = m_fields[0]->GetTotPoints();
340 // Allocation of geometry memory
341 m_GeometricInfo = Array<OneD, Array<OneD, NekDouble>>(5);
342 for (int i = 0; i < m_GeometricInfo.size(); i++)
343 {
344 m_GeometricInfo[i] = Array<OneD, NekDouble>(phystot, 0.0);
345 }
346
347 bool waveSpace = m_fields[0]->GetWaveSpace();
348 m_fields[0]->SetWaveSpace(false);
349
350 // Calculate derivatives of x transformation --> m_GeometricInfo 0-1
352 m_GeometricInfo[0]);
354 m_GeometricInfo[1]);
355
356 // Calculate derivatives of y transformation m_GeometricInfo 2-3
358 m_GeometricInfo[2]);
360 m_GeometricInfo[3]);
361
362 // Calculate fx*gy-gx*fy --> m_GeometricInfo4
363 Vmath::Vmul(phystot, m_GeometricInfo[1], 1, m_GeometricInfo[2], 1,
364 m_GeometricInfo[4], 1);
365 Vmath::Vvtvm(phystot, m_GeometricInfo[0], 1, m_GeometricInfo[3], 1,
366 m_GeometricInfo[4], 1, m_GeometricInfo[4], 1);
367 //
370
371 m_fields[0]->SetWaveSpace(waveSpace);
372}
Array< OneD, Array< OneD, NekDouble > > m_coords
Array with the Cartesian coordinates.
Definition: Mapping.h:406

References CalculateChristoffel(), CalculateMetricTensor(), Nektar::MultiRegions::DirCartesianMap, Nektar::GlobalMapping::Mapping::m_coords, Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Vmul(), and Vmath::Vvtvm().

Friends And Related Function Documentation

◆ MemoryManager< MappingXYofXY >

friend class MemoryManager< MappingXYofXY >
friend

Definition at line 1 of file MappingXYofXY.h.

Member Data Documentation

◆ className

std::string Nektar::GlobalMapping::MappingXYofXY::className
static
Initial value:
=
"X = X(x,y), Y = Y(x,y)")
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: MappingXYofXY.h:56
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
MappingFactory & GetMappingFactory()
Declaration of the mapping factory singleton.
Definition: Mapping.cpp:49

Name of the class.

Definition at line 68 of file MappingXYofXY.h.

◆ m_Christoffel

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingXYofXY::m_Christoffel
protected

◆ m_metricTensor

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingXYofXY::m_metricTensor
protected

Definition at line 77 of file MappingXYofXY.h.

Referenced by CalculateMetricTensor(), v_GetInvMetricTensor(), and v_GetMetricTensor().