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 66 of file Mapping.h.

Constructor & Destructor Documentation

◆ ~Mapping()

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

Destructor.

Definition at line 70 of file Mapping.h.

71 {
72 }

◆ Mapping()

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

Constructor.

Definition at line 55 of file Mapping.cpp.

57 : m_session(pSession), m_fields(pFields)
58{
59 switch (m_fields[0]->GetExpType())
60 {
62 {
64 }
65 break;
66
68 {
70 }
71 break;
72
76 {
78 }
79 break;
80
81 default:
82 ASSERTL0(0, "Dimension not supported");
83 break;
84 }
85
87}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:412
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Mapping.h:400
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:404
LibUtilities::FieldIOSharedPtr m_fld
Definition: Mapping.h:402
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:195

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 209 of file Mapping.h.

212 {
213 v_ApplyChristoffelContravar(inarray, outarray);
214 }
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 227 of file Mapping.h.

230 {
231 v_ApplyChristoffelCovar(inarray, outarray);
232 }
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 541 of file Mapping.cpp.

544{
545 if (inarray == outarray)
546 {
547 int physTot = m_fields[0]->GetTotPoints();
548 int nvel = m_nConvectiveFields;
549
550 for (int i = 0; i < nvel; i++)
551 {
552 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
553 }
555 }
556 else
557 {
558 v_ContravarFromCartesian(inarray, outarray);
559 }
560}
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:435
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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 479 of file Mapping.cpp.

482{
483 if (inarray == outarray)
484 {
485 int physTot = m_fields[0]->GetTotPoints();
486 int nvel = m_nConvectiveFields;
487
488 for (int i = 0; i < nvel; i++)
489 {
490 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
491 }
492 v_ContravarToCartesian(m_tmp, outarray);
493 }
494 else
495 {
496 v_ContravarToCartesian(inarray, outarray);
497 }
498}
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 572 of file Mapping.cpp.

575{
576 if (inarray == outarray)
577 {
578 int physTot = m_fields[0]->GetTotPoints();
579 int nvel = m_nConvectiveFields;
580
581 for (int i = 0; i < nvel; i++)
582 {
583 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
584 }
585 v_CovarFromCartesian(m_tmp, outarray);
586 }
587 else
588 {
589 v_CovarFromCartesian(inarray, outarray);
590 }
591}
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 510 of file Mapping.cpp.

513{
514 if (inarray == outarray)
515 {
516 int physTot = m_fields[0]->GetTotPoints();
517 int nvel = m_nConvectiveFields;
518
519 for (int i = 0; i < nvel; i++)
520 {
521 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
522 }
523 v_CovarToCartesian(m_tmp, outarray);
524 }
525 else
526 {
527 v_CovarToCartesian(inarray, outarray);
528 }
529}
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 322 of file Mapping.h.

325 {
326 v_CurlCurlField(inarray, outarray, generalized);
327 }
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:987

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 264 of file Mapping.h.

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

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 167 of file Mapping.h.

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

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 400 of file Mapping.cpp.

