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

#include <MappingGeneral.h>

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

Static Public Member Functions

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

Static Public Attributes

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

Protected Member Functions

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

Protected Attributes

Array< OneD, Array< OneD, NekDouble > > m_metricTensor
 
Array< OneD, Array< OneD, NekDouble > > m_invMetricTensor
 
Array< OneD, Array< OneD, NekDouble > > m_deriv
 
Array< OneD, Array< OneD, NekDouble > > m_invDeriv
 
Array< OneD, Array< OneD, NekDouble > > m_Christoffel
 
Array< OneD, NekDoublem_jac
 
- Protected Attributes inherited from Nektar::GlobalMapping::Mapping
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 
Array< OneD, Array< OneD, NekDouble > > m_coords
 Array with the Cartesian coordinates. More...
 
Array< OneD, Array< OneD, NekDouble > > m_coordsVel
 Array with the velocity of the coordinates. More...
 
Array< OneD, Array< OneD, NekDouble > > m_GeometricInfo
 Array with metric terms of the mapping. More...
 
int m_nConvectiveFields
 Number of velocity components. More...
 
std::string m_funcName
 Name of the function containing the coordinates. More...
 
std::string m_velFuncName
 Name of the function containing the velocity of the coordinates. More...
 
bool m_constantJacobian
 Flag defining if the Jacobian is constant. More...
 
bool m_timeDependent
 Flag defining if the Mapping is time-dependent. More...
 
bool m_fromFunction
 Flag defining if the Mapping is defined by a function. More...
 
Array< OneD, Array< OneD, NekDouble > > m_wk1
 
Array< OneD, Array< OneD, NekDouble > > m_wk2
 
Array< OneD, Array< OneD, NekDouble > > m_tmp
 

Friends

class MemoryManager< MappingGeneral >
 

Additional Inherited Members

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

Detailed Description

This class implements the most general mapping, defined by the transformation

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

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

\[ \bar{z} = \bar{z}(x,y,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 51 of file MappingGeneral.h.

Constructor & Destructor Documentation

◆ MappingGeneral()

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

Definition at line 59 of file MappingGeneral.cpp.

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

Member Function Documentation

◆ CalculateChristoffel()

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

Definition at line 381 of file MappingGeneral.cpp.

382 {
383  int physTot = m_fields[0]->GetTotPoints();
384  int nvel = m_nConvectiveFields;
385 
386  Array<OneD, Array<OneD, NekDouble>> gradG(nvel * nvel * nvel);
387  Array<OneD, Array<OneD, NekDouble>> tmp(nvel * nvel * nvel);
388  m_Christoffel = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel * nvel);
389  // Allocate memory
390  for (int i = 0; i < gradG.size(); i++)
391  {
392  gradG[i] = Array<OneD, NekDouble>(physTot, 0.0);
393  tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
394  m_Christoffel[i] = Array<OneD, NekDouble>(physTot, 0.0);
395  }
396 
397  // Set wavespace to false and store current value
398  bool waveSpace = m_fields[0]->GetWaveSpace();
399  m_fields[0]->SetWaveSpace(false);
400 
401  // Calculate gradients of g_ij
402  for (int i = 0; i < nvel; i++)
403  {
404  for (int j = 0; j < nvel; j++)
405  {
406  for (int k = 0; k < nvel; k++)
407  {
408  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[k],
409  m_metricTensor[i * nvel + j],
410  gradG[i * nvel * nvel + j * nvel + k]);
411  }
412  }
413  }
414 
415  // Calculate tmp[p,j,k] = 1/2( gradG[pj,k]+ gradG[pk,j]-gradG[jk,p])
416  for (int p = 0; p < nvel; p++)
417  {
418  for (int j = 0; j < nvel; j++)
419  {
420  for (int k = 0; k < nvel; k++)
421  {
422  Vmath::Vadd(physTot, gradG[p * nvel * nvel + j * nvel + k], 1,
423  gradG[p * nvel * nvel + k * nvel + j], 1,
424  tmp[p * nvel * nvel + j * nvel + k], 1);
425  Vmath::Vsub(physTot, tmp[p * nvel * nvel + j * nvel + k], 1,
426  gradG[j * nvel * nvel + k * nvel + p], 1,
427  tmp[p * nvel * nvel + j * nvel + k], 1);
428  Vmath::Smul(physTot, 0.5, tmp[p * nvel * nvel + j * nvel + k],
429  1, tmp[p * nvel * nvel + j * nvel + k], 1);
430  }
431  }
432  }
433 
434  // Calculate Christoffel symbols = g^ip tmp[p,j,k]
435  for (int i = 0; i < nvel; i++)
436  {
437  for (int j = 0; j < nvel; j++)
438  {
439  for (int k = 0; k < nvel; k++)
440  {
441  for (int p = 0; p < nvel; p++)
442  {
443  Vmath::Vvtvp(
444  physTot, m_invMetricTensor[i * nvel + p], 1,
445  tmp[p * nvel * nvel + j * nvel + k], 1,
446  m_Christoffel[i * nvel * nvel + j * nvel + k], 1,
447  m_Christoffel[i * nvel * nvel + j * nvel + k], 1);
448  }
449  }
450  }
451  }
452  // Restore wavespace
453  m_fields[0]->SetWaveSpace(waveSpace);
454 }
Array< OneD, Array< OneD, NekDouble > > m_metricTensor
Array< OneD, Array< OneD, NekDouble > > m_invMetricTensor
Array< OneD, Array< OneD, NekDouble > > m_Christoffel
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:414
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:406
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:91
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:574
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:359
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:248
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:419

