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)
 
GLOBAL_MAPPING_EXPORT void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping) override
 
GLOBAL_MAPPING_EXPORT void v_ContravarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_CovarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_GetJacobian (Array< OneD, NekDouble > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_GetMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelContravar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
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_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pMapping)
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarToCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_ContravarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_CovarFromCartesian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_GetCartesianCoordinates (Array< OneD, NekDouble > &out0, Array< OneD, NekDouble > &out1, Array< OneD, NekDouble > &out2)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetCoordVelocity (Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetJacobian (Array< OneD, NekDouble > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_DotGradJacobian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_GetMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_GetInvMetricTensor (Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_LowerIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_RaiseIndex (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelContravar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_ApplyChristoffelCovar (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)=0
 
virtual GLOBAL_MAPPING_EXPORT void v_Divergence (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_VelocityLaplacian (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble alpha)
 
virtual GLOBAL_MAPPING_EXPORT void v_gradgradU (const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual GLOBAL_MAPPING_EXPORT void v_CurlCurlField (Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const bool generalized)
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateMapping (const NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &coords=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &coordsVel=NullNekDoubleArrayOfArray)
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateGeomInfo ()=0
 
virtual GLOBAL_MAPPING_EXPORT void v_UpdateBCs (const NekDouble time)
 

Protected Attributes

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 49 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 57 of file MappingGeneral.cpp.

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

Member Function Documentation

◆ CalculateChristoffel()

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

Definition at line 379 of file MappingGeneral.cpp.

380{
381 int physTot = m_fields[0]->GetTotPoints();
382 int nvel = m_nConvectiveFields;
383
384 Array<OneD, Array<OneD, NekDouble>> gradG(nvel * nvel * nvel);
385 Array<OneD, Array<OneD, NekDouble>> tmp(nvel * nvel * nvel);
386 m_Christoffel = Array<OneD, Array<OneD, NekDouble>>(nvel * nvel * nvel);
387 // Allocate memory
388 for (int i = 0; i < gradG.size(); i++)
389 {
390 gradG[i] = Array<OneD, NekDouble>(physTot, 0.0);
391 tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
392 m_Christoffel[i] = Array<OneD, NekDouble>(physTot, 0.0);
393 }
394
395 // Set wavespace to false and store current value
396 bool waveSpace = m_fields[0]->GetWaveSpace();
397 m_fields[0]->SetWaveSpace(false);
398
399 // Calculate gradients of g_ij
400 for (int i = 0; i < nvel; i++)
401 {
402 for (int j = 0; j < nvel; j++)
403 {
404 for (int k = 0; k < nvel; k++)
405 {
407 m_metricTensor[i * nvel + j],
408 gradG[i * nvel * nvel + j * nvel + k]);
409 }
410 }
411 }
412
413 // Calculate tmp[p,j,k] = 1/2( gradG[pj,k]+ gradG[pk,j]-gradG[jk,p])
414 for (int p = 0; p < nvel; p++)
415 {
416 for (int j = 0; j < nvel; j++)
417 {
418 for (int k = 0; k < nvel; k++)
419 {
420 Vmath::Vadd(physTot, gradG[p * nvel * nvel + j * nvel + k], 1,
421 gradG[p * nvel * nvel + k * nvel + j], 1,
422 tmp[p * nvel * nvel + j * nvel + k], 1);
423 Vmath::Vsub(physTot, tmp[p * nvel * nvel + j * nvel + k], 1,
424 gradG[j * nvel * nvel + k * nvel + p], 1,
425 tmp[p * nvel * nvel + j * nvel + k], 1);
426 Vmath::Smul(physTot, 0.5, tmp[p * nvel * nvel + j * nvel + k],
427 1, tmp[p * nvel * nvel + j * nvel + k], 1);
428 }
429 }
430 }
431
432 // Calculate Christoffel symbols = g^ip tmp[p,j,k]
433 for (int i = 0; i < nvel; i++)
434 {
435 for (int j = 0; j < nvel; j++)
436 {
437 for (int k = 0; k < nvel; k++)
438 {
439 for (int p = 0; p < nvel; p++)
440 {
442 physTot, m_invMetricTensor[i * nvel + p], 1,
443 tmp[p * nvel * nvel + j * nvel + k], 1,
444 m_Christoffel[i * nvel * nvel + j * nvel + k], 1,
445 m_Christoffel[i * nvel * nvel + j * nvel + k], 1);
446 }
447 }
448 }
449 }
450 // Restore wavespace
451 m_fields[0]->SetWaveSpace(waveSpace);
452}
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:412
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Definition: Mapping.h:404
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

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

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

References Nektar::UnitTests::d(), 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 56 of file MappingGeneral.h.

60 {
63 p->InitObject(pFields, pMapping);
64 return p;
65 }
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:51

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 191 of file MappingGeneral.cpp.

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

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 215 of file MappingGeneral.cpp.

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

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 115 of file MappingGeneral.cpp.

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

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 79 of file MappingGeneral.cpp.

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

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 133 of file MappingGeneral.cpp.

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

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

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

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 174 of file MappingGeneral.cpp.

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

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 151 of file MappingGeneral.cpp.

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

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 157 of file MappingGeneral.cpp.

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

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 67 of file MappingGeneral.cpp.

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

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.
MappingFactory & GetMappingFactory()
Declaration of the mapping factory singleton.
Definition: Mapping.cpp:49

Name of the class.

Definition at line 68 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 82 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 77 of file MappingGeneral.h.

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