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.

72 {}

◆ Mapping()

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

Constructor.

Definition at line 58 of file Mapping.cpp.

60  : m_session(pSession), m_fields(pFields)
61 {
62  switch (m_fields[0]->GetExpType())
63  {
64  case MultiRegions::e1D:
65  {
67  }
68  break;
69 
70  case MultiRegions::e2D:
71  {
73  }
74  break;
75 
76  case MultiRegions::e3D:
79  {
81  }
82  break;
83 
84  default:
85  ASSERTL0(0,"Dimension not supported");
86  break;
87  }
88 
90 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:416
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Mapping.h:404
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:408
LibUtilities::FieldIOSharedPtr m_fld
Definition: Mapping.h:406
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 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 549 of file Mapping.cpp.

552 {
553  if(inarray == outarray)
554  {
555  int physTot = m_fields[0]->GetTotPoints();
556  int nvel = m_nConvectiveFields;
557 
558  for (int i=0; i< nvel; i++)
559  {
560  Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
561  }
562  v_ContravarFromCartesian( m_tmp, outarray);
563  }
564  else
565  {
566  v_ContravarFromCartesian( inarray, outarray);
567  }
568 }
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:440
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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

492 {
493  if(inarray == outarray)
494  {
495  int physTot = m_fields[0]->GetTotPoints();
496  int nvel = m_nConvectiveFields;
497 
498  for (int i=0; i< nvel; i++)
499  {
500  Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
501  }
502  v_ContravarToCartesian( m_tmp, outarray);
503  }
504  else
505  {
506  v_ContravarToCartesian( inarray, outarray);
507  }
508 }
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 579 of file Mapping.cpp.

582 {
583  if(inarray == outarray)
584  {
585  int physTot = m_fields[0]->GetTotPoints();
586  int nvel = m_nConvectiveFields;
587 
588  for (int i=0; i< nvel; i++)
589  {
590  Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
591  }
592  v_CovarFromCartesian( m_tmp, outarray);
593  }
594  else
595  {
596  v_CovarFromCartesian( inarray, outarray);
597  }
598 }
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 519 of file Mapping.cpp.

522 {
523  if(inarray == outarray)
524  {
525  int physTot = m_fields[0]->GetTotPoints();
526  int nvel = m_nConvectiveFields;
527 
528  for (int i=0; i< nvel; i++)
529  {
530  Vmath::Vcopy(physTot, inarray[i], 1, m_tmp[i], 1);
531  }
532  v_CovarToCartesian( m_tmp, outarray);
533  }
534  else
535  {
536  v_CovarToCartesian( inarray, outarray);
537  }
538 }
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.

328  {
329  v_CurlCurlField( inarray, outarray, generalized);
330  }
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:1012

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:766

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:690

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.

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

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

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

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.

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

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:678

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 154 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:902

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

365  {
366  return m_constantJacobian;
367  }
bool m_constantJacobian
Flag defining if the Jacobian is constant.
Definition: Mapping.h:426

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:100

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

371  {
372  return m_isDefined;
373  }

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

353  {
354  return m_fromFunction;
355  }
bool m_fromFunction
Flag defining if the Mapping is defined by a function.
Definition: Mapping.h:430

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

341  {
342  return m_timeDependent;
343  }
bool m_timeDependent
Flag defining if the Mapping is time-dependent.
Definition: Mapping.h:428

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

271 {
272  if (!m_init)
273  {
274  TiXmlElement* vMapping = NULL;
275  string vType;
276  if (pSession->DefinesElement("Nektar/Mapping"))
277  {
278  vMapping = pSession->GetElement("Nektar/Mapping");
279  vType = vMapping->Attribute("TYPE");
280  m_isDefined = true;
281  }
282  else
283  {
284  vType = "Translation";
285  }
286 
288  vType, pSession, pFields,
289  vMapping);
290 
291  m_init = true;
292  }
293 
294  return m_mappingPtr;
295 }
static MappingSharedPtr m_mappingPtr
Definition: Mapping.h:433
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
MappingFactory & GetMappingFactory()
Declaration of the mapping factory singleton.
Definition: Mapping.cpp:52

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::GlobalMapping::GetMappingFactory(), 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 610 of file Mapping.cpp.

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

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

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

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

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

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

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

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