References Nektar::MultiRegions::DirCartesianMap, m_Christoffel, Nektar::GlobalMapping::Mapping::m_fields, m_invMetricTensor, m_metricTensor, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, CellMLToNektar.cellml_metadata::p, Vmath::Smul(), Vmath::Vadd(), Vmath::Vsub(), and Vmath::Vvtvp().

Referenced by v_UpdateGeomInfo().

◆ CalculateMetricTerms()

void Nektar::GlobalMapping::MappingGeneral::CalculateMetricTerms ( )
protected

Definition at line 247 of file MappingGeneral.cpp.

248 {
249  int physTot = m_fields[0]->GetTotPoints();
250  int nvel = m_nConvectiveFields;
251 
252  // Set wavespace to false and store current value
253  bool wavespace = m_fields[0]->GetWaveSpace();
254  m_fields[0]->SetWaveSpace(false);
255 
256  // Allocate memory
257  m_metricTensor = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel);
258  m_invMetricTensor = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel);
259  m_deriv = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel);
260  m_invDeriv = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel);
261  for (int i = 0; i < m_metricTensor.size(); i++)
262  {
263  m_metricTensor[i] = Array<OneD, NekDouble>(physTot, 0.0);
264  m_invMetricTensor[i] = Array<OneD, NekDouble>(physTot, 0.0);
265  m_deriv[i] = Array<OneD, NekDouble>(physTot, 0.0);
266  m_invDeriv[i] = Array<OneD, NekDouble>(physTot, 0.0);
267  }
268  m_jac = Array<OneD, NekDouble>(physTot, 0.0);
269 
270  // First, calculate derivatives of the mapping -> dX^i/dx^j = c^i_j
271  for (int i = 0; i < nvel; i++)
272  {
273  for (int j = 0; j < nvel; j++)
274  {
275  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j],
276  m_coords[i], m_deriv[i * nvel + j]);
277  }
278  }
279  // In Homogeneous case, m_deriv(2,2) needs to be set to 1
280  // because differentiation in wavespace is not valid for non-periodic field
281  if (m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
282  {
283  Vmath::Fill(physTot, 1.0, m_deriv[2 * nvel + 2], 1);
284  }
285 
286  // Now calculate the metric tensor --> g_ij = sum_k { c^k_i c^k_j }
287  for (int i = 0; i < nvel; i++)
288  {
289  for (int j = 0; j < nvel; j++)
290  {
291  for (int k = 0; k < nvel; k++)
292  {
293  Vmath::Vvtvp(physTot, m_deriv[k * nvel + i], 1,
294  m_deriv[k * nvel + j], 1,
295  m_metricTensor[i * nvel + j], 1,
296  m_metricTensor[i * nvel + j], 1);
297  }
298  }
299  }
300 
301  // Put the adjoint of g in m_invMetricTensor
302  switch (nvel)
303  {
304  case 1:
305  Vmath::Fill(physTot, 1.0, m_invMetricTensor[0], 1);
306  break;
307  case 2:
308  Vmath::Vcopy(physTot, m_metricTensor[1 * nvel + 1], 1,
309  m_invMetricTensor[0 * nvel + 0], 1);
310  Vmath::Smul(physTot, -1.0, m_metricTensor[0 * nvel + 1], 1,
311  m_invMetricTensor[1 * nvel + 0], 1);
312  Vmath::Smul(physTot, -1.0, m_metricTensor[1 * nvel + 0], 1,
313  m_invMetricTensor[0 * nvel + 1], 1);
314  Vmath::Vcopy(physTot, m_metricTensor[0 * nvel + 0], 1,
315  m_invMetricTensor[1 * nvel + 1], 1);
316  break;
317  case 3:
318  {
319  int a, b, c, d, e, i, j;
320 
321  // Compute g^{ij} by computing Cofactors(g_ij)^T
322  for (i = 0; i < nvel; ++i)
323  {
324  for (j = 0; j < nvel; ++j)
325  {
326  a = ((i + 1) % nvel) * nvel + ((j + 1) % nvel);
327  b = ((i + 1) % nvel) * nvel + ((j + 2) % nvel);
328  c = ((i + 2) % nvel) * nvel + ((j + 1) % nvel);
329  d = ((i + 2) % nvel) * nvel + ((j + 2) % nvel);
330  e = i * nvel + j;
331  // a*d - b*c
332  Vmath::Vmul(physTot, m_metricTensor[b], 1,
333  m_metricTensor[c], 1, m_invMetricTensor[e], 1);
334  Vmath::Vvtvm(physTot, m_metricTensor[a], 1,
335  m_metricTensor[d], 1, m_invMetricTensor[e], 1,
336  m_invMetricTensor[e], 1);
337  }
338  }
339  break;
340  }
341  }
342 
343  // Compute g = det(g_{ij}) (= Jacobian squared) and store
344  // temporarily in m_jac.
345  for (int i = 0; i < nvel; ++i)
346  {
347  Vmath::Vvtvp(physTot, m_metricTensor[i], 1, m_invMetricTensor[i * nvel],
348  1, m_jac, 1, m_jac, 1);
349  }
350 
351  // Calculate g^ij (the inverse of g_ij) by dividing by jac
352  for (int i = 0; i < nvel * nvel; ++i)
353  {
354  Vmath::Vdiv(physTot, m_invMetricTensor[i], 1, m_jac, 1,
355  m_invMetricTensor[i], 1);
356  }
357 
358  // Compute the Jacobian = sqrt(g)
359  Vmath::Vsqrt(physTot, m_jac, 1, m_jac, 1);
360 
361  // Calculate the derivatives of the inverse transformation
362  // c'^j_i = dx^j/dX^i = sum_k {g^jk c^i_k}
363  for (int i = 0; i < nvel; ++i)
364  {
365  for (int j = 0; j < nvel; ++j)
366  {
367  for (int k = 0; k < nvel; ++k)
368  {
369  Vmath::Vvtvp(physTot, m_deriv[i * nvel + k], 1,
370  m_invMetricTensor[j * nvel + k], 1,
371  m_invDeriv[i * nvel + j], 1,
372  m_invDeriv[i * nvel + j], 1);
373  }
374  }
375  }
376 
377  // Restore value of wavespace
378  m_fields[0]->SetWaveSpace(wavespace);
379 }
Array< OneD, Array< OneD, NekDouble > > m_invDeriv
Array< OneD, Array< OneD, NekDouble > > m_deriv
Array< OneD, NekDouble > m_jac
Array< OneD, Array< OneD, NekDouble > > m_coords
Array with the Cartesian coordinates.
Definition: Mapping.h:408
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:534
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:209
void Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.cpp:598
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:284
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::e3DH1D, Vmath::Fill(), Nektar::GlobalMapping::Mapping::m_coords, m_deriv, Nektar::GlobalMapping::Mapping::m_fields, m_invDeriv, m_invMetricTensor, m_jac, m_metricTensor, Nektar::GlobalMapping::Mapping::m_nConvectiveFields, Vmath::Smul(), Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vsqrt(), Vmath::Vvtvm(), and Vmath::Vvtvp().

