Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public 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>

Collaboration diagram for Nektar::SpatialDomains::GeomFactors:
Collaboration graph
[legend]

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.
 GeomFactors (const GeomFactors &S)
 Copy constructor.
 ~GeomFactors ()
 Destructor.
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}$.
const Array< OneD, const
NekDouble
GetJac (const LibUtilities::PointsKeyVector &keyTgt)
 Return the Jacobian of the mapping and cache the result.
const Array< TwoD, const
NekDouble
GetGmat (const LibUtilities::PointsKeyVector &keyTgt)
 Return the Laplacian coefficients $g_{ij}$.
const Array< TwoD, const
NekDouble
GetDerivFactors (const LibUtilities::PointsKeyVector &keyTgt)
 Return the derivative of the reference coordinates with respect to the mapping, $\frac{\partial \xi_i}{\partial \chi_j}$.
GeomType GetGtype ()
 Returns whether the geometry is regular or deformed.
bool IsValid () const
 Determine if element is valid and not self-intersecting.
int GetCoordim () const
 Return the number of dimensions of the coordinate system.
size_t GetHash ()
 Computes a hash of this GeomFactors element.

Protected Attributes

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

Private Member Functions

void CheckIfValid ()
 Tests if the element is valid and not self-intersecting.
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.
Array< TwoD, NekDoubleComputeGmat (const LibUtilities::PointsKeyVector &keyTgt) const
 Computes the Laplacian coefficients $g_{ij}$.
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}$.
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.
void Adjoint (const Array< TwoD, const NekDouble > &src, Array< TwoD, NekDouble > &tgt) const
 Compute the transpose of the cofactors matrix.

Friends

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

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

Constructor & Destructor Documentation

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

References CheckIfValid().

:
m_type(gtype),
m_expDim(xmap->GetShapeDimension()),
m_coordDim(coordim),
m_valid(true),
m_xmap(xmap),
m_coords(coords)
{
}
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 109 of file GeomFactors.cpp.

:
m_type(S.m_type),
m_expDim(S.m_expDim),
m_coordDim(S.m_coordDim),
m_valid(S.m_valid),
m_xmap(S.m_xmap),
m_coords(S.m_coords)
{
}
Nektar::SpatialDomains::GeomFactors::~GeomFactors ( )

Destructor.

Definition at line 123 of file GeomFactors.cpp.

{
}

Member Function Documentation

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

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

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

{
ASSERTL1(src.num_elements() == tgt.num_elements(),
"Source matrix is of different size to destination"
"matrix for computing adjoint.");
int n = src[0].num_elements();
switch (m_expDim)
{
case 1:
Vmath::Fill (n, 1.0, &tgt[0][0], 1);
break;
case 2:
Vmath::Vcopy(n, &src[3][0], 1, &tgt[0][0], 1);
Vmath::Smul (n, -1.0, &src[1][0], 1, &tgt[1][0], 1);
Vmath::Smul (n, -1.0, &src[2][0], 1, &tgt[2][0], 1);
Vmath::Vcopy(n, &src[0][0], 1, &tgt[3][0], 1);
break;
case 3:
{
int a, b, c, d, e, i, j;
// Compute g^{ij} by computing Cofactors(g_ij)^T
for (i = 0; i < m_expDim; ++i)
{
for (j = 0; j < m_expDim; ++j)
{
a = ((i+1)%m_expDim) * m_expDim + ((j+1)%m_expDim);
b = ((i+1)%m_expDim) * m_expDim + ((j+2)%m_expDim);
c = ((i+2)%m_expDim) * m_expDim + ((j+1)%m_expDim);
d = ((i+2)%m_expDim) * m_expDim + ((j+2)%m_expDim);
e = j*m_expDim + i;
Vmath::Vvtvvtm(n, &src[a][0], 1, &src[d][0], 1,
&src[b][0], 1, &src[c][0], 1,
&tgt[e][0], 1);
}
}
break;
}
}
}
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 481 of file GeomFactors.cpp.

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

Referenced by GeomFactors().

