Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Friends | List of all members
Nektar::GlobalMapping::MappingXofXZ Class Reference

#include <MappingXofXZ.h>

Inheritance diagram for Nektar::GlobalMapping::MappingXofXZ:
Inheritance graph
[legend]
Collaboration diagram for Nektar::GlobalMapping::MappingXofXZ:
Collaboration graph
[legend]

Static Public Member Functions

static GLOBAL_MAPPING_EXPORT
MappingSharedPtr 
create (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::GlobalMapping::Mapping
static GLOBAL_MAPPING_EXPORT
MappingSharedPtr 
Load (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Return a pointer to the mapping, creating it on first call. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

 MappingXofXZ (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 
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)
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetJacobian (Array< OneD, NekDouble > &outarray)
 
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)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)
 
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)
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateGeomInfo ()
 
- Protected Member Functions inherited from Nektar::GlobalMapping::Mapping
GLOBAL_MAPPING_EXPORT Mapping (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Constructor. More...
 
GLOBAL_MAPPING_EXPORT void EvaluateFunction (Array< OneD, MultiRegions::ExpListSharedPtr > pFields, LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
 
GLOBAL_MAPPING_EXPORT void EvaluateTimeFunction (LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
 
virtual GLOBAL_MAPPING_EXPORT void v_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_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 (const 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_UpdateBCs (const NekDouble time)
 

Friends

class MemoryManager< MappingXofXZ >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::GlobalMapping::Mapping
virtual GLOBAL_MAPPING_EXPORT ~Mapping ()
 Destructor. More...
 
GLOBAL_MAPPING_EXPORT void InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
 Initialise the mapping object. More...
 
GLOBAL_MAPPING_EXPORT void ReplaceField (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Replace the Expansion List used by the mapping. More...
 
GLOBAL_MAPPING_EXPORT void Output (LibUtilities::FieldMetaDataMap &fieldMetaDataMap, const std::string &outname)
 Output function called when a chk or fld file is written. More...
 
GLOBAL_MAPPING_EXPORT void ContravarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Convert a contravariant vector to the Cartesian system. More...
 
GLOBAL_MAPPING_EXPORT void CovarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Convert a covariant vector to the Cartesian system. More...
 
GLOBAL_MAPPING_EXPORT void ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Convert a contravariant vector to the transformed system. More...
 
GLOBAL_MAPPING_EXPORT void CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Convert a covariant vector to the transformed system. More...
 
GLOBAL_MAPPING_EXPORT void GetCartesianCoordinates (Array< OneD, NekDouble > &out0, Array< OneD, NekDouble > &out1, Array< OneD, NekDouble > &out2)
 Get the Cartesian coordinates in the field. More...
 
GLOBAL_MAPPING_EXPORT void GetJacobian (Array< OneD, NekDouble > &outarray)
 Get the Jacobian of the transformation. More...
 
GLOBAL_MAPPING_EXPORT void DotGradJacobian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the dot product with the gradient of the Jacobian. More...
 
GLOBAL_MAPPING_EXPORT void GetMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)
 Get the metric tensor $g_{ij}$. More...
 
GLOBAL_MAPPING_EXPORT void GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)
 Get the inverse of metric tensor $g^{ij}$. More...
 
GLOBAL_MAPPING_EXPORT void LowerIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Lower index of vector: $v_{i} = g_{ij}*v^{j}$. More...
 
GLOBAL_MAPPING_EXPORT void RaiseIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Raise index of vector: $v^{i} = g^{ij}*v_{j}$. More...
 
GLOBAL_MAPPING_EXPORT void ApplyChristoffelContravar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Apply the Christoffel symbols to a contravariant vector. More...
 
GLOBAL_MAPPING_EXPORT void ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Apply the Christoffel symbols to a covariant vector. More...
 
GLOBAL_MAPPING_EXPORT void GetCoordVelocity (Array< OneD, Array< OneD, NekDouble > > &outarray)
 Obtain the velocity of the coordinates. More...
 
GLOBAL_MAPPING_EXPORT void Divergence (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the generalised divergence operator. More...
 
GLOBAL_MAPPING_EXPORT void VelocityLaplacian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble alpha=0.0)
 Generalised (correction to the) velocity Laplacian operator. More...
 
GLOBAL_MAPPING_EXPORT void gradgradU (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Second order covariant derivatives of a contravariant vector. More...
 
GLOBAL_MAPPING_EXPORT void CurlCurlField (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const bool generalized)
 CurlCurl calculated on the whole field. More...
 
GLOBAL_MAPPING_EXPORT bool IsTimeDependent ()
 Get flag defining if mapping is time-dependent. More...
 
GLOBAL_MAPPING_EXPORT void SetTimeDependent (const bool value)
 Set flag defining if mapping is time-dependent. More...
 
GLOBAL_MAPPING_EXPORT bool IsFromFunction ()
 Get flag defining if mapping is defined by a function. More...
 
GLOBAL_MAPPING_EXPORT void SetFromFunction (const bool value)
 Set flag defining if mapping is defined by a function. More...
 
GLOBAL_MAPPING_EXPORT bool HasConstantJacobian ()
 Get flag defining if mapping has constant Jacobian. More...
 
GLOBAL_MAPPING_EXPORT bool IsDefined ()
 Get flag determining if the mapping was defined or is trivial. More...
 
GLOBAL_MAPPING_EXPORT void UpdateBCs (const NekDouble time)
 Update the Dirichlet Boundary Conditions when using Mappings. More...
 
GLOBAL_MAPPING_EXPORT void UpdateMapping (const NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &coords=NullNekDoubleArrayofArray, const Array< OneD, Array< OneD, NekDouble > > &coordsVel=NullNekDoubleArrayofArray)
 Update the Mapping with new coordinates. More...
 
GLOBAL_MAPPING_EXPORT void UpdateGeomInfo ()
 Recompute the metric terms of the Mapping. More...
 
- Protected Attributes inherited from Nektar::GlobalMapping::Mapping
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_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...
 
string m_funcName
 Name of the function containing the coordinates. More...
 
string m_velFuncName
 Name of the function containing the velocity of the coordinates. More...
 
bool m_constantJacobian
 Flag defining if the Jacobian is constant. More...
 
bool m_timeDependent
 Flag defining if the Mapping is time-dependent. More...
 
bool m_fromFunction
 Flag defining if the Mapping is defined by a function. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_wk1
 
Array< OneD, Array< OneD,
NekDouble > > 
m_wk2
 
Array< OneD, Array< OneD,
NekDouble > > 
m_tmp
 
- Static Protected Attributes inherited from Nektar::GlobalMapping::Mapping
static MappingSharedPtr m_mappingPtr = MappingSharedPtr()
 
static bool m_init = false
 
static bool m_isDefined = false
 

Detailed Description

This class implements a mapping defined by a transformation of the type

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

\[ \bar{y} = y \]

\[ \bar{z} = z \]

where $(\bar{x},\bar{y},\bar{z})$ are the Cartesian (physical) coordinates and $(x,y,z)$ are the transformed (computational) coordinates.

Definition at line 52 of file MappingXofXZ.h.

Constructor & Destructor Documentation

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

Definition at line 58 of file MappingXofXZ.cpp.

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

Member Function Documentation

static GLOBAL_MAPPING_EXPORT MappingSharedPtr Nektar::GlobalMapping::MappingXofXZ::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const TiXmlElement *  pMapping 
)
inlinestatic

Creates an instance of this class.

Definition at line 60 of file MappingXofXZ.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::GlobalMapping::MappingSharedPtr.

64  {
67  pFields);
68  p->InitObject(pFields, pMapping);
69  return p;
70  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
GLOBAL_MAPPING_EXPORT typedef boost::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:51
void Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelContravar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 314 of file MappingXofXZ.cpp.

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

317 {
318  int physTot = m_fields[0]->GetTotPoints();
319  int nvel = m_nConvectiveFields;
320  Array<OneD, NekDouble> wk(physTot, 0.0);
321 
322  for (int i = 0; i< nvel; i++)
323  {
324  for (int j = 0; j< nvel; j++)
325  {
326  outarray[i*nvel+j] = Array<OneD, NekDouble>(physTot,0.0);
327  }
328  }
329 
330  // Calculate non-zero terms
331 
332  // outarray(0,0) = U1 * fxx/fx + U3 * fxz/fx
333  Vmath::Vdiv(physTot,m_GeometricInfo[2],1,m_GeometricInfo[0],1,wk,1);
334  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+0],1);
335  Vmath::Vdiv(physTot,m_GeometricInfo[3],1,m_GeometricInfo[0],1,wk,1);
336  Vmath::Vvtvp(physTot,wk,1,inarray[2],1,outarray[0*nvel+0],1,
337  outarray[0*nvel+0],1);
338 
339  // outarray(0,2) = U1 * fxz/fx + U3 * fzz/fx
340  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+2],1);
341  Vmath::Vdiv(physTot,m_GeometricInfo[4],1,m_GeometricInfo[0],1,wk,1);
342  Vmath::Vvtvp(physTot,wk,1,inarray[2],1,outarray[0*nvel+2],1,
343  outarray[0*nvel+2],1);
344 
345 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:428
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:227
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_ApplyChristoffelCovar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 347 of file MappingXofXZ.cpp.

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

