Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::SpatialDomains::GeomFactors Class Reference

Calculation and storage of geometric factors associated with the mapping from StdRegions reference elements to a given LocalRegions physical element in the mesh. More...

#include <GeomFactors.h>

Public Member Functions

 GeomFactors (const GeomType gtype, const int coordim, const StdRegions::StdExpansionSharedPtr &xmap, const Array< OneD, Array< OneD, NekDouble > > &coords)
 Constructor for GeomFactors class. More...
 
 GeomFactors (const GeomFactors &S)
 Copy constructor. More...
 
 ~GeomFactors ()
 Destructor. More...
 
DerivStorage GetDeriv (const LibUtilities::PointsKeyVector &keyTgt)
 Return the derivative of the mapping with respect to the reference coordinates, \(\frac{\partial \chi_i}{\partial \xi_j}\). More...
 
const Array< OneD, const NekDoubleGetJac (const LibUtilities::PointsKeyVector &keyTgt)
 Return the Jacobian of the mapping and cache the result. More...
 
const Array< TwoD, const NekDoubleGetGmat (const LibUtilities::PointsKeyVector &keyTgt)
 Return the Laplacian coefficients \(g_{ij}\). More...
 
const Array< TwoD, const NekDoubleGetDerivFactors (const LibUtilities::PointsKeyVector &keyTgt)
 Return the derivative of the reference coordinates with respect to the mapping, \(\frac{\partial \xi_i}{\partial \chi_j}\). More...
 
