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 ()=default
 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 71 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 86 of file GeomFactors.cpp.

89 : m_type(gtype), m_expDim(xmap->GetShapeDimension()), m_coordDim(coordim),
90 m_valid(true), m_xmap(xmap), m_coords(coords)
91{
93}
int m_coordDim
Dimension of coordinate system.
Definition: GeomFactors.h:131
void CheckIfValid()
Tests if the element is valid and not self-intersecting.
int m_expDim
Dimension of expansion.
Definition: GeomFactors.h:129
bool m_valid
Validity of element (Jacobian positive)
Definition: GeomFactors.h:133
StdRegions::StdExpansionSharedPtr m_xmap
Stores information about the expansion.
Definition: GeomFactors.h:139
GeomType m_type
Type of geometry (e.g. eRegular, eDeformed, eMovingRegular).
Definition: GeomFactors.h:127
Array< OneD, Array< OneD, NekDouble > > m_coords
Stores coordinates of the geometry.
Definition: GeomFactors.h:141

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 99 of file GeomFactors.cpp.

100 : m_type(S.m_type), m_expDim(S.m_expDim), m_coordDim(S.m_coordDim),
101 m_valid(S.m_valid), m_xmap(S.m_xmap), m_coords(S.m_coords)
102{
103}

◆ ~GeomFactors()

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

Destructor.

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 668 of file GeomFactors.cpp.

670{
671 ASSERTL1(src.size() == tgt.size(),
672 "Source matrix is of different size to destination"
673 "matrix for computing adjoint.");
674
675 int n = src[0].size();
676 switch (m_expDim)
677 {
678 case 1:
679 Vmath::Fill(n, 1.0, &tgt[0][0], 1);
680 break;
681 case 2:
682 Vmath::Vcopy(n, &src[3][0], 1, &tgt[0][0], 1);
683 Vmath::Smul(n, -1.0, &src[1][0], 1, &tgt[1][0], 1);
684 Vmath::Smul(n, -1.0, &src[2][0], 1, &tgt[2][0], 1);
685 Vmath::Vcopy(n, &src[0][0], 1, &tgt[3][0], 1);
686 break;
687 case 3:
688 {
689 int a, b, c, d, e, i, j;
690
691 // Compute g^{ij} by computing Cofactors(g_ij)^T
692 for (i = 0; i < m_expDim; ++i)
693 {
694 for (j = 0; j < m_expDim; ++j)
695 {
696 a = ((i + 1) % m_expDim) * m_expDim + ((j + 1) % m_expDim);
697 b = ((i + 1) % m_expDim) * m_expDim + ((j + 2) % m_expDim);
698 c = ((i + 2) % m_expDim) * m_expDim + ((j + 1) % m_expDim);
699 d = ((i + 2) % m_expDim) * m_expDim + ((j + 2) % m_expDim);
700 e = j * m_expDim + i;
701 Vmath::Vvtvvtm(n, &src[a][0], 1, &src[d][0], 1, &src[b][0],
702 1, &src[c][0], 1, &tgt[e][0], 1);
703 }
704 }
705 break;
706 }
707 }
708}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
std::vector< double > d(NPUPPER *NPUPPER)
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.hpp:456
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 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 ASSERTL1, Nektar::UnitTests::d(), 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 570 of file GeomFactors.cpp.