Referenced by v_UpdateGeomInfo().

◆ create()

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

Creates an instance of this class.

Definition at line 58 of file MappingGeneral.h.

62  {
65  p->InitObject(pFields, pMapping);
66  return p;
67  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:50

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

◆ v_ApplyChristoffelContravar()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 193 of file MappingGeneral.cpp.

196 {
197  int physTot = m_fields[0]->GetTotPoints();
198  int nvel = m_nConvectiveFields;
199 
200  // Calculate {i,jk} u^j
201  for (int i = 0; i < nvel; i++)
202  {
203  for (int k = 0; k < nvel; k++)
204  {
205  outarray[i * nvel + k] = Array<OneD, NekDouble>(physTot, 0.0);
206  for (int j = 0; j < nvel; j++)
207  {
208  Vmath::Vvtvp(physTot, inarray[j], 1,
209  m_Christoffel[i * nvel * nvel + j * nvel + k], 1,
210  outarray[i * nvel + k], 1, outarray[i * nvel + k],
211  1);
212  }
213  }
214  }
215 }

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

◆ v_ApplyChristoffelCovar()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 217 of file MappingGeneral.cpp.

220 {
221  int physTot = m_fields[0]->GetTotPoints();
222  int nvel = m_nConvectiveFields;
223 
224  // Calculate {i,jk} u_i
225  for (int j = 0; j < nvel; j++)
226  {
227  for (int k = 0; k < nvel; k++)
228  {
229  outarray[j * nvel + k] = Array<OneD, NekDouble>(physTot, 0.0);
230  for (int i = 0; i < nvel; i++)
231  {
232  Vmath::Vvtvp(physTot, inarray[i], 1,
233  m_Christoffel[i * nvel * nvel + j * nvel + k], 1,
234  outarray[j * nvel + k], 1, outarray[j * nvel + k],
235  1);
236  }
237  }
238  }
239 }

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

◆ v_ContravarFromCartesian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 117 of file MappingGeneral.cpp.

120 {
121  int physTot = m_fields[0]->GetTotPoints();
122  int nvel = m_nConvectiveFields;
123 
124  for (int i = 0; i < nvel; i++)
125  {
126  outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
127  for (int j = 0; j < nvel; j++)
128  {
129  Vmath::Vvtvp(physTot, inarray[j], 1, m_invDeriv[j * nvel + i], 1,
130  outarray[i], 1, outarray[i], 1);
131  }
132  }
133 }

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

◆ v_ContravarToCartesian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 81 of file MappingGeneral.cpp.

84 {
85  int physTot = m_fields[0]->GetTotPoints();
86  int nvel = m_nConvectiveFields;
87 
88  for (int i = 0; i < nvel; i++)
89  {
90  outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
91  for (int j = 0; j < nvel; j++)
92  {
93  Vmath::Vvtvp(physTot, inarray[j], 1, m_deriv[i * nvel + j], 1,
94  outarray[i], 1, outarray[i], 1);
95  }
96  }
97 }

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

◆ v_CovarFromCartesian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 135 of file MappingGeneral.cpp.

138 {
139  int physTot = m_fields[0]->GetTotPoints();
140  int nvel = m_nConvectiveFields;
141 
142  for (int i = 0; i < nvel; i++)
143  {
144  outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
145  for (int j = 0; j < nvel; j++)
146  {
147  Vmath::Vvtvp(physTot, inarray[j], 1, m_deriv[j * nvel + i], 1,
148  outarray[i], 1, outarray[i], 1);
149  }
150  }
151 }

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

◆ v_CovarToCartesian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 99 of file MappingGeneral.cpp.

102 {
103  int physTot = m_fields[0]->GetTotPoints();
104  int nvel = m_nConvectiveFields;
105 
106  for (int i = 0; i < nvel; i++)
107  {
108  outarray[i] = Array<OneD, NekDouble>(physTot, 0.0);
109  for (int j = 0; j < nvel; j++)
110  {
111  Vmath::Vvtvp(physTot, inarray[j], 1, m_invDeriv[i * nvel + j], 1,
112  outarray[i], 1, outarray[i], 1);
113  }
114  }
115 }

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

◆ v_GetInvMetricTensor()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 176 of file MappingGeneral.cpp.

178 {
179  int physTot = m_fields[0]->GetTotPoints();
180  int nvel = m_nConvectiveFields;
181 
182  for (int i = 0; i < nvel; i++)
183  {
184  for (int j = 0; j < nvel; j++)
185  {
186  outarray[i * nvel + j] = Array<OneD, NekDouble>(physTot, 0.0);
187  Vmath::Vcopy(physTot, m_invMetricTensor[i * nvel + j], 1,
188  outarray[i * nvel + j], 1);
189  }
190  }
191 }

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

◆ v_GetJacobian()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 153 of file MappingGeneral.cpp.

154 {
155  int physTot = m_fields[0]->GetTotPoints();
156  Vmath::Vcopy(physTot, m_jac, 1, outarray, 1);
157 }

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

◆ v_GetMetricTensor()

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

Implements Nektar::GlobalMapping::Mapping.

Definition at line 159 of file MappingGeneral.cpp.

161 {
162  int physTot = m_fields[0]->GetTotPoints();
163  int nvel = m_nConvectiveFields;
164 
165  for (int i = 0; i < nvel; i++)
166  {
167  for (int j = 0; j < nvel; j++)
168  {
169  outarray[i * nvel + j] = Array<OneD, NekDouble>(physTot, 0.0);
170  Vmath::Vcopy(physTot, m_metricTensor[i * nvel + j], 1,
171  outarray[i * nvel + j], 1);
172  }
173  }
174 }

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

◆ v_InitObject()

void Nektar::GlobalMapping::MappingGeneral::v_InitObject ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const TiXmlElement *  pMapping 
)
overrideprotectedvirtual