359  {
360  m_fromFunction = value;
361  }

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

347  {
348  m_timeDependent = value;
349  }

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

381  {
382  v_UpdateBCs(time);
383  }
virtual GLOBAL_MAPPING_EXPORT void v_UpdateBCs(const NekDouble time)
Definition: Mapping.cpp:1075

References v_UpdateBCs().

◆ UpdateGeomInfo()

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

Recompute the metric terms of the Mapping.

Definition at line 396 of file Mapping.h.

397  {
399  }
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 386 of file Mapping.h.

391  {
392  v_UpdateMapping( time, coords, coordsVel);
393  }
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:1230

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

1016 {
1017  int physTot = m_fields[0]->GetTotPoints();
1018  int nvel = m_nConvectiveFields;
1019 
1020  // Set wavespace to false and store current value
1021  bool wavespace = m_fields[0]->GetWaveSpace();
1022  m_fields[0]->SetWaveSpace(false);
1023 
1024  // For implicit treatment of viscous terms, we want the generalized
1025  // curlcurl and for explicit treatment, we want the cartesian one.
1026  if (generalized)
1027  {
1028  // Get the second derivatives u^i_{,jk}
1029  Array<OneD, Array<OneD, NekDouble> > tmp(nvel);
1030  Array<OneD, Array<OneD, NekDouble> > ddU(nvel*nvel*nvel);
1031  gradgradU(inarray, ddU);
1032 
1033  // Raise index to obtain A^{ip}_{k} = g^pj u^i_{,jk}
1034  for (int i = 0; i < nvel; ++i)
1035  {
1036  for (int k = 0; k < nvel; ++k)
1037  {
1038  // Copy to wk
1039  for (int j = 0; j < nvel; ++j)
1040  {
1041  tmp[j] = ddU[i*nvel*nvel+j*nvel+k];
1042  }
1043  RaiseIndex(tmp, m_tmp);
1044  for (int p=0; p<nvel; ++p)
1045  {
1046  Vmath::Vcopy(physTot, m_tmp[p], 1,
1047  ddU[i*nvel*nvel+p*nvel+k], 1);
1048  }
1049  }
1050  }
1051  // The curlcurl is g^ji u^k_{kj} - g^jk u^i_kj = A^{ki}_k - A^{ik}_k
1052  for (int i = 0; i < nvel; ++i)
1053  {
1054  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
1055  for (int k = 0; k < nvel; ++k)
1056  {
1057  Vmath::Vadd(physTot, outarray[i], 1,
1058  ddU[k*nvel*nvel+i*nvel+k], 1,
1059  outarray[i], 1);
1060  Vmath::Vsub(physTot, outarray[i], 1,
1061  ddU[i*nvel*nvel+k*nvel+k], 1,
1062  outarray[i], 1);
1063  }
1064  }
1065  }
1066  else
1067  {
1068  m_fields[0]->CurlCurl(inarray, outarray);
1069  }
1070 
1071  // Restore value of wavespace
1072  m_fields[0]->SetWaveSpace(wavespace);
1073 }
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:641
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:322
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:372

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

769 {
770  int physTot = m_fields[0]->GetTotPoints();
771  Array<OneD, NekDouble> wk(physTot, 0.0);
772 
773  Vmath::Zero(physTot, outarray, 1);
774 
775  // Set wavespace to false and store current value
776  bool wavespace = m_fields[0]->GetWaveSpace();
777  m_fields[0]->SetWaveSpace(false);
778 
779  // Get Mapping Jacobian
780  Array<OneD, NekDouble> Jac(physTot, 0.0);
781  GetJacobian(Jac);
782 
783  for(int i = 0; i < m_nConvectiveFields; ++i)
784  {
785  Vmath::Vmul(physTot,Jac, 1, inarray[i], 1, wk, 1); // J*Ui
786  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i],
787  wk, wk); // (J*Ui)_i
788  Vmath::Vadd(physTot, wk, 1, outarray, 1, outarray, 1);
789  }
790  Vmath::Vdiv(physTot,outarray,1,Jac,1,outarray,1); //1/J*(J*Ui)_i
791 
792  m_fields[0]->SetWaveSpace(wavespace);
793 }
GLOBAL_MAPPING_EXPORT void GetJacobian(Array< OneD, NekDouble > &outarray)
Get the Jacobian of the transformation.
Definition: Mapping.h:154
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:192
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:257

