Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::GlobalMapping::MappingGeneral:
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

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)
 
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_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_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_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 (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)
 

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

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 (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...
 
- 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 52 of file MappingGeneral.h.

Constructor & Destructor Documentation

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

Member Function Documentation

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

Definition at line 402 of file MappingGeneral.cpp.

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

Referenced by v_UpdateGeomInfo().

403 {
404  int physTot = m_fields[0]->GetTotPoints();
405  int nvel = m_nConvectiveFields;
406 
407  Array<OneD, Array<OneD, NekDouble> > gradG(nvel*nvel*nvel);
408  Array<OneD, Array<OneD, NekDouble> > tmp(nvel*nvel*nvel);
410  // Allocate memory
411  for (int i = 0; i < gradG.num_elements(); i++)
412  {
413  gradG[i] = Array<OneD, NekDouble>(physTot, 0.0);
414  tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
415  m_Christoffel[i] = Array<OneD, NekDouble>(physTot, 0.0);
416  }
417 
418  // Set wavespace to false and store current value
419  bool waveSpace = m_fields[0]->GetWaveSpace();
420  m_fields[0]->SetWaveSpace(false);
421 
422  //Calculate gradients of g_ij
423  for (int i = 0; i <nvel; i++)
424  {
425  for(int j=0; j<nvel; j++)
426  {
427  for(int k=0; k<nvel; k++)
428  {
429  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[k],
430  m_metricTensor[i*nvel+j],
431  gradG[i*nvel*nvel + j*nvel + k]);
432  }
433  }
434  }
435 
436  // Calculate tmp[p,j,k] = 1/2( gradG[pj,k]+ gradG[pk,j]-gradG[jk,p])
437  for (int p = 0; p <nvel; p++)
438  {
439  for (int j = 0; j < nvel; j++)
440  {
441  for (int k = 0; k < nvel; k++)
442  {
443  Vmath::Vadd(physTot, gradG[p*nvel*nvel + j*nvel + k], 1,
444  gradG[p*nvel*nvel + k*nvel + j], 1,
445  tmp[p*nvel*nvel + j*nvel + k], 1);
446  Vmath::Vsub(physTot, tmp[p*nvel*nvel + j*nvel + k], 1,
447  gradG[j*nvel*nvel + k*nvel + p], 1,
448  tmp[p*nvel*nvel + j*nvel + k], 1);
449  Vmath::Smul(physTot, 0.5, tmp[p*nvel*nvel + j*nvel + k], 1,
450  tmp[p*nvel*nvel + j*nvel + k], 1);
451  }
452  }
453  }
454 
455  // Calculate Christoffel symbols = g^ip tmp[p,j,k]
456  for (int i = 0; i <nvel; i++)
457  {
458  for (int j = 0; j < nvel; j++)
459  {
460  for (int k = 0; k < nvel; k++)
461  {
462  for (int p = 0; p < nvel; p++)
463  {
464  Vmath::Vvtvp(physTot, m_invMetricTensor[i*nvel+p], 1,
465  tmp[p*nvel*nvel + j*nvel + k], 1,
466  m_Christoffel[i*nvel*nvel+j*nvel+k], 1,
467  m_Christoffel[i*nvel*nvel+j*nvel+k], 1);
468  }
469  }
470  }
471  }
472  // Restore wavespace
473  m_fields[0]->SetWaveSpace(waveSpace);
474 
475 }
Array< OneD, Array< OneD, NekDouble > > m_invMetricTensor
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
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, Array< OneD, NekDouble > > m_Christoffel
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
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
Array< OneD, Array< OneD, NekDouble > > m_metricTensor
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 Nektar::GlobalMapping::MappingGeneral::CalculateMetricTerms ( )
protected

Definition at line 264 of file MappingGeneral.cpp.

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