405{
406 ASSERTL0(pSession->DefinesFunction(pFunctionName),
407 "Function '" + pFunctionName + "' does not exist.");
408
409 unsigned int nq = pFields[0]->GetNpoints();
410 if (pArray.size() != nq)
411 {
412 pArray = Array<OneD, NekDouble>(nq);
413 }
414
416 vType = pSession->GetFunctionType(pFunctionName, pFieldName);
418 {
419 Array<OneD, NekDouble> x0(nq);
420 Array<OneD, NekDouble> x1(nq);
421 Array<OneD, NekDouble> x2(nq);
422
423 pFields[0]->GetCoords(x0, x1, x2);
425 pSession->GetFunction(pFunctionName, pFieldName);
426
427 ffunc->Evaluate(x0, x1, x2, pTime, pArray);
428 }
429 else if (vType == LibUtilities::eFunctionTypeFile)
430 {
431 std::string filename =
432 pSession->GetFunctionFilename(pFunctionName, pFieldName);
433
434 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
435 std::vector<std::vector<NekDouble>> FieldData;
436 Array<OneD, NekDouble> vCoeffs(pFields[0]->GetNcoeffs());
437 Vmath::Zero(vCoeffs.size(), vCoeffs, 1);
438
440 LibUtilities::FieldIO::CreateForFile(pSession, filename);
441 fld->Import(filename, FieldDef, FieldData);
442
443 int idx = -1;
444 for (int i = 0; i < FieldDef.size(); ++i)
445 {
446 for (int j = 0; j < FieldDef[i]->m_fields.size(); ++j)
447 {
448 if (FieldDef[i]->m_fields[j] == pFieldName)
449 {
450 idx = j;
451 }
452 }
453
454 if (idx >= 0)
455 {
456 pFields[0]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
457 FieldDef[i]->m_fields[idx],
458 vCoeffs);
459 }
460 else
461 {
462 cout << "Field " + pFieldName + " not found." << endl;
463 }
464 }
465 pFields[0]->BwdTrans(vCoeffs, pArray);
466 }
467}
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:224
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

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 382 of file Mapping.cpp.

386{
387 ASSERTL0(pSession->DefinesFunction(pFunctionName),
388 "Function '" + pFunctionName + "' does not exist.");
389
391 pSession->GetFunction(pFunctionName, pFieldName);
392
393 Array<OneD, NekDouble> x0(1, 0.0);
394 Array<OneD, NekDouble> x1(1, 0.0);
395 Array<OneD, NekDouble> x2(1, 0.0);
396
397 ffunc->Evaluate(x0, x1, x2, pTime, pArray);
398}

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 133 of file Mapping.h.

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

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 242 of file Mapping.h.

244 {
245 v_GetCoordVelocity(outarray);
246 }
virtual GLOBAL_MAPPING_EXPORT void v_GetCoordVelocity(Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Mapping.cpp:668

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 182 of file Mapping.h.

184 {
185 v_GetInvMetricTensor(outarray);
186 }
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 153 of file Mapping.h.

154 {
155 v_GetJacobian(outarray);
156 }
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 175 of file Mapping.h.

177 {
178 v_GetMetricTensor(outarray);
179 }
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 304 of file Mapping.h.

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

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 360 of file Mapping.h.

361 {
362 return m_constantJacobian;
363 }
bool m_constantJacobian
Flag defining if the Jacobian is constant.
Definition: Mapping.h:421

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 75 of file Mapping.h.

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

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 366 of file Mapping.h.

367 {
368 return m_isDefined;
369 }

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 348 of file Mapping.h.

349 {
350 return m_fromFunction;
351 }
bool m_fromFunction
Flag defining if the Mapping is defined by a function.
Definition: Mapping.h:425

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 336 of file Mapping.h.

337 {
338 return m_timeDependent;
339 }
bool m_timeDependent
Flag defining if the Mapping is time-dependent.
Definition: Mapping.h:423

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 266 of file Mapping.cpp.

269{
270 if (!m_init)
271 {
272 TiXmlElement *vMapping = nullptr;
273 string vType;
274 if (pSession->DefinesElement("Nektar/Mapping"))
275 {
276 vMapping = pSession->GetElement("Nektar/Mapping");
277 vType = vMapping->Attribute("TYPE");
278 m_isDefined = true;
280 vMapping, pSession->GetTimeLevel());
281 }
282 else
283 {
284 vType = "Translation";
285 }
286
288 pFields, vMapping);
289
290 m_init = true;
291 }
292
293 return m_mappingPtr;
294}
static MappingSharedPtr m_mappingPtr
Definition: Mapping.h:428
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
Get XML elment time level (Parallel-in-Time)
MappingFactory & GetMappingFactory()
Declaration of the mapping factory singleton.
Definition: Mapping.cpp:49

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 603 of file Mapping.cpp.