{
// Jacobian test only makes sense when expdim = coorddim
// If one-dimensional then element is valid.
if (m_coordDim != m_expDim || m_expDim == 1)
{
m_valid = true;
return;
}
int nqtot = 1;
for (int i = 0; i < m_expDim; ++i)
{
p[i] = m_xmap->GetBasis(i)->GetPointsKey();
nqtot *= p[i].GetNumPoints();
}
int pts = (m_type == eRegular || m_type == eMovingRegular)
? 1 : nqtot;
Array<OneD, NekDouble> jac(pts, 0.0);
DerivStorage deriv = GetDeriv(p);
switch (m_expDim)
{
case 2:
{
Vmath::Vvtvvtm(pts, &deriv[0][0][0], 1, &deriv[1][1][0], 1,
&deriv[1][0][0], 1, &deriv[0][1][0], 1,
&jac[0], 1);
break;
}
case 3:
{
Array<OneD, NekDouble> tmp(pts, 0.0);
Vmath::Vvtvvtm(pts, &deriv[1][1][0], 1, &deriv[2][2][0], 1,
&deriv[2][1][0], 1, &deriv[1][2][0], 1,
&tmp[0], 1);
Vmath::Vvtvp (pts, &deriv[0][0][0], 1, &tmp[0], 1,
&jac[0], 1, &jac[0], 1);
Vmath::Vvtvvtm(pts, &deriv[2][1][0], 1, &deriv[0][2][0], 1,
&deriv[0][1][0], 1, &deriv[2][2][0], 1,
&tmp[0], 1);
Vmath::Vvtvp (pts, &deriv[1][0][0], 1, &tmp[0], 1,
&jac[0], 1, &jac[0], 1);
Vmath::Vvtvvtm(pts, &deriv[0][1][0], 1, &deriv[1][2][0], 1,
&deriv[1][1][0], 1, &deriv[0][2][0], 1,
&tmp[0], 1);
Vmath::Vvtvp (pts, &deriv[2][0][0], 1, &tmp[0], 1,
&jac[0], 1, &jac[0], 1);
break;
}
}
if (Vmath::Vmin(pts, &jac[0], 1) < 0)
{
m_valid = false;
}
}
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 170 of file GeomFactors.cpp.

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

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

{
ASSERTL1(keyTgt.size() == m_expDim,
"Dimension of target point distribution does not match "
"expansion dimension.");
int i = 0, j = 0;
int nqtot_map = 1;
int nqtot_tbasis = 1;
// Allocate storage and compute number of points
for (i = 0; i < m_expDim; ++i)
{
map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
nqtot_map *= map_points[i].GetNumPoints();
nqtot_tbasis *= keyTgt[i].GetNumPoints();
deriv[i] = Array<OneD, Array<OneD,NekDouble> >(m_coordDim);
d_map[i] = Array<OneD, Array<OneD,NekDouble> >(m_coordDim);
}
// Calculate local derivatives
for(i = 0; i < m_coordDim; ++i)
{
Array<OneD, NekDouble> tmp(nqtot_map);
// Transform from coefficient space to physical space
m_xmap->BwdTrans(m_coords[i], tmp);
// Allocate storage and take the derivative (calculated at the
// points as specified in 'Coords')
for (j = 0; j < m_expDim; ++j)
{
d_map[j][i] = Array<OneD,NekDouble>(nqtot_map);
deriv[j][i] = Array<OneD,NekDouble>(nqtot_tbasis);
m_xmap->StdPhysDeriv(j, tmp, d_map[j][i]);
}
}
for (i = 0; i < m_coordDim; ++i)
{
// Interpolate the derivatives:
// - from the points as defined in the mapping ('Coords')
// - to the points we at which we want to know the metrics
// ('tbasis')
bool same = true;
for (j = 0; j < m_expDim; ++j)
{
same = same && (map_points[j] == keyTgt[j]);
}
if( same )
{
for (j = 0; j < m_expDim; ++j)
{
deriv[j][i] = d_map[j][i];
}
}
else
{
for (j = 0; j < m_expDim; ++j)
{
Interp(map_points, d_map[j][i], keyTgt, deriv[j][i]);
}
}
}
return deriv;
}
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.

Definition at line 398 of file GeomFactors.cpp.

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

Referenced by GetDerivFactors().

