Nektar++
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
Nektar::GlobalMapping::Mapping Class Referenceabstract

Base class for mapping to be applied to the coordinate system. More...

#include <Mapping.h>

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

Public Member Functions

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 Public Member Functions

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

Protected Member Functions

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

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

static MappingSharedPtr m_mappingPtr = MappingSharedPtr()
 
static bool m_init = false
 
static bool m_isDefined = false
 

Detailed Description

Base class for mapping to be applied to the coordinate system.

Definition at line 68 of file Mapping.h.

Constructor & Destructor Documentation

◆ ~Mapping()

virtual GLOBAL_MAPPING_EXPORT Nektar::GlobalMapping::Mapping::~Mapping ( )
inlinevirtual

Destructor.

Definition at line 72 of file Mapping.h.

73 {
74 }

◆ Mapping()

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

Constructor.

Definition at line 59 of file Mapping.cpp.

61 : m_session(pSession), m_fields(pFields)
62{
63 switch (m_fields[0]->GetExpType())
64 {
66 {
68 }
69 break;
70
72 {
74 }
75 break;
76
80 {
82 }
83 break;
84
85 default:
86 ASSERTL0(0, "Dimension not supported");
87 break;
88 }
89
91}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:414
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Mapping.h:402
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:406
LibUtilities::FieldIOSharedPtr m_fld
Definition: Mapping.h:404
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:197

References ASSERTL0, Nektar::LibUtilities::FieldIO::CreateDefault(), Nektar::MultiRegions::e1D, Nektar::MultiRegions::e2D, Nektar::MultiRegions::e3D, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_fields, m_fld, and m_nConvectiveFields.

Member Function Documentation

◆ ApplyChristoffelContravar()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::ApplyChristoffelContravar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
inline

Apply the Christoffel symbols to a contravariant vector.

This function is used apply the Christoffel symbols \( \left( i,pk\right)\) to a contravariant vector \(u^p\), following the relation

\[ (out)^{i}_{k} = \left( i,pk\right)u^p \]

Parameters
inarrayContravariant vector \(u^p\)
outarrayResult of applying Christoffel symbols to \(u^p\)

Definition at line 211 of file Mapping.h.