void GetMovingFrames (const LibUtilities::PointsKeyVector &keyTgt, const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Returns moving frames. More...
 
GeomType GetGtype ()
 Returns whether the geometry is regular or deformed. More...
 
bool IsValid () const
 Determine if element is valid and not self-intersecting. More...
 
int GetCoordim () const
 Return the number of dimensions of the coordinate system. More...
 
size_t GetHash ()
 Computes a hash of this GeomFactors element. More...
 

Protected Member Functions

StdRegions::StdExpansionSharedPtrGetXmap ()
 Return the Xmap;. More...
 

Protected Attributes

GeomType m_type
 Type of geometry (e.g. eRegular, eDeformed, eMovingRegular). More...
 
int m_expDim
 Dimension of expansion. More...
 
int m_coordDim
 Dimension of coordinate system. More...
 
bool m_valid
 Validity of element (Jacobian positive) More...
 
enum GeomMMF m_MMFDir
 Principle tangent direction for MMF. More...
 
StdRegions::StdExpansionSharedPtr m_xmap
 Stores information about the expansion. More...
 
Array< OneD, Array< OneD, NekDouble > > m_coords
 Stores coordinates of the geometry. More...
 
std::map< LibUtilities::PointsKeyVector, Array< OneD, NekDouble > > m_jacCache
 Jacobian vector cache. More...
 
std::map< LibUtilities::PointsKeyVector, Array< TwoD, NekDouble > > m_derivFactorCache
 DerivFactors vector cache. More...
 

Private Member Functions

void CheckIfValid ()
 Tests if the element is valid and not self-intersecting. More...
 
DerivStorage ComputeDeriv (const LibUtilities::PointsKeyVector &keyTgt) const
 
Array< OneD, NekDoubleComputeJac (const LibUtilities::PointsKeyVector &keyTgt) const
 Return the Jacobian of the mapping and cache the result. More...
 
Array< TwoD, NekDoubleComputeGmat (const LibUtilities::PointsKeyVector &keyTgt) const
 Computes the Laplacian coefficients \(g_{ij}\). More...
 
Array< TwoD, NekDoubleComputeDerivFactors (const LibUtilities::PointsKeyVector &keyTgt) const
 Return the derivative of the reference coordinates with respect to the mapping, \(\frac{\partial \xi_i}{\partial \chi_j}\). More...
 
void ComputeMovingFrames (const LibUtilities::PointsKeyVector &keyTgt, const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &movingframes)
 
void Interp (const LibUtilities::PointsKeyVector &src_points, const Array< OneD, const NekDouble > &src, const LibUtilities::PointsKeyVector &tgt_points, Array< OneD, NekDouble > &tgt) const
 Perform interpolation of data between two point distributions. More...
 
void Adjoint (const Array< TwoD, const NekDouble > &src, Array< TwoD, NekDouble > &tgt) const
 Compute the transpose of the cofactors matrix. More...
 
void ComputePrincipleDirection (const LibUtilities::PointsKeyVector &keyTgt, const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &output)
 
void VectorNormalise (Array< OneD, Array< OneD, NekDouble > > &array)
 
void VectorCrossProd (const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &v3)
 Computes the vector cross-product in 3D of v1 and v2, storing the result in v3. More...
 

Friends

bool operator== (const GeomFactors &lhs, const GeomFactors &rhs)
 Tests if two GeomFactors classes are equal. More...
 

Detailed Description

Calculation and storage of geometric factors associated with the mapping from StdRegions reference elements to a given LocalRegions physical element in the mesh.

This class stores the various geometric factors associated with a specific element, necessary for fundamental integration and differentiation operations as well as edge and surface normals.

Initially, these algorithms are provided with a mapping from the reference region element to the physical element. Practically, this is represented using a corresponding reference region element for each coordinate component. Note that for straight-sided elements, these elements will be of linear order. Curved elements are represented using higher-order coordinate mappings. This geometric order is in contrast to the order of the spectral/hp expansion order on the element.

For application of the chain rule during differentiation we require the partial derivatives

\[\frac{\partial \xi_i}{\partial \chi_j}\]

evaluated at the physical points of the expansion basis. We also construct the inverse metric tensor \(g^{ij}\) which, in the case of a domain embedded in a higher-dimensional space, supports the conversion of covariant quantities to contravariant quantities. When the expansion dimension is equal to the coordinate dimension the Jacobian of the mapping \(\chi_j\) is a square matrix and consequently the required terms are the entries of the inverse of the Jacobian. However, in general this is not the case, so we therefore implement the construction of these terms following the derivation in Cantwell, et. al. [CaYaKiPeSh13]. Given the coordinate maps \(\chi_i\), this comprises of five steps

  1. Compute the terms of the Jacobian \(\frac{\partial \chi_i}{\partial \xi_j}\).
  2. Compute the metric tensor \(g_{ij}=\mathbf{t}_{(i)}\cdot\mathbf{t}_{(j)}\).
  3. Compute the square of the Jacobian determinant \(g=|\mathbf{g}|\).
  4. Compute the inverse metric tensor \(g^{ij}\).
  5. Compute the terms \(\frac{\partial \xi_i}{\partial \chi_j}\).

Definition at line 75 of file GeomFactors.h.

Constructor & Destructor Documentation

◆ GeomFactors() [1/2]

Nektar::SpatialDomains::GeomFactors::GeomFactors ( const GeomType  gtype,
const int  coordim,
const StdRegions::StdExpansionSharedPtr xmap,
const Array< OneD, Array< OneD, NekDouble > > &  coords 
)

Constructor for GeomFactors class.

Parameters
gtypeSpecified whether the geometry is regular or deformed.
coordimSpecifies the dimension of the coordinate system.
CoordsCoordinate maps of the element.

Definition at line 88 of file GeomFactors.cpp.

92  :
93  m_type(gtype),
94  m_expDim(xmap->GetShapeDimension()),
95  m_coordDim(coordim),
96  m_valid(true),
97  m_xmap(xmap),
98  m_coords(coords)
99  {
100  CheckIfValid();
101  }
int m_coordDim
Dimension of coordinate system.
Definition: GeomFactors.h:139
void CheckIfValid()
Tests if the element is valid and not self-intersecting.
int m_expDim
Dimension of expansion.
Definition: GeomFactors.h:137
bool m_valid
Validity of element (Jacobian positive)
Definition: GeomFactors.h:141
StdRegions::StdExpansionSharedPtr m_xmap
Stores information about the expansion.
Definition: GeomFactors.h:147
GeomType m_type
Type of geometry (e.g. eRegular, eDeformed, eMovingRegular).
Definition: GeomFactors.h:135
Array< OneD, Array< OneD, NekDouble > > m_coords
Stores coordinates of the geometry.
Definition: GeomFactors.h:149

References CheckIfValid().

◆ GeomFactors() [2/2]

Nektar::SpatialDomains::GeomFactors::GeomFactors ( const GeomFactors S)

Copy constructor.

Parameters
SAn instance of a GeomFactors class from which to construct a new instance.

Definition at line 108 of file GeomFactors.cpp.

108  :
109  m_type(S.m_type),
110  m_expDim(S.m_expDim),
111  m_coordDim(S.m_coordDim),
112  m_valid(S.m_valid),
113  m_xmap(S.m_xmap),
114  m_coords(S.m_coords)
115  {
116  }

◆ ~GeomFactors()

Nektar::SpatialDomains::GeomFactors::~GeomFactors ( )

Destructor.

Definition at line 122 of file GeomFactors.cpp.

123  {
124  }

Member Function Documentation

◆ Adjoint()

void Nektar::SpatialDomains::GeomFactors::Adjoint ( const Array< TwoD, const NekDouble > &  src,
Array< TwoD, NekDouble > &  tgt 
) const
private

Compute the transpose of the cofactors matrix.

Input and output arrays are of dimension (m_expDim*m_expDim) x num_points. The first index of the input and output arrays are ordered row-by-row.

Parameters
srcInput data array.
tgtStorage for adjoint matrix data.

Definition at line 718 of file GeomFactors.cpp.

721  {
722  ASSERTL1(src.size() == tgt.size(),
723  "Source matrix is of different size to destination"
724  "matrix for computing adjoint.");
725 
726  int n = src[0].size();
727  switch (m_expDim)
728  {
729  case 1:
730  Vmath::Fill (n, 1.0, &tgt[0][0], 1);
731  break;
732  case 2:
733  Vmath::Vcopy(n, &src[3][0], 1, &tgt[0][0], 1);
734  Vmath::Smul (n, -1.0, &src[1][0], 1, &tgt[1][0], 1);
735  Vmath::Smul (n, -1.0, &src[2][0], 1, &tgt[2][0], 1);
736  Vmath::Vcopy(n, &src[0][0], 1, &tgt[3][0], 1);
737  break;
738  case 3:
739  {
740  int a, b, c, d, e, i, j;
741 
742  // Compute g^{ij} by computing Cofactors(g_ij)^T
743  for (i = 0; i < m_expDim; ++i)
744  {
745  for (j = 0; j < m_expDim; ++j)
746  {
747  a = ((i+1)%m_expDim) * m_expDim + ((j+1)%m_expDim);
748  b = ((i+1)%m_expDim) * m_expDim + ((j+2)%m_expDim);
749  c = ((i+2)%m_expDim) * m_expDim + ((j+1)%m_expDim);
750  d = ((i+2)%m_expDim) * m_expDim + ((j+2)%m_expDim);
751  e = j*m_expDim + i;
752  Vmath::Vvtvvtm(n, &src[a][0], 1, &src[d][0], 1,
753  &src[b][0], 1, &src[c][0], 1,
754  &tgt[e][0], 1);
755  }
756  }
757  break;
758  }
759  }
760  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
void Vvtvvtm(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtm (vector times vector minus vector times vector):
Definition: Vmath.cpp:658
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
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:1199

References ASSERTL1, Vmath::Fill(), m_expDim, Vmath::Smul(), Vmath::Vcopy(), and Vmath::Vvtvvtm().

Referenced by ComputeDerivFactors(), ComputeGmat(), and ComputeJac().

◆ CheckIfValid()

void Nektar::SpatialDomains::GeomFactors::CheckIfValid ( )
private

Tests if the element is valid and not self-intersecting.

Constructs the Jacobian as per Spencer's book p158 and tests if negative.

Definition at line 610 of file GeomFactors.cpp.

611  {
612  // Jacobian test only makes sense when expdim = coorddim
613  // If one-dimensional then element is valid.
614  if (m_coordDim != m_expDim || m_expDim == 1)
615  {
616  m_valid = true;
617  return;
618  }
619 
621  int nqtot = 1;
622  for (int i = 0; i < m_expDim; ++i)
623  {
624  p[i] = m_xmap->GetBasis(i)->GetPointsKey();
625  nqtot *= p[i].GetNumPoints();
626  }
627  int pts = (m_type == eRegular || m_type == eMovingRegular)
628  ? 1 : nqtot;
629  Array<OneD, NekDouble> jac(pts, 0.0);
630 
631  DerivStorage deriv = GetDeriv(p);
632 
633  switch (m_expDim)
634  {
635  case 2:
636  {
637  Vmath::Vvtvvtm(pts, &deriv[0][0][0], 1, &deriv[1][1][0], 1,
638  &deriv[1][0][0], 1, &deriv[0][1][0], 1,
639  &jac[0], 1);
640  break;
641  }
642  case 3:
643  {
644  Array<OneD, NekDouble> tmp(pts, 0.0);
645 
646  Vmath::Vvtvvtm(pts, &deriv[1][1][0], 1, &deriv[2][2][0], 1,
647  &deriv[2][1][0], 1, &deriv[1][2][0], 1,
648  &tmp[0], 1);
649  Vmath::Vvtvp (pts, &deriv[0][0][0], 1, &tmp[0], 1,
650  &jac[0], 1, &jac[0], 1);
651 
652  Vmath::Vvtvvtm(pts, &deriv[2][1][0], 1, &deriv[0][2][0], 1,
653  &deriv[0][1][0], 1, &deriv[2][2][0], 1,
654  &tmp[0], 1);
655  Vmath::Vvtvp (pts, &deriv[1][0][0], 1, &tmp[0], 1,
656  &jac[0], 1, &jac[0], 1);
657 
658  Vmath::Vvtvvtm(pts, &deriv[0][1][0], 1, &deriv[1][2][0], 1,
659  &deriv[1][1][0], 1, &deriv[0][2][0], 1,
660  &tmp[0], 1);
661  Vmath::Vvtvp (pts, &deriv[2][0][0], 1, &tmp[0], 1,
662  &jac[0], 1, &jac[0], 1);
663 
664  break;
665  }
666  }
667 
668  if (Vmath::Vmin(pts, &jac[0], 1) < 0)
669  {
670  m_valid = false;
671  }
672  }
DerivStorage GetDeriv(const LibUtilities::PointsKeyVector &keyTgt)
Return the derivative of the mapping with respect to the reference coordinates, .
Definition: GeomFactors.h:232
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > DerivStorage
Storage type for derivative of mapping.
Definition: GeomFactors.h:70
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.cpp:992
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:513

References Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, GetDeriv(), m_coordDim, m_expDim, m_type, m_valid, m_xmap, CellMLToNektar.cellml_metadata::p, Vmath::Vmin(), Vmath::Vvtvp(), and Vmath::Vvtvvtm().

Referenced by GeomFactors().

◆ ComputeDeriv()

DerivStorage Nektar::SpatialDomains::GeomFactors::ComputeDeriv ( const LibUtilities::PointsKeyVector keyTgt) const
private

Derivatives are computed at the geometry point distributions and interpolated to the target point distributions.

Parameters
tpointsTarget point distributions.
Returns
Derivative of coordinate map evaluated at target point distributions.

Definition at line 169 of file GeomFactors.cpp.

171  {
172  ASSERTL1(keyTgt.size() == m_expDim,
173  "Dimension of target point distribution does not match "
174  "expansion dimension.");
175 
176  int i = 0, j = 0;
177  int nqtot_map = 1;
178  int nqtot_tbasis = 1;
182 
183  // Allocate storage and compute number of points
184  for (i = 0; i < m_expDim; ++i)
185  {
186  map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
187  nqtot_map *= map_points[i].GetNumPoints();
188  nqtot_tbasis *= keyTgt[i].GetNumPoints();
189  deriv[i] = Array<OneD, Array<OneD,NekDouble> >(m_coordDim);
190  d_map[i] = Array<OneD, Array<OneD,NekDouble> >(m_coordDim);
191  }
192 
193  // Calculate local derivatives
194  for(i = 0; i < m_coordDim; ++i)
195  {
196  Array<OneD, NekDouble> tmp(nqtot_map);
197  // Transform from coefficient space to physical space
198  m_xmap->BwdTrans(m_coords[i], tmp);
199 
200  // Allocate storage and take the derivative (calculated at the
201  // points as specified in 'Coords')
202  for (j = 0; j < m_expDim; ++j)
203  {
204  d_map[j][i] = Array<OneD,NekDouble>(nqtot_map);
205  deriv[j][i] = Array<OneD,NekDouble>(nqtot_tbasis);
206  m_xmap->StdPhysDeriv(j, tmp, d_map[j][i]);
207  }
208  }
209 
210  for (i = 0; i < m_coordDim; ++i)
211  {
212  // Interpolate the derivatives:
213  // - from the points as defined in the mapping ('Coords')
214  // - to the points at which we want to know the metrics
215  // ('tbasis')
216  bool same = true;
217  for (j = 0; j < m_expDim; ++j)
218  {
219  same = same && (map_points[j] == keyTgt[j]);
220  }
221  if( same )
222  {
223  for (j = 0; j < m_expDim; ++j)
224  {
225  deriv[j][i] = d_map[j][i];
226  }
227  }
228  else
229  {
230  for (j = 0; j < m_expDim; ++j)
231  {
232  Interp(map_points, d_map[j][i], keyTgt, deriv[j][i]);
233  }
234  }
235  }
236 
237  return deriv;
238  }
void Interp(const LibUtilities::PointsKeyVector &src_points, const Array< OneD, const NekDouble > &src, const LibUtilities::PointsKeyVector &tgt_points, Array< OneD, NekDouble > &tgt) const
Perform interpolation of data between two point distributions.

References ASSERTL1, Interp(), m_coordDim, m_coords, m_expDim, and m_xmap.

Referenced by ComputeDerivFactors(), ComputeGmat(), ComputeJac(), ComputeMovingFrames(), and GetDeriv().

◆ ComputeDerivFactors()

Array< TwoD, NekDouble > Nektar::SpatialDomains::GeomFactors::ComputeDerivFactors ( const LibUtilities::PointsKeyVector keyTgt) const
private

Return the derivative of the reference coordinates with respect to the mapping, \(\frac{\partial \xi_i}{\partial \chi_j}\).

Parameters
keyTgtTarget point distributions.
Returns
Derivative factors evaluated at the target point distributions. A 1D example: /f$ Jac =(\partial x/ \partial \xi) /f$ ; /f$ factor = 1/Jac = (\partial \xi/ \partial x) /f$

Definition at line 399 of file GeomFactors.cpp.

401  {
402  ASSERTL1(keyTgt.size() == m_expDim,
403  "Dimension of target point distribution does not match "
404  "expansion dimension.");
405 
406  int i = 0, j = 0, k = 0, l = 0;
407  int ptsTgt = 1;
408 
409  if (m_type == eDeformed)
410  {
411  // Allocate storage and compute number of points
412  for (i = 0; i < m_expDim; ++i)
413  {
414  ptsTgt *= keyTgt[i].GetNumPoints();
415  }
416  }
417 
418  // Get derivative at geometry points
419  DerivStorage deriv = ComputeDeriv(keyTgt);
420 
421  Array<TwoD, NekDouble> tmp (m_expDim*m_expDim, ptsTgt, 0.0);
422  Array<TwoD, NekDouble> gmat(m_expDim*m_expDim, ptsTgt, 0.0);
423  Array<OneD, NekDouble> jac (ptsTgt, 0.0);
424  Array<TwoD, NekDouble> factors(m_expDim*m_coordDim, ptsTgt, 0.0);
425 
426  // Compute g_{ij} as t_i \cdot t_j and store in tmp
427  for (i = 0, l = 0; i < m_expDim; ++i)
428  {
429  for (j = 0; j < m_expDim; ++j, ++l)
430  {
431  for (k = 0; k < m_coordDim; ++k)
432  {
433  Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
434  &deriv[j][k][0], 1,
435  &tmp[l][0], 1,
436  &tmp[l][0], 1);
437  }
438  }
439  }
440 
441  Adjoint(tmp, gmat);
442 
443  // Compute g = det(g_{ij}) (= Jacobian squared) and store
444  // temporarily in m_jac.
445  for (i = 0; i < m_expDim; ++i)
446  {
447  Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i*m_expDim][0], 1,
448  &jac[0], 1, &jac[0], 1);
449  }
450 
451  for (i = 0; i < m_expDim*m_expDim; ++i)
452  {
453  Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
454  }
455 
456  // Compute the Jacobian = sqrt(g)
457  Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
458 
459  // Compute the derivative factors
460  for (k = 0, l = 0; k < m_coordDim; ++k)
461  {
462  for (j = 0; j < m_expDim; ++j, ++l)
463  {
464  for (i = 0; i < m_expDim; ++i)
465  {
466  Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
467  &gmat[m_expDim*i+j][0], 1,
468  &factors[l][0], 1,
469  &factors[l][0], 1);
470  }
471  }
472  }
473 
474  return factors;
475  }
DerivStorage ComputeDeriv(const LibUtilities::PointsKeyVector &keyTgt) const
void Adjoint(const Array< TwoD, const NekDouble > &src, Array< TwoD, NekDouble > &tgt) const
Compute the transpose of the cofactors matrix.
@ eDeformed
Geometry is curved or has non-constant factors.
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:475
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:257