{
ASSERTL1(keyTgt.size() == m_expDim,
"Dimension of target point distribution does not match "
"expansion dimension.");
int i = 0, j = 0, k = 0, l = 0;
int ptsTgt = 1;
if (m_type == eDeformed)
{
// Allocate storage and compute number of points
for (i = 0; i < m_expDim; ++i)
{
ptsTgt *= keyTgt[i].GetNumPoints();
}
}
// Get derivative at geometry points
DerivStorage deriv = ComputeDeriv(keyTgt);
Array<TwoD, NekDouble> tmp (m_expDim*m_expDim, ptsTgt, 0.0);
Array<TwoD, NekDouble> gmat(m_expDim*m_expDim, ptsTgt, 0.0);
Array<OneD, NekDouble> jac (ptsTgt, 0.0);
Array<TwoD, NekDouble> factors(m_expDim*m_coordDim, ptsTgt, 0.0);
// Compute g_{ij} as t_i \cdot t_j and store in tmp
for (i = 0, l = 0; i < m_expDim; ++i)
{
for (j = 0; j < m_expDim; ++j, ++l)
{
for (k = 0; k < m_coordDim; ++k)
{
Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
&deriv[j][k][0], 1,
&tmp[l][0], 1,
&tmp[l][0], 1);
}
}
}
Adjoint(tmp, gmat);
// Compute g = det(g_{ij}) (= Jacobian squared) and store
// temporarily in m_jac.
for (i = 0; i < m_expDim; ++i)
{
Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i*m_expDim][0], 1,
&jac[0], 1, &jac[0], 1);
}
for (i = 0; i < m_expDim*m_expDim; ++i)
{
Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
}
// Compute the Jacobian = sqrt(g)
Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
// Compute the derivative factors
for (k = 0, l = 0; k < m_coordDim; ++k)
{
for (j = 0; j < m_expDim; ++j, ++l)
{
for (i = 0; i < m_expDim; ++i)
{
Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
&gmat[m_expDim*i+j][0], 1,
&factors[l][0], 1,
&factors[l][0], 1);
}
}
}
return factors;
}
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 333 of file GeomFactors.cpp.

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

Referenced by GetGmat().

{
ASSERTL1(keyTgt.size() == m_expDim,
"Dimension of target point distribution does not match "
"expansion dimension.");
int i = 0, j = 0, k = 0, l = 0;
int ptsTgt = 1;
if (m_type == eDeformed)
{
// Allocate storage and compute number of points
for (i = 0; i < m_expDim; ++i)
{
ptsTgt *= keyTgt[i].GetNumPoints();
}
}
// Get derivative at geometry points
DerivStorage deriv = ComputeDeriv(keyTgt);
Array<TwoD, NekDouble> tmp (m_expDim*m_expDim, ptsTgt, 0.0);
Array<TwoD, NekDouble> gmat(m_expDim*m_expDim, ptsTgt, 0.0);
Array<OneD, NekDouble> jac (ptsTgt, 0.0);
// Compute g_{ij} as t_i \cdot t_j and store in tmp
for (i = 0, l = 0; i < m_expDim; ++i)
{
for (j = 0; j < m_expDim; ++j, ++l)
{
for (k = 0; k < m_coordDim; ++k)
{
Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
&deriv[j][k][0], 1,
&tmp[l][0], 1,
&tmp[l][0], 1);
}
}
}
Adjoint(tmp, gmat);
// Compute g = det(g_{ij}) (= Jacobian squared) and store
// temporarily in m_jac.
for (i = 0; i < m_expDim; ++i)
{
Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i*m_expDim][0], 1,
&jac[0], 1, &jac[0], 1);
}
for (i = 0; i < m_expDim*m_expDim; ++i)
{
Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
}
return gmat;
}
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 253 of file GeomFactors.cpp.

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

Referenced by GetJac(), and Nektar::SpatialDomains::operator==().

{
ASSERTL1(keyTgt.size() == m_expDim,
"Dimension of target point distribution does not match "
"expansion dimension.");
int i = 0, j = 0, k = 0, l = 0;
int ptsTgt = 1;
if (m_type == eDeformed)
{
// Allocate storage and compute number of points
for (i = 0; i < m_expDim; ++i)
{
ptsTgt *= keyTgt[i].GetNumPoints();
}
}
// Get derivative at geometry points
DerivStorage deriv = ComputeDeriv(keyTgt);
Array<TwoD, NekDouble> tmp (m_expDim*m_expDim, ptsTgt, 0.0);
Array<TwoD, NekDouble> gmat(m_expDim*m_expDim, ptsTgt, 0.0);
Array<OneD, NekDouble> jac (ptsTgt, 0.0);
// Compute g_{ij} as t_i \cdot t_j and store in tmp
for (i = 0, l = 0; i < m_expDim; ++i)
{
for (j = 0; j < m_expDim; ++j, ++l)
{
for (k = 0; k < m_coordDim; ++k)
{
Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
&deriv[j][k][0], 1,
&tmp[l][0], 1,
&tmp[l][0], 1);
}
}
}
Adjoint(tmp, gmat);
// Compute g = det(g_{ij}) (= Jacobian squared) and store
// temporarily in m_jac.
for (i = 0; i < m_expDim; ++i)
{
Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i*m_expDim][0], 1,
&jac[0], 1, &jac[0], 1);
}
// Compute the Jacobian = sqrt(g)
Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
return jac;
}
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 295 of file GeomFactors.h.

References m_coordDim.