214 {
215 v_ApplyChristoffelContravar(inarray, outarray);
216 }
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelContravar(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0

References v_ApplyChristoffelContravar().

Referenced by v_gradgradU(), and v_VelocityLaplacian().

◆ ApplyChristoffelCovar()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::ApplyChristoffelCovar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
inline

Apply the Christoffel symbols to a covariant vector.

This function is used apply the Christoffel symbols \( \left( p,ik\right)\) to a covariant vector \(u_p\), following the relation

\[ (out)_{ik} = \left( p,ik\right)u_p \]

Parameters
inarrayContravariant vector \(u_p\)
outarrayResult of applying Christoffel symbols to \(u_p\)

Definition at line 229 of file Mapping.h.

232 {
233 v_ApplyChristoffelCovar(inarray, outarray);
234 }
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0

References v_ApplyChristoffelCovar().

Referenced by v_gradgradU().

◆ ContravarFromCartesian()

void Nektar::GlobalMapping::Mapping::ContravarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Convert a contravariant vector to the transformed system.

This function converts a contravariant vector in Cartesian space \(\bar{v}^{i}\) to the corresponding \(v^{j}\) in transformed space using the relation

\[v^{j} = \frac{\partial x^j}{\partial \bar{x}^i}\bar{v}^{i}\]

Parameters
inarrayComponents of the vector in Cartesian space ( \(\bar{v}^{i}\))
outarrayComponents of the vector in transformed space ( \(v^{j}\))

Definition at line 547 of file Mapping.cpp.

550{
551 if (inarray == outarray)
552 {
553 int physTot = m_fields[0]->GetTotPoints();
554 int nvel = m_nConvectiveFields;
555
556 for (int i = 0; i < nvel; i++)
557 {
558 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
559 }
561 }
562 else
563 {
564 v_ContravarFromCartesian(inarray, outarray);
565 }
566}
virtual GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
Array< OneD, Array< OneD, NekDouble > > m_tmp
Definition: Mapping.h:437
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

References m_fields, m_nConvectiveFields, m_tmp, v_ContravarFromCartesian(), and Vmath::Vcopy().

Referenced by v_UpdateBCs().

◆ ContravarToCartesian()

void Nektar::GlobalMapping::Mapping::ContravarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Convert a contravariant vector to the Cartesian system.

This function converts a contravariant vector in transformed space \(v^{j}\) to the corresponding \(\bar{v}^{i}\) in Cartesian (physical) space using the relation

\[\bar{v}^{i} = \frac{\partial \bar{x}^i}{\partial x^j}v^{j}\]

Parameters
inarrayComponents of the vector in transformed space ( \(v^{j}\))
outarrayComponents of the vector in Cartesian space ( \(\bar{v}^{i}\))

Definition at line 485 of file Mapping.cpp.

488{
489 if (inarray == outarray)
490 {
491 int physTot = m_fields[0]->GetTotPoints();
492 int nvel = m_nConvectiveFields;
493
494 for (int i = 0; i < nvel; i++)
495 {
496 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
497 }
498 v_ContravarToCartesian(m_tmp, outarray);
499 }
500 else
501 {
502 v_ContravarToCartesian(inarray, outarray);
503 }
504}
virtual GLOBAL_MAPPING_EXPORT void v_ContravarToCartesian(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0

References m_fields, m_nConvectiveFields, m_tmp, v_ContravarToCartesian(), and Vmath::Vcopy().

◆ CovarFromCartesian()

void Nektar::GlobalMapping::Mapping::CovarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Convert a covariant vector to the transformed system.

This function converts a covariant vector in Cartesian space \(\bar{v}_{i}\) to the corresponding \(v_{j}\) in transformed space using the relation

\[\bar{v}_{j} = \frac{\partial \bar{x}^i}{\partial x^j}v_{i}\]

Parameters
inarrayComponents of the vector in Cartesian space ( \(\bar{v}_{i}\))
outarrayComponents of the vector in transformed space ( \(v_{j}\))

Definition at line 578 of file Mapping.cpp.

581{
582 if (inarray == outarray)
583 {
584 int physTot = m_fields[0]->GetTotPoints();
585 int nvel = m_nConvectiveFields;
586
587 for (int i = 0; i < nvel; i++)
588 {
589 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
590 }
591 v_CovarFromCartesian(m_tmp, outarray);
592 }
593 else
594 {
595 v_CovarFromCartesian(inarray, outarray);
596 }
597}
virtual GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0

References m_fields, m_nConvectiveFields, m_tmp, v_CovarFromCartesian(), and Vmath::Vcopy().

◆ CovarToCartesian()

void Nektar::GlobalMapping::Mapping::CovarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Convert a covariant vector to the Cartesian system.

This function converts a covariant vector in transformed space \(v_{j}\) to the corresponding \(\bar{v}_{i}\) in Cartesian (physical) space using the relation

\[\bar{v}_{i} = \frac{\partial x^j}{\partial \bar{x}^i}v_{j}\]

Parameters
inarrayComponents of the vector in transformed space ( \(v_{j}\))
outarrayComponents of the vector in Cartesian space ( \(\bar{v}_{i}\))

Definition at line 516 of file Mapping.cpp.

519{
520 if (inarray == outarray)
521 {
522 int physTot = m_fields[0]->GetTotPoints();
523 int nvel = m_nConvectiveFields;
524
525 for (int i = 0; i < nvel; i++)
526 {
527 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
528 }
529 v_CovarToCartesian(m_tmp, outarray);
530 }
531 else
532 {
533 v_CovarToCartesian(inarray, outarray);
534 }
535}
virtual GLOBAL_MAPPING_EXPORT void v_CovarToCartesian(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0

References m_fields, m_nConvectiveFields, m_tmp, v_CovarToCartesian(), and Vmath::Vcopy().

◆ CurlCurlField()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::CurlCurlField ( Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const bool  generalized 
)
inline

CurlCurl calculated on the whole field.

This function can be used to compute both the generalised CurlCurl or the typical (Cartesian) CurlCurl, depending on the flag generalized

Parameters
inarrayContravariant vector \(u^i\)
outarrayCurlCurl of \(u\)
generalizedFlag defining if generalised or typical CurlCurl

Definition at line 324 of file Mapping.h.

327 {
328 v_CurlCurlField(inarray, outarray, generalized);
329 }
virtual GLOBAL_MAPPING_EXPORT void v_CurlCurlField(Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const bool generalized)
Definition: Mapping.cpp:993

References v_CurlCurlField().

◆ Divergence()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::Divergence ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inline

Calculate the generalised divergence operator.

This function is used to calculate the generalised divergence of a contravariant vector, defined as

\[ D = u^i_{,i} = \frac{1}{J}\frac{\partial(Ju^i)}{\partial x^i}\]

Parameters
inarrayContravariant vector \(u^i\)
outarrayDivergence of \(u^i\)

Definition at line 266 of file Mapping.h.

269 {
270 v_Divergence(inarray, outarray);
271 }
virtual GLOBAL_MAPPING_EXPORT void v_Divergence(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Mapping.cpp:756

References v_Divergence().

◆ DotGradJacobian()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::DotGradJacobian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inline

Calculate the dot product with the gradient of the Jacobian.

This function calculates the dot product of an array against the gradient of the Jacobian of the Mapping.

Parameters
inarrayInput array
outarrayOutput array

Definition at line 169 of file Mapping.h.

172 {
173 v_DotGradJacobian(inarray, outarray);
174 }
virtual GLOBAL_MAPPING_EXPORT void v_DotGradJacobian(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Mapping.cpp:685

References v_DotGradJacobian().

◆ EvaluateFunction()

void Nektar::GlobalMapping::Mapping::EvaluateFunction ( Array< OneD, MultiRegions::ExpListSharedPtr pFields,
LibUtilities::SessionReaderSharedPtr  pSession,
std::string  pFieldName,
Array< OneD, NekDouble > &  pArray,
const std::string &  pFunctionName,
NekDouble  pTime = NekDouble(0) 
)
protected

Definition at line 406 of file Mapping.cpp.

411{
412 ASSERTL0(pSession->DefinesFunction(pFunctionName),
413 "Function '" + pFunctionName + "' does not exist.");
414
415 unsigned int nq = pFields[0]->GetNpoints();
416 if (pArray.size() != nq)
417 {
418 pArray = Array<OneD, NekDouble>(nq);
419 }
420
422 vType = pSession->GetFunctionType(pFunctionName, pFieldName);
424 {
425 Array<OneD, NekDouble> x0(nq);
426 Array<OneD, NekDouble> x1(nq);
427 Array<OneD, NekDouble> x2(nq);
428
429 pFields[0]->GetCoords(x0, x1, x2);
431 pSession->GetFunction(pFunctionName, pFieldName);
432
433 ffunc->Evaluate(x0, x1, x2, pTime, pArray);
434 }
435 else if (vType == LibUtilities::eFunctionTypeFile)
436 {
437 std::string filename =
438 pSession->GetFunctionFilename(pFunctionName, pFieldName);
439
440 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
441 std::vector<std::vector<NekDouble>> FieldData;
442 Array<OneD, NekDouble> vCoeffs(pFields[0]->GetNcoeffs());
443 Vmath::Zero(vCoeffs.size(), vCoeffs, 1);
444
446 LibUtilities::FieldIO::CreateForFile(pSession, filename);
447 fld->Import(filename, FieldDef, FieldData);
448
449 int idx = -1;
450 for (int i = 0; i < FieldDef.size(); ++i)
451 {
452 for (int j = 0; j < FieldDef[i]->m_fields.size(); ++j)
453 {
454 if (FieldDef[i]->m_fields[j] == pFieldName)
455 {
456 idx = j;
457 }
458 }
459
460 if (idx >= 0)
461 {
462 pFields[0]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
463 FieldDef[i]->m_fields[idx],
464 vCoeffs);
465 }
466 else
467 {
468 cout << "Field " + pFieldName + " not found." << endl;
469 }
470 }
471 pFields[0]->BwdTrans(vCoeffs, pArray);
472 }
473}
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:226
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:327
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

References ASSERTL0, Nektar::LibUtilities::FieldIO::CreateForFile(), Nektar::LibUtilities::eFunctionTypeExpression, Nektar::LibUtilities::eFunctionTypeFile, m_fields, and Vmath::Zero().

Referenced by v_InitObject(), and v_UpdateMapping().

◆ EvaluateTimeFunction()

void Nektar::GlobalMapping::Mapping::EvaluateTimeFunction ( LibUtilities::SessionReaderSharedPtr  pSession,
std::string  pFieldName,
Array< OneD, NekDouble > &  pArray,
const std::string &  pFunctionName,
NekDouble  pTime = NekDouble(0) 
)
protected

Definition at line 388 of file Mapping.cpp.

392{
393 ASSERTL0(pSession->DefinesFunction(pFunctionName),
394 "Function '" + pFunctionName + "' does not exist.");
395
397 pSession->GetFunction(pFunctionName, pFieldName);
398
399 Array<OneD, NekDouble> x0(1, 0.0);
400 Array<OneD, NekDouble> x1(1, 0.0);
401 Array<OneD, NekDouble> x2(1, 0.0);
402
403 ffunc->Evaluate(x0, x1, x2, pTime, pArray);
404}

References ASSERTL0.

◆ GetCartesianCoordinates()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::GetCartesianCoordinates ( Array< OneD, NekDouble > &  out0,
Array< OneD, NekDouble > &  out1,
Array< OneD, NekDouble > &  out2 
)
inline

Get the Cartesian coordinates in the field.

This function is used to obtain the Cartesian coordinates associated withthe Mapping

Parameters
out0Coordinates in the x-direction
out1Coordinates in the y-direction
out2Coordinates in the z-direction

Definition at line 135 of file Mapping.h.

138 {
139 v_GetCartesianCoordinates(out0, out1, out2);
140 }
virtual GLOBAL_MAPPING_EXPORT void v_GetCartesianCoordinates(Array< OneD, NekDouble > &out0, Array< OneD, NekDouble > &out1, Array< OneD, NekDouble > &out2)
Definition: Mapping.cpp:659

References v_GetCartesianCoordinates().

Referenced by v_UpdateBCs().

◆ GetCoordVelocity()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::GetCoordVelocity ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
inline

Obtain the velocity of the coordinates.

This function is used to obtain the velocity of the coordinates associated with the Mapping

Parameters
outarrayVelocity of the coordinates

Definition at line 244 of file Mapping.h.

246 {
247 v_GetCoordVelocity(outarray);
248 }
virtual GLOBAL_MAPPING_EXPORT void v_GetCoordVelocity(Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Mapping.cpp:674

References v_GetCoordVelocity().

Referenced by v_UpdateBCs().

◆ GetInvMetricTensor()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::GetInvMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
inline

Get the inverse of metric tensor \(g^{ij}\).

Definition at line 184 of file Mapping.h.

186 {
187 v_GetInvMetricTensor(outarray);
188 }
virtual GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)=0

References v_GetInvMetricTensor().

Referenced by Nektar::GlobalMapping::MappingXYofXY::CalculateChristoffel(), and v_RaiseIndex().

◆ GetJacobian()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::GetJacobian ( Array< OneD, NekDouble > &  outarray)
inline

Get the Jacobian of the transformation.

This function is used to obtain the Jacobian of the Mapping.

Parameters
outarrayArray containing the Jacobian

Definition at line 155 of file Mapping.h.

156 {
157 v_GetJacobian(outarray);
158 }
virtual GLOBAL_MAPPING_EXPORT void v_GetJacobian(Array< OneD, NekDouble > &outarray)=0

References v_GetJacobian().

Referenced by v_Divergence(), v_DotGradJacobian(), and Nektar::GlobalMapping::MappingXYofXY::v_GetInvMetricTensor().

◆ GetMetricTensor()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::GetMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
inline

Get the metric tensor \(g_{ij}\).

Definition at line 177 of file Mapping.h.

179 {
180 v_GetMetricTensor(outarray);
181 }
virtual GLOBAL_MAPPING_EXPORT void v_GetMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)=0

References v_GetMetricTensor().

Referenced by Nektar::GlobalMapping::MappingXYofXY::CalculateChristoffel(), and v_LowerIndex().

◆ gradgradU()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::gradgradU ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
inline

Second order covariant derivatives of a contravariant vector.

This function computes the second order covariant derivatives of a contravariant vector, resulting in \(u^{i}_{,jk}\)

Parameters
inarrayContravariant vector \(u^i\)
outarraySecond order derivatives \(u^{i}_{,jk}\)

Definition at line 306 of file Mapping.h.

309 {
310 v_gradgradU(inarray, outarray);
311 }
virtual GLOBAL_MAPPING_EXPORT void v_gradgradU(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Mapping.cpp:885

References v_gradgradU().

Referenced by v_CurlCurlField().

◆ HasConstantJacobian()

GLOBAL_MAPPING_EXPORT bool Nektar::GlobalMapping::Mapping::HasConstantJacobian ( )
inline

Get flag defining if mapping has constant Jacobian.

Definition at line 362 of file Mapping.h.

363 {
364 return m_constantJacobian;
365 }
bool m_constantJacobian
Flag defining if the Jacobian is constant.
Definition: Mapping.h:423

References m_constantJacobian.

Referenced by v_DotGradJacobian().

◆ InitObject()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::InitObject ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const TiXmlElement *  pMapping 
)
inline

Initialise the mapping object.

Definition at line 77 of file Mapping.h.

80 {
81 v_InitObject(pFields, pMapping);
82 }
virtual GLOBAL_MAPPING_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
Definition: Mapping.cpp:101

References v_InitObject().

Referenced by ReplaceField().

◆ IsDefined()

GLOBAL_MAPPING_EXPORT bool Nektar::GlobalMapping::Mapping::IsDefined ( )
inline

Get flag determining if the mapping was defined or is trivial.

Definition at line 368 of file Mapping.h.

369 {
370 return m_isDefined;
371 }

References m_isDefined.

◆ IsFromFunction()

GLOBAL_MAPPING_EXPORT bool Nektar::GlobalMapping::Mapping::IsFromFunction ( )
inline

Get flag defining if mapping is defined by a function.

Definition at line 350 of file Mapping.h.

351 {
352 return m_fromFunction;
353 }
bool m_fromFunction
Flag defining if the Mapping is defined by a function.
Definition: Mapping.h:427

References m_fromFunction.

◆ IsTimeDependent()

GLOBAL_MAPPING_EXPORT bool Nektar::GlobalMapping::Mapping::IsTimeDependent ( void  )
inline

Get flag defining if mapping is time-dependent.

Definition at line 338 of file Mapping.h.

339 {
340 return m_timeDependent;
341 }
bool m_timeDependent
Flag defining if the Mapping is time-dependent.
Definition: Mapping.h:425

References m_timeDependent.

Referenced by v_UpdateBCs().

◆ Load()

MappingSharedPtr Nektar::GlobalMapping::Mapping::Load ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields 
)
static

Return a pointer to the mapping, creating it on first call.

This function is responsible for loading the Mapping, guaranteeing that a single instance of this class exists. When it is first called, it creates a Mapping and returns a pointer to it. On subsequent calls, it just returns the pointer.

Parameters
pSessionSession reader object
pFieldsFields which will be used by the Mapping
Returns
Pointer to the Mapping

Definition at line 272 of file Mapping.cpp.

275{
276 if (!m_init)
277 {
278 TiXmlElement *vMapping = NULL;
279 string vType;
280 if (pSession->DefinesElement("Nektar/Mapping"))
281 {
282 vMapping = pSession->GetElement("Nektar/Mapping");
283 vType = vMapping->Attribute("TYPE");
284 m_isDefined = true;
286 vMapping, pSession->GetTimeLevel());
287 }
288 else
289 {
290 vType = "Translation";
291 }
292
294 pFields, vMapping);
295
296 m_init = true;
297 }
298
299 return m_mappingPtr;
300}
static MappingSharedPtr m_mappingPtr
Definition: Mapping.h:430
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool disableCheck=true)
Get XML elment time level (Parallel-in-Time)
MappingFactory & GetMappingFactory()
Declaration of the mapping factory singleton.
Definition: Mapping.cpp:53

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::GlobalMapping::GetMappingFactory(), Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel(), m_init, m_isDefined, and m_mappingPtr.