References Adjoint(), ASSERTL1, ComputeDeriv(), Nektar::SpatialDomains::eDeformed, m_coordDim, m_expDim, m_type, Vmath::Vdiv(), Vmath::Vsqrt(), and Vmath::Vvtvp().

Referenced by GetDerivFactors().

◆ ComputeGmat()

Array< TwoD, NekDouble > Nektar::SpatialDomains::GeomFactors::ComputeGmat ( const LibUtilities::PointsKeyVector keyTgt) const
private

Computes the Laplacian coefficients \(g_{ij}\).

This routine returns a two-dimensional array of values specifying the inverse metric terms associated with the coordinate mapping of the corresponding reference region to the physical element. These terms correspond to the \(g^{ij}\) terms in [CaYaKiPeSh13] and, in the case of an embedded manifold, map covariant quantities to contravariant quantities. The leading index of the array is the index of the term in the tensor numbered as

\[\left(\begin{array}{ccc} 0 & 1 & 2 \\ 1 & 3 & 4 \\ 2 & 4 & 5 \end{array}\right)\]

. The second dimension is either of size 1 in the case of elements having GeomType eRegular, or of size equal to the number of quadrature points for eDeformed elements.

See also
Wikipedia "Covariance and Contravariance of Vectors"
Returns
Two-dimensional array containing the inverse metric tensor of the coordinate mapping.