This function initialises the Mapping object. It computes the coordinates and velocity coordinates, initialises the workspace variables, and calls UpdateGeomInfo, which will perform the calculations specific for each type of Mapping.

Parameters
pFieldsExpList array used in the mapping
pMappingxml element describing the mapping

Reimplemented from Nektar::GlobalMapping::Mapping.

Definition at line 69 of file MappingGeneral.cpp.

72 {
73  Mapping::v_InitObject(pFields, pMapping);
74 
75  m_constantJacobian = false;
76 
78  "General Mapping needs at least 2 velocity components.");
79 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
virtual GLOBAL_MAPPING_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
Definition: Mapping.cpp:101
bool m_constantJacobian
Flag defining if the Jacobian is constant.
Definition: Mapping.h:423

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

◆ v_UpdateGeomInfo()

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

Friends And Related Function Documentation

◆ MemoryManager< MappingGeneral >

friend class MemoryManager< MappingGeneral >
friend

Definition at line 1 of file MappingGeneral.h.

Member Data Documentation

◆ className

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

Name of the class.

Definition at line 70 of file MappingGeneral.h.

◆ m_Christoffel

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

◆ m_deriv

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingGeneral::m_deriv
protected

◆ m_invDeriv

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingGeneral::m_invDeriv
protected

◆ m_invMetricTensor

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingGeneral::m_invMetricTensor
protected

◆ m_jac

Array<OneD, NekDouble> Nektar::GlobalMapping::MappingGeneral::m_jac
protected

Definition at line 84 of file MappingGeneral.h.

Referenced by CalculateMetricTerms(), and v_GetJacobian().

◆ m_metricTensor

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

Definition at line 79 of file MappingGeneral.h.

Referenced by CalculateChristoffel(), CalculateMetricTerms(), and v_GetMetricTensor().