605{
606 if (inarray == outarray)
607 {
608 int physTot = m_fields[0]->GetTotPoints();
609 int nvel = m_nConvectiveFields;
610
611 for (int i = 0; i < nvel; i++)
612 {
613 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
614 }
615 v_LowerIndex(m_tmp, outarray);
616 }
617 else
618 {
619 v_LowerIndex(inarray, outarray);
620 }
621}
virtual GLOBAL_MAPPING_EXPORT void v_LowerIndex(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Mapping.cpp:708

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 304 of file Mapping.cpp.

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

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 633 of file Mapping.cpp.

635{
636 if (inarray == outarray)
637 {
638 int physTot = m_fields[0]->GetTotPoints();
639 int nvel = m_nConvectiveFields;
640
641 for (int i = 0; i < nvel; i++)
642 {
643 Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
644 }
645 v_RaiseIndex(m_tmp, outarray);
646 }
647 else
648 {
649 v_RaiseIndex(inarray, outarray);
650 }
651}
virtual GLOBAL_MAPPING_EXPORT void v_RaiseIndex(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: Mapping.cpp:729

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 241 of file Mapping.cpp.

243{
244 m_fields = pFields;
245
246 TiXmlElement *vMapping = nullptr;
247
248 if (m_session->DefinesElement("Nektar/Mapping"))
249 {
250 vMapping = m_session->GetElement("Nektar/Mapping");
252 vMapping, m_session->GetTimeLevel());
253 }
254 InitObject(pFields, vMapping);
255}
GLOBAL_MAPPING_EXPORT void InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
Initialise the mapping object.
Definition: Mapping.h:75

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 354 of file Mapping.h.

355 {
356 m_fromFunction = value;
357 }

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 342 of file Mapping.h.

343 {
344 m_timeDependent = value;
345 }

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 376 of file Mapping.h.

377 {
378 v_UpdateBCs(time);
379 }
virtual GLOBAL_MAPPING_EXPORT void v_UpdateBCs(const NekDouble time)
Definition: Mapping.cpp:1049

References v_UpdateBCs().

◆ UpdateGeomInfo()

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

Recompute the metric terms of the Mapping.

Definition at line 393 of file Mapping.h.

394 {
396 }
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 382 of file Mapping.h.

388 {
389 v_UpdateMapping(time, coords, coordsVel);
390 }
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:1207

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 987 of file Mapping.cpp.

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

752{
753 int physTot = m_fields[0]->GetTotPoints();
754 Array<OneD, NekDouble> wk(physTot, 0.0);
755
756 Vmath::Zero(physTot, outarray, 1);
757
758 // Set wavespace to false and store current value
759 bool wavespace = m_fields[0]->GetWaveSpace();
760 m_fields[0]->SetWaveSpace(false);
761
762 // Get Mapping Jacobian
763 Array<OneD, NekDouble> Jac(physTot, 0.0);
764 GetJacobian(Jac);
765
766 for (int i = 0; i < m_nConvectiveFields; ++i)
767 {
768 Vmath::Vmul(physTot, Jac, 1, inarray[i], 1, wk, 1); // J*Ui
769 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i], wk,
770 wk); // (J*Ui)_i
771 Vmath::Vadd(physTot, wk, 1, outarray, 1, outarray, 1);
772 }
773 Vmath::Vdiv(physTot, outarray, 1, Jac, 1, outarray, 1); // 1/J*(J*Ui)_i
774
775 m_fields[0]->SetWaveSpace(wavespace);
776}
GLOBAL_MAPPING_EXPORT void GetJacobian(Array< OneD, NekDouble > &outarray)
Get the Jacobian of the transformation.
Definition: Mapping.h:153
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.hpp:126

References Nektar::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 679 of file Mapping.cpp.