571{
572 // Jacobian test only makes sense when expdim = coorddim
573 // If one-dimensional then element is valid.
574 if (m_coordDim != m_expDim || m_expDim == 1)
575 {
576 m_valid = true;
577 return;
578 }
579
581 int nqtot = 1;
582 for (int i = 0; i < m_expDim; ++i)
583 {
584 p[i] = m_xmap->GetBasis(i)->GetPointsKey();
585 nqtot *= p[i].GetNumPoints();
586 }
587 int pts = (m_type == eRegular || m_type == eMovingRegular) ? 1 : nqtot;
588 Array<OneD, NekDouble> jac(pts, 0.0);
589
590 DerivStorage deriv = GetDeriv(p);
591
592 switch (m_expDim)
593 {
594 case 2:
595 {
596 Vmath::Vvtvvtm(pts, &deriv[0][0][0], 1, &deriv[1][1][0], 1,
597 &deriv[1][0][0], 1, &deriv[0][1][0], 1, &jac[0], 1);
598 break;
599 }
600 case 3:
601 {
602 Array<OneD, NekDouble> tmp(pts, 0.0);
603
604 Vmath::Vvtvvtm(pts, &deriv[1][1][0], 1, &deriv[2][2][0], 1,
605 &deriv[2][1][0], 1, &deriv[1][2][0], 1, &tmp[0], 1);
606 Vmath::Vvtvp(pts, &deriv[0][0][0], 1, &tmp[0], 1, &jac[0], 1,
607 &jac[0], 1);
608
609 Vmath::Vvtvvtm(pts, &deriv[2][1][0], 1, &deriv[0][2][0], 1,
610 &deriv[0][1][0], 1, &deriv[2][2][0], 1, &tmp[0], 1);
611 Vmath::Vvtvp(pts, &deriv[1][0][0], 1, &tmp[0], 1, &jac[0], 1,
612 &jac[0], 1);
613
614 Vmath::Vvtvvtm(pts, &deriv[0][1][0], 1, &deriv[1][2][0], 1,
615 &deriv[1][1][0], 1, &deriv[0][2][0], 1, &tmp[0], 1);
616 Vmath::Vvtvp(pts, &deriv[2][0][0], 1, &tmp[0], 1, &jac[0], 1,
617 &jac[0], 1);
618
619 break;
620 }
621 }
622
623 if (Vmath::Vmin(pts, &jac[0], 1) < 0)
624 {
625 m_valid = false;
626 }
627}
DerivStorage GetDeriv(const LibUtilities::PointsKeyVector &keyTgt)
Return the derivative of the mapping with respect to the reference coordinates, .
Definition: GeomFactors.h:216
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > DerivStorage
Storage type for derivative of mapping.
Definition: GeomFactors.h:66
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.hpp:725
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

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 146 of file GeomFactors.cpp.

148{
149 ASSERTL1(keyTgt.size() == m_expDim,
150 "Dimension of target point distribution does not match "
151 "expansion dimension.");
152
153 int i = 0, j = 0;
154 int nqtot_map = 1;
155 int nqtot_tbasis = 1;
159
160 // Allocate storage and compute number of points
161 for (i = 0; i < m_expDim; ++i)
162 {
163 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
164 nqtot_map *= map_points[i].GetNumPoints();
165 nqtot_tbasis *= keyTgt[i].GetNumPoints();
166 deriv[i] = Array<OneD, Array<OneD, NekDouble>>(m_coordDim);
167 d_map[i] = Array<OneD, Array<OneD, NekDouble>>(m_coordDim);
168 }
169
170 // Calculate local derivatives
171 for (i = 0; i < m_coordDim; ++i)
172 {
173 Array<OneD, NekDouble> tmp(nqtot_map);
174 // Transform from coefficient space to physical space
175 m_xmap->BwdTrans(m_coords[i], tmp);
176
177 // Allocate storage and take the derivative (calculated at the
178 // points as specified in 'Coords')
179 for (j = 0; j < m_expDim; ++j)
180 {
181 d_map[j][i] = Array<OneD, NekDouble>(nqtot_map);
182 deriv[j][i] = Array<OneD, NekDouble>(nqtot_tbasis);
183 m_xmap->StdPhysDeriv(j, tmp, d_map[j][i]);
184 }
185 }
186
187 for (i = 0; i < m_coordDim; ++i)
188 {
189 // Interpolate the derivatives:
190 // - from the points as defined in the mapping ('Coords')
191 // - to the points at which we want to know the metrics
192 // ('tbasis')
193 bool same = true;
194 for (j = 0; j < m_expDim; ++j)
195 {
196 same = same && (map_points[j] == keyTgt[j]);
197 }
198 if (same)
199 {
200 for (j = 0; j < m_expDim; ++j)
201 {
202 deriv[j][i] = d_map[j][i];
203 }
204 }
205 else
206 {
207 for (j = 0; j < m_expDim; ++j)
208 {
209 Interp(map_points, d_map[j][i], keyTgt, deriv[j][i]);
210 }
211 }
212 }
213
214 return deriv;
215}
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 369 of file GeomFactors.cpp.