265 {
266  int physTot = m_fields[0]->GetTotPoints();
267  int nvel = m_nConvectiveFields;
268 
269  // Set wavespace to false and store current value
270  bool wavespace = m_fields[0]->GetWaveSpace();
271  m_fields[0]->SetWaveSpace(false);
272 
273  // Allocate memory
278  for (int i = 0; i < m_metricTensor.num_elements(); i++)
279  {
280  m_metricTensor[i] = Array<OneD, NekDouble>(physTot, 0.0);
281  m_invMetricTensor[i] = Array<OneD, NekDouble>(physTot, 0.0);
282  m_deriv[i] = Array<OneD, NekDouble>(physTot, 0.0);
283  m_invDeriv[i] = Array<OneD, NekDouble>(physTot, 0.0);
284  }
285  m_jac = Array<OneD, NekDouble>(physTot, 0.0);
286 
287  // First, calculate derivatives of the mapping -> dX^i/dx^j = c^i_j
288  for( int i = 0; i<nvel; i++)
289  {
290  for( int j = 0; j<nvel; j++)
291  {
292  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j],
293  m_coords[i],
294  m_deriv[i*nvel+j]);
295  }
296  }
297  // In Homogeneous case, m_deriv(2,2) needs to be set to 1
298  // because differentiation in wavespace is not valid for non-periodic field
299  if (m_fields[0]->GetExpType() == MultiRegions::e3DH1D)
300  {
301  Vmath::Fill(physTot, 1.0, m_deriv[2*nvel+2], 1);
302  }
303 
304  // Now calculate the metric tensor --> g_ij = sum_k { c^k_i c^k_j }
305  for( int i = 0; i<nvel; i++)
306  {
307  for( int j = 0; j<nvel; j++)
308  {
309  for( int k = 0; k<nvel; k++)
310  {
311  Vmath::Vvtvp(physTot, m_deriv[k*nvel+i], 1,
312  m_deriv[k*nvel+j], 1,
313  m_metricTensor[i*nvel+j], 1,
314  m_metricTensor[i*nvel+j], 1);
315  }
316  }
317  }
318 
319  // Put the adjoint of g in m_invMetricTensor
320  switch (nvel)
321  {
322  case 1:
323  Vmath::Fill (physTot, 1.0, m_invMetricTensor[0], 1);
324  break;
325  case 2:
326  Vmath::Vcopy(physTot, m_metricTensor[1*nvel+1], 1,
327  m_invMetricTensor[0*nvel+0], 1);
328  Vmath::Smul (physTot, -1.0, m_metricTensor[0*nvel+1], 1,
329  m_invMetricTensor[1*nvel+0], 1);
330  Vmath::Smul (physTot, -1.0, m_metricTensor[1*nvel+0], 1,
331  m_invMetricTensor[0*nvel+1], 1);
332  Vmath::Vcopy(physTot, m_metricTensor[0*nvel+0], 1,
333  m_invMetricTensor[1*nvel+1], 1);
334  break;
335  case 3:
336  {
337  int a, b, c, d, e, i, j;
338 
339  // Compute g^{ij} by computing Cofactors(g_ij)^T
340  for (i = 0; i < nvel; ++i)
341  {
342  for (j = 0; j < nvel; ++j)
343  {
344  a = ((i+1)%nvel) * nvel + ((j+1)%nvel);
345  b = ((i+1)%nvel) * nvel + ((j+2)%nvel);
346  c = ((i+2)%nvel) * nvel + ((j+1)%nvel);
347  d = ((i+2)%nvel) * nvel + ((j+2)%nvel);
348  e = i*nvel + j;
349  // a*d - b*c
350  Vmath::Vmul(physTot, m_metricTensor[b], 1,
351  m_metricTensor[c], 1,
352  m_invMetricTensor[e], 1);
353  Vmath::Vvtvm(physTot, m_metricTensor[a], 1,
354  m_metricTensor[d], 1,
355  m_invMetricTensor[e], 1,
356  m_invMetricTensor[e], 1);
357  }
358  }
359  break;
360  }
361  }
362 
363  // Compute g = det(g_{ij}) (= Jacobian squared) and store
364  // temporarily in m_jac.
365  for (int i = 0; i < nvel; ++i)
366  {
367  Vmath::Vvtvp(physTot, m_metricTensor[i], 1,
368  m_invMetricTensor[i*nvel], 1,
369  m_jac, 1, m_jac, 1);
370  }
371 
372  // Calculate g^ij (the inverse of g_ij) by dividing by jac
373  for (int i = 0; i < nvel*nvel; ++i)
374  {
375  Vmath::Vdiv(physTot, m_invMetricTensor[i], 1, m_jac, 1,
376  m_invMetricTensor[i], 1);
377  }
378 
379  // Compute the Jacobian = sqrt(g)
380  Vmath::Vsqrt(physTot, m_jac, 1, m_jac, 1);
381 
382  // Calculate the derivatives of the inverse transformation
383  // c'^j_i = dx^j/dX^i = sum_k {g^jk c^i_k}
384  for (int i = 0; i < nvel; ++i)
385  {
386  for (int j = 0; j < nvel; ++j)
387  {
388  for (int k = 0; k < nvel; ++k)
389  {
390  Vmath::Vvtvp(physTot, m_deriv[i*nvel+k], 1,
391  m_invMetricTensor[j*nvel+k], 1,
392  m_invDeriv[i*nvel+j], 1,
393  m_invDeriv[i*nvel+j], 1);
394  }
395  }
396  }
397 
398  // Restore value of wavespace
399  m_fields[0]->SetWaveSpace(wavespace);
400 }
Array< OneD, Array< OneD, NekDouble > > m_invMetricTensor
Array< OneD, Array< OneD, NekDouble > > m_coords
Array with the Cartesian coordinates.
Definition: Mapping.h:411
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:394
Array< OneD, Array< OneD, NekDouble > > m_invDeriv
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, Array< OneD, NekDouble > > m_deriv
Array< OneD, NekDouble > m_jac
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
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 plus vector): z = w*x - y
Definition: Vmath.cpp:451
Array< OneD, Array< OneD, NekDouble > > m_metricTensor
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
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 60 of file MappingGeneral.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::MappingGeneral::v_ApplyChristoffelContravar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 209 of file MappingGeneral.cpp.

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