Referenced by Nektar::FieldUtils::ProcessMapping::GetMapping(), Nektar::MappingExtrapolate::MappingExtrapolate(), Nektar::SolverUtils::DriverAdaptive::v_Execute(), Nektar::SolverUtils::FilterAeroForces::v_Initialise(), Nektar::VCSMapping::v_InitObject(), Nektar::ForcingMovingBody::v_InitObject(), and Nektar::SolverUtils::EquationSystem::WriteFld().

◆ LowerIndex()

void Nektar::GlobalMapping::Mapping::LowerIndex ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Lower index of vector: \(v_{i} = g_{ij}*v^{j}\).

This function lowers the index of the contravariant vector \(v^{i}\), transforming it into its associated covariant vector \(v_{j}\) according to the relation

\[v_{j} = g_{ij}v^{i}\]

where \(g_{ij}\) is the metric tensor.

Parameters
inarrayComponents of the contravariant vector \(v^{i}\)
outarrayComponents of the covariant vector \(v_{j}\)

Definition at line 609 of file Mapping.cpp.

611{
612 if (inarray == outarray)
613 {
614 int physTot = m_fields[0]->GetTotPoints();
615 int nvel = m_nConvectiveFields;
616
617 for (int i = 0; i < nvel; i++)
618 {
619 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
620 }
621 v_LowerIndex(m_tmp, outarray);
622 }
623 else
624 {
625 v_LowerIndex(inarray, outarray);
626 }
627}
virtual GLOBAL_MAPPING_EXPORT void v_LowerIndex(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Mapping.cpp:714

References m_fields, m_nConvectiveFields, m_tmp, v_LowerIndex(), and Vmath::Vcopy().

◆ Output()

void Nektar::GlobalMapping::Mapping::Output ( LibUtilities::FieldMetaDataMap fieldMetaDataMap,
const std::string &  outname 
)

Output function called when a chk or fld file is written.

This function should be called before writing a fld or chk file. It updates the metadata with information from the Mapping. Also, if the mapping is not defined by a function, it writes a .map file containing the coordinates and velocity of the coordinates.

Parameters
fieldMetaDataMapMetadata of the output file
outnameName of the output file

Definition at line 310 of file Mapping.cpp.

312{
313 // Only do anything if mapping exists
314 if (m_isDefined)
315 {
316 fieldMetaDataMap["MappingCartesianVel"] = std::string("False");
317 if (m_fromFunction)
318 {
319 // Add metadata
320 fieldMetaDataMap["MappingType"] = std::string("Expression");
321 fieldMetaDataMap["MappingExpression"] = m_funcName;
322 if (m_timeDependent)
323 {
324 fieldMetaDataMap["MappingVelExpression"] = m_velFuncName;
325 }
326 }
327 else
328 {
329 int expdim = m_fields[0]->GetGraph()->GetMeshDimension();
330 string fieldNames[3] = {"x", "y", "z"};
331 string velFieldNames[3] = {"vx", "vy", "vz"};
332
333 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
334 m_fields[0]->GetFieldDefinitions();
335 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
336
337 int ncoeffs = m_fields[0]->GetNcoeffs();
338 Array<OneD, NekDouble> fieldcoeffs(ncoeffs);
339
340 bool wavespace = m_fields[0]->GetWaveSpace();
341 m_fields[0]->SetWaveSpace(false);
342 // copy coordinates Data into FieldData and set variable
343 for (int j = 0; j < expdim; ++j)
344 {
345 m_fields[0]->FwdTransLocalElmt(m_coords[j], fieldcoeffs);
346
347 for (int i = 0; i < FieldDef.size(); ++i)
348 {
349 // Could do a search here to find correct variable
350 FieldDef[i]->m_fields.push_back(fieldNames[j]);
351 m_fields[0]->AppendFieldData(FieldDef[i], FieldData[i],
352 fieldcoeffs);
353 }
354 }
355 if (m_timeDependent)
356 {
357 // copy coordinates velocity Data into FieldData and set
358 // variable
359 for (int j = 0; j < expdim; ++j)
360 {
361 m_fields[0]->FwdTransLocalElmt(m_coordsVel[j], fieldcoeffs);
362
363 for (int i = 0; i < FieldDef.size(); ++i)
364 {
365 // Could do a search here to find correct variable
366 FieldDef[i]->m_fields.push_back(velFieldNames[j]);
367 m_fields[0]->AppendFieldData(FieldDef[i], FieldData[i],
368 fieldcoeffs);
369 }
370 }
371 }
372
373 std::string outfile = outname;
374 outfile.erase(outfile.end() - 4, outfile.end());
375 outfile += ".map";
376
377 m_fld->Write(outfile, FieldDef, FieldData, fieldMetaDataMap);
378
379 // Write metadata to orginal output
380 fieldMetaDataMap["MappingType"] = std::string("File");
381 fieldMetaDataMap["FileName"] = outfile;
382
383 m_fields[0]->SetWaveSpace(wavespace);
384 }
385 }
386}
Array< OneD, Array< OneD, NekDouble > > m_coords
Array with the Cartesian coordinates.
Definition: Mapping.h:408
Array< OneD, Array< OneD, NekDouble > > m_coordsVel
Array with the velocity of the coordinates.
Definition: Mapping.h:410
std::string m_velFuncName
Name of the function containing the velocity of the coordinates.
Definition: Mapping.h:419
std::string m_funcName
Name of the function containing the coordinates.
Definition: Mapping.h:417

References m_coords, m_coordsVel, m_fields, m_fld, m_fromFunction, m_funcName, m_isDefined, m_timeDependent, and m_velFuncName.

◆ RaiseIndex()

void Nektar::GlobalMapping::Mapping::RaiseIndex ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Raise index of vector: \(v^{i} = g^{ij}*v_{j}\).

This function raises the index of the covariant vector \(v_{j}\), transforming it into its associated contravariant vector \(v^{i}\) according to the relation

\[v^{i} = g^{ij}v_{j}\]

where \(g^{ij}\) is the inverse of the metric tensor.

Parameters
inarrayComponents of the contravariant vector \(v^{i}\)
outarrayComponents of the covariant vector \(v_{j}\)

Definition at line 639 of file Mapping.cpp.

641{
642 if (inarray == outarray)
643 {
644 int physTot = m_fields[0]->GetTotPoints();
645 int nvel = m_nConvectiveFields;
646
647 for (int i = 0; i < nvel; i++)
648 {
649 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
650 }
651 v_RaiseIndex(m_tmp, outarray);
652 }
653 else
654 {
655 v_RaiseIndex(inarray, outarray);
656 }
657}
virtual GLOBAL_MAPPING_EXPORT void v_RaiseIndex(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Mapping.cpp:735

References m_fields, m_nConvectiveFields, m_tmp, v_RaiseIndex(), and Vmath::Vcopy().

Referenced by v_CurlCurlField(), and v_VelocityLaplacian().

◆ ReplaceField()

void Nektar::GlobalMapping::Mapping::ReplaceField ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields)

Replace the Expansion List used by the mapping.

This function replaces the expansion list contained in m_fields, and then proceeds reinitialising the Mapping with this new field.

Parameters
pFieldsNew field to be used by the Mapping

Definition at line 247 of file Mapping.cpp.

249{
250 m_fields = pFields;
251
252 TiXmlElement *vMapping = NULL;
253
254 if (m_session->DefinesElement("Nektar/Mapping"))
255 {
256 vMapping = m_session->GetElement("Nektar/Mapping");
258 vMapping, m_session->GetTimeLevel());
259 }
260 InitObject(pFields, vMapping);
261}
GLOBAL_MAPPING_EXPORT void InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
Initialise the mapping object.
Definition: Mapping.h:77

References Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel(), InitObject(), m_fields, and m_session.

◆ SetFromFunction()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::SetFromFunction ( const bool  value)
inline

Set flag defining if mapping is defined by a function.

Definition at line 356 of file Mapping.h.

357 {
358 m_fromFunction = value;
359 }

References m_fromFunction.

◆ SetTimeDependent()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::SetTimeDependent ( const bool  value)
inline

Set flag defining if mapping is time-dependent.

Definition at line 344 of file Mapping.h.

345 {
346 m_timeDependent = value;
347 }

References m_timeDependent.

◆ UpdateBCs()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::UpdateBCs ( const NekDouble  time)
inline

Update the Dirichlet Boundary Conditions when using Mappings.

Definition at line 378 of file Mapping.h.

379 {
380 v_UpdateBCs(time);
381 }
virtual GLOBAL_MAPPING_EXPORT void v_UpdateBCs(const NekDouble time)
Definition: Mapping.cpp:1055