References Nektar::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::MappingXYofZ, Nektar::GlobalMapping::MappingXofZ, Nektar::GlobalMapping::MappingXofXZ, and Nektar::GlobalMapping::MappingTranslation.

Definition at line 690 of file Mapping.cpp.

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

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

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

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

680 {
681  int physTot = m_fields[0]->GetTotPoints();
682 
683  for(int i = 0; i < m_nConvectiveFields; ++i)
684  {
685  outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
686  Vmath::Vcopy(physTot, m_coordsVel[i], 1, outarray[i], 1);
687  }
688 }

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

905 {
906  int physTot = m_fields[0]->GetTotPoints();
907  int nvel = m_nConvectiveFields;
908 
909  // Declare variables
910  outarray = Array<OneD, Array<OneD, NekDouble> > (nvel*nvel*nvel);
911  Array<OneD, Array<OneD, NekDouble> > tmp(nvel);
912  for (int i=0; i< nvel*nvel*nvel; i++)
913  {
914  outarray[i] = Array<OneD, NekDouble>(physTot,0.0);
915  }
916 
917  // Set wavespace to false and store current value
918  bool wavespace = m_fields[0]->GetWaveSpace();
919  m_fields[0]->SetWaveSpace(false);
920 
921  // Calculate vector gradient u^i_(,j) = du^i/dx^j + {i,pj}*u^p
923  for (int i=0; i< nvel; i++)
924  {
925  if (nvel == 2)
926  {
927  m_fields[0]->PhysDeriv(inarray[i],
928  m_wk2[i*nvel+0],
929  m_wk2[i*nvel+1]);
930  }
931  else
932  {
933  m_fields[0]->PhysDeriv(inarray[i],
934  m_wk2[i*nvel+0],
935  m_wk2[i*nvel+1],
936  m_wk2[i*nvel+2]);
937  }
938  for (int j=0; j< nvel; j++)
939  {
940  Vmath::Vadd(physTot,m_wk1[i*nvel+j],1,m_wk2[i*nvel+j],1,
941  m_wk1[i*nvel+j], 1);
942  }
943  }
944 
945  //
946  // Calculate (u^i_,j),k
947  //
948 
949  // Step 1 : d(u^i_,j))/d(x^k)
950  for (int i=0; i< nvel; i++)
951  {
952  for (int j=0; j< nvel; j++)
953  {
954  if (nvel == 2)
955  {
956  m_fields[0]->PhysDeriv(m_wk1[i*nvel+j],
957  outarray[i*nvel*nvel+j*nvel+0],
958  outarray[i*nvel*nvel+j*nvel+1]);
959  }
960  else
961  {
962  m_fields[0]->PhysDeriv(m_wk1[i*nvel+j],
963  outarray[i*nvel*nvel+j*nvel+0],
964  outarray[i*nvel*nvel+j*nvel+1],
965  outarray[i*nvel*nvel+j*nvel+2]);
966  }
967  }
968  }
969 
970  // Step 2: d(u^i_,j)/d(x^k) - {p,jk}*u^i_,p
971  for (int i=0; i< nvel; i++)
972  {
973  for (int p=0; p< nvel; p++)
974  {
975  tmp[p] = m_wk1[i*nvel+p];
976  }
978  for (int j=0; j< nvel; j++)
979  {
980  for (int k=0; k< nvel; k++)
981  {
982  Vmath::Vsub(physTot,outarray[i*nvel*nvel+j*nvel+k],1,
983  m_wk2[j*nvel+k],1,
984  outarray[i*nvel*nvel+j*nvel+k], 1);
985  }
986  }
987  }
988 
989  // Step 3: d(u^i_,j)/d(x^k) - {p,jk}*u^i_,p + {i,pk} u^p_,j
990  for (int j=0; j< nvel; j++)
991  {
992  for (int p=0; p< nvel; p++)
993  {
994  tmp[p] = m_wk1[p*nvel+j];
995  }
997  for (int i=0; i< nvel; i++)
998  {
999  for (int k=0; k< nvel; k++)
1000  {
1001  Vmath::Vadd(physTot,outarray[i*nvel*nvel+j*nvel+k],1,
1002  m_wk2[i*nvel+k],1,
1003  outarray[i*nvel*nvel+j*nvel+k], 1);
1004  }
1005  }
1006  }
1007 
1008  // Restore value of wavespace
1009  m_fields[0]->SetWaveSpace(wavespace);
1010 }
Array< OneD, Array< OneD, NekDouble > > m_wk2
Definition: Mapping.h:439
Array< OneD, Array< OneD, NekDouble > > m_wk1
Definition: Mapping.h:438
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::MappingXYofZ, Nektar::GlobalMapping::MappingXYofXY, Nektar::GlobalMapping::MappingXofZ, Nektar::GlobalMapping::MappingXofXZ, Nektar::GlobalMapping::MappingTranslation, and Nektar::GlobalMapping::MappingGeneral.