Definition at line 332 of file GeomFactors.cpp.

334  {
335  ASSERTL1(keyTgt.size() == m_expDim,
336  "Dimension of target point distribution does not match "
337  "expansion dimension.");
338 
339  int i = 0, j = 0, k = 0, l = 0;
340  int ptsTgt = 1;
341 
342  if (m_type == eDeformed)
343  {
344  // Allocate storage and compute number of points
345  for (i = 0; i < m_expDim; ++i)
346  {
347  ptsTgt *= keyTgt[i].GetNumPoints();
348  }
349  }
350 
351  // Get derivative at geometry points
352  DerivStorage deriv = ComputeDeriv(keyTgt);
353 
354  Array<TwoD, NekDouble> tmp (m_expDim*m_expDim, ptsTgt, 0.0);
355  Array<TwoD, NekDouble> gmat(m_expDim*m_expDim, ptsTgt, 0.0);
356  Array<OneD, NekDouble> jac (ptsTgt, 0.0);
357 
358  // Compute g_{ij} as t_i \cdot t_j and store in tmp
359  for (i = 0, l = 0; i < m_expDim; ++i)
360  {
361  for (j = 0; j < m_expDim; ++j, ++l)
362  {
363  for (k = 0; k < m_coordDim; ++k)
364  {
365  Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
366  &deriv[j][k][0], 1,
367  &tmp[l][0], 1,
368  &tmp[l][0], 1);
369  }
370  }
371  }
372 
373  Adjoint(tmp, gmat);
374 
375  // Compute g = det(g_{ij}) (= Jacobian squared) and store
376  // temporarily in m_jac.
377  for (i = 0; i < m_expDim; ++i)
378  {
379  Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i*m_expDim][0], 1,
380  &jac[0], 1, &jac[0], 1);
381  }
382 
383  for (i = 0; i < m_expDim*m_expDim; ++i)
384  {
385  Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
386  }
387 
388  return gmat;
389  }

References Adjoint(), ASSERTL1, ComputeDeriv(), Nektar::SpatialDomains::eDeformed, m_coordDim, m_expDim, m_type, Vmath::Vdiv(), and Vmath::Vvtvp().

Referenced by GetGmat().

◆ ComputeJac()

Array< OneD, NekDouble > Nektar::SpatialDomains::GeomFactors::ComputeJac ( const LibUtilities::PointsKeyVector keyTgt) const
private

Return the Jacobian of the mapping and cache the result.

This routine returns an array of values specifying the Jacobian of the mapping at quadrature points in the element. The array is either of size 1 in the case of elements having GeomType eRegular, or of size equal to the number of quadrature points for eDeformed elements.

Returns
Array containing the Jacobian of the coordinate mapping at the quadrature points of the element.
See also
GeomType

Definition at line 252 of file GeomFactors.cpp.

254  {
255  ASSERTL1(keyTgt.size() == m_expDim,
256  "Dimension of target point distribution does not match "
257  "expansion dimension.");
258 
259  int i = 0, j = 0, k = 0, l = 0;
260  int ptsTgt = 1;
261 
262  if (m_type == eDeformed)
263  {
264  // Allocate storage and compute number of points
265  for (i = 0; i < m_expDim; ++i)
266  {
267  ptsTgt *= keyTgt[i].GetNumPoints();
268  }
269  }
270 
271  // Get derivative at geometry points
272  DerivStorage deriv = ComputeDeriv(keyTgt);
273 
274  Array<TwoD, NekDouble> tmp (m_expDim*m_expDim, ptsTgt, 0.0);
275  Array<TwoD, NekDouble> gmat(m_expDim*m_expDim, ptsTgt, 0.0);
276  Array<OneD, NekDouble> jac (ptsTgt, 0.0);
277 
278  // Compute g_{ij} as t_i \cdot t_j and store in tmp
279  for (i = 0, l = 0; i < m_expDim; ++i)
280  {
281  for (j = 0; j < m_expDim; ++j, ++l)
282  {
283  for (k = 0; k < m_coordDim; ++k)
284  {
285  Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
286  &deriv[j][k][0], 1,
287  &tmp[l][0], 1,
288  &tmp[l][0], 1);
289  }
290  }
291  }
292 
293  Adjoint(tmp, gmat);
294 
295  // Compute g = det(g_{ij}) (= Jacobian squared) and store
296  // temporarily in m_jac.
297  for (i = 0; i < m_expDim; ++i)
298  {
299  Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i*m_expDim][0], 1,
300  &jac[0], 1, &jac[0], 1);
301  }
302 
303  // Compute the Jacobian = sqrt(g)
304  Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
305 
306  return jac;
307  }