350 {
351  int physTot = m_fields[0]->GetTotPoints();
352  int nvel = m_nConvectiveFields;
353  Array<OneD, NekDouble> wk(physTot, 0.0);
354 
355  for (int i = 0; i< nvel; i++)
356  {
357  for (int j = 0; j< nvel; j++)
358  {
359  outarray[i*nvel+j] = Array<OneD, NekDouble>(physTot,0.0);
360  }
361  }
362 
363  // Calculate non-zero terms
364 
365  // outarray(0,0) = U1 * fxx/fx
366  Vmath::Vdiv(physTot,m_GeometricInfo[2],1,m_GeometricInfo[0],1,wk,1);
367  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+0],1);
368 
369  //outarray(0,2) = outarray(2,0) = U1 * fxz/fx
370  Vmath::Vdiv(physTot,m_GeometricInfo[3],1,m_GeometricInfo[0],1,wk,1);
371  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[0*nvel+2],1);
372  Vmath::Vcopy(physTot,outarray[0*nvel+2],1,outarray[2*nvel+0],1);
373 
374  // outarray(2,2) = U1 * fzz/fx
375  Vmath::Vdiv(physTot,m_GeometricInfo[4],1,m_GeometricInfo[0],1,wk,1);
376  Vmath::Vmul(physTot,wk,1,inarray[0],1,outarray[2*nvel+2],1);
377 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:227
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_ContravarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 121 of file MappingXofXZ.cpp.

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