Definition at line 100 of file Mapping.cpp.

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

Definition at line 720 of file Mapping.cpp.

723 {
724  int physTot = m_fields[0]->GetTotPoints();
725  int nvel = m_nConvectiveFields;
726 
727  Array<OneD, Array<OneD, NekDouble> > g(nvel*nvel);
728 
729  GetMetricTensor(g);
730 
731  for (int i=0; i< nvel; i++)
732  {
733  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
734  for (int j=0; j< nvel; j++)
735  {
736  Vmath::Vvtvp(physTot, g[i*nvel+j], 1, inarray[j], 1,
737  outarray[i], 1,
738  outarray[i], 1);
739  }
740  }
741 }
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::MappingXofZ, Nektar::GlobalMapping::MappingXofXZ, and Nektar::GlobalMapping::MappingTranslation.

Definition at line 743 of file Mapping.cpp.

746 {
747  int physTot = m_fields[0]->GetTotPoints();
748  int nvel = m_nConvectiveFields;
749 
750  Array<OneD, Array<OneD, NekDouble> > g(nvel*nvel);
751 
753 
754  for (int i=0; i< nvel; i++)
755  {
756  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
757  for (int j=0; j< nvel; j++)
758  {
759  Vmath::Vvtvp(physTot, g[i*nvel+j], 1, inarray[j], 1,
760  outarray[i], 1,
761  outarray[i], 1);
762  }
763  }
764 }
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 1075 of file Mapping.cpp.

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