References v_UpdateBCs().

◆ UpdateGeomInfo()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::UpdateGeomInfo ( )
inline

Recompute the metric terms of the Mapping.

Definition at line 395 of file Mapping.h.

396 {
398 }
virtual GLOBAL_MAPPING_EXPORT void v_UpdateGeomInfo()=0

References v_UpdateGeomInfo().

Referenced by v_InitObject(), and v_UpdateMapping().

◆ UpdateMapping()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::UpdateMapping ( const NekDouble  time,
const Array< OneD, Array< OneD, NekDouble > > &  coords = NullNekDoubleArrayOfArray,
const Array< OneD, Array< OneD, NekDouble > > &  coordsVel = NullNekDoubleArrayOfArray 
)
inline

Update the Mapping with new coordinates.

Definition at line 384 of file Mapping.h.

390 {
391 v_UpdateMapping(time, coords, coordsVel);
392 }
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)
Definition: Mapping.cpp:1213

References v_UpdateMapping().

◆ v_ApplyChristoffelContravar()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_ApplyChristoffelContravar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedpure virtual

◆ v_ApplyChristoffelCovar()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_ApplyChristoffelCovar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedpure virtual

◆ v_ContravarFromCartesian()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_ContravarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedpure virtual

◆ v_ContravarToCartesian()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_ContravarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedpure virtual

◆ v_CovarFromCartesian()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_CovarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedpure virtual

◆ v_CovarToCartesian()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_CovarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedpure virtual

◆ v_CurlCurlField()

void Nektar::GlobalMapping::Mapping::v_CurlCurlField ( Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const bool  generalized 
)
protectedvirtual

Definition at line 993 of file Mapping.cpp.

996{
997 int physTot = m_fields[0]->GetTotPoints();
998 int nvel = m_nConvectiveFields;
999
1000 // Set wavespace to false and store current value
1001 bool wavespace = m_fields[0]->GetWaveSpace();
1002 m_fields[0]->SetWaveSpace(false);
1003
1004 // For implicit treatment of viscous terms, we want the generalized
1005 // curlcurl and for explicit treatment, we want the cartesian one.
1006 if (generalized)
1007 {
1008 // Get the second derivatives u^i_{,jk}
1009 Array<OneD, Array<OneD, NekDouble>> tmp(nvel);
1010 Array<OneD, Array<OneD, NekDouble>> ddU(nvel * nvel * nvel);
1011 gradgradU(inarray, ddU);
1012
1013 // Raise index to obtain A^{ip}_{k} = g^pj u^i_{,jk}
1014 for (int i = 0; i < nvel; ++i)
1015 {
1016 for (int k = 0; k < nvel; ++k)
1017 {
1018 // Copy to wk
1019 for (int j = 0; j < nvel; ++j)
1020 {
1021 tmp[j] = ddU[i * nvel * nvel + j * nvel + k];
1022 }
1023 RaiseIndex(tmp, m_tmp);
1024 for (int p = 0; p < nvel; ++p)
1025 {
1026 Vmath::Vcopy(physTot, m_tmp[p], 1,
1027 ddU[i * nvel * nvel + p * nvel + k], 1);
1028 }
1029 }
1030 }
1031 // The curlcurl is g^ji u^k_{kj} - g^jk u^i_kj = A^{ki}_k - A^{ik}_k
1032 for (int i = 0; i < nvel; ++i)
1033 {
1034 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
1035 for (int k = 0; k < nvel; ++k)
1036 {
1037 Vmath::Vadd(physTot, outarray[i], 1,
1038 ddU[k * nvel * nvel + i * nvel + k], 1, outarray[i],
1039 1);
1040 Vmath::Vsub(physTot, outarray[i], 1,
1041 ddU[i * nvel * nvel + k * nvel + k], 1, outarray[i],
1042 1);
1043 }
1044 }
1045 }
1046 else
1047 {
1048 m_fields[0]->CurlCurl(inarray, outarray);
1049 }
1050
1051 // Restore value of wavespace
1052 m_fields[0]->SetWaveSpace(wavespace);
1053}
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.
Definition: Mapping.h:306
GLOBAL_MAPPING_EXPORT void RaiseIndex(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Raise index of vector: .
Definition: Mapping.cpp:639
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:354
void 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:414

References gradgradU(), m_fields, m_nConvectiveFields, m_tmp, CellMLToNektar.cellml_metadata::p, RaiseIndex(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vsub().

Referenced by CurlCurlField().

◆ v_Divergence()

void Nektar::GlobalMapping::Mapping::v_Divergence ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 756 of file Mapping.cpp.

758{
759 int physTot = m_fields[0]->GetTotPoints();
760 Array<OneD, NekDouble> wk(physTot, 0.0);
761
762 Vmath::Zero(physTot, outarray, 1);
763
764 // Set wavespace to false and store current value
765 bool wavespace = m_fields[0]->GetWaveSpace();
766 m_fields[0]->SetWaveSpace(false);
767
768 // Get Mapping Jacobian
769 Array<OneD, NekDouble> Jac(physTot, 0.0);
770 GetJacobian(Jac);
771
772 for (int i = 0; i < m_nConvectiveFields; ++i)
773 {
774 Vmath::Vmul(physTot, Jac, 1, inarray[i], 1, wk, 1); // J*Ui
775 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i], wk,
776 wk); // (J*Ui)_i
777 Vmath::Vadd(physTot, wk, 1, outarray, 1, outarray, 1);
778 }
779 Vmath::Vdiv(physTot, outarray, 1, Jac, 1, outarray, 1); // 1/J*(J*Ui)_i
780
781 m_fields[0]->SetWaveSpace(wavespace);
782}
GLOBAL_MAPPING_EXPORT void GetJacobian(Array< OneD, NekDouble > &outarray)
Get the Jacobian of the transformation.
Definition: Mapping.h:155
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
void 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:280

References Nektar::MultiRegions::DirCartesianMap, GetJacobian(), m_fields, m_nConvectiveFields, Vmath::Vadd(), Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Zero().

Referenced by Divergence().

◆ v_DotGradJacobian()

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

Reimplemented in Nektar::GlobalMapping::MappingTranslation, Nektar::GlobalMapping::MappingXofXZ, Nektar::GlobalMapping::MappingXofZ, and Nektar::GlobalMapping::MappingXYofZ.

Definition at line 685 of file Mapping.cpp.

688{
689 int physTot = m_fields[0]->GetTotPoints();
690
691 outarray = Array<OneD, NekDouble>(physTot, 0.0);
692 if (!HasConstantJacobian())
693 {
694 // Set wavespace to false and store current value
695 bool wavespace = m_fields[0]->GetWaveSpace();
696 m_fields[0]->SetWaveSpace(false);
697
698 // Get Mapping Jacobian
699 Array<OneD, NekDouble> Jac(physTot, 0.0);
700 GetJacobian(Jac);
701
702 // Calculate inarray . grad(Jac)
703 Array<OneD, NekDouble> wk(physTot, 0.0);
704 for (int i = 0; i < m_nConvectiveFields; ++i)
705 {
706 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i], Jac, wk);
707 Vmath::Vvtvp(physTot, inarray[i], 1, wk, 1, outarray, 1, outarray,
708 1);
709 }
710 m_fields[0]->SetWaveSpace(wavespace);
711 }
712}
GLOBAL_MAPPING_EXPORT bool HasConstantJacobian()
Get flag defining if mapping has constant Jacobian.
Definition: Mapping.h:362
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569

References Nektar::MultiRegions::DirCartesianMap, GetJacobian(), HasConstantJacobian(), m_fields, m_nConvectiveFields, and Vmath::Vvtvp().

Referenced by DotGradJacobian().

◆ v_GetCartesianCoordinates()

void Nektar::GlobalMapping::Mapping::v_GetCartesianCoordinates ( Array< OneD, NekDouble > &  out0,
Array< OneD, NekDouble > &  out1,
Array< OneD, NekDouble > &  out2 
)
protectedvirtual

Definition at line 659 of file Mapping.cpp.

662{
663 int physTot = m_fields[0]->GetTotPoints();
664
665 out0 = Array<OneD, NekDouble>(physTot, 0.0);
666 out1 = Array<OneD, NekDouble>(physTot, 0.0);
667 out2 = Array<OneD, NekDouble>(physTot, 0.0);
668
669 Vmath::Vcopy(physTot, m_coords[0], 1, out0, 1);
670 Vmath::Vcopy(physTot, m_coords[1], 1, out1, 1);
671 Vmath::Vcopy(physTot, m_coords[2], 1, out2, 1);
672}

References m_coords, m_fields, and Vmath::Vcopy().

Referenced by GetCartesianCoordinates().

◆ v_GetCoordVelocity()

void Nektar::GlobalMapping::Mapping::v_GetCoordVelocity ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
protectedvirtual

Definition at line 674 of file Mapping.cpp.

675{
676 int physTot = m_fields[0]->GetTotPoints();
677
678 for (int i = 0; i < m_nConvectiveFields; ++i)
679 {
680 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
681 Vmath::Vcopy(physTot, m_coordsVel[i], 1, outarray[i], 1);
682 }
683}

References m_coordsVel, m_fields, m_nConvectiveFields, and Vmath::Vcopy().

Referenced by GetCoordVelocity().

◆ v_GetInvMetricTensor()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_GetInvMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
protectedpure virtual

◆ v_GetJacobian()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_GetJacobian ( Array< OneD, NekDouble > &  outarray)
protectedpure virtual

◆ v_GetMetricTensor()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_GetMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
protectedpure virtual

◆ v_gradgradU()

void Nektar::GlobalMapping::Mapping::v_gradgradU ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Definition at line 885 of file Mapping.cpp.