212 {
213  int physTot = m_fields[0]->GetTotPoints();
214  int nvel = m_nConvectiveFields;
215 
216  // Calculate {i,jk} u^j
217  for (int i = 0; i <nvel; i++)
218  {
219  for (int k = 0; k < nvel; k++)
220  {
221  outarray[i*nvel+k] = Array<OneD, NekDouble>(physTot,0.0);
222  for (int j = 0; j < nvel; j++)
223  {
224  Vmath::Vvtvp(physTot, inarray[j], 1,
225  m_Christoffel[i*nvel*nvel+j*nvel+k], 1,
226  outarray[i*nvel+k], 1,
227  outarray[i*nvel+k], 1);
228  }
229  }
230  }
231 
232 }
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
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
Array< OneD, Array< OneD, NekDouble > > m_Christoffel
void Nektar::GlobalMapping::MappingGeneral::v_ApplyChristoffelCovar ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 234 of file MappingGeneral.cpp.

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

237 {
238  int physTot = m_fields[0]->GetTotPoints();
239  int nvel = m_nConvectiveFields;
240 
241  // Calculate {i,jk} u_i
242  for (int j = 0; j < nvel; j++)
243  {
244  for (int k = 0; k < nvel; k++)
245  {
246  outarray[j*nvel+k] = Array<OneD, NekDouble>(physTot,0.0);
247  for (int i = 0; i <nvel; i++)
248  {
249  Vmath::Vvtvp(physTot, inarray[i], 1,
250  m_Christoffel[i*nvel*nvel+j*nvel+k], 1,
251  outarray[j*nvel+k], 1,
252  outarray[j*nvel+k], 1);
253  }
254  }
255  }
256 }
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
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
Array< OneD, Array< OneD, NekDouble > > m_Christoffel
void Nektar::GlobalMapping::MappingGeneral::v_ContravarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 124 of file MappingGeneral.cpp.

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

127 {
128  int physTot = m_fields[0]->GetTotPoints();
129  int nvel = m_nConvectiveFields;
130 
131  for (int i=0; i<nvel; i++)
132  {
133  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
134  for (int j=0; j<nvel; j++)
135  {
136  Vmath::Vvtvp(physTot, inarray[j], 1,
137  m_invDeriv[j*nvel+i], 1,
138  outarray[i], 1,
139  outarray[i], 1);
140  }
141 
142  }
143 }
Array< OneD, Array< OneD, NekDouble > > m_invDeriv
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
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Nektar::GlobalMapping::MappingGeneral::v_ContravarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 82 of file MappingGeneral.cpp.

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

85 {
86  int physTot = m_fields[0]->GetTotPoints();
87  int nvel = m_nConvectiveFields;
88 
89  for (int i=0; i<nvel; i++)
90  {
91  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
92  for (int j=0; j<nvel; j++)
93  {
94  Vmath::Vvtvp(physTot, inarray[j], 1,
95  m_deriv[i*nvel+j], 1,
96  outarray[i], 1,
97  outarray[i], 1);
98  }
99 
100  }
101 }
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
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
Array< OneD, Array< OneD, NekDouble > > m_deriv
void Nektar::GlobalMapping::MappingGeneral::v_CovarFromCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 145 of file MappingGeneral.cpp.

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