References Adjoint(), ASSERTL1, ComputeDeriv(), Nektar::SpatialDomains::eDeformed, m_coordDim, m_expDim, m_type, Vmath::Vsqrt(), and Vmath::Vvtvp().

Referenced by GetJac().

◆ ComputeMovingFrames()

void Nektar::SpatialDomains::GeomFactors::ComputeMovingFrames ( const LibUtilities::PointsKeyVector keyTgt,
const SpatialDomains::GeomMMF  MMFdir,
const Array< OneD, const NekDouble > &  CircCentre,
Array< OneD, Array< OneD, NekDouble > > &  movingframes 
)
private

Definition at line 477 of file GeomFactors.cpp.

482  {
483  ASSERTL1(keyTgt.size() == m_expDim,
484  "Dimension of target point distribution does not match "
485  "expansion dimension.");
486 
487  int i = 0, k = 0;
488  int ptsTgt = 1;
489  int nq = 1;
490 
491  for (i = 0; i < m_expDim; ++i)
492  {
493  nq *= keyTgt[i].GetNumPoints();
494  }
495 
496  if (m_type == eDeformed)
497  {
498  // Allocate storage and compute number of points
499  for (i = 0; i < m_expDim; ++i)
500  {
501  ptsTgt *= keyTgt[i].GetNumPoints();
502  }
503  }
504 
505  // Get derivative at geometry points
506  DerivStorage deriv = ComputeDeriv(keyTgt);
507 
508  // number of moving frames is requited to be 3, even for surfaces
509  int MFdim = 3;
510 
511  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > MFtmp(MFdim);
512 
513  // Compute g_{ij} as t_i \cdot t_j and store in tmp
514  for (i = 0; i < MFdim; ++i)
515  {
516  MFtmp[i] = Array<OneD, Array<OneD, NekDouble> >(m_coordDim);
517  for (k = 0; k < m_coordDim; ++k)
518  {
519  MFtmp[i][k] = Array<OneD, NekDouble>(nq);
520  }
521  }
522 
523  // Compute g_{ij} as t_i \cdot t_j and store in tmp
524  for (i = 0; i < MFdim-1; ++i)
525  {
526  for (k = 0; k < m_coordDim; ++k)
527  {
528  if (m_type == eDeformed)
529  {
530  Vmath::Vcopy(ptsTgt, &deriv[i][k][0], 1,
531  &MFtmp[i][k][0], 1);
532  }
533  else
534  {
535  Vmath::Fill(nq, deriv[i][k][0], MFtmp[i][k], 1);
536  }
537  }
538  }
539 
540  // Direction of MF1 is preserved: MF2 is considered in the same
541  // tangent plane as MF1. MF3 is computed by cross product of MF1
542  // and MF2. MF2 is consequently computed as the cross product of
543  // MF3 and MF1.
544  Array<OneD, Array<OneD,NekDouble> > PrincipleDir(m_coordDim);
545  for (int k=0; k<m_coordDim; k++)
546  {
547  PrincipleDir[k] = Array<OneD, NekDouble>(nq);
548  }
549 
550  if( !(MMFdir == eLOCAL) )
551  {
552  ComputePrincipleDirection(keyTgt, MMFdir,
553  factors, PrincipleDir);
554  }
555 
556  // MF3 = MF1 \times MF2
557  VectorCrossProd(MFtmp[0], MFtmp[1], MFtmp[2]);
558 
559  // Normalizing MF3
560  VectorNormalise(MFtmp[2]);
561 
562  if ( !(MMFdir == eLOCAL) )
563  {
564  Array<OneD, NekDouble> temp(nq, 0.0);
565 
566  // Reorient MF1 along the PrincipleDir
567  for (i = 0; i < m_coordDim; ++i)
568  {
569  Vmath::Vvtvp(nq, MFtmp[2][i], 1,
570  PrincipleDir[i], 1,
571  temp, 1,
572  temp, 1);
573  }
574  Vmath::Neg(nq, temp, 1);
575 
576  // u2 = v2 - < u1 , v2 > ( u1 / < u1, u1 > )
577  for (i = 0; i < m_coordDim; ++i)
578  {
579  Vmath::Vvtvp(nq, temp, 1,
580  MFtmp[2][i], 1,
581  PrincipleDir[i], 1,
582  MFtmp[0][i], 1);
583  }
584  }
585 
586  // Normalizing MF1
587  VectorNormalise(MFtmp[0]);
588 
589  // MF2 = MF3 \times MF1
590  VectorCrossProd(MFtmp[2], MFtmp[0], MFtmp[1]);
591 
592  // Normalizing MF2
593  VectorNormalise(MFtmp[1]);
594 
595  for (i = 0; i < MFdim; ++i)
596  {
597  for (k = 0; k < m_coordDim; ++k)
598  {
599  Vmath::Vcopy(nq, &MFtmp[i][k][0], 1,
600  &movingframes[i*m_coordDim+k][0], 1);
601  }
602  }
603  }
void ComputePrincipleDirection(const LibUtilities::PointsKeyVector &keyTgt, const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &output)
void VectorNormalise(Array< OneD, Array< OneD, NekDouble > > &array)
void VectorCrossProd(const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &v3)
Computes the vector cross-product in 3D of v1 and v2, storing the result in v3.
@ eLOCAL
No Principal direction.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461