1234 {
1235  if (m_fromFunction)
1236  {
1237  std::string s_FieldStr;
1238  string fieldNames[3] = {"x", "y", "z"};
1239  string velFieldNames[3] = {"vx", "vy", "vz"};
1240  // Check if function from session file defines each component
1241  // and evaluate them, otherwise there is no need to update
1242  // coords
1243  for(int i = 0; i < 3; i++)
1244  {
1245  s_FieldStr = fieldNames[i];
1246  if ( m_session->DefinesFunction(m_funcName, s_FieldStr))
1247  {
1248  EvaluateFunction(m_fields, m_session, s_FieldStr, m_coords[i],
1249  m_funcName, time);
1250  if ( i==2 && m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
1251  {
1252  ASSERTL0 (false,
1253  "3DH1D does not support mapping in the z-direction.");
1254  }
1255  }
1256  s_FieldStr = velFieldNames[i];
1257  if ( m_session->DefinesFunction(m_velFuncName, s_FieldStr))
1258  {
1259  EvaluateFunction(m_fields, m_session, s_FieldStr,
1260  m_coordsVel[i], m_velFuncName, time);
1261  if ( i==2 && m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
1262  {
1263  ASSERTL0 (false,
1264  "3DH1D does not support mapping in the z-direction.");
1265  }
1266  }
1267  }
1268  }
1269  else
1270  {
1271  int physTot = m_fields[0]->GetTotPoints();
1272  int nvel = m_nConvectiveFields;
1273  // Copy coordinates
1274  for(int i = 0; i < 3; i++)
1275  {
1276  Vmath::Vcopy(physTot, coords[i], 1, m_coords[i], 1);
1277  }
1278 
1279  for(int i = 0; i < nvel; i++)
1280  {
1281  Vmath::Vcopy(physTot, coordsVel[i], 1, m_coordsVel[i], 1);
1282  }
1283  }
1284 
1285  // Update the information required by the specific mapping
1286  UpdateGeomInfo();
1287 }

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

799 {
800  int physTot = m_fields[0]->GetTotPoints();
801  int nvel = m_nConvectiveFields;
802 
803  Array<OneD, Array<OneD, NekDouble> > tmp(nvel);
804 
805  // Set wavespace to false and store current value
806  bool wavespace = m_fields[0]->GetWaveSpace();
807  m_fields[0]->SetWaveSpace(false);
808 
809  // Calculate vector gradient wk2 = u^i_(,k) = du^i/dx^k + {i,jk}*u^j
811  for (int i=0; i< nvel; i++)
812  {
813  if(nvel == 2)
814  {
815  m_fields[0]->PhysDeriv(inarray[i],
816  m_wk2[i*nvel+0],
817  m_wk2[i*nvel+1]);
818  }
819  else
820  {
821  m_fields[0]->PhysDeriv(inarray[i],
822  m_wk2[i*nvel+0],
823  m_wk2[i*nvel+1],
824  m_wk2[i*nvel+2]);
825  }
826  for (int k=0; k< nvel; k++)
827  {
828  Vmath::Vadd(physTot,m_wk1[i*nvel+k],1,m_wk2[i*nvel+k],1,
829  m_wk1[i*nvel+k], 1);
830  }
831  }
832  // Calculate wk1 = A^(ij) = g^(jk)*u^i_(,k)
833  for (int i=0; i< nvel; i++)
834  {
835  for (int k=0; k< nvel; k++)
836  {
837  tmp[k] = m_wk1[i*nvel+k];
838  }
839  RaiseIndex(tmp, m_tmp);
840  for (int j=0; j<nvel; j++)
841  {
842  Vmath::Vcopy(physTot, m_tmp[j], 1, m_wk1[i*nvel+j], 1);
843  }
844  }
845  //
846  // Calculate L(U)^i = (A^(ij))_(,j) - alpha*d^2(u^i)/dx^jdx^j
847  //
848 
849  // Step 1 :
850  // d(A^(ij) - alpha*du^i/dx^j)/d(x^j)
851  for (int i=0; i< nvel; i++)
852  {
853  outarray[i] = Array<OneD, NekDouble>(physTot,0.0);
854  for (int j=0; j< nvel; j++)
855  {
856  Vmath::Smul(physTot, alpha, m_wk2[i*nvel+j], 1,
857  m_tmp[0], 1);
858  Vmath::Vsub(physTot, m_wk1[i*nvel+j], 1, m_tmp[0], 1,
859  m_tmp[0], 1);
860 
861  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j],
862  m_tmp[0],
863  m_tmp[1]);
864  Vmath::Vadd(physTot,outarray[i],1,m_tmp[1],1,outarray[i], 1);
865  }
866  }
867 
868  // Step 2: d(A^(ij))/d(x^j) + {j,pj}*A^(ip)
869  for (int i=0; i< nvel; i++)
870  {
871  for (int p=0; p< nvel; p++)
872  {
873  tmp[p] = m_wk1[i*nvel+p];
874  }
876  for (int j=0; j< nvel; j++)
877  {
878  Vmath::Vadd(physTot,outarray[i],1,m_wk2[j*nvel+j],1,
879  outarray[i], 1);
880  }
881  }
882 
883  // Step 3: d(A^(ij))/d(x^j) + {j,pj}*A^(ip) + {i,pj} A^(pj)
884  for (int j=0; j< nvel; j++)
885  {
886  for (int p=0; p< nvel; p++)
887  {
888  tmp[p] = m_wk1[p*nvel+j];
889  }
891  for (int i=0; i< nvel; i++)
892  {
893  Vmath::Vadd(physTot,outarray[i], 1, m_wk2[i*nvel+j], 1,
894  outarray[i], 1);
895  }
896  }
897 
898  // Restore value of wavespace
899  m_fields[0]->SetWaveSpace(wavespace);
900 }
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:225

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:795

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 412 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 408 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 406 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 430 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 419 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 414 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 434 of file Mapping.h.

Referenced by Load().

◆ m_isDefined

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

Definition at line 435 of file Mapping.h.

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

◆ m_mappingPtr

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

Definition at line 433 of file Mapping.h.

Referenced by Load().

◆ m_nConvectiveFields

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

Number of velocity components.

Definition at line 416 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 428 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 421 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 438 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 439 of file Mapping.h.

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