371{
372 ASSERTL1(keyTgt.size() == m_expDim,
373 "Dimension of target point distribution does not match "
374 "expansion dimension.");
375
376 int i = 0, j = 0, k = 0, l = 0;
377 int ptsTgt = 1;
378
379 if (m_type == eDeformed)
380 {
381 // Allocate storage and compute number of points
382 for (i = 0; i < m_expDim; ++i)
383 {
384 ptsTgt *= keyTgt[i].GetNumPoints();
385 }
386 }
387
388 // Get derivative at geometry points
389 DerivStorage deriv = ComputeDeriv(keyTgt);
390
391 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
392 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
393 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
394 Array<TwoD, NekDouble> factors(m_expDim * m_coordDim, ptsTgt, 0.0);
395
396 // Compute g_{ij} as t_i \cdot t_j and store in tmp
397 for (i = 0, l = 0; i < m_expDim; ++i)
398 {
399 for (j = 0; j < m_expDim; ++j, ++l)
400 {
401 for (k = 0; k < m_coordDim; ++k)
402 {
403 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
404 &tmp[l][0], 1, &tmp[l][0], 1);
405 }
406 }
407 }
408
409 Adjoint(tmp, gmat);
410
411 // Compute g = det(g_{ij}) (= Jacobian squared) and store
412 // temporarily in m_jac.
413 for (i = 0; i < m_expDim; ++i)
414 {
415 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
416 1, &jac[0], 1);
417 }
418
419 for (i = 0; i < m_expDim * m_expDim; ++i)
420 {
421 Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
422 }
423
424 // Compute the Jacobian = sqrt(g)
425 Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
426
427 // Compute the derivative factors
428 for (k = 0, l = 0; k < m_coordDim; ++k)
429 {
430 for (j = 0; j < m_expDim; ++j, ++l)
431 {
432 for (i = 0; i < m_expDim; ++i)
433 {
434 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
435 &gmat[m_expDim * i + j][0], 1, &factors[l][0], 1,
436 &factors[l][0], 1);
437 }
438 }
439 }
440
441 return factors;
442}
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.
StdRegions::ConstFactorMap factors
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.hpp:340
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

References Adjoint(), ASSERTL1, ComputeDeriv(), Nektar::SpatialDomains::eDeformed, Nektar::VarcoeffHashingTest::factors, 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"] (http://en.wikipedia.org/wiki/Covariance_and_contravariance_of_vectors)
Returns
Two-dimensional array containing the inverse metric tensor of the coordinate mapping.

Definition at line 305 of file GeomFactors.cpp.

307{
308 ASSERTL1(keyTgt.size() == m_expDim,
309 "Dimension of target point distribution does not match "
310 "expansion dimension.");
311
312 int i = 0, j = 0, k = 0, l = 0;
313 int ptsTgt = 1;
314
315 if (m_type == eDeformed)
316 {
317 // Allocate storage and compute number of points
318 for (i = 0; i < m_expDim; ++i)
319 {
320 ptsTgt *= keyTgt[i].GetNumPoints();
321 }
322 }
323
324 // Get derivative at geometry points
325 DerivStorage deriv = ComputeDeriv(keyTgt);
326
327 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
328 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
329 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
330
331 // Compute g_{ij} as t_i \cdot t_j and store in tmp
332 for (i = 0, l = 0; i < m_expDim; ++i)
333 {
334 for (j = 0; j < m_expDim; ++j, ++l)
335 {
336 for (k = 0; k < m_coordDim; ++k)
337 {
338 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
339 &tmp[l][0], 1, &tmp[l][0], 1);
340 }
341 }
342 }
343
344 Adjoint(tmp, gmat);
345
346 // Compute g = det(g_{ij}) (= Jacobian squared) and store
347 // temporarily in m_jac.
348 for (i = 0; i < m_expDim; ++i)
349 {
350 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
351 1, &jac[0], 1);
352 }
353
354 for (i = 0; i < m_expDim * m_expDim; ++i)
355 {
356 Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
357 }
358
359 return gmat;
360}

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 228 of file GeomFactors.cpp.