148 {
149  int physTot = m_fields[0]->GetTotPoints();
150  int nvel = m_nConvectiveFields;
151 
152  for (int i=0; i<nvel; i++)
153  {
154  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
155  for (int j=0; j<nvel; j++)
156  {
157  Vmath::Vvtvp(physTot, inarray[j], 1,
158  m_deriv[j*nvel+i], 1,
159  outarray[i], 1,
160  outarray[i], 1);
161  }
162 
163  }
164 }
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
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
Array< OneD, Array< OneD, NekDouble > > m_deriv
void Nektar::GlobalMapping::MappingGeneral::v_CovarToCartesian ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 103 of file MappingGeneral.cpp.

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

106 {
107  int physTot = m_fields[0]->GetTotPoints();
108  int nvel = m_nConvectiveFields;
109 
110  for (int i=0; i<nvel; i++)
111  {
112  outarray[i] = Array<OneD, NekDouble> (physTot, 0.0);
113  for (int j=0; j<nvel; j++)
114  {
115  Vmath::Vvtvp(physTot, inarray[j], 1,
116  m_invDeriv[i*nvel+j], 1,
117  outarray[i], 1,
118  outarray[i], 1);
119  }
120 
121  }
122 }
Array< OneD, Array< OneD, NekDouble > > m_invDeriv
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
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
void Nektar::GlobalMapping::MappingGeneral::v_GetInvMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 191 of file MappingGeneral.cpp.

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

193 {
194  int physTot = m_fields[0]->GetTotPoints();
195  int nvel = m_nConvectiveFields;
196 
197  for (int i=0; i<nvel; i++)
198  {
199  for (int j=0; j<nvel; j++)
200  {
201  outarray[i*nvel+j] = Array<OneD, NekDouble> (physTot, 0.0);
202  Vmath::Vcopy(physTot, m_invMetricTensor[i*nvel+j], 1,
203  outarray[i*nvel+j], 1);
204  }
205 
206  }
207 }
Array< OneD, Array< OneD, NekDouble > > m_invMetricTensor
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 Nektar::GlobalMapping::MappingGeneral::v_GetJacobian ( Array< OneD, NekDouble > &  outarray)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 166 of file MappingGeneral.cpp.

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

168 {
169  int physTot = m_fields[0]->GetTotPoints();
170  Vmath::Vcopy(physTot, m_jac, 1, outarray, 1);
171 }
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
Array< OneD, NekDouble > m_jac
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::GlobalMapping::MappingGeneral::v_GetMetricTensor ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
protectedvirtual

Implements Nektar::GlobalMapping::Mapping.

Definition at line 173 of file MappingGeneral.cpp.

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

175 {
176  int physTot = m_fields[0]->GetTotPoints();
177  int nvel = m_nConvectiveFields;
178 
179  for (int i=0; i<nvel; i++)
180  {
181  for (int j=0; j<nvel; j++)
182  {
183  outarray[i*nvel+j] = Array<OneD, NekDouble> (physTot, 0.0);
184  Vmath::Vcopy(physTot, m_metricTensor[i*nvel+j], 1,
185  outarray[i*nvel+j], 1);
186  }
187 
188  }
189 }
int m_nConvectiveFields
Number of velocity components.
Definition: Mapping.h:417
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:409
Array< OneD, Array< OneD, NekDouble > > m_metricTensor
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::GlobalMapping::MappingGeneral::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 70 of file MappingGeneral.cpp.

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

73 {
74  Mapping::v_InitObject(pFields, pMapping);
75 
76  m_constantJacobian = false;
77 
79  "General Mapping needs at least 2 velocity components.");
80 }
#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::MappingGeneral::v_UpdateGeomInfo ( )
protectedvirtual

Friends And Related Function Documentation

friend class MemoryManager< MappingGeneral >
friend

Definition at line 56 of file MappingGeneral.h.

Member Data Documentation

std::string Nektar::GlobalMapping::MappingGeneral::className
static
Initial value:
=
MappingGeneral::create, "X = X(x,y,z), Y = Y(x,y,z), Z=Z(x,y,z)")

Name of the class.

Definition at line 73 of file MappingGeneral.h.

Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingGeneral::m_Christoffel
protected
Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingGeneral::m_deriv
protected
Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingGeneral::m_invDeriv
protected
Array<OneD, Array<OneD, NekDouble> > Nektar::GlobalMapping::MappingGeneral::m_invMetricTensor
protected
Array<OneD, NekDouble> Nektar::GlobalMapping::MappingGeneral::m_jac
protected

Definition at line 87 of file MappingGeneral.h.

Referenced by CalculateMetricTerms(), and v_GetJacobian().

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

Definition at line 82 of file MappingGeneral.h.

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