682{
683 int physTot = m_fields[0]->GetTotPoints();
684
685 outarray = Array<OneD, NekDouble>(physTot, 0.0);
686 if (!HasConstantJacobian())
687 {
688 // Set wavespace to false and store current value
689 bool wavespace = m_fields[0]->GetWaveSpace();
690 m_fields[0]->SetWaveSpace(false);
691
692 // Get Mapping Jacobian
693 Array<OneD, NekDouble> Jac(physTot, 0.0);
694 GetJacobian(Jac);
695
696 // Calculate inarray . grad(Jac)
697 Array<OneD, NekDouble> wk(physTot, 0.0);
698 for (int i = 0; i < m_nConvectiveFields; ++i)
699 {
700 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i], Jac, wk);
701 Vmath::Vvtvp(physTot, inarray[i], 1, wk, 1, outarray, 1, outarray,
702 1);
703 }
704 m_fields[0]->SetWaveSpace(wavespace);
705 }
706}
GLOBAL_MAPPING_EXPORT bool HasConstantJacobian()
Get flag defining if mapping has constant Jacobian.
Definition: Mapping.h:360
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

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 653 of file Mapping.cpp.

656{
657 int physTot = m_fields[0]->GetTotPoints();
658
659 out0 = Array<OneD, NekDouble>(physTot, 0.0);
660 out1 = Array<OneD, NekDouble>(physTot, 0.0);
661 out2 = Array<OneD, NekDouble>(physTot, 0.0);
662
663 Vmath::Vcopy(physTot, m_coords[0], 1, out0, 1);
664 Vmath::Vcopy(physTot, m_coords[1], 1, out1, 1);
665 Vmath::Vcopy(physTot, m_coords[2], 1, out2, 1);
666}

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 668 of file Mapping.cpp.

669{
670 int physTot = m_fields[0]->GetTotPoints();
671
672 for (int i = 0; i < m_nConvectiveFields; ++i)
673 {
674 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
675 Vmath::Vcopy(physTot, m_coordsVel[i], 1, outarray[i], 1);
676 }
677}

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 879 of file Mapping.cpp.

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

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 97 of file Mapping.cpp.

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

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 708 of file Mapping.cpp.

710{
711 int physTot = m_fields[0]->GetTotPoints();
712 int nvel = m_nConvectiveFields;
713
714 Array<OneD, Array<OneD, NekDouble>> g(nvel * nvel);
715
717
718 for (int i = 0; i < nvel; i++)
719 {
720 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
721 for (int j = 0; j < nvel; j++)
722 {
723 Vmath::Vvtvp(physTot, g[i * nvel + j], 1, inarray[j], 1,
724 outarray[i], 1, outarray[i], 1);
725 }
726 }
727}
GLOBAL_MAPPING_EXPORT void GetMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)
Get the metric tensor .
Definition: Mapping.h:175

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 729 of file Mapping.cpp.

731{
732 int physTot = m_fields[0]->GetTotPoints();
733 int nvel = m_nConvectiveFields;
734
735 Array<OneD, Array<OneD, NekDouble>> g(nvel * nvel);
736
738
739 for (int i = 0; i < nvel; i++)
740 {
741 outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
742 for (int j = 0; j < nvel; j++)
743 {
744 Vmath::Vvtvp(physTot, g[i * nvel + j], 1, inarray[j], 1,
745 outarray[i], 1, outarray[i], 1);
746 }
747 }
748}
GLOBAL_MAPPING_EXPORT void GetInvMetricTensor(Array< OneD, Array< OneD, NekDouble > > &outarray)
Get the inverse of metric tensor .
Definition: Mapping.h:182

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 1049 of file Mapping.cpp.

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

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 1207 of file Mapping.cpp.

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

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 778 of file Mapping.cpp.

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

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 287 of file Mapping.h.

291 {
292 v_VelocityLaplacian(inarray, outarray, alpha);
293 }
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:778

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 408 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 404 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 402 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 425 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 415 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 410 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 429 of file Mapping.h.

Referenced by Load().

◆ m_isDefined

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

Definition at line 430 of file Mapping.h.

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

◆ m_mappingPtr

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

Definition at line 428 of file Mapping.h.

Referenced by Load().

◆ m_nConvectiveFields

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

Number of velocity components.

Definition at line 412 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 423 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 417 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 433 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 434 of file Mapping.h.

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