230{
231 ASSERTL1(keyTgt.size() == m_expDim,
232 "Dimension of target point distribution does not match "
233 "expansion dimension.");
234
235 int i = 0, j = 0, k = 0, l = 0;
236 int ptsTgt = 1;
237
238 if (m_type == eDeformed)
239 {
240 // Allocate storage and compute number of points
241 for (i = 0; i < m_expDim; ++i)
242 {
243 ptsTgt *= keyTgt[i].GetNumPoints();
244 }
245 }
246
247 // Get derivative at geometry points
248 DerivStorage deriv = ComputeDeriv(keyTgt);
249
250 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
251 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
252 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
253
254 // Compute g_{ij} as t_i \cdot t_j and store in tmp
255 for (i = 0, l = 0; i < m_expDim; ++i)
256 {
257 for (j = 0; j < m_expDim; ++j, ++l)
258 {
259 for (k = 0; k < m_coordDim; ++k)
260 {
261 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
262 &tmp[l][0], 1, &tmp[l][0], 1);
263 }
264 }
265 }
266
267 Adjoint(tmp, gmat);
268
269 // Compute g = det(g_{ij}) (= Jacobian squared) and store
270 // temporarily in m_jac.
271 for (i = 0; i < m_expDim; ++i)
272 {
273 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
274 1, &jac[0], 1);
275 }
276
277 // Compute the Jacobian = sqrt(g)
278 Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
279
280 return jac;
281}

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 444 of file GeomFactors.cpp.

449{
450 ASSERTL1(keyTgt.size() == m_expDim,
451 "Dimension of target point distribution does not match "
452 "expansion dimension.");
453
454 int i = 0, k = 0;
455 int ptsTgt = 1;
456 int nq = 1;
457
458 for (i = 0; i < m_expDim; ++i)
459 {
460 nq *= keyTgt[i].GetNumPoints();
461 }
462
463 if (m_type == eDeformed)
464 {
465 // Allocate storage and compute number of points
466 for (i = 0; i < m_expDim; ++i)
467 {
468 ptsTgt *= keyTgt[i].GetNumPoints();
469 }
470 }
471
472 // Get derivative at geometry points
473 DerivStorage deriv = ComputeDeriv(keyTgt);
474
475 // number of moving frames is requited to be 3, even for surfaces
476 int MFdim = 3;
477
478 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> MFtmp(MFdim);
479
480 // Compute g_{ij} as t_i \cdot t_j and store in tmp
481 for (i = 0; i < MFdim; ++i)
482 {
483 MFtmp[i] = Array<OneD, Array<OneD, NekDouble>>(m_coordDim);
484 for (k = 0; k < m_coordDim; ++k)
485 {
486 MFtmp[i][k] = Array<OneD, NekDouble>(nq);
487 }
488 }
489
490 // Compute g_{ij} as t_i \cdot t_j and store in tmp
491 for (i = 0; i < MFdim - 1; ++i)
492 {
493 for (k = 0; k < m_coordDim; ++k)
494 {
495 if (m_type == eDeformed)
496 {
497 Vmath::Vcopy(ptsTgt, &deriv[i][k][0], 1, &MFtmp[i][k][0], 1);
498 }
499 else
500 {
501 Vmath::Fill(nq, deriv[i][k][0], MFtmp[i][k], 1);
502 }
503 }
504 }
505
506 // Direction of MF1 is preserved: MF2 is considered in the same
507 // tangent plane as MF1. MF3 is computed by cross product of MF1
508 // and MF2. MF2 is consequently computed as the cross product of
509 // MF3 and MF1.
510 Array<OneD, Array<OneD, NekDouble>> PrincipleDir(m_coordDim);
511 for (k = 0; k < m_coordDim; k++)
512 {
513 PrincipleDir[k] = Array<OneD, NekDouble>(nq);
514 }
515
516 if (!(MMFdir == eLOCAL))
517 {
518 ComputePrincipleDirection(keyTgt, MMFdir, factors, PrincipleDir);
519 }
520
521 // MF3 = MF1 \times MF2
522 VectorCrossProd(MFtmp[0], MFtmp[1], MFtmp[2]);
523
524 // Normalizing MF3
525 VectorNormalise(MFtmp[2]);
526
527 if (!(MMFdir == eLOCAL))
528 {
529 Array<OneD, NekDouble> temp(nq, 0.0);
530
531 // Reorient MF1 along the PrincipleDir
532 for (i = 0; i < m_coordDim; ++i)
533 {
534 Vmath::Vvtvp(nq, MFtmp[2][i], 1, PrincipleDir[i], 1, temp, 1, temp,
535 1);
536 }
537 Vmath::Neg(nq, temp, 1);
538
539 // u2 = v2 - < u1 , v2 > ( u1 / < u1, u1 > )
540 for (i = 0; i < m_coordDim; ++i)
541 {
542 Vmath::Vvtvp(nq, temp, 1, MFtmp[2][i], 1, PrincipleDir[i], 1,
543 MFtmp[0][i], 1);
544 }
545 }
546
547 // Normalizing MF1
548 VectorNormalise(MFtmp[0]);
549
550 // MF2 = MF3 \times MF1
551 VectorCrossProd(MFtmp[2], MFtmp[0], MFtmp[1]);
552
553 // Normalizing MF2
554 VectorNormalise(MFtmp[1]);
555
556 for (i = 0; i < MFdim; ++i)
557 {
558 for (k = 0; k < m_coordDim; ++k)
559 {
560 Vmath::Vcopy(nq, &MFtmp[i][k][0], 1,
561 &movingframes[i * m_coordDim + k][0], 1);
562 }
563 }
564}
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.hpp:292