887{
888 int physTot = m_fields[0]->GetTotPoints();
889 int nvel = m_nConvectiveFields;
890
891 // Declare variables
892 outarray = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel * nvel);
893 Array<OneD, Array<OneD, NekDouble>> tmp(nvel);
894 for (int i = 0; i < nvel * nvel * nvel; i++)
895 {
896 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
897 }
898
899 // Set wavespace to false and store current value
900 bool wavespace = m_fields[0]->GetWaveSpace();
901 m_fields[0]->SetWaveSpace(false);
902
903 // Calculate vector gradient u^i_(,j) = du^i/dx^j + {i,pj}*u^p
905 for (int i = 0; i < nvel; i++)
906 {
907 if (nvel == 2)
908 {
909 m_fields[0]->PhysDeriv(inarray[i], m_wk2[i * nvel + 0],
910 m_wk2[i * nvel + 1]);
911 }
912 else
913 {
914 m_fields[0]->PhysDeriv(inarray[i], m_wk2[i * nvel + 0],
915 m_wk2[i * nvel + 1], m_wk2[i * nvel + 2]);
916 }
917 for (int j = 0; j < nvel; j++)
918 {
919 Vmath::Vadd(physTot, m_wk1[i * nvel + j], 1, m_wk2[i * nvel + j], 1,
920 m_wk1[i * nvel + j], 1);
921 }
922 }
923
924 //
925 // Calculate (u^i_,j),k
926 //
927
928 // Step 1 : d(u^i_,j))/d(x^k)
929 for (int i = 0; i < nvel; i++)
930 {
931 for (int j = 0; j < nvel; j++)
932 {
933 if (nvel == 2)
934 {
935 m_fields[0]->PhysDeriv(
936 m_wk1[i * nvel + j],
937 outarray[i * nvel * nvel + j * nvel + 0],
938 outarray[i * nvel * nvel + j * nvel + 1]);
939 }
940 else
941 {
942 m_fields[0]->PhysDeriv(
943 m_wk1[i * nvel + j],
944 outarray[i * nvel * nvel + j * nvel + 0],
945 outarray[i * nvel * nvel + j * nvel + 1],
946 outarray[i * nvel * nvel + j * nvel + 2]);
947 }
948 }
949 }
950
951 // Step 2: d(u^i_,j)/d(x^k) - {p,jk}*u^i_,p
952 for (int i = 0; i < nvel; i++)
953 {
954 for (int p = 0; p < nvel; p++)
955 {
956 tmp[p] = m_wk1[i * nvel + p];
957 }
959 for (int j = 0; j < nvel; j++)
960 {
961 for (int k = 0; k < nvel; k++)
962 {
963 Vmath::Vsub(physTot, outarray[i * nvel * nvel + j * nvel + k],
964 1, m_wk2[j * nvel + k], 1,
965 outarray[i * nvel * nvel + j * nvel + k], 1);
966 }
967 }
968 }
969
970 // Step 3: d(u^i_,j)/d(x^k) - {p,jk}*u^i_,p + {i,pk} u^p_,j
971 for (int j = 0; j < nvel; j++)
972 {
973 for (int p = 0; p < nvel; p++)
974 {
975 tmp[p] = m_wk1[p * nvel + j];
976 }
978 for (int i = 0; i < nvel; i++)
979 {
980 for (int k = 0; k < nvel; k++)
981 {
982 Vmath::Vadd(physTot, outarray[i * nvel * nvel + j * nvel + k],
983 1, m_wk2[i * nvel + k], 1,
984 outarray[i * nvel * nvel + j * nvel + k], 1);
985 }
986 }
987 }
988
989 // Restore value of wavespace
990 m_fields[0]->SetWaveSpace(wavespace);
991}
Array< OneD, Array< OneD, NekDouble > > m_wk2
Definition: Mapping.h:436
Array< OneD, Array< OneD, NekDouble > > m_wk1
Definition: Mapping.h:435
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.
Definition: Mapping.h:211
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.
Definition: Mapping.h:229

References ApplyChristoffelContravar(), ApplyChristoffelCovar(), m_fields, m_nConvectiveFields, m_wk1, m_wk2, CellMLToNektar.cellml_metadata::p, Vmath::Vadd(), and Vmath::Vsub().

Referenced by gradgradU().

◆ v_InitObject()

void Nektar::GlobalMapping::Mapping::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 in Nektar::GlobalMapping::MappingGeneral, Nektar::GlobalMapping::MappingTranslation, Nektar::GlobalMapping::MappingXofXZ, Nektar::GlobalMapping::MappingXofZ, Nektar::GlobalMapping::MappingXYofXY, and Nektar::GlobalMapping::MappingXYofZ.

Definition at line 101 of file Mapping.cpp.

104{
105 boost::ignore_unused(pFields);
106
107 int phystot = m_fields[0]->GetTotPoints();
108 m_fromFunction = true;
109 // Initialise variables
110 m_coords = Array<OneD, Array<OneD, NekDouble>>(3);
111 m_coordsVel = Array<OneD, Array<OneD, NekDouble>>(3);
112 Array<OneD, Array<OneD, NekDouble>> coords(3);
113 for (int i = 0; i < 3; i++)
114 {
115 m_coords[i] = Array<OneD, NekDouble>(phystot);
116 m_coordsVel[i] = Array<OneD, NekDouble>(phystot);
117 coords[i] = Array<OneD, NekDouble>(phystot);
118 }
119
120 // Check if mapping is defined as time-dependent
121 const TiXmlElement *timeDep = pMapping->FirstChildElement("TIMEDEPENDENT");
122 if (timeDep)
123 {
124 string sTimeDep = timeDep->GetText();
125 m_timeDependent = (boost::iequals(sTimeDep, "true")) ||
126 (boost::iequals(sTimeDep, "yes"));
127 }
128 else
129 {
130 m_timeDependent = false;
131 }
132
133 // Load coordinates
134 string fieldNames[3] = {"x", "y", "z"};
135 const TiXmlElement *funcNameElmt = pMapping->FirstChildElement("COORDS");
136 if (funcNameElmt)
137 {
138 m_funcName = funcNameElmt->GetText();
139 ASSERTL0(m_session->DefinesFunction(m_funcName),
140 "Function '" + m_funcName + "' not defined.");
141
142 // Get coordinates in the domain
143 m_fields[0]->GetCoords(coords[0], coords[1], coords[2]);
144
145 std::string s_FieldStr;
146 // Check if function from session file defines each component
147 // and evaluate them, otherwise use trivial transformation
148 for (int i = 0; i < 3; i++)
149 {
150 s_FieldStr = fieldNames[i];
151 if (m_session->DefinesFunction(m_funcName, s_FieldStr))
152 {
154 m_funcName);
155 if (i == 2 && m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
156 {
157 ASSERTL0(
158 false,
159 "3DH1D does not support mapping in the z-direction.");
160 }
161 }
162 else
163 {
164 // This coordinate is not defined, so use (x^i)' = x^i
165 Vmath::Vcopy(phystot, coords[i], 1, m_coords[i], 1);
166 }
167 }
168 }
169 else
170 {
171 m_fields[0]->GetCoords(coords[0], coords[1], coords[2]);
172 for (int i = 0; i < 3; i++)
173 {
174 // Use (x^i)' = x^i as default. This can be useful if we
175 // have a time-dependent mapping, and then only the
176 // initial mapping will be trivial
177 Vmath::Vcopy(phystot, coords[i], 1, m_coords[i], 1);
178 }
179 }
180
181 // Load coordinate velocity if they are defined,
182 // otherwise use zero to make it general
183 string velFieldNames[3] = {"vx", "vy", "vz"};
184 const TiXmlElement *velFuncNameElmt = pMapping->FirstChildElement("VEL");
185 if (velFuncNameElmt)
186 {
187 m_velFuncName = velFuncNameElmt->GetText();
188 ASSERTL0(m_session->DefinesFunction(m_velFuncName),
189 "Function '" + m_velFuncName + "' not defined.");
190
191 std::string s_FieldStr;
192 // Check if function from session file defines each component
193 // and evaluate them, otherwise use 0
194 for (int i = 0; i < 3; i++)
195 {
196 s_FieldStr = velFieldNames[i];
197 if (m_session->DefinesFunction(m_velFuncName, s_FieldStr))
198 {
201 if (i == 2 && m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
202 {
203 ASSERTL0(
204 false,
205 "3DH1D does not support mapping in the z-direction.");
206 }
207 }
208 else
209 {
210 // This coordinate velocity is not defined, so use 0
211 Vmath::Zero(phystot, m_coordsVel[i], 1);
212 }
213 }
214 }
215 else
216 {
217 for (int i = 0; i < 3; i++)
218 {
219 Vmath::Zero(phystot, m_coordsVel[i], 1);
220 }
221 }
222
223 // Initialise workspace variables
224 int nvel = m_nConvectiveFields;
225 m_wk1 = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel);
226 m_wk2 = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel);
227 m_tmp = Array<OneD, Array<OneD, NekDouble>>(nvel);
228 for (int i = 0; i < nvel; i++)
229 {
230 m_tmp[i] = Array<OneD, NekDouble>(phystot, 0.0);
231 for (int j = 0; j < nvel; j++)
232 {
233 m_wk1[i * nvel + j] = Array<OneD, NekDouble>(phystot, 0.0);
234 m_wk2[i * nvel + j] = Array<OneD, NekDouble>(phystot, 0.0);
235 }
236 }
237
238 // Calculate information required by the particular mapping
240}
GLOBAL_MAPPING_EXPORT void UpdateGeomInfo()
Recompute the metric terms of the Mapping.
Definition: Mapping.h:395
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))
Definition: Mapping.cpp:406