124 {
125  int physTot = m_fields[0]->GetTotPoints();
126  Array<OneD, NekDouble> wk(physTot, 0.0);
127 
128  // U1 = u1/fx - fz/fx * u3
129  Vmath::Vdiv(physTot, inarray[0], 1,
130  m_GeometricInfo[0], 1, outarray[0], 1);
131  Vmath::Vdiv(physTot, m_GeometricInfo[1], 1,
132  m_GeometricInfo[0], 1, wk, 1);
133  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, wk, 1);
134  Vmath::Vsub(physTot, outarray[0], 1, wk, 1, outarray[0], 1);
135 
136  // U2 = u2
137  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
138 
139  // U3 = u3
140  Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
141 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:227
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
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:329
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_ContravarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 81 of file MappingXofXZ.cpp.

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

84 {
85  int physTot = m_fields[0]->GetTotPoints();
86 
87  // U1 = fx*u1 + fz*u3
88  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, inarray[0], 1,
89  outarray[0], 1);
90  Vmath::Vvtvp(physTot, m_GeometricInfo[1], 1, inarray[2], 1,
91  outarray[0], 1, outarray[0],1);
92 
93  // U2 = u2
94  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
95 
96  // U3 = u3
97  Vmath::Vcopy(physTot, inarray[2], 1, outarray[2], 1);
98 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:428
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_CovarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 143 of file MappingXofXZ.cpp.

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