References ASSERTL1, ComputeDeriv(), ComputePrincipleDirection(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::eLOCAL, Vmath::Fill(), m_coordDim, m_expDim, m_type, Vmath::Neg(), Vmath::Vcopy(), VectorCrossProd(), VectorNormalise(), and Vmath::Vvtvp().

Referenced by GetMovingFrames().

◆ ComputePrincipleDirection()

void Nektar::SpatialDomains::GeomFactors::ComputePrincipleDirection ( const LibUtilities::PointsKeyVector keyTgt,
const SpatialDomains::GeomMMF  MMFdir,
const Array< OneD, const NekDouble > &  CircCentre,
Array< OneD, Array< OneD, NekDouble > > &  output 
)
private

Definition at line 766 of file GeomFactors.cpp.

771  {
772  int nq = output[0].size();
773 
774  output = Array<OneD,Array<OneD,NekDouble> >(m_coordDim);
775  for (int i = 0; i < m_coordDim; ++i)
776  {
777  output[i] = Array<OneD, NekDouble> (nq, 0.0);
778  }
779 
780  // Construction of Connection
781  switch(MMFdir)
782  {
783  // projection to x-axis
784  case eTangentX:
785  {
786  Vmath::Fill(nq, 1.0, output[0], 1);
787  break;
788  }
789  case eTangentY:
790  {
791  Vmath::Fill(nq, 1.0, output[1], 1);
792  break;
793  }
794  case eTangentXY:
795  {
796  Vmath::Fill(nq, sqrt(2.0), output[0], 1);
797  Vmath::Fill(nq, sqrt(2.0), output[1], 1);
798  break;
799  }
800  case eTangentZ:
801  {
802  Vmath::Fill(nq, 1.0, output[2], 1);
803  break;
804  }
805  case eTangentCircular:
806  {
807  // Tangent direction depends on spatial location.
808  Array<OneD, Array<OneD, NekDouble> > x(m_coordDim);
809  for (int k = 0; k < m_coordDim; k++)
810  {
811  x[k] = Array<OneD, NekDouble>(nq);
812  }
813 
814  // m_coords are StdExpansions which store the mapping
815  // between the std element and the local element. Bwd
816  // transforming the std element minimum basis gives a
817  // minimum physical basis for geometry. Need to then
818  // interpolate this up to the quadrature basis.
819  int nqtot_map = 1;
821  for (int i = 0; i < m_expDim; ++i)
822  {
823  map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
824  nqtot_map *= map_points[i].GetNumPoints();
825  }
826  Array<OneD, NekDouble> tmp(nqtot_map);
827  for (int k = 0; k < m_coordDim; k++)
828  {
829  m_xmap->BwdTrans(m_coords[k], tmp);
830  Interp(map_points, tmp, keyTgt, x[k]);
831  }
832 
833  // circular around the center of the domain
834  NekDouble radius, xc=0.0, yc=0.0, xdis, ydis;
835  NekDouble la, lb;
836 
837  ASSERTL1(factors.size() >= 4,
838  "factors is too short.");
839 
840  la = factors[0];
841  lb = factors[1];
842  xc = factors[2];
843  yc = factors[3];
844 
845  for (int i = 0; i < nq; i++)
846  {
847  xdis = x[0][i]-xc;
848  ydis = x[1][i]-yc;
849  radius = sqrt(xdis*xdis/la/la+ydis*ydis/lb/lb);
850  output[0][i] = ydis/radius;
851  output[1][i] = -1.0*xdis/radius;
852  }
853  break;
854  }
855  case eTangentIrregular:
856  {
857  // Tangent direction depends on spatial location.
858  Array<OneD, Array<OneD, NekDouble> > x(m_coordDim);
859  for (int k = 0; k < m_coordDim; k++)
860  {
861  x[k] = Array<OneD, NekDouble>(nq);
862  }
863 
864  int nqtot_map = 1;
866  for (int i = 0; i < m_expDim; ++i)
867  {
868  map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
869  nqtot_map *= map_points[i].GetNumPoints();
870  }
871  Array<OneD, NekDouble> tmp(nqtot_map);
872  for (int k = 0; k < m_coordDim; k++)
873  {
874  m_xmap->BwdTrans(m_coords[k], tmp);
875  Interp(map_points, tmp, keyTgt, x[k]);
876  }
877 
878  // circular around the center of the domain
879  NekDouble xtan, ytan, mag;
880  for (int i = 0; i < nq; i++)
881  {
882  xtan = -1.0*(x[1][i]*x[1][i]*x[1][i] + x[1][i]);
883  ytan = 2.0*x[0][i];
884  mag = sqrt(xtan*xtan + ytan*ytan);
885  output[0][i] = xtan/mag;
886  output[1][i] = ytan/mag;
887  }
888  break;
889  }
890  case eTangentNonconvex:
891  {
892  // Tangent direction depends on spatial location.
893  Array<OneD, Array<OneD, NekDouble> > x(m_coordDim);
894  for (int k = 0; k < m_coordDim; k++)
895  {
896  x[k] = Array<OneD, NekDouble>(nq);
897  }
898 
899  int nqtot_map = 1;
901  for (int i = 0; i < m_expDim; ++i)
902  {
903  map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
904  nqtot_map *= map_points[i].GetNumPoints();
905  }
906  Array<OneD, NekDouble> tmp(nqtot_map);
907  for (int k = 0; k < m_coordDim; k++)
908  {
909  m_xmap->BwdTrans(m_coords[k], tmp);
910  Interp(map_points, tmp, keyTgt, x[k]);
911  }
912 
913  // circular around the center of the domain
914  NekDouble xtan, ytan, mag;
915  for (int i = 0; i < nq; i++)
916  {
917  xtan = -2.0*x[1][i]*x[1][i]*x[1][i] + x[1][i];
918  ytan = sqrt(3.0)*x[0][i];
919  mag = sqrt(xtan*xtan + ytan*ytan);
920  output[0][i] = xtan/mag;
921  output[1][i] = ytan/mag;
922  }
923  break;
924  }
925  default:
926  {
927  break;
928  }
929  }
930  }
@ eTangentIrregular
Circular around the centre of domain.
@ eTangentX
X coordinate direction.
@ eTangentCircular
Circular around the centre of domain.
@ eTangentNonconvex
Circular around the centre of domain.
@ eTangentXY
XY direction.
@ eTangentZ
Z coordinate direction.
@ eTangentY
Y coordinate direction.
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

References ASSERTL1, Nektar::SpatialDomains::eTangentCircular, Nektar::SpatialDomains::eTangentIrregular, Nektar::SpatialDomains::eTangentNonconvex, Nektar::SpatialDomains::eTangentX, Nektar::SpatialDomains::eTangentXY, Nektar::SpatialDomains::eTangentY, Nektar::SpatialDomains::eTangentZ, Vmath::Fill(), Interp(), m_coordDim, m_coords, m_expDim, m_xmap, and tinysimd::sqrt().

Referenced by ComputeMovingFrames().

◆ GetCoordim()

int Nektar::SpatialDomains::GeomFactors::GetCoordim ( ) const
inline

Return the number of dimensions of the coordinate system.

This is greater than or equal to the expansion dimension.

Returns
The dimension of the coordinate system.

Definition at line 328 of file GeomFactors.h.

329  {
330  return m_coordDim;
331  }

References m_coordDim.

◆ GetDeriv()

DerivStorage Nektar::SpatialDomains::GeomFactors::GetDeriv ( const LibUtilities::PointsKeyVector tpoints)
inline

Return the derivative of the mapping with respect to the reference coordinates, \(\frac{\partial \chi_i}{\partial \xi_j}\).

Parameters
tpointsTarget point distributions.
Returns
Derivative evaluated at target point distributions.
See also
GeomFactors::ComputeDeriv

Definition at line 232 of file GeomFactors.h.

234  {
235  return ComputeDeriv(tpoints);
236  }

References ComputeDeriv().

Referenced by CheckIfValid().

◆ GetDerivFactors()

const Array< TwoD, const NekDouble > Nektar::SpatialDomains::GeomFactors::GetDerivFactors ( const LibUtilities::PointsKeyVector keyTgt)
inline

Return the derivative of the reference coordinates with respect to the mapping, \(\frac{\partial \xi_i}{\partial \chi_j}\).

Returns cached value if available, otherwise computes derivative factors and stores result in cache.

Parameters
keyTgtTarget point distributions.
Returns
Derivative factors evaluated at target point distributions.
See also
GeomFactors::ComputeDerivFactors

Definition at line 286 of file GeomFactors.h.

288  {
289  auto x = m_derivFactorCache.find(keyTgt);
290 
291  if (x != m_derivFactorCache.end())
292  {
293  return x->second;
294  }
295 
296  m_derivFactorCache[keyTgt] = ComputeDerivFactors(keyTgt);
297 
298  return m_derivFactorCache[keyTgt];
299 
300  }
Array< TwoD, NekDouble > ComputeDerivFactors(const LibUtilities::PointsKeyVector &keyTgt) const
Return the derivative of the reference coordinates with respect to the mapping, .
std::map< LibUtilities::PointsKeyVector, Array< TwoD, NekDouble > > m_derivFactorCache
DerivFactors vector cache.
Definition: GeomFactors.h:155

References ComputeDerivFactors(), and m_derivFactorCache.

◆ GetGmat()

const Array< TwoD, const NekDouble > Nektar::SpatialDomains::GeomFactors::GetGmat ( const LibUtilities::PointsKeyVector keyTgt)
inline

Return the Laplacian coefficients \(g_{ij}\).

Parameters
keyTgtTarget point distributions.
Returns
Inverse metric tensor evaluated at target point distributions.
See also
GeomFactors::ComputeGmat

Definition at line 270 of file GeomFactors.h.

272  {
273  return ComputeGmat(keyTgt);
274  }
Array< TwoD, NekDouble > ComputeGmat(const LibUtilities::PointsKeyVector &keyTgt) const
Computes the Laplacian coefficients .

References ComputeGmat().

◆ GetGtype()

GeomType Nektar::SpatialDomains::GeomFactors::GetGtype ( )
inline

Returns whether the geometry is regular or deformed.

A geometric shape is considered regular if it has constant geometric information, and deformed if this information changes throughout the shape.

Returns
The type of geometry.
See also
GeomType

Definition at line 318 of file GeomFactors.h.

319  {
320  return m_type;
321  }

References m_type.

◆ GetHash()

size_t Nektar::SpatialDomains::GeomFactors::GetHash ( )
inline

Computes a hash of this GeomFactors element.

The hash is computed from the geometry type, expansion dimension, coordinate dimension and Jacobian.

Returns
Hash of this GeomFactors object.

Definition at line 348 of file GeomFactors.h.

349  {
350  LibUtilities::PointsKeyVector ptsKeys = m_xmap->GetPointsKeys();
351  const Array<OneD, const NekDouble> jac = GetJac(ptsKeys);
352 
353  size_t hash = 0;
354  hash_combine(hash, (int)m_type, m_expDim, m_coordDim);
355  if (m_type == eDeformed)
356  {
357  hash_range(hash, jac.begin(), jac.end());
358  }
359  else
360  {
361  hash_combine(hash, jac[0]);
362  }
363  return hash;
364  }
const Array< OneD, const NekDouble > GetJac(const LibUtilities::PointsKeyVector &keyTgt)
Return the Jacobian of the mapping and cache the result.
Definition: GeomFactors.h:248
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:69
void hash_combine(std::size_t &seed)
Definition: HashUtils.hpp:46

References Nektar::SpatialDomains::eDeformed, GetJac(), Nektar::hash_combine(), Nektar::hash_range(), m_coordDim, m_expDim, m_type, and m_xmap.

◆ GetJac()

const Array< OneD, const NekDouble > Nektar::SpatialDomains::GeomFactors::GetJac ( const LibUtilities::PointsKeyVector keyTgt)
inline

Return the Jacobian of the mapping and cache the result.

Returns cached value if available, otherwise computes Jacobian and stores result in cache.

Parameters
keyTgtTarget point distributions.
Returns
Jacobian evaluated at target point distributions.
See also
GeomFactors::ComputeJac

Definition at line 248 of file GeomFactors.h.

250  {
251  auto x = m_jacCache.find(keyTgt);
252 
253  if (x != m_jacCache.end())
254  {
255  return x->second;
256  }
257 
258  m_jacCache[keyTgt] = ComputeJac(keyTgt);
259 
260  return m_jacCache[keyTgt];
261 
262  }
Array< OneD, NekDouble > ComputeJac(const LibUtilities::PointsKeyVector &keyTgt) const
Return the Jacobian of the mapping and cache the result.
std::map< LibUtilities::PointsKeyVector, Array< OneD, NekDouble > > m_jacCache
Jacobian vector cache.
Definition: GeomFactors.h:152

References ComputeJac(), and m_jacCache.

Referenced by GetHash().

◆ GetMovingFrames()

void Nektar::SpatialDomains::GeomFactors::GetMovingFrames ( const LibUtilities::PointsKeyVector keyTgt,
const SpatialDomains::GeomMMF  MMFdir,
const Array< OneD, const NekDouble > &  CircCentre,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
inline

Returns moving frames.

Definition at line 302 of file GeomFactors.h.

306  {
307  ComputeMovingFrames(keyTgt,MMFdir,CircCentre,outarray);
308  }
void ComputeMovingFrames(const LibUtilities::PointsKeyVector &keyTgt, const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &movingframes)

References ComputeMovingFrames().

◆ GetXmap()

StdRegions::StdExpansionSharedPtr & Nektar::SpatialDomains::GeomFactors::GetXmap ( void  )
inlineprotected

Return the Xmap;.

Definition at line 366 of file GeomFactors.h.

367  {
368  return m_xmap;
369  }

References m_xmap.

◆ Interp()

void Nektar::SpatialDomains::GeomFactors::Interp ( const LibUtilities::PointsKeyVector src_points,
const Array< OneD, const NekDouble > &  src,
const LibUtilities::PointsKeyVector tgt_points,
Array< OneD, NekDouble > &  tgt 
) const
private

Perform interpolation of data between two point distributions.

Parameters
map_pointsSource data point distribution.
srcSource data to be interpolated.
tpointsTarget data point distribution.
tgtTarget data storage.

Definition at line 681 of file GeomFactors.cpp.

686  {
687  ASSERTL1(src_points.size() == tgt_points.size(),
688  "Dimension of target point distribution does not match "
689  "expansion dimension.");
690 
691  switch (m_expDim)
692  {
693  case 1:
694  LibUtilities::Interp1D(src_points[0], src,
695  tgt_points[0], tgt);
696  break;
697  case 2:
698  LibUtilities::Interp2D(src_points[0], src_points[1], src,
699  tgt_points[0], tgt_points[1], tgt);
700  break;
701  case 3:
702  LibUtilities::Interp3D(src_points[0], src_points[1],
703  src_points[2], src,
704  tgt_points[0], tgt_points[1],
705  tgt_points[2], tgt);
706  break;
707  }
708  }
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:53
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
Definition: Interp.cpp:185
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:115

References ASSERTL1, Nektar::LibUtilities::Interp1D(), Nektar::LibUtilities::Interp2D(), Nektar::LibUtilities::Interp3D(), and m_expDim.

Referenced by ComputeDeriv(), and ComputePrincipleDirection().

◆ IsValid()

bool Nektar::SpatialDomains::GeomFactors::IsValid ( ) const
inline

Determine if element is valid and not self-intersecting.

The validity test is performed by testing if the Jacobian is negative at any point in the shape.

Returns
True if the element is not self-intersecting.

Definition at line 338 of file GeomFactors.h.

339  {
340  return m_valid;
341  }

References m_valid.

◆ VectorCrossProd()

void Nektar::SpatialDomains::GeomFactors::VectorCrossProd ( const Array< OneD, const Array< OneD, NekDouble > > &  v1,
const Array< OneD, const Array< OneD, NekDouble > > &  v2,
Array< OneD, Array< OneD, NekDouble > > &  v3 
)
private

Computes the vector cross-product in 3D of v1 and v2, storing the result in v3.

Parameters
v1First input vector.
v2Second input vector.
v3Output vector computed to be orthogonal to both v1 and v2.

Definition at line 978 of file GeomFactors.cpp.

982  {
983  ASSERTL0(v1.size() == 3,
984  "Input 1 has dimension not equal to 3.");
985  ASSERTL0(v2.size() == 3,
986  "Input 2 has dimension not equal to 3.");
987  ASSERTL0(v3.size() == 3,
988  "Output vector has dimension not equal to 3.");
989 
990  int nq = v1[0].size();
991  Array<OneD, NekDouble> temp(nq);
992 
993  Vmath::Vmul (nq, v1[2], 1, v2[1], 1, temp, 1);
994  Vmath::Vvtvm(nq, v1[1], 1, v2[2], 1, temp, 1, v3[0], 1);
995 
996  Vmath::Vmul (nq, v1[0], 1, v2[2], 1, temp, 1);
997  Vmath::Vvtvm(nq, v1[2], 1, v2[0], 1, temp, 1, v3[1], 1);
998 
999  Vmath::Vmul (nq, v1[1], 1, v2[0], 1, temp, 1);
1000  Vmath::Vvtvm(nq, v1[0], 1, v2[1], 1, temp, 1, v3[2], 1);
1001  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:192
void 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:541

References ASSERTL0, Vmath::Vmul(), and Vmath::Vvtvm().

Referenced by ComputeMovingFrames().

◆ VectorNormalise()

void Nektar::SpatialDomains::GeomFactors::VectorNormalise ( Array< OneD, Array< OneD, NekDouble > > &  array)
private

Definition at line 936 of file GeomFactors.cpp.

938  {
939  int ndim = array.size();
940  ASSERTL0(ndim > 0, "Number of components must be > 0.");
941  for (int i = 1; i < ndim; ++i)
942  {
943  ASSERTL0(array[i].size() == array[0].size(),
944  "Array size mismatch in coordinates.");
945  }
946 
947  int nq = array[0].size();
948  Array<OneD, NekDouble> norm (nq, 0.0);
949 
950  // Compute the norm of each vector.
951  for (int i = 0; i < ndim; ++i)
952  {
953  Vmath::Vvtvp(nq, array[i], 1,
954  array[i], 1,
955  norm, 1,
956  norm, 1);
957  }
958 
959  Vmath::Vsqrt(nq, norm, 1, norm, 1);
960 
961  // Normalise the vectors by the norm
962  for (int i = 0; i < ndim; ++i)
963  {
964  Vmath::Vdiv(nq, array[i], 1, norm, 1, array[i], 1);
965  }
966  }

References ASSERTL0, Vmath::Vdiv(), Vmath::Vsqrt(), and Vmath::Vvtvp().

Referenced by ComputeMovingFrames().

Friends And Related Function Documentation

◆ operator==

bool operator== ( const GeomFactors lhs,
const GeomFactors rhs 
)
friend

Tests if two GeomFactors classes are equal.

Member data equivalence is tested in the following order: shape type, expansion dimension, coordinate dimension and coordinates.

Definition at line 131 of file GeomFactors.cpp.

132  {
133  if(!(lhs.m_type == rhs.m_type))
134  {
135  return false;
136  }
137 
138  if(!(lhs.m_expDim == rhs.m_expDim))
139  {
140  return false;
141  }
142 
143  if(!(lhs.m_coordDim == rhs.m_coordDim))
144  {
145  return false;
146  }
147 
148  const Array<OneD, const NekDouble> jac_lhs =
149  lhs.ComputeJac(lhs.m_xmap->GetPointsKeys());
150  const Array<OneD, const NekDouble> jac_rhs =
151  rhs.ComputeJac(rhs.m_xmap->GetPointsKeys());
152  if(!(jac_lhs == jac_rhs))
153  {
154  return false;
155  }
156 
157  return true;
158  }

Member Data Documentation

◆ m_coordDim

int Nektar::SpatialDomains::GeomFactors::m_coordDim
protected

◆ m_coords

Array<OneD, Array<OneD, NekDouble> > Nektar::SpatialDomains::GeomFactors::m_coords
protected

Stores coordinates of the geometry.

Definition at line 149 of file GeomFactors.h.

Referenced by ComputeDeriv(), and ComputePrincipleDirection().

◆ m_derivFactorCache

std::map<LibUtilities::PointsKeyVector, Array<TwoD, NekDouble> > Nektar::SpatialDomains::GeomFactors::m_derivFactorCache
protected

DerivFactors vector cache.

Definition at line 155 of file GeomFactors.h.

Referenced by GetDerivFactors().

◆ m_expDim

int Nektar::SpatialDomains::GeomFactors::m_expDim
protected

◆ m_jacCache

std::map<LibUtilities::PointsKeyVector, Array<OneD, NekDouble> > Nektar::SpatialDomains::GeomFactors::m_jacCache
protected

Jacobian vector cache.

Definition at line 152 of file GeomFactors.h.

Referenced by GetJac().

◆ m_MMFDir

enum GeomMMF Nektar::SpatialDomains::GeomFactors::m_MMFDir
protected

Principle tangent direction for MMF.

Definition at line 141 of file GeomFactors.h.

◆ m_type

GeomType Nektar::SpatialDomains::GeomFactors::m_type
protected

Type of geometry (e.g. eRegular, eDeformed, eMovingRegular).

Definition at line 135 of file GeomFactors.h.

Referenced by CheckIfValid(), ComputeDerivFactors(), ComputeGmat(), ComputeJac(), ComputeMovingFrames(), GetGtype(), and GetHash().

◆ m_valid

bool Nektar::SpatialDomains::GeomFactors::m_valid
protected

Validity of element (Jacobian positive)

Definition at line 141 of file GeomFactors.h.

Referenced by CheckIfValid(), and IsValid().

◆ m_xmap

StdRegions::StdExpansionSharedPtr Nektar::SpatialDomains::GeomFactors::m_xmap
protected

Stores information about the expansion.

Definition at line 147 of file GeomFactors.h.

Referenced by CheckIfValid(), ComputeDeriv(), ComputePrincipleDirection(), GetHash(), and GetXmap().