References ASSERTL0, Nektar::MultiRegions::e3DH1D, EvaluateFunction(), m_coords, m_coordsVel, m_fields, m_fromFunction, m_funcName, m_nConvectiveFields, m_session, m_timeDependent, m_tmp, m_velFuncName, m_wk1, m_wk2, UpdateGeomInfo(), Vmath::Vcopy(), and Vmath::Zero().

Referenced by InitObject(), Nektar::GlobalMapping::MappingGeneral::v_InitObject(), Nektar::GlobalMapping::MappingTranslation::v_InitObject(), Nektar::GlobalMapping::MappingXofXZ::v_InitObject(), Nektar::GlobalMapping::MappingXofZ::v_InitObject(), Nektar::GlobalMapping::MappingXYofXY::v_InitObject(), and Nektar::GlobalMapping::MappingXYofZ::v_InitObject().

◆ v_LowerIndex()

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

Reimplemented in Nektar::GlobalMapping::MappingTranslation, Nektar::GlobalMapping::MappingXofXZ, and Nektar::GlobalMapping::MappingXofZ.

Definition at line 714 of file Mapping.cpp.

716{
717 int physTot = m_fields[0]->GetTotPoints();
718 int nvel = m_nConvectiveFields;
719
720 Array<OneD, Array<OneD, NekDouble>> g(nvel * nvel);
721
723
724 for (int i = 0; i < nvel; i++)
725 {
726 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
727 for (int j = 0; j < nvel; j++)
728 {
729 Vmath::Vvtvp(physTot, g[i * nvel + j], 1, inarray[j], 1,
730 outarray[i], 1, outarray[i], 1);
731 }
732 }
733}
GLOBAL_MAPPING_EXPORT void GetMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)
Get the metric tensor .
Definition: Mapping.h:177

References GetMetricTensor(), m_fields, m_nConvectiveFields, and Vmath::Vvtvp().

Referenced by LowerIndex().

◆ v_RaiseIndex()

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

Reimplemented in Nektar::GlobalMapping::MappingTranslation, Nektar::GlobalMapping::MappingXofXZ, and Nektar::GlobalMapping::MappingXofZ.

Definition at line 735 of file Mapping.cpp.

737{
738 int physTot = m_fields[0]->GetTotPoints();
739 int nvel = m_nConvectiveFields;
740
741 Array<OneD, Array<OneD, NekDouble>> g(nvel * nvel);
742
744
745 for (int i = 0; i < nvel; i++)
746 {
747 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
748 for (int j = 0; j < nvel; j++)
749 {
750 Vmath::Vvtvp(physTot, g[i * nvel + j], 1, inarray[j], 1,
751 outarray[i], 1, outarray[i], 1);
752 }
753 }
754}
GLOBAL_MAPPING_EXPORT void GetInvMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)
Get the inverse of metric tensor .
Definition: Mapping.h:184

References GetInvMetricTensor(), m_fields, m_nConvectiveFields, and Vmath::Vvtvp().

Referenced by RaiseIndex().

◆ v_UpdateBCs()

void Nektar::GlobalMapping::Mapping::v_UpdateBCs ( const NekDouble  time)
protectedvirtual

Definition at line 1055 of file Mapping.cpp.