146 {
147  int physTot = m_fields[0]->GetTotPoints();
148 
149  // U1 = u1*fx
150  Vmath::Vmul(physTot, inarray[0], 1, m_GeometricInfo[0], 1,
151  outarray[0], 1);
152 
153  // U2 = u2
154  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
155 
156  // U3 = u3 + fz*u1
157  Vmath::Vmul(physTot, m_GeometricInfo[1], 1,
158  inarray[0], 1, outarray[2], 1);
159  Vmath::Vadd(physTot, inarray[2], 1, outarray[2], 1, outarray[2], 1);
160 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_CovarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 100 of file MappingXofXZ.cpp.

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

103 {
104  int physTot = m_fields[0]->GetTotPoints();
105  Array<OneD, NekDouble> wk(physTot, 0.0);
106 
107  // U1 = u1/fx
108  Vmath::Vdiv(physTot, inarray[0], 1, m_GeometricInfo[0], 1,
109  outarray[0], 1);
110 
111  // U2 = u2
112  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
113 
114  // U3 = u3 - fz/fx*u1
115  Vmath::Vdiv(physTot, m_GeometricInfo[1], 1,
116  m_GeometricInfo[0], 1, wk, 1);
117  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, wk, 1);
118  Vmath::Vsub(physTot, inarray[2], 1, wk, 1, outarray[2], 1);
119 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:227
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
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:329
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_DotGradJacobian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 169 of file MappingXofXZ.cpp.

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

172 {
173  int physTot = m_fields[0]->GetTotPoints();
174 
175  Vmath::Vmul(physTot, m_GeometricInfo[2], 1, inarray[0], 1, outarray, 1);
176  Vmath::Vvtvp(physTot, m_GeometricInfo[3], 1, inarray[2], 1,
177  outarray, 1, outarray,1);
178 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:428
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_GetInvMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 215 of file MappingXofXZ.cpp.

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

217 {
218  int physTot = m_fields[0]->GetTotPoints();
219  int nvel = m_nConvectiveFields;
220  Array<OneD, NekDouble> wk(physTot, 0.0);
221 
222  for (int i=0; i<nvel*nvel; i++)
223  {
224  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
225  }
226  // Fill diagonal with 1.0
227  for (int i=0; i<nvel; i++)
228  {
229  Vmath::Sadd(physTot, 1.0, outarray[i+nvel*i], 1,
230  outarray[i+nvel*i], 1);
231  }
232 
233  // G^{13} and G^{31} = -fz/fx
234  Vmath::Vdiv(physTot,m_GeometricInfo[1],1,
235  m_GeometricInfo[0],1,wk,1); // fz/fx
236  Vmath::Neg(physTot, wk, 1);
237  Vmath::Vcopy(physTot, wk, 1, outarray[0*nvel+2], 1);
238  Vmath::Vcopy(physTot, wk, 1, outarray[2*nvel+0], 1);
239 
240  // G^{11} = (1+fz^2)/(fx^2)
241  Vmath::Vmul(physTot, m_GeometricInfo[1], 1,
242  m_GeometricInfo[1], 1, wk, 1); // fz^2
243  Vmath::Vadd(physTot, wk, 1, outarray[0*nvel+0], 1,
244  outarray[0*nvel+0], 1);
245 
246  Vmath::Vmul(physTot, m_GeometricInfo[0], 1,
247  m_GeometricInfo[0], 1, wk, 1); // fx^2
248  Vmath::Vdiv(physTot, outarray[0*nvel+0], 1, wk,1,
249  outarray[0*nvel+0], 1);
250 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:227
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_GetJacobian ( Array< OneD, NekDouble > &  outarray)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 162 of file MappingXofXZ.cpp.

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

164 {
165  int physTot = m_fields[0]->GetTotPoints();
166  Vmath::Vcopy(physTot, m_GeometricInfo[0], 1, outarray, 1);
167 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::GlobalMapping::MappingXofXZ::v_GetMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 180 of file MappingXofXZ.cpp.

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

182 {
183  int physTot = m_fields[0]->GetTotPoints();
184  int nvel = m_nConvectiveFields;
185  Array<OneD, NekDouble> wk(physTot, 0.0);
186 
187  for (int i=0; i<nvel*nvel; i++)
188  {
189  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
190  }
191  // Fill G^{22} and G^{33} with 1.0
192  for (int i=1; i<nvel; i++)
193  {
194  Vmath::Sadd(physTot, 1.0, outarray[i+nvel*i], 1,
195  outarray[i+nvel*i], 1);
196  }
197 
198  // G_{13} and G_{31} = fz*fx
199  Vmath::Vmul(physTot,m_GeometricInfo[1],1,
200  m_GeometricInfo[0],1,wk,1); // fz*fx
201  Vmath::Vcopy(physTot, wk, 1, outarray[0*nvel+2], 1);
202  Vmath::Vcopy(physTot, wk, 1, outarray[2*nvel+0], 1);
203 
204  // G^{11} = (fx^2)
205  Vmath::Vmul(physTot, m_GeometricInfo[0], 1,
206  m_GeometricInfo[0], 1, outarray[0*nvel+0], 1);
207 
208  // G^{33} = (1+fz^2)
209  Vmath::Vmul(physTot, m_GeometricInfo[1], 1,
210  m_GeometricInfo[1], 1, wk, 1); // fz^2
211  Vmath::Vadd(physTot, wk, 1, outarray[2*nvel+2], 1,
212  outarray[2*nvel+2], 1);
213 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
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:169
void Nektar::GlobalMapping::MappingXofXZ::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 from Nektar::GlobalMapping::Mapping.

Definition at line 69 of file MappingXofXZ.cpp.

References ASSERTL0, Nektar::GlobalMapping::Mapping::m_constantJacobian, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, and Nektar::GlobalMapping::Mapping::v_InitObject().

72 {
73  Mapping::v_InitObject(pFields, pMapping);
74 
75  m_constantJacobian = false;
76 
78  "Mapping X = X(x,z) needs 3 velocity components.");
79 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
virtual GLOBAL_MAPPING_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
Definition: Mapping.cpp:98
bool m_constantJacobian
Flag defining if the Jacobian is constant.
Definition: Mapping.h:427
void Nektar::GlobalMapping::MappingXofXZ::v_LowerIndex ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 252 of file MappingXofXZ.cpp.

References Nektar::GlobalMapping::Mapping::m_fields, Nektar::GlobalMapping::Mapping::m_GeometricInfo, Vmath::Sadd(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vmul().

255 {
256  int physTot = m_fields[0]->GetTotPoints();
257  Array<OneD, NekDouble> wk(physTot, 0.0);
258 
259  // out[0] = in[0]*fx^2 + in[2] * fz*fx
260  Vmath::Vmul(physTot,m_GeometricInfo[1],1,m_GeometricInfo[0],1,
261  wk,1); // fz*fx
262  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, outarray[0], 1); //in[2]*fz*fx
263  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[2], 1); //in[0]*fz*fx
264 
265  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1,
266  wk, 1); //fx^2
267  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, wk, 1); //in[0]*fx^2
268 
269  Vmath::Vadd(physTot, outarray[0], 1, wk, 1, outarray[0], 1);
270 
271  // out[1] = in[1]
272  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
273 
274  // out[2] = fx*fz*in[0] + (1+fz^2)*in[2]
275  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1,
276  wk, 1); // fz^2
277  Vmath::Sadd(physTot, 1.0, wk, 1, wk, 1); // 1+fz^2
278  Vmath::Vmul(physTot, wk, 1, inarray[2],1, wk, 1); // (1+fz^2)*in[2]
279 
280  Vmath::Vadd(physTot, wk, 1, outarray[2],1, outarray[2], 1);
281 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_RaiseIndex ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 283 of file MappingXofXZ.cpp.

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

286 {
287  int physTot = m_fields[0]->GetTotPoints();
288  Array<OneD, NekDouble> wk(physTot, 0.0);
289  Array<OneD, NekDouble> wk_2(physTot, 0.0);
290 
291  // out[2] = in[2] - in[0] * fz/fx
292  Vmath::Vdiv(physTot,m_GeometricInfo[1],1,m_GeometricInfo[0],1,
293  wk,1);
294  Vmath::Vmul(physTot, wk, 1, inarray[0], 1, outarray[2], 1);
295  Vmath::Vsub(physTot, inarray[2], 1, outarray[2], 1, outarray[2], 1);
296 
297  // out[0] = in[0]*(1+fz^2)/(fx^2) - in[2] * fz/fx
298  Vmath::Vmul(physTot, wk, 1, inarray[2], 1, outarray[0], 1);
299  Vmath::Vmul(physTot, m_GeometricInfo[1], 1, m_GeometricInfo[1], 1,
300  wk, 1);
301  Vmath::Sadd(physTot, 1.0, wk, 1, wk, 1);
302  Vmath::Vmul(physTot, m_GeometricInfo[0], 1, m_GeometricInfo[0], 1,
303  wk_2, 1);
304  Vmath::Vdiv(physTot, wk, 1, wk_2,1, wk, 1);
305  Vmath::Vmul(physTot, wk, 1, inarray[0],1, wk, 1);
306  Vmath::Vsub(physTot, wk, 1, outarray[0], 1, outarray[0], 1);
307 
308  // out[1] = in[1]
309  Vmath::Vcopy(physTot, inarray[1], 1, outarray[1], 1);
310 
311 
312 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
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:227
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
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:329
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:169
void Nektar::GlobalMapping::MappingXofXZ::v_UpdateGeomInfo ( )
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 379 of file MappingXofXZ.cpp.

References Nektar::MultiRegions::DirCartesianMap, Nektar::GlobalMapping::Mapping::m_coords, Nektar::GlobalMapping::Mapping::m_fields, and Nektar::GlobalMapping::Mapping::m_GeometricInfo.

380 {
381  int phystot = m_fields[0]->GetTotPoints();
382  // Allocation of geometry memory
384  for (int i = 0; i < m_GeometricInfo.num_elements(); i++)
385  {
386  m_GeometricInfo[i] = Array<OneD, NekDouble>(phystot, 0.0);
387  }
388 
389  bool waveSpace = m_fields[0]->GetWaveSpace();
390  m_fields[0]->SetWaveSpace(false);
391 
392  // Calculate derivatives of transformation
393  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
394  m_coords[0], m_GeometricInfo[0]); //f_x
395  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
396  m_coords[0], m_GeometricInfo[1]); //f_z
397 
398  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
399  m_GeometricInfo[0], m_GeometricInfo[2]); //f_xx
400  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
401  m_GeometricInfo[0], m_GeometricInfo[3]); //f_xz
402  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
403  m_GeometricInfo[1], m_GeometricInfo[4]); //f_zz
404 
405  m_fields[0]->SetWaveSpace(waveSpace);
406 }
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
Array with metric terms of the mapping.
Definition: Mapping.h:415
Array< OneD, Array< OneD, NekDouble > > m_coords
Array with the Cartesian coordinates.
Definition: Mapping.h:411
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86

Friends And Related Function Documentation

friend class MemoryManager< MappingXofXZ >
friend

Definition at line 56 of file MappingXofXZ.h.

Member Data Documentation

std::string Nektar::GlobalMapping::MappingXofXZ::className
static
Initial value:

Name of the class.

Definition at line 73 of file MappingXofXZ.h.