References ASSERTL1, ComputeDeriv(), ComputePrincipleDirection(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::eLOCAL, Nektar::VarcoeffHashingTest::factors, 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 713 of file GeomFactors.cpp.

718{
719 int nq = output[0].size();
720
721 output = Array<OneD, Array<OneD, NekDouble>>(m_coordDim);
722 for (int i = 0; i < m_coordDim; ++i)
723 {
724 output[i] = Array<OneD, NekDouble>(nq, 0.0);
725 }
726
727 // Construction of Connection
728 switch (MMFdir)
729 {
730 // projection to x-axis
731 case eTangentX:
732 {
733 Vmath::Fill(nq, 1.0, output[0], 1);
734 break;
735 }
736 case eTangentY:
737 {
738 Vmath::Fill(nq, 1.0, output[1], 1);
739 break;
740 }
741 case eTangentXY:
742 {
743 Vmath::Fill(nq, sqrt(2.0), output[0], 1);
744 Vmath::Fill(nq, sqrt(2.0), output[1], 1);
745 break;
746 }
747 case eTangentZ:
748 {
749 Vmath::Fill(nq, 1.0, output[2], 1);
750 break;
751 }
752 case eTangentCircular:
753 {
754 // Tangent direction depends on spatial location.
755 Array<OneD, Array<OneD, NekDouble>> x(m_coordDim);
756 for (int k = 0; k < m_coordDim; k++)
757 {
758 x[k] = Array<OneD, NekDouble>(nq);
759 }
760
761 // m_coords are StdExpansions which store the mapping
762 // between the std element and the local element. Bwd
763 // transforming the std element minimum basis gives a
764 // minimum physical basis for geometry. Need to then
765 // interpolate this up to the quadrature basis.
766 int nqtot_map = 1;
768 for (int i = 0; i < m_expDim; ++i)
769 {
770 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
771 nqtot_map *= map_points[i].GetNumPoints();
772 }
773 Array<OneD, NekDouble> tmp(nqtot_map);
774 for (int k = 0; k < m_coordDim; k++)
775 {
776 m_xmap->BwdTrans(m_coords[k], tmp);
777 Interp(map_points, tmp, keyTgt, x[k]);
778 }
779
780 // circular around the center of the domain
781 NekDouble radius, xc = 0.0, yc = 0.0, xdis, ydis;
782 NekDouble la, lb;
783
784 ASSERTL1(factors.size() >= 4, "factors is too short.");
785
786 la = factors[0];
787 lb = factors[1];
788 xc = factors[2];
789 yc = factors[3];
790
791 for (int i = 0; i < nq; i++)
792 {
793 xdis = x[0][i] - xc;
794 ydis = x[1][i] - yc;
795 radius = sqrt(xdis * xdis / la / la + ydis * ydis / lb / lb);
796 output[0][i] = ydis / radius;
797 output[1][i] = -1.0 * xdis / radius;
798 }
799 break;
800 }
802 {
803 // Tangent direction depends on spatial location.
804 Array<OneD, Array<OneD, NekDouble>> x(m_coordDim);
805 for (int k = 0; k < m_coordDim; k++)
806 {
807 x[k] = Array<OneD, NekDouble>(nq);
808 }
809
810 int nqtot_map = 1;
812 for (int i = 0; i < m_expDim; ++i)
813 {
814 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
815 nqtot_map *= map_points[i].GetNumPoints();
816 }
817 Array<OneD, NekDouble> tmp(nqtot_map);
818 for (int k = 0; k < m_coordDim; k++)
819 {
820 m_xmap->BwdTrans(m_coords[k], tmp);
821 Interp(map_points, tmp, keyTgt, x[k]);
822 }
823
824 // circular around the center of the domain
825 NekDouble xtan, ytan, mag;
826 for (int i = 0; i < nq; i++)
827 {
828 xtan = -1.0 * (x[1][i] * x[1][i] * x[1][i] + x[1][i]);
829 ytan = 2.0 * x[0][i];
830 mag = sqrt(xtan * xtan + ytan * ytan);
831 output[0][i] = xtan / mag;
832 output[1][i] = ytan / mag;
833 }
834 break;
835 }
837 {
838 // Tangent direction depends on spatial location.
839 Array<OneD, Array<OneD, NekDouble>> x(m_coordDim);
840 for (int k = 0; k < m_coordDim; k++)
841 {
842 x[k] = Array<OneD, NekDouble>(nq);
843 }
844
845 int nqtot_map = 1;
847 for (int i = 0; i < m_expDim; ++i)
848 {
849 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
850 nqtot_map *= map_points[i].GetNumPoints();
851 }
852 Array<OneD, NekDouble> tmp(nqtot_map);
853 for (int k = 0; k < m_coordDim; k++)
854 {
855 m_xmap->BwdTrans(m_coords[k], tmp);
856 Interp(map_points, tmp, keyTgt, x[k]);
857 }
858
859 // circular around the center of the domain
860 NekDouble xtan, ytan, mag;
861 for (int i = 0; i < nq; i++)
862 {
863 xtan = -2.0 * x[1][i] * x[1][i] * x[1][i] + x[1][i];
864 ytan = sqrt(3.0) * x[0][i];
865 mag = sqrt(xtan * xtan + ytan * ytan);
866 output[0][i] = xtan / mag;
867 output[1][i] = ytan / mag;
868 }
869 break;
870 }
871 default:
872 {
873 break;
874 }
875 }
876}
@ 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:285

References ASSERTL1, Nektar::SpatialDomains::eTangentCircular, Nektar::SpatialDomains::eTangentIrregular, Nektar::SpatialDomains::eTangentNonconvex, Nektar::SpatialDomains::eTangentX, Nektar::SpatialDomains::eTangentXY, Nektar::SpatialDomains::eTangentY, Nektar::SpatialDomains::eTangentZ, Nektar::VarcoeffHashingTest::factors, 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 307 of file GeomFactors.h.

308{
309 return m_coordDim;
310}

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 216 of file GeomFactors.h.

218{
219 return ComputeDeriv(tpoints);
220}

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 267 of file GeomFactors.h.

269{
270 auto x = m_derivFactorCache.find(keyTgt);
271
272 if (x != m_derivFactorCache.end())
273 {
274 return x->second;
275 }
276
277 m_derivFactorCache[keyTgt] = ComputeDerivFactors(keyTgt);
278
279 return m_derivFactorCache[keyTgt];
280}
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:146

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 252 of file GeomFactors.h.

254{
255 return ComputeGmat(keyTgt);
256}
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 298 of file GeomFactors.h.

299{
300 return m_type;
301}

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 327 of file GeomFactors.h.

328{
329 LibUtilities::PointsKeyVector ptsKeys = m_xmap->GetPointsKeys();
330 const Array<OneD, const NekDouble> jac = GetJac(ptsKeys);
331
332 size_t hash = 0;
334 if (m_type == eDeformed)
335 {
336 hash_range(hash, jac.begin(), jac.end());
337 }
338 else
339 {
340 hash_combine(hash, jac[0]);
341 }
342 return hash;
343}
const Array< OneD, const NekDouble > GetJac(const LibUtilities::PointsKeyVector &keyTgt)
Return the Jacobian of the mapping and cache the result.
Definition: GeomFactors.h:231
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:64
void hash_combine(std::size_t &seed)
Definition: HashUtils.hpp:44

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 231 of file GeomFactors.h.

233{
234 auto x = m_jacCache.find(keyTgt);
235
236 if (x != m_jacCache.end())
237 {
238 return x->second;
239 }
240
241 m_jacCache[keyTgt] = ComputeJac(keyTgt);
242
243 return m_jacCache[keyTgt];
244}
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:143

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 282 of file GeomFactors.h.

287{
288 ComputeMovingFrames(keyTgt, MMFdir, CircCentre, outarray);
289}
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 345 of file GeomFactors.h.

346{
347 return m_xmap;
348}

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 635 of file GeomFactors.cpp.

639{
640 ASSERTL1(src_points.size() == tgt_points.size(),
641 "Dimension of target point distribution does not match "
642 "expansion dimension.");
643
644 switch (m_expDim)
645 {
646 case 1:
647 LibUtilities::Interp1D(src_points[0], src, tgt_points[0], tgt);
648 break;
649 case 2:
650 LibUtilities::Interp2D(src_points[0], src_points[1], src,
651 tgt_points[0], tgt_points[1], tgt);
652 break;
653 case 3:
654 LibUtilities::Interp3D(src_points[0], src_points[1], src_points[2],
655 src, tgt_points[0], tgt_points[1],
656 tgt_points[2], tgt);
657 break;
658 }
659}
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:47
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:162
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:101

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 317 of file GeomFactors.h.

318{
319 return m_valid;
320}

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 918 of file GeomFactors.cpp.

922{
923 ASSERTL0(v1.size() == 3, "Input 1 has dimension not equal to 3.");
924 ASSERTL0(v2.size() == 3, "Input 2 has dimension not equal to 3.");
925 ASSERTL0(v3.size() == 3, "Output vector has dimension not equal to 3.");
926
927 int nq = v1[0].size();
928 Array<OneD, NekDouble> temp(nq);
929
930 Vmath::Vmul(nq, v1[2], 1, v2[1], 1, temp, 1);
931 Vmath::Vvtvm(nq, v1[1], 1, v2[2], 1, temp, 1, v3[0], 1);
932
933 Vmath::Vmul(nq, v1[0], 1, v2[2], 1, temp, 1);
934 Vmath::Vvtvm(nq, v1[2], 1, v2[0], 1, temp, 1, v3[1], 1);
935
936 Vmath::Vmul(nq, v1[1], 1, v2[0], 1, temp, 1);
937 Vmath::Vvtvm(nq, v1[0], 1, v2[1], 1, temp, 1, v3[2], 1);
938}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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

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 881 of file GeomFactors.cpp.

882{
883 int ndim = array.size();
884 ASSERTL0(ndim > 0, "Number of components must be > 0.");
885 for (int i = 1; i < ndim; ++i)
886 {
887 ASSERTL0(array[i].size() == array[0].size(),
888 "Array size mismatch in coordinates.");
889 }
890
891 int nq = array[0].size();
892 Array<OneD, NekDouble> norm(nq, 0.0);
893
894 // Compute the norm of each vector.
895 for (int i = 0; i < ndim; ++i)
896 {
897 Vmath::Vvtvp(nq, array[i], 1, array[i], 1, norm, 1, norm, 1);
898 }
899
900 Vmath::Vsqrt(nq, norm, 1, norm, 1);
901
902 // Normalise the vectors by the norm
903 for (int i = 0; i < ndim; ++i)
904 {
905 Vmath::Vdiv(nq, array[i], 1, norm, 1, array[i], 1);
906 }
907}

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 109 of file GeomFactors.cpp.

110{
111 if (!(lhs.m_type == rhs.m_type))
112 {
113 return false;
114 }
115
116 if (!(lhs.m_expDim == rhs.m_expDim))
117 {
118 return false;
119 }
120
121 if (!(lhs.m_coordDim == rhs.m_coordDim))
122 {
123 return false;
124 }
125
126 const Array<OneD, const NekDouble> jac_lhs =
127 lhs.ComputeJac(lhs.m_xmap->GetPointsKeys());
128 const Array<OneD, const NekDouble> jac_rhs =
129 rhs.ComputeJac(rhs.m_xmap->GetPointsKeys());
130 if (!(jac_lhs == jac_rhs))
131 {
132 return false;
133 }
134
135 return true;
136}

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 141 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 146 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 143 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 136 of file GeomFactors.h.

◆ m_type

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

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

Definition at line 127 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 133 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 139 of file GeomFactors.h.

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