1056{
1057 int physTot = m_fields[0]->GetTotPoints();
1058 int nvel = m_nConvectiveFields;
1059 int nfields = m_fields.size();
1060 int nbnds = m_fields[0]->GetBndConditions().size();
1061
1062 // Declare variables
1063 Array<OneD, const SpatialDomains::BoundaryConditionShPtr> BndConds;
1064 Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
1065
1066 Array<OneD, bool> isDirichlet(nfields);
1067
1068 Array<OneD, Array<OneD, NekDouble>> values(nfields);
1069 for (int i = 0; i < nfields; i++)
1070 {
1071 values[i] = Array<OneD, NekDouble>(physTot, 0.0);
1072 }
1073
1074 Array<OneD, Array<OneD, NekDouble>> tmp(nvel);
1075 Array<OneD, Array<OneD, NekDouble>> tmp2(nvel);
1076 Array<OneD, Array<OneD, NekDouble>> coordVel(nvel);
1077 for (int i = 0; i < nvel; i++)
1078 {
1079 tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
1080 tmp2[i] = Array<OneD, NekDouble>(physTot, 0.0);
1081 coordVel[i] = Array<OneD, NekDouble>(physTot, 0.0);
1082 }
1083
1084 // Get coordinates velocity in transformed system (for MovingBody regions)
1085 GetCoordVelocity(tmp);
1086 ContravarFromCartesian(tmp, coordVel);
1087
1088 // Get Cartesian coordinates for evaluating boundary conditions
1089 Array<OneD, Array<OneD, NekDouble>> coords(3);
1090 for (int dir = 0; dir < 3; dir++)
1091 {
1092 coords[dir] = Array<OneD, NekDouble>(physTot, 0.0);
1093 }
1094 GetCartesianCoordinates(coords[0], coords[1], coords[2]);
1095
1096 // Loop boundary conditions looking for Dirichlet bc's
1097 for (int n = 0; n < nbnds; ++n)
1098 {
1099 // Evaluate original Dirichlet boundary conditions in whole domain
1100 for (int i = 0; i < nfields; ++i)
1101 {
1102 BndConds = m_fields[i]->GetBndConditions();
1103 BndExp = m_fields[i]->GetBndCondExpansions();
1104 if (BndConds[n]->GetBoundaryConditionType() ==
1106 {
1107 isDirichlet[i] = true;
1108 // If we have the a velocity component
1109 // check if all vel bc's are also Dirichlet
1110 if (i < nvel)
1111 {
1112 for (int j = 0; j < nvel; ++j)
1113 {
1115 ->GetBndConditions()[n]
1116 ->GetBoundaryConditionType() ==
1118 "Mapping only supported when all velocity "
1119 "components have the same type of boundary "
1120 "conditions");
1121 }
1122 }
1123 // Check if bc is time-dependent
1124 ASSERTL0(!BndConds[n]->IsTimeDependent(),
1125 "Time-dependent Dirichlet boundary conditions not "
1126 "supported with mapping yet.");
1127
1128 // Get boundary condition
1129 LibUtilities::Equation condition =
1130 std::static_pointer_cast<
1131 SpatialDomains::DirichletBoundaryCondition>(BndConds[n])
1132 ->m_dirichletCondition;
1133 // Evaluate
1134 condition.Evaluate(coords[0], coords[1], coords[2], time,
1135 values[i]);
1136 }
1137 else
1138 {
1139 isDirichlet[i] = false;
1140 }
1141 }
1142 // Convert velocity vector to transformed system
1143 if (isDirichlet[0])
1144 {
1145 for (int i = 0; i < nvel; ++i)
1146 {
1147 Vmath::Vcopy(physTot, values[i], 1, tmp[i], 1);
1148 }
1149 ContravarFromCartesian(tmp, tmp2);
1150 for (int i = 0; i < nvel; ++i)
1151 {
1152 Vmath::Vcopy(physTot, tmp2[i], 1, values[i], 1);
1153 }
1154 }
1155
1156 // Now, project result to boundary
1157 for (int i = 0; i < nfields; ++i)
1158 {
1159 BndConds = m_fields[i]->GetBndConditions();
1160 BndExp = m_fields[i]->GetBndCondExpansions();
1161
1162 if (BndConds[n]->GetUserDefined() == "" ||
1163 BndConds[n]->GetUserDefined() == "MovingBody")
1164 {
1165 m_fields[i]->ExtractPhysToBnd(n, values[i],
1166 BndExp[n]->UpdatePhys());
1167
1168 // Apply MovingBody correction
1169 if ((i < nvel) && BndConds[n]->GetUserDefined() == "MovingBody")
1170 {
1171 // Get coordinate velocity on boundary
1172 Array<OneD, NekDouble> coordVelBnd(
1173 BndExp[n]->GetTotPoints());
1174 m_fields[i]->ExtractPhysToBnd(n, coordVel[i], coordVelBnd);
1175
1176 // Apply correction
1177 Vmath::Vadd(BndExp[n]->GetTotPoints(), coordVelBnd, 1,
1178 BndExp[n]->UpdatePhys(), 1,
1179 BndExp[n]->UpdatePhys(), 1);
1180 }
1181 }
1182 }
1183 }
1184
1185 // Finally, perform FwdTrans in all fields
1186 for (int i = 0; i < m_fields.size(); ++i)
1187 {
1188 // Get boundary condition information
1189 BndConds = m_fields[i]->GetBndConditions();
1190 BndExp = m_fields[i]->GetBndCondExpansions();
1191
1192 for (int n = 0; n < BndConds.size(); ++n)
1193 {
1194 if (BndConds[n]->GetBoundaryConditionType() ==
1196 {
1197 if (BndConds[n]->GetUserDefined() == "" ||
1198 BndConds[n]->GetUserDefined() == "MovingBody")
1199 {
1200 if (m_fields[i]->GetExpType() == MultiRegions::e3DH1D)
1201 {
1202 BndExp[n]->SetWaveSpace(false);
1203 }
1204
1205 BndExp[n]->FwdTransBndConstrained(
1206 BndExp[n]->GetPhys(), BndExp[n]->UpdateCoeffs());
1207 }
1208 }
1209 }
1210 }
1211}
GLOBAL_MAPPING_EXPORT void GetCartesianCoordinates(Array< OneD, NekDouble > &out0, Array< OneD, NekDouble > &out1, Array< OneD, NekDouble > &out2)
Get the Cartesian coordinates in the field.
Definition: Mapping.h:135
GLOBAL_MAPPING_EXPORT bool IsTimeDependent()
Get flag defining if mapping is time-dependent.
Definition: Mapping.h:338
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.
Definition: Mapping.cpp:547
GLOBAL_MAPPING_EXPORT void GetCoordVelocity(Array< OneD, Array< OneD, NekDouble > > &outarray)
Obtain the velocity of the coordinates.
Definition: Mapping.h:244

References ASSERTL0, ContravarFromCartesian(), Nektar::MultiRegions::e3DH1D, Nektar::SpatialDomains::eDirichlet, Nektar::LibUtilities::Equation::Evaluate(), GetCartesianCoordinates(), GetCoordVelocity(), IsTimeDependent(), m_fields, m_nConvectiveFields, Vmath::Vadd(), and Vmath::Vcopy().

Referenced by UpdateBCs().

◆ v_UpdateGeomInfo()

virtual GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::v_UpdateGeomInfo ( )
protectedpure virtual

◆ v_UpdateMapping()

void Nektar::GlobalMapping::Mapping::v_UpdateMapping ( const NekDouble  time,
const Array< OneD, Array< OneD, NekDouble > > &  coords = NullNekDoubleArrayOfArray,
const Array< OneD, Array< OneD, NekDouble > > &  coordsVel = NullNekDoubleArrayOfArray 
)
protectedvirtual

Definition at line 1213 of file Mapping.cpp.

1216{
1217 if (m_fromFunction)
1218 {
1219 std::string s_FieldStr;
1220 string fieldNames[3] = {"x", "y", "z"};
1221 string velFieldNames[3] = {"vx", "vy", "vz"};
1222 // Check if function from session file defines each component
1223 // and evaluate them, otherwise there is no need to update
1224 // coords
1225 for (int i = 0; i < 3; i++)
1226 {
1227 s_FieldStr = fieldNames[i];
1228 if (m_session->DefinesFunction(m_funcName, s_FieldStr))
1229 {
1231 m_funcName, time);
1232 if (i == 2 && m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
1233 {
1234 ASSERTL0(
1235 false,
1236 "3DH1D does not support mapping in the z-direction.");
1237 }
1238 }
1239 s_FieldStr = velFieldNames[i];
1240 if (m_session->DefinesFunction(m_velFuncName, s_FieldStr))
1241 {
1242 EvaluateFunction(m_fields, m_session, s_FieldStr,
1243 m_coordsVel[i], m_velFuncName, time);
1244 if (i == 2 && m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
1245 {
1246 ASSERTL0(
1247 false,
1248 "3DH1D does not support mapping in the z-direction.");
1249 }
1250 }
1251 }
1252 }
1253 else
1254 {
1255 int physTot = m_fields[0]->GetTotPoints();
1256 int nvel = m_nConvectiveFields;
1257 // Copy coordinates
1258 for (int i = 0; i < 3; i++)
1259 {
1260 Vmath::Vcopy(physTot, coords[i], 1, m_coords[i], 1);
1261 }
1262
1263 for (int i = 0; i < nvel; i++)
1264 {
1265 Vmath::Vcopy(physTot, coordsVel[i], 1, m_coordsVel[i], 1);
1266 }
1267 }
1268
1269 // Update the information required by the specific mapping
1271}

References ASSERTL0, Nektar::MultiRegions::e3DH1D, EvaluateFunction(), m_coords, m_coordsVel, m_fields, m_fromFunction, m_funcName, m_nConvectiveFields, m_session, m_velFuncName, UpdateGeomInfo(), and Vmath::Vcopy().

Referenced by UpdateMapping().

◆ v_VelocityLaplacian()

void Nektar::GlobalMapping::Mapping::v_VelocityLaplacian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  alpha 
)
protectedvirtual

Definition at line 784 of file Mapping.cpp.

787{
788 int physTot = m_fields[0]->GetTotPoints();
789 int nvel = m_nConvectiveFields;
790
791 Array<OneD, Array<OneD, NekDouble>> tmp(nvel);
792
793 // Set wavespace to false and store current value
794 bool wavespace = m_fields[0]->GetWaveSpace();
795 m_fields[0]->SetWaveSpace(false);
796
797 // Calculate vector gradient wk2 = u^i_(,k) = du^i/dx^k + {i,jk}*u^j
799 for (int i = 0; i < nvel; i++)
800 {
801 if (nvel == 2)
802 {
803 m_fields[0]->PhysDeriv(inarray[i], m_wk2[i * nvel + 0],
804 m_wk2[i * nvel + 1]);
805 }
806 else
807 {
808 m_fields[0]->PhysDeriv(inarray[i], m_wk2[i * nvel + 0],
809 m_wk2[i * nvel + 1], m_wk2[i * nvel + 2]);
810 }
811 for (int k = 0; k < nvel; k++)
812 {
813 Vmath::Vadd(physTot, m_wk1[i * nvel + k], 1, m_wk2[i * nvel + k], 1,
814 m_wk1[i * nvel + k], 1);
815 }
816 }
817 // Calculate wk1 = A^(ij) = g^(jk)*u^i_(,k)
818 for (int i = 0; i < nvel; i++)
819 {
820 for (int k = 0; k < nvel; k++)
821 {
822 tmp[k] = m_wk1[i * nvel + k];
823 }
824 RaiseIndex(tmp, m_tmp);
825 for (int j = 0; j < nvel; j++)
826 {
827 Vmath::Vcopy(physTot, m_tmp[j], 1, m_wk1[i * nvel + j], 1);
828 }
829 }
830 //
831 // Calculate L(U)^i = (A^(ij))_(,j) - alpha*d^2(u^i)/dx^jdx^j
832 //
833
834 // Step 1 :
835 // d(A^(ij) - alpha*du^i/dx^j)/d(x^j)
836 for (int i = 0; i < nvel; i++)
837 {
838 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
839 for (int j = 0; j < nvel; j++)
840 {
841 Vmath::Smul(physTot, alpha, m_wk2[i * nvel + j], 1, m_tmp[0], 1);
842 Vmath::Vsub(physTot, m_wk1[i * nvel + j], 1, m_tmp[0], 1, m_tmp[0],
843 1);
844
845 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], m_tmp[0],
846 m_tmp[1]);
847 Vmath::Vadd(physTot, outarray[i], 1, m_tmp[1], 1, outarray[i], 1);
848 }
849 }
850
851 // Step 2: d(A^(ij))/d(x^j) + {j,pj}*A^(ip)
852 for (int i = 0; i < nvel; i++)
853 {
854 for (int p = 0; p < nvel; p++)
855 {
856 tmp[p] = m_wk1[i * nvel + p];
857 }
859 for (int j = 0; j < nvel; j++)
860 {
861 Vmath::Vadd(physTot, outarray[i], 1, m_wk2[j * nvel + j], 1,
862 outarray[i], 1);
863 }
864 }
865
866 // Step 3: d(A^(ij))/d(x^j) + {j,pj}*A^(ip) + {i,pj} A^(pj)
867 for (int j = 0; j < nvel; j++)
868 {
869 for (int p = 0; p < nvel; p++)
870 {
871 tmp[p] = m_wk1[p * nvel + j];
872 }
874 for (int i = 0; i < nvel; i++)
875 {
876 Vmath::Vadd(physTot, outarray[i], 1, m_wk2[i * nvel + j], 1,
877 outarray[i], 1);
878 }
879 }
880
881 // Restore value of wavespace
882 m_fields[0]->SetWaveSpace(wavespace);
883}
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245

References ApplyChristoffelContravar(), Nektar::MultiRegions::DirCartesianMap, m_fields, m_nConvectiveFields, m_tmp, m_wk1, m_wk2, CellMLToNektar.cellml_metadata::p, RaiseIndex(), Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vsub().

Referenced by VelocityLaplacian().

◆ VelocityLaplacian()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::Mapping::VelocityLaplacian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  alpha = 0.0 
)
inline

Generalised (correction to the) velocity Laplacian operator.

This function is used to calculate a correction defined as the difference between the generalised Laplacian and the original Laplacian multiplied by a constant \(\alpha\), resulting in

\[ L^i = g^{jk}u^{i}_{,jk} - \alpha \frac{\partial^2 x^i} {\partial x^j \partial x^j}\]

By default, \(\alpha\) is zero, resulting in the generalised Laplacian.

Parameters
inarrayContravariant vector \(u^i\)
outarrayResult of the operation
alphaThe constant \(\alpha\)

Definition at line 289 of file Mapping.h.

293 {
294 v_VelocityLaplacian(inarray, outarray, alpha);
295 }
virtual GLOBAL_MAPPING_EXPORT void v_VelocityLaplacian(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble alpha)
Definition: Mapping.cpp:784

References v_VelocityLaplacian().

Member Data Documentation

◆ m_constantJacobian

bool Nektar::GlobalMapping::Mapping::m_constantJacobian
protected

◆ m_coords

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::Mapping::m_coords
protected

◆ m_coordsVel

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::Mapping::m_coordsVel
protected

Array with the velocity of the coordinates.

Definition at line 410 of file Mapping.h.

Referenced by Output(), v_GetCoordVelocity(), v_InitObject(), Nektar::GlobalMapping::MappingTranslation::v_InitObject(), and v_UpdateMapping().

◆ m_fields

Array<OneD, MultiRegions::ExpListSharedPtr> Nektar::GlobalMapping::Mapping::m_fields
protected

Definition at line 406 of file Mapping.h.