{
return m_coordDim;
}
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 204 of file GeomFactors.h.

References ComputeDeriv().

Referenced by CheckIfValid().

{
return ComputeDeriv(tpoints);
}
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 260 of file GeomFactors.h.

References ComputeDerivFactors(), and m_derivFactorCache.

{
Array<TwoD, NekDouble> >::const_iterator x;
if ((x = m_derivFactorCache.find(keyTgt)) != m_derivFactorCache.end())
{
return x->second;
}
return m_derivFactorCache[keyTgt];
}
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 244 of file GeomFactors.h.

References ComputeGmat().

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

References m_type.

{
return m_type;
}
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 315 of file GeomFactors.h.

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

{
LibUtilities::PointsKeyVector ptsKeys = m_xmap->GetPointsKeys();
const Array<OneD, const NekDouble> jac = GetJac(ptsKeys);
size_t hash = 0;
boost::hash_combine(hash, (int)m_type);
boost::hash_combine(hash, m_expDim);
boost::hash_combine(hash, m_coordDim);
if (m_type == eDeformed)
{
boost::hash_range(hash, jac.begin(), jac.end());
}
else
{
boost::hash_combine(hash, jac[0]);
}
return hash;
}
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 220 of file GeomFactors.h.

References ComputeJac(), and m_jacCache.

Referenced by GetHash().

{
Array<OneD, NekDouble> >::const_iterator x;
if ((x = m_jacCache.find(keyTgt)) != m_jacCache.end())
{
return x->second;
}
m_jacCache[keyTgt] = ComputeJac(keyTgt);
return m_jacCache[keyTgt];
}
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 552 of file GeomFactors.cpp.

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

Referenced by ComputeDeriv().

{
ASSERTL1(src_points.size() == tgt_points.size(),
"Dimension of target point distribution does not match "
"expansion dimension.");
switch (m_expDim)
{
case 1:
LibUtilities::Interp1D(src_points[0], src,
tgt_points[0], tgt);
break;
case 2:
LibUtilities::Interp2D(src_points[0], src_points[1], src,
tgt_points[0], tgt_points[1], tgt);
break;
case 3:
LibUtilities::Interp3D(src_points[0], src_points[1],
src_points[2], src,
tgt_points[0], tgt_points[1],
tgt_points[2], tgt);
break;
}
}
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 305 of file GeomFactors.h.

References m_valid.

{
return m_valid;
}

Friends And Related Function Documentation

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

{
if(!(lhs.m_type == rhs.m_type))
{
return false;
}
if(!(lhs.m_expDim == rhs.m_expDim))
{
return false;
}
if(!(lhs.m_coordDim == rhs.m_coordDim))
{
return false;
}
const Array<OneD, const NekDouble> jac_lhs =
lhs.ComputeJac(lhs.m_xmap->GetPointsKeys());
const Array<OneD, const NekDouble> jac_rhs =
rhs.ComputeJac(rhs.m_xmap->GetPointsKeys());
if(!(jac_lhs == jac_rhs))
{
return false;
}
return true;
}

Member Data Documentation

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

Dimension of coordinate system.

Definition at line 137 of file GeomFactors.h.

Referenced by CheckIfValid(), ComputeDeriv(), ComputeDerivFactors(), ComputeGmat(), ComputeJac(), GetCoordim(), GetHash(), and Nektar::SpatialDomains::operator==().

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

Stores coordinates of the geometry.

Definition at line 143 of file GeomFactors.h.

Referenced by ComputeDeriv().

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

DerivFactors vector cache.

Definition at line 149 of file GeomFactors.h.

Referenced by GetDerivFactors().

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

Dimension of expansion.

Definition at line 135 of file GeomFactors.h.

Referenced by Adjoint(), CheckIfValid(), ComputeDeriv(), ComputeDerivFactors(), ComputeGmat(), ComputeJac(), GetHash(), Interp(), and Nektar::SpatialDomains::operator==().

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

Jacobian vector cache.

Definition at line 146 of file GeomFactors.h.

Referenced by GetJac().

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

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

Definition at line 133 of file GeomFactors.h.

Referenced by CheckIfValid(), ComputeDerivFactors(), ComputeGmat(), ComputeJac(), GetGtype(), GetHash(), and Nektar::SpatialDomains::operator==().

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

Validity of element (Jacobian positive)

Definition at line 139 of file GeomFactors.h.

Referenced by CheckIfValid(), and IsValid().

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

Stores information about the expansion.

Definition at line 141 of file GeomFactors.h.

Referenced by CheckIfValid(), ComputeDeriv(), GetHash(), and Nektar::SpatialDomains::operator==().