Referenced by Nektar::GlobalMapping::MappingGeneral::CalculateChristoffel(), Nektar::GlobalMapping::MappingXYofXY::CalculateChristoffel(), Nektar::GlobalMapping::MappingXYofXY::CalculateMetricTensor(), Nektar::GlobalMapping::MappingGeneral::CalculateMetricTerms(), ContravarFromCartesian(), ContravarToCartesian(), CovarFromCartesian(), CovarToCartesian(), EvaluateFunction(), LowerIndex(), Mapping(), Output(), RaiseIndex(), ReplaceField(), Nektar::GlobalMapping::MappingGeneral::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingTranslation::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXofZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXYofXY::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXYofZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingGeneral::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingTranslation::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXofZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXYofXY::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXYofZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingGeneral::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingTranslation::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXofZ::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingGeneral::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingTranslation::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXofZ::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingGeneral::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingTranslation::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXofZ::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingGeneral::v_CovarToCartesian(), Nektar::GlobalMapping::MappingTranslation::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXofZ::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_CovarToCartesian(), v_CurlCurlField(), v_Divergence(), v_DotGradJacobian(), Nektar::GlobalMapping::MappingTranslation::v_DotGradJacobian(), Nektar::GlobalMapping::MappingXofXZ::v_DotGradJacobian(), Nektar::GlobalMapping::MappingXofZ::v_DotGradJacobian(), Nektar::GlobalMapping::MappingXYofZ::v_DotGradJacobian(), v_GetCartesianCoordinates(), v_GetCoordVelocity(), Nektar::GlobalMapping::MappingGeneral::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingTranslation::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXofXZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXofZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXYofXY::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXYofZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingGeneral::v_GetJacobian(), Nektar::GlobalMapping::MappingTranslation::v_GetJacobian(), Nektar::GlobalMapping::MappingXofXZ::v_GetJacobian(), Nektar::GlobalMapping::MappingXofZ::v_GetJacobian(), Nektar::GlobalMapping::MappingXYofXY::v_GetJacobian(), Nektar::GlobalMapping::MappingXYofZ::v_GetJacobian(), Nektar::GlobalMapping::MappingGeneral::v_GetMetricTensor(), Nektar::GlobalMapping::MappingTranslation::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXofXZ::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXofZ::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXYofXY::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXYofZ::v_GetMetricTensor(), v_gradgradU(), v_InitObject(), Nektar::GlobalMapping::MappingTranslation::v_InitObject(), v_LowerIndex(), Nektar::GlobalMapping::MappingTranslation::v_LowerIndex(), Nektar::GlobalMapping::MappingXofXZ::v_LowerIndex(), Nektar::GlobalMapping::MappingXofZ::v_LowerIndex(), v_RaiseIndex(), Nektar::GlobalMapping::MappingTranslation::v_RaiseIndex(), Nektar::GlobalMapping::MappingXofXZ::v_RaiseIndex(), Nektar::GlobalMapping::MappingXofZ::v_RaiseIndex(), v_UpdateBCs(), Nektar::GlobalMapping::MappingXofXZ::v_UpdateGeomInfo(), Nektar::GlobalMapping::MappingXofZ::v_UpdateGeomInfo(), Nektar::GlobalMapping::MappingXYofXY::v_UpdateGeomInfo(), Nektar::GlobalMapping::MappingXYofZ::v_UpdateGeomInfo(), v_UpdateMapping(), and v_VelocityLaplacian().

◆ m_fld

LibUtilities::FieldIOSharedPtr Nektar::GlobalMapping::Mapping::m_fld
protected

Definition at line 404 of file Mapping.h.

Referenced by Mapping(), and Output().

◆ m_fromFunction

bool Nektar::GlobalMapping::Mapping::m_fromFunction
protected

Flag defining if the Mapping is defined by a function.

Definition at line 427 of file Mapping.h.

Referenced by IsFromFunction(), Output(), SetFromFunction(), v_InitObject(), and v_UpdateMapping().

◆ m_funcName

std::string Nektar::GlobalMapping::Mapping::m_funcName
protected

Name of the function containing the coordinates.

Definition at line 417 of file Mapping.h.

Referenced by Output(), v_InitObject(), and v_UpdateMapping().

◆ m_GeometricInfo

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::Mapping::m_GeometricInfo
protected

Array with metric terms of the mapping.

Definition at line 412 of file Mapping.h.

Referenced by Nektar::GlobalMapping::MappingXYofXY::CalculateMetricTensor(), Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXofZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXYofZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXofZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXYofZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXofXZ::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXofZ::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXofZ::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXofZ::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXofZ::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXYofZ::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXofXZ::v_DotGradJacobian(), Nektar::GlobalMapping::MappingXofXZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXofZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXYofZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXofXZ::v_GetJacobian(), Nektar::GlobalMapping::MappingXYofXY::v_GetJacobian(), Nektar::GlobalMapping::MappingXofXZ::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXofZ::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXYofZ::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXofXZ::v_LowerIndex(), Nektar::GlobalMapping::MappingXofZ::v_LowerIndex(), Nektar::GlobalMapping::MappingXofXZ::v_RaiseIndex(), Nektar::GlobalMapping::MappingXofZ::v_RaiseIndex(), Nektar::GlobalMapping::MappingXofXZ::v_UpdateGeomInfo(), Nektar::GlobalMapping::MappingXofZ::v_UpdateGeomInfo(), Nektar::GlobalMapping::MappingXYofXY::v_UpdateGeomInfo(), and Nektar::GlobalMapping::MappingXYofZ::v_UpdateGeomInfo().

◆ m_init

bool Nektar::GlobalMapping::Mapping::m_init = false
staticprotected

Definition at line 431 of file Mapping.h.

Referenced by Load().

◆ m_isDefined

bool Nektar::GlobalMapping::Mapping::m_isDefined = false
staticprotected

Definition at line 432 of file Mapping.h.

Referenced by IsDefined(), Load(), and Output().

◆ m_mappingPtr

MappingSharedPtr Nektar::GlobalMapping::Mapping::m_mappingPtr = MappingSharedPtr()
staticprotected

Definition at line 430 of file Mapping.h.

Referenced by Load().

◆ m_nConvectiveFields

int Nektar::GlobalMapping::Mapping::m_nConvectiveFields
protected

Number of velocity components.

Definition at line 414 of file Mapping.h.

Referenced by Nektar::GlobalMapping::MappingGeneral::CalculateChristoffel(), Nektar::GlobalMapping::MappingXYofXY::CalculateChristoffel(), Nektar::GlobalMapping::MappingGeneral::CalculateMetricTerms(), ContravarFromCartesian(), ContravarToCartesian(), CovarFromCartesian(), CovarToCartesian(), LowerIndex(), Mapping(), RaiseIndex(), Nektar::GlobalMapping::MappingGeneral::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingTranslation::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXofZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXYofXY::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingXYofZ::v_ApplyChristoffelContravar(), Nektar::GlobalMapping::MappingGeneral::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingTranslation::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXofZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXYofXY::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingXYofZ::v_ApplyChristoffelCovar(), Nektar::GlobalMapping::MappingGeneral::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingTranslation::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_ContravarFromCartesian(), Nektar::GlobalMapping::MappingGeneral::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingTranslation::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_ContravarToCartesian(), Nektar::GlobalMapping::MappingGeneral::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingTranslation::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_CovarFromCartesian(), Nektar::GlobalMapping::MappingGeneral::v_CovarToCartesian(), Nektar::GlobalMapping::MappingTranslation::v_CovarToCartesian(), Nektar::GlobalMapping::MappingXYofXY::v_CovarToCartesian(), v_CurlCurlField(), v_Divergence(), v_DotGradJacobian(), v_GetCoordVelocity(), Nektar::GlobalMapping::MappingGeneral::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingTranslation::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXofXZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXofZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXYofXY::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingXYofZ::v_GetInvMetricTensor(), Nektar::GlobalMapping::MappingGeneral::v_GetMetricTensor(), Nektar::GlobalMapping::MappingTranslation::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXofXZ::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXofZ::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXYofXY::v_GetMetricTensor(), Nektar::GlobalMapping::MappingXYofZ::v_GetMetricTensor(), v_gradgradU(), v_InitObject(), Nektar::GlobalMapping::MappingGeneral::v_InitObject(), Nektar::GlobalMapping::MappingTranslation::v_InitObject(), Nektar::GlobalMapping::MappingXofXZ::v_InitObject(), Nektar::GlobalMapping::MappingXofZ::v_InitObject(), Nektar::GlobalMapping::MappingXYofXY::v_InitObject(), Nektar::GlobalMapping::MappingXYofZ::v_InitObject(), v_LowerIndex(), Nektar::GlobalMapping::MappingTranslation::v_LowerIndex(), v_RaiseIndex(), Nektar::GlobalMapping::MappingTranslation::v_RaiseIndex(), v_UpdateBCs(), v_UpdateMapping(), and v_VelocityLaplacian().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::GlobalMapping::Mapping::m_session
protected

◆ m_timeDependent

bool Nektar::GlobalMapping::Mapping::m_timeDependent
protected

Flag defining if the Mapping is time-dependent.

Definition at line 425 of file Mapping.h.

Referenced by IsTimeDependent(), Output(), SetTimeDependent(), v_InitObject(), and Nektar::GlobalMapping::MappingTranslation::v_InitObject().

◆ m_tmp

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::Mapping::m_tmp
protected

◆ m_velFuncName

std::string Nektar::GlobalMapping::Mapping::m_velFuncName
protected

Name of the function containing the velocity of the coordinates.

Definition at line 419 of file Mapping.h.

Referenced by Output(), v_InitObject(), and v_UpdateMapping().

◆ m_wk1

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::Mapping::m_wk1
protected

Definition at line 435 of file Mapping.h.

Referenced by v_gradgradU(), v_InitObject(), and v_VelocityLaplacian().

◆ m_wk2

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::Mapping::m_wk2
protected

Definition at line 436 of file Mapping.h.

Referenced by v_gradgradU(), v_InitObject(), and v_VelocityLaplacian().