Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Typedefs | Enumerations | Functions | Variables
Nektar::StdRegions Namespace Reference

The namespace associated with the the StdRegions library (StdRegions introduction) More...

Classes

struct  IndexValue
class  IndexMapKey
class  StdExpansion
 The base class for all shapes. More...
class  StdExpansion0D
class  StdExpansion1D
class  StdExpansion2D
class  StdExpansion3D
class  StdHexExp
 Class representing a hexehedral element in reference space. More...
class  StdLinSysKey
class  StdMatrixKey
class  StdNodalPrismExp
class  StdNodalTetExp
class  StdNodalTriExp
class  StdPointExp
class  StdPrismExp
 Class representing a prismatic element in reference space. More...
struct  cmpop
class  StdPyrExp
class  StdQuadExp
class  StdSegExp
 Class representing a segment element in reference space. More...
class  StdTetExp
class  StdTriExp

Typedefs

typedef Array< OneD, IndexValueIndexMapValues
typedef boost::shared_ptr
< IndexMapKey
IndexMapKeySharedPtr
typedef boost::shared_ptr
< IndexMapValues
IndexMapValuesSharedPtr
typedef Array< OneD, Array
< OneD, NekDouble > > 
NormalVector
typedef boost::shared_ptr
< StdExpansion
StdExpansionSharedPtr
typedef std::vector
< StdExpansionSharedPtr
StdExpansionVector
typedef std::vector
< StdExpansionSharedPtr >
::iterator 
StdExpansionVectorIter
typedef boost::shared_ptr
< StdExpansion0D
StdExpansion0DSharedPtr
typedef boost::shared_ptr
< StdExpansion1D
StdExpansion1DSharedPtr
typedef boost::shared_ptr
< StdExpansion2D
StdExpansion2DSharedPtr
typedef boost::shared_ptr
< StdExpansion3D
StdExpansion3DSharedPtr
typedef boost::shared_ptr
< StdHexExp
StdHexExpSharedPtr
typedef boost::shared_ptr
< StdMatrixKey
StdMatrixKeySharedPtr
typedef boost::shared_ptr
< StdNodalPrismExp
StdNodalPrismExpSharedPtr
typedef boost::shared_ptr
< StdNodalTetExp
StdNodalTetExpSharedPtr
typedef boost::shared_ptr
< StdNodalTriExp
StdNodalTriExpSharedPtr
typedef boost::shared_ptr
< StdPointExp
StdPointExpSharedPtr
typedef boost::shared_ptr
< StdPrismExp
StdPrismExpSharedPtr
typedef boost::tuple< unsigned
int, unsigned int, unsigned
int, unsigned int > 
Mode
typedef boost::shared_ptr
< StdPyrExp
StdPyrExpSharedPtr
typedef boost::shared_ptr
< StdQuadExp
StdQuadExpSharedPtr
typedef std::map
< StdRegions::VarCoeffType,
Array< OneD, NekDouble > > 
VarCoeffMap
typedef std::map
< ConstFactorType, NekDouble
ConstFactorMap
typedef boost::shared_ptr
< StdSegExp
StdSegExpSharedPtr
typedef boost::shared_ptr
< StdTetExp
StdTetExpSharedPtr
typedef boost::shared_ptr
< StdTriExp
StdTriExpSharedPtr

Enumerations

enum  ElementType {
  eStdSegExp, eSegExp, eStdQuadExp, eStdTriExp,
  eStdNodalTriExp, eQuadExp, eTriExp, eNodalTriExp,
  eStdHexExp, eStdPrismExp, eStdPyrExp, eStdTetExp,
  eStdNodalTetExp, eHexExp, ePrismExp, ePyrExp,
  eTetExp, eNodalTetExp, SIZE_ElementType
}
enum  MatrixType {
  eMass, eInvMass, eLaplacian, eLaplacian00,
  eLaplacian01, eLaplacian02, eLaplacian10, eLaplacian11,
  eLaplacian12, eLaplacian20, eLaplacian21, eLaplacian22,
  eInvLaplacianWithUnityMean, eWeakDeriv0, eWeakDeriv1, eWeakDeriv2,
  eWeakDirectionalDeriv, eMassLevelCurvature, eLinearAdvectionReaction, eLinearAdvectionDiffusionReaction,
  eNBasisTrans, eInvNBasisTrans, eBwdTrans, eIProductWRTBase,
  eIProductWRTDerivBase0, eIProductWRTDerivBase1, eIProductWRTDerivBase2, eHelmholtz,
  eHybridDGHelmholtz, eInvHybridDGHelmholtz, eHybridDGHelmBndLam, eHybridDGLamToQ0,
  eHybridDGLamToQ1, eHybridDGLamToQ2, eHybridDGLamToU, eFwdTrans,
  ePreconR, ePreconRT, ePreconRMass, ePreconRTMass,
  ePreconLinearSpace, ePreconLinearSpaceMass, eInterpGauss, eGaussDG,
  SIZE_MatrixType
}
enum  VarCoeffType {
  eVarCoeffMass, eVarCoeffLaplacian, eVarCoeffWeakDeriv, eVarCoeffD00,
  eVarCoeffD11, eVarCoeffD22, eVarCoeffD01, eVarCoeffD02,
  eVarCoeffD12, eVarCoeffVelX, eVarCoeffVelY
}
enum  ConstFactorType {
  eFactorLambda, eFactorTau, eFactorTime, eFactorSVVCutoffRatio,
  eFactorSVVDiffCoeff, eFactorGaussVertex, eFactorGaussEdge
}
enum  IndexMapType {
  eEdgeToElement, eFaceToElement, eEdgeInterior, eFaceInterior,
  eBoundary, eVertex
}
enum  Orientation {
  eNoOrientation, eFwd, eBwd, eForwards,
  eBackwards, eDir1FwdDir1_Dir2FwdDir2, eDir1FwdDir1_Dir2BwdDir2, eDir1BwdDir1_Dir2FwdDir2,
  eDir1BwdDir1_Dir2BwdDir2, eDir1FwdDir2_Dir2FwdDir1, eDir1FwdDir2_Dir2BwdDir1, eDir1BwdDir2_Dir2FwdDir1,
  eDir1BwdDir2_Dir2BwdDir1, SIZE_Orientation
}

Functions

bool operator< (const IndexMapKey &lhs, const IndexMapKey &rhs)
bool operator== (const IndexMapKey &lhs, const IndexMapKey &rhs)
std::ostream & operator<< (std::ostream &os, const IndexMapKey &rhs)
template<typename T >
NekVector< T > GetColumn (const NekMatrix< T > &matA, int n)
NekMatrix< NekDouble > & SetColumn (NekMatrix< NekDouble > &matA, int n, const NekVector< NekDouble > &x)
NekVector< NekDoubleGetE (int rows, int n)
NekMatrix< NekDoubleInvert (const NekMatrix< NekDouble > &matA)
NekMatrix< NekDoubleGetTranspose (const NekMatrix< NekDouble > &matA)
int GetSize (const ConstArray< OneD, NekDouble > &x)
int GetSize (const NekVector< NekDouble > &x)
NekVector< NekDoubleToVector (const ConstArray< OneD, NekDouble > &x)
Array< OneD, NekDoubleToArray (const NekVector< NekDouble > &x)
NekVector< NekDoubleHadamard (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y)
NekVector< NekDoubleVectorPower (const NekVector< NekDouble > &x, NekDouble p)
std::string MatrixToString (const NekMatrix< NekDouble > &A, int precision, NekDouble expSigFigs)
std::string VectorToString (const NekVector< NekDouble > &v, int precision, NekDouble expSigFigs)
int GetTriNumPoints (int degree)
int GetDegree (int nBasisFunctions)
int GetTetNumPoints (int degree)
int GetTetDegree (int nBasisFunc)
NekDouble MakeRound (NekDouble x)
NekMatrix< NekDoubleGetMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, int degree)
NekMatrix< NekDoubleGetMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z, int degree)
NekMatrix< NekDoubleGetMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z)
NekMatrix< NekDoubleGetMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y)
NekMatrix< NekDoubleGetXDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, int degree)
NekMatrix< NekDoubleGetXDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y)
NekMatrix< NekDoubleGetTetXDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z, int degree)
NekMatrix< NekDoubleGetTetXDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z)
NekMatrix< NekDoubleGetTetYDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z, int degree)
NekMatrix< NekDoubleGetTetYDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z)
NekMatrix< NekDoubleGetTetZDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z, int degree)
NekMatrix< NekDoubleGetTetZDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, const NekVector< NekDouble > &z)
NekMatrix< NekDoubleGetYDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y, int degree)
NekMatrix< NekDoubleGetYDerivativeOfMonomialVandermonde (const NekVector< NekDouble > &x, const NekVector< NekDouble > &y)
NekVector< NekDoubleGetIntegralOfMonomialVandermonde (int degree)
int GetSize (const Array< OneD, const NekDouble > &x)
NekVector< NekDoubleToVector (const Array< OneD, const NekDouble > &x)
bool operator< (const StdMatrixKey &lhs, const StdMatrixKey &rhs)
bool operator== (const StdMatrixKey &lhs, const StdMatrixKey &rhs)
std::ostream & operator<< (std::ostream &os, const StdMatrixKey &rhs)
template<class InputIterator , class EqualityComparable >
InputIterator find (InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)

Variables

const int g_shapenverts [LibUtilities::SIZE_ShapeType] = {0,2,3,4,4,5,6,8}
const int g_shapenedges [LibUtilities::SIZE_ShapeType] = {0,1,3,4,6,8,9,12}
const int g_shapenfaces [LibUtilities::SIZE_ShapeType] = {0,0,0,0,4,5,5,6}
const char *const ElementTypeMap []
const char *const MatrixTypeMap []
const char *const VarCoeffTypeMap []
static VarCoeffMap NullVarCoeffMap
const char *const ConstFactorTypeMap []
static ConstFactorMap NullConstFactorMap
const char *const IndexMapTypeMap []
const char *const OrientationMap []

Detailed Description

The namespace associated with the the StdRegions library (StdRegions introduction)

Typedef Documentation

Definition at line 246 of file StdRegions.hpp.

Definition at line 125 of file IndexMapKey.h.

Definition at line 55 of file IndexMapKey.h.

Definition at line 126 of file IndexMapKey.h.

typedef boost::tuple< unsigned int, unsigned int, unsigned int, unsigned int> Nektar::StdRegions::Mode

Definition at line 49 of file StdPyrExp.h.

Definition at line 59 of file StdExpansion.h.

Definition at line 96 of file StdExpansion0D.h.

Definition at line 113 of file StdExpansion1D.h.

Definition at line 198 of file StdExpansion2D.h.

Definition at line 49 of file StdExpansion3D.h.

Definition at line 1767 of file StdExpansion.h.

Definition at line 1768 of file StdExpansion.h.

Definition at line 1769 of file StdExpansion.h.

Definition at line 285 of file StdHexExp.h.

Definition at line 196 of file StdMatrixKey.h.

Definition at line 141 of file StdNodalPrismExp.h.

Definition at line 142 of file StdNodalTetExp.h.

Definition at line 184 of file StdNodalTriExp.h.

Definition at line 123 of file StdPointExp.h.

Definition at line 244 of file StdPrismExp.h.

Definition at line 254 of file StdPyrExp.h.

Definition at line 253 of file StdQuadExp.h.

Definition at line 47 of file StdSegExp.h.

Definition at line 250 of file StdTetExp.h.

Definition at line 258 of file StdTriExp.h.

Definition at line 223 of file StdRegions.hpp.

Enumeration Type Documentation

Enumerator:
eFactorLambda 
eFactorTau 
eFactorTime 
eFactorSVVCutoffRatio 
eFactorSVVDiffCoeff 
eFactorGaussVertex 
eFactorGaussEdge 

Definition at line 226 of file StdRegions.hpp.

Enumerator:
eStdSegExp 
eSegExp 
eStdQuadExp 
eStdTriExp 
eStdNodalTriExp 
eQuadExp 
eTriExp 
eNodalTriExp 
eStdHexExp 
eStdPrismExp 
eStdPyrExp 
eStdTetExp 
eStdNodalTetExp 
eHexExp 
ePrismExp 
ePyrExp 
eTetExp 
eNodalTetExp 
SIZE_ElementType 

Definition at line 52 of file StdRegions.hpp.

Enumerator:
eEdgeToElement 
eFaceToElement 
eEdgeInterior 
eFaceInterior 
eBoundary 
eVertex 

Definition at line 249 of file StdRegions.hpp.

Todo:
we need to tidy up matrix construction approach probably using a factory type approach
Enumerator:
eMass 
eInvMass 
eLaplacian 
eLaplacian00 
eLaplacian01 
eLaplacian02 
eLaplacian10 
eLaplacian11 
eLaplacian12 
eLaplacian20 
eLaplacian21 
eLaplacian22 
eInvLaplacianWithUnityMean 
eWeakDeriv0 
eWeakDeriv1 
eWeakDeriv2 
eWeakDirectionalDeriv 
eMassLevelCurvature 
eLinearAdvectionReaction 
eLinearAdvectionDiffusionReaction 
eNBasisTrans 
eInvNBasisTrans 
eBwdTrans 
eIProductWRTBase 
eIProductWRTDerivBase0 
eIProductWRTDerivBase1 
eIProductWRTDerivBase2 
eHelmholtz 
eHybridDGHelmholtz 
eInvHybridDGHelmholtz 
eHybridDGHelmBndLam 
eHybridDGLamToQ0 
eHybridDGLamToQ1 
eHybridDGLamToQ2 
eHybridDGLamToU 
eFwdTrans 
ePreconR 
ePreconRT 
ePreconRMass 
ePreconRTMass 
ePreconLinearSpace 
ePreconLinearSpaceMass 
eInterpGauss 
eGaussDG 
SIZE_MatrixType 

Definition at line 101 of file StdRegions.hpp.

Enumerator:
eNoOrientation 
eFwd 
eBwd 
eForwards 
eBackwards 
eDir1FwdDir1_Dir2FwdDir2 
eDir1FwdDir1_Dir2BwdDir2 
eDir1BwdDir1_Dir2FwdDir2 
eDir1BwdDir1_Dir2BwdDir2 
eDir1FwdDir2_Dir2FwdDir1 
eDir1FwdDir2_Dir2BwdDir1 
eDir1BwdDir2_Dir2FwdDir1 
eDir1BwdDir2_Dir2BwdDir1 
SIZE_Orientation 

Definition at line 269 of file StdRegions.hpp.

{
eDir1FwdDir1_Dir2FwdDir2, // These flags are interpreted as
eDir1FwdDir1_Dir2BwdDir2, // taking the second direction to the
eDir1BwdDir1_Dir2FwdDir2, // first direction. So Dir1FwdDir2 takes
eDir1BwdDir1_Dir2BwdDir2, // direction 2 and makes it backward
eDir1FwdDir2_Dir2FwdDir1, // to direction 1 in the mapped face.
eDir1FwdDir2_Dir2BwdDir1, // Note be careful not to flip this
eDir1BwdDir2_Dir2FwdDir1, // convention especially when using
eDir1BwdDir2_Dir2BwdDir1, // transposed mappings.
};
Enumerator:
eVarCoeffMass 
eVarCoeffLaplacian 
eVarCoeffWeakDeriv 
eVarCoeffD00 
eVarCoeffD11 
eVarCoeffD22 
eVarCoeffD01 
eVarCoeffD02 
eVarCoeffD12 
eVarCoeffVelX 
eVarCoeffVelY 

Definition at line 195 of file StdRegions.hpp.

Function Documentation

template<class InputIterator , class EqualityComparable >
InputIterator Nektar::StdRegions::find ( InputIterator  first,
InputIterator  last,
InputIterator  startingpoint,
const EqualityComparable &  value 
)

Definition at line 310 of file StdRegions.hpp.

Referenced by Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), Nektar::SpatialDomains::MeshGraph::GetCompositeList(), Nektar::LibUtilities::AnalyticExpressionEvaluator::GetConstant(), Nektar::LibUtilities::SessionReader::GetGlobalSysSolnInfo(), Nektar::LibUtilities::SessionReader::GetSolverInfoAsEnum(), Nektar::LibUtilities::SessionReader::GetValueAsEnum(), Nektar::SpatialDomains::PointGeom::IsElmtConnected(), Nektar::Utilities::operator==(), Nektar::Utilities::Face::operator==(), Nektar::Utilities::ProcessCyl::Process(), Nektar::SpatialDomains::Domain::Read(), Nektar::SpatialDomains::BoundaryConditions::ReadBoundaryConditions(), Nektar::SpatialDomains::MeshGraph::ReadCurves(), Nektar::SpatialDomains::MeshGraph::ReadExpansions(), Nektar::LibUtilities::SessionReader::ReadGlobalSysSolnInfo(), Nektar::LibUtilities::SessionReader::ReadVariables(), Nektar::LibUtilities::SessionReader::RegisterEnumValue(), Nektar::StdRegions::StdPyrExp::StdPyrExp(), Nektar::SolverUtils::FilterAeroForces::v_Initialise(), Nektar::SpatialDomains::QuadGeom::v_IsElmtConnected(), Nektar::SpatialDomains::TriGeom::v_IsElmtConnected(), Nektar::SpatialDomains::SegGeom::v_IsElmtConnected(), Nektar::SpatialDomains::Geometry3D::v_IsElmtConnected(), and Nektar::Utilities::OutputNekpp::WriteXmlComposites().

{
InputIterator val;
if(startingpoint == first)
{
val = find(first,last,value);
}
else
{
val = find(startingpoint,last,value);
if(val == last)
{
val = find(first,startingpoint,value);
if(val == startingpoint)
{
val = last;
}
}
}
return val;
}
template<typename T >
NekVector< T > Nektar::StdRegions::GetColumn ( const NekMatrix< T > &  matA,
int  n 
)

Definition at line 52 of file StdExpUtil.cpp.

{
NekVector<T> v(matA.GetRows());
for( int i=0; i<matA.GetRows(); ++i )
{
v[i] = matA(i,n);
}
return v;
}
int Nektar::StdRegions::GetDegree ( int  nBasisFunctions)

Definition at line 207 of file StdExpUtil.cpp.

References ASSERTL1, GetTriNumPoints(), and MakeRound().

Referenced by GetMonomialVandermonde(), GetXDerivativeOfMonomialVandermonde(), and GetYDerivativeOfMonomialVandermonde().

{
int degree = int(MakeRound((-3.0 + sqrt(1.0 + 8.0*nBasisFunctions))/2.0));
// TODO: Find out why ASSERTL0 and ASSERTL1 don't work
ASSERTL1( GetTriNumPoints(degree) == nBasisFunctions, "The number of points defines an expansion of fractional degree, which is not supported." );
return degree;
}
NekVector< NekDouble > Nektar::StdRegions::GetE ( int  rows,
int  n 
)

Definition at line 72 of file StdExpUtil.cpp.

Referenced by Invert().

{
NekVector<NekDouble> e(rows, 0.0);
e(n) = 1;
return e;
}
NekVector< NekDouble > Nektar::StdRegions::GetIntegralOfMonomialVandermonde ( int  degree)

Definition at line 497 of file StdExpUtil.cpp.

{
int cols = (degree + 1) * (degree + 2) / 2;
NekVector<NekDouble>integralMVandermonde(cols, 0.0);
for(int d=0, n=0; d <= degree; ++d)
{
for(int p=0; p <= d; ++p, ++n)
{
int q = d - p;
int sp = 1 - 2 * (p % 2);
int sq = 1 - 2 * (q % 2);
integralMVandermonde(n) = NekDouble(sp * (p+1) + sq * (q+1) + sp * sq * (p+q+2)) / ((p+1) * (q+1) * (p+q+2));
}
}
return integralMVandermonde;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
int  degree 
)

Definition at line 236 of file StdExpUtil.cpp.

References GetSize().

Referenced by GetMonomialVandermonde().

{
int rows = GetSize(x);
int cols = (degree + 1)*(degree + 2)/2;
NekMatrix<NekDouble>matV(rows,cols, 0.0);
for(int d=0, n=0; d <= degree; ++d)
{
for(int p=0; p <= d; ++p, ++n)
{
int q = d - p;
for(int i=0; i<rows; ++i)
{
matV(i, n) = pow(x[i], p) * pow(y[i], q);
}
}
}
return matV;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z,
int  degree 
)

Definition at line 257 of file StdExpUtil.cpp.

References GetSize().

{
int rows = GetSize(x);
int cols = (degree + 1) * (degree + 2) * (degree + 3)/6;
NekMatrix<NekDouble> matV(rows, cols, 0.0);
for(int d=0, n=0; d<=degree; ++d)
{
for(int p=0; p <= d; ++p)
{
for(int q=0; q <= d - p; ++q, ++n)
{
int r = d - p - q;
// Set n-th column of V to the monomial vector
for(int i=0; i<rows; ++i)
{
matV(i,n) = pow(x[i],p) * pow(y[i],q) * pow(z[i],r);
}
}
}
}
return matV;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z 
)

Definition at line 284 of file StdExpUtil.cpp.

References GetMonomialVandermonde(), GetSize(), and GetTetDegree().

{
int rows = GetSize(x);
int degree = GetTetDegree(rows);
return GetMonomialVandermonde(x, y, z, degree);
}
NekMatrix< NekDouble > Nektar::StdRegions::GetMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y 
)

Definition at line 291 of file StdExpUtil.cpp.

References GetDegree(), GetMonomialVandermonde(), and GetSize().

{
int rows = GetSize(x);
int degree = GetDegree(rows);
return GetMonomialVandermonde(x, y, degree);
}
int Nektar::StdRegions::GetSize ( const Array< OneD, const NekDouble > &  x)

Definition at line 111 of file NodalUtil.cpp.

Referenced by Nektar::NekMatrix< DataType, StandardMatrixTag >::CalculateIndex(), Nektar::LibUtilities::DubinerPoly(), Nektar::LibUtilities::DubinerPolyXDerivative(), Nektar::LibUtilities::DubinerPolyYDerivative(), Nektar::LibUtilities::GetInterpolationMatrix(), Nektar::LibUtilities::GetMonomialVandermonde(), Nektar::LibUtilities::GetTetInterpolationMatrix(), Nektar::LibUtilities::GetTetVandermonde(), Nektar::LibUtilities::GetTetXDerivativeMatrix(), Nektar::LibUtilities::GetTetXDerivativeOfMonomialVandermonde(), Nektar::LibUtilities::GetTetYDerivativeOfMonomialVandermonde(), Nektar::LibUtilities::GetTetZDerivativeMatrix(), Nektar::LibUtilities::GetTetZDerivativeOfMonomialVandermonde(), Nektar::LibUtilities::GetVandermonde(), Nektar::LibUtilities::GetVandermondeForTetXDerivative(), Nektar::LibUtilities::GetVandermondeForTetYDerivative(), Nektar::LibUtilities::GetVandermondeForTetZDerivative(), Nektar::LibUtilities::GetVandermondeForXDerivative(), Nektar::LibUtilities::GetVandermondeForYDerivative(), Nektar::LibUtilities::GetXDerivativeMatrix(), Nektar::LibUtilities::GetXDerivativeOfMonomialVandermonde(), Nektar::LibUtilities::GetYDerivativeMatrix(), Nektar::LibUtilities::GetYDerivativeOfMonomialVandermonde(), Nektar::LibUtilities::Hadamard(), Nektar::LibUtilities::JacobiPoly(), Nektar::LibUtilities::JacobiPolyDerivative(), Nektar::LibUtilities::LegendrePoly(), Nektar::LibUtilities::LegendrePolyDerivative(), Nektar::LibUtilities::TetrahedralBasis(), Nektar::LibUtilities::TetXDerivative(), Nektar::LibUtilities::TetYDerivative(), Nektar::LibUtilities::TetZDerivative(), Nektar::LibUtilities::ToArray(), ToVector(), and Nektar::LibUtilities::VectorPower().

{
return x.num_elements();
}
int Nektar::StdRegions::GetSize ( const ConstArray< OneD, NekDouble > &  x)
int Nektar::StdRegions::GetSize ( const NekVector< NekDouble > &  x)

Definition at line 116 of file StdExpUtil.cpp.

References Nektar::NekVector< DataType >::GetRows().

{
return x.GetRows();
}
int Nektar::StdRegions::GetTetDegree ( int  nBasisFunc)

Definition at line 222 of file StdExpUtil.cpp.

References ASSERTL1, GetTetNumPoints(), and MakeRound().

Referenced by GetMonomialVandermonde(), GetTetXDerivativeOfMonomialVandermonde(), GetTetYDerivativeOfMonomialVandermonde(), and GetTetZDerivativeOfMonomialVandermonde().

{
NekDouble eq = pow( 81.0 * nBasisFunc + 3.0 * sqrt(-3.0 + 729.0 * nBasisFunc * nBasisFunc), 1.0/3.0);
int degree = int(MakeRound(eq/3.0 + 1.0/eq - 1.0)) - 1;
ASSERTL1( GetTetNumPoints(degree) == nBasisFunc, "The number of points defines an expansion of fractional degree, which is not supported." );
return degree;
}
int Nektar::StdRegions::GetTetNumPoints ( int  degree)

Definition at line 216 of file StdExpUtil.cpp.

Referenced by GetTetDegree().

{
return (degree+1) * (degree+2) * (degree+3) / 6;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetTetXDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z,
int  degree 
)

Definition at line 336 of file StdExpUtil.cpp.

References GetSize().

Referenced by GetTetXDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int cols = (degree + 1) * (degree + 2) * (degree + 3) / 6;
NekMatrix<NekDouble> matVx(rows, cols, 0.0);
for(int d=0, n=0; d<=degree; ++d)
{
for(int p=0; p <= d; ++p)
{
for(int q=0; q <= d - p; ++q, ++n)
{
int r = d - p - q;
if(p > 0)
{
// Set n-th column of V to the monomial vector
for(int i=0; i<rows; ++i)
{
matVx(i,n) = p * pow(x[i],p-1.0) * pow(y[i],q) * pow(z[i],r);
}
}
else{
for(int j=0; j<rows; ++j)
{
matVx(j, n) = 0.0;
}
}
}
}
}
return matVx;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetTetXDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z 
)

Definition at line 369 of file StdExpUtil.cpp.

References GetSize(), GetTetDegree(), and GetTetXDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int degree = GetTetDegree(rows);
return GetTetXDerivativeOfMonomialVandermonde(x, y, z, degree);
}
NekMatrix< NekDouble > Nektar::StdRegions::GetTetYDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z,
int  degree 
)

Definition at line 377 of file StdExpUtil.cpp.

References GetSize().

Referenced by GetTetYDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int cols = (degree + 1) * (degree + 2) * (degree + 3) / 6;
NekMatrix<NekDouble> matVy(rows, cols, 0.0);
for(int d=0, n=0; d<=degree; ++d)
{
for(int p=0; p <= d; ++p)
{
for(int q=0; q <= d - p; ++q, ++n)
{
int r = d - p - q;
if(q > 0)
{
// Set n-th column of V to the monomial vector
for(int i=0; i<rows; ++i)
{
matVy(i,n) = q * pow(x[i],p) * pow(y[i],q-1.0) * pow(z[i],r);
}
}
else
{
for(int j=0; j<rows; ++j)
{
matVy(j, n) = 0.0;
}
}
}
}
}
return matVy;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetTetYDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z 
)

Definition at line 411 of file StdExpUtil.cpp.

References GetSize(), GetTetDegree(), and GetTetYDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int degree = GetTetDegree(rows);
return GetTetYDerivativeOfMonomialVandermonde(x, y, z, degree);
}
NekMatrix< NekDouble > Nektar::StdRegions::GetTetZDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z,
int  degree 
)

Definition at line 419 of file StdExpUtil.cpp.

References GetSize().

Referenced by GetTetZDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int cols = (degree + 1) * (degree + 2) * (degree + 3) / 6;
NekMatrix<NekDouble> matVz(rows, cols, 0.0);
for(int d=0, n=0; d<=degree; ++d)
{
for(int p=0; p <= d; ++p)
{
for(int q=0; q <= d - p; ++q, ++n)
{
int r = d - p - q;
if(r > 0)
{
// Set n-th column of V to the monomial vector
for(int i=0; i<rows; ++i)
{
matVz(i,n) = r * pow(x[i],p) * pow(y[i],q) * pow(z[i],r-1.0);
}
}
else{
for(int j=0; j<rows; ++j)
{
matVz(j, n) = 0.0;
}
}
}
}
}
return matVz;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetTetZDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
const NekVector< NekDouble > &  z 
)

Definition at line 452 of file StdExpUtil.cpp.

References GetSize(), GetTetDegree(), and GetTetZDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int degree = GetTetDegree(rows);
return GetTetZDerivativeOfMonomialVandermonde(x, y, z, degree);
}
NekMatrix< NekDouble > Nektar::StdRegions::GetTranspose ( const NekMatrix< NekDouble > &  matA)

Definition at line 96 of file StdExpUtil.cpp.

{
int rows = matA.GetRows(), columns = matA.GetColumns();
NekMatrix<NekDouble> matX(rows,columns);
for( int i=0; i<rows; ++i )
{
for( int j=0; j<columns; ++j )
{
matX(j,i) = matA(i,j);
}
}
return matX;
}
int Nektar::StdRegions::GetTriNumPoints ( int  degree)

Definition at line 202 of file StdExpUtil.cpp.

Referenced by GetDegree().

{
return (degree+1) * (degree+2) / 2;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetXDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
int  degree 
)

Definition at line 298 of file StdExpUtil.cpp.

References GetSize().

Referenced by GetXDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int cols = (degree + 1) * (degree + 2) / 2;
NekMatrix<NekDouble> matVx(rows, cols, 0.0);
for(int d=0, n=0; d <= degree; ++d)
{
for(int p=0; p <= d; ++p, ++n)
{
int q = d - p;
if(p > 0)
{
for(int i=0; i<rows; ++i)
{
matVx(i, n) = p * pow(x[i], p-1.0) * pow(y[i],q);
}
}
else
{
for(int j=0; j<rows; ++j)
{
matVx(j, n) = 0.0;
}
}
}
}
return matVx;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetXDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y 
)

Definition at line 329 of file StdExpUtil.cpp.

References GetDegree(), GetSize(), and GetXDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int degree = GetDegree(rows);
}
NekMatrix< NekDouble > Nektar::StdRegions::GetYDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y,
int  degree 
)

Definition at line 460 of file StdExpUtil.cpp.

References GetSize().

Referenced by GetYDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int cols = (degree + 1) * (degree + 2) / 2;
NekMatrix<NekDouble> matVy(rows, cols, 0.0);
for(int d=0, n=0; d <= degree; ++d)
{
for(int p=0; p <= d; ++p, ++n)
{
int q = d - p;
if(q > 0)
{
for(int i=0; i<rows; ++i)
{
matVy(i, n) = q * pow(x[i], p) * pow(y[i], q-1.0);
}
}
else
{
for(int j=0; j<rows; ++j)
{
matVy(j, n) = 0.0;
}
}
}
}
return matVy;
}
NekMatrix< NekDouble > Nektar::StdRegions::GetYDerivativeOfMonomialVandermonde ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y 
)

Definition at line 490 of file StdExpUtil.cpp.

References GetDegree(), GetSize(), and GetYDerivativeOfMonomialVandermonde().

{
int rows = GetSize(x);
int degree = GetDegree(rows);
}
NekVector< NekDouble > Nektar::StdRegions::Hadamard ( const NekVector< NekDouble > &  x,
const NekVector< NekDouble > &  y 
)

Definition at line 131 of file StdExpUtil.cpp.

References GetSize().

{
int size = GetSize(x);
NekVector<NekDouble> z(size);
for( int i=0; i<size; ++i )
{
z[i] = x[i] * y[i];
}
return z;
}
NekMatrix< NekDouble > Nektar::StdRegions::Invert ( const NekMatrix< NekDouble > &  matA)

Definition at line 79 of file StdExpUtil.cpp.

References GetE(), SetColumn(), and Nektar::LinearSystem::Solve().

{
int rows = matA.GetRows(), columns = matA.GetColumns();
NekMatrix<NekDouble> matX(rows,columns);
// The linear system solver
LinearSystem<NekMatrix<NekDouble> > matL( SharedNekMatrixPtr(new NekMatrix<NekDouble>(matA)) );
// Solve each column for the identity matrix
for( int j=0; j<columns; ++j )
{
SetColumn( matX, j, matL.Solve( GetE(rows,j) ) );
}
return matX;
}
NekDouble Nektar::StdRegions::MakeRound ( NekDouble  x)

Definition at line 231 of file StdExpUtil.cpp.

Referenced by GetDegree(), GetTetDegree(), MatrixToString(), and VectorToString().

{
return floor(x + 0.5);
}
std::string Nektar::StdRegions::MatrixToString ( const NekMatrix< NekDouble > &  A,
int  precision,
NekDouble  expSigFigs 
)

Definition at line 157 of file StdExpUtil.cpp.

References MakeRound().

{
stringstream s;
s << setprecision(precision);
int M = int(A.GetRows()), N = int(A.GetColumns());
for(int i=0; i<M; ++i )
{
for(int j=0; j<N; ++j)
{
NekDouble a = MakeRound(expSigFigs * A(i, j)) / expSigFigs;
s << setw(7) << right << a;
if( j < N-1 )
{
s << ", ";
}
}
if( i < M-1 )
{
s << "\n";
}
}
return s.str();
}
bool Nektar::StdRegions::operator< ( const IndexMapKey &  lhs,
const IndexMapKey &  rhs 
)

Definition at line 90 of file IndexMapKey.cpp.

References Nektar::StdRegions::IndexMapKey::m_entityID, Nektar::StdRegions::IndexMapKey::m_indexMapType, Nektar::StdRegions::IndexMapKey::m_orientation, Nektar::StdRegions::IndexMapKey::m_p, Nektar::StdRegions::IndexMapKey::m_q, Nektar::StdRegions::IndexMapKey::m_r, and Nektar::StdRegions::IndexMapKey::m_shapeType.

{
if(lhs.m_indexMapType < rhs.m_indexMapType)
{
return true;
}
if(lhs.m_indexMapType > rhs.m_indexMapType)
{
return false;
}
if(lhs.m_shapeType < rhs.m_shapeType)
{
return true;
}
if(lhs.m_shapeType > rhs.m_shapeType)
{
return false;
}
if(lhs.m_p < rhs.m_p)
{
return true;
}
if(lhs.m_p > rhs.m_p)
{
return false;
}
if(lhs.m_q < rhs.m_q)
{
return true;
}
if(lhs.m_q > rhs.m_q)
{
return false;
}
if(lhs.m_r < rhs.m_r)
{
return true;
}
if(lhs.m_r > rhs.m_r)
{
return false;
}
if(lhs.m_entityID < rhs.m_entityID)
{
return true;
}
if(lhs.m_entityID > rhs.m_entityID)
{
return false;
}
if(lhs.m_orientation < rhs.m_orientation)
{
return true;
}
if(lhs.m_orientation > rhs.m_orientation)
{
return false;
}
return false;
}
bool Nektar::StdRegions::operator< ( const StdMatrixKey &  lhs,
const StdMatrixKey &  rhs 
)

Definition at line 100 of file StdMatrixKey.cpp.

References Nektar::StdRegions::StdMatrixKey::m_base, Nektar::StdRegions::StdMatrixKey::m_factors, Nektar::StdRegions::StdMatrixKey::m_matrixType, Nektar::StdRegions::StdMatrixKey::m_ncoeffs, Nektar::StdRegions::StdMatrixKey::m_nodalPointsType, Nektar::StdRegions::StdMatrixKey::m_shapeType, Nektar::StdRegions::StdMatrixKey::m_varcoeff_hashes, Nektar::StdRegions::StdMatrixKey::m_varcoeffs, and Nektar::LibUtilities::ShapeTypeDimMap.

{
if(lhs.m_matrixType < rhs.m_matrixType)
{
return true;
}
if(lhs.m_matrixType > rhs.m_matrixType)
{
return false;
}
if(lhs.m_ncoeffs < rhs.m_ncoeffs) // probably do not need this check since checking the m_base below?
{
return true;
}
if(lhs.m_ncoeffs > rhs.m_ncoeffs)
{
return false;
}
for(unsigned int i = 0; i < LibUtilities::ShapeTypeDimMap[lhs.m_shapeType]; ++i)
{
if(lhs.m_base[i].get() < rhs.m_base[i].get())
{
return true;
}
if(lhs.m_base[i].get() > rhs.m_base[i].get())
{
return false;
}
}
if(lhs.m_factors.size() < rhs.m_factors.size())
{
return true;
}
else if(lhs.m_factors.size() > rhs.m_factors.size())
{
return false;
}
else
{
ConstFactorMap::const_iterator x, y;
for(x = lhs.m_factors.begin(), y = rhs.m_factors.begin();
x != lhs.m_factors.end(); ++x, ++y)
{
if (x->second < y->second)
{
return true;
}
if (x->second > y->second)
{
return false;
}
}
}
if(lhs.m_varcoeffs.size() < rhs.m_varcoeffs.size())
{
return true;
}
if(lhs.m_varcoeffs.size() > rhs.m_varcoeffs.size())
{
return false;
}
for (unsigned int i = 0; i < lhs.m_varcoeff_hashes.size(); ++i)
{
if(lhs.m_varcoeff_hashes[i] < rhs.m_varcoeff_hashes[i])
{
return true;
}
if(lhs.m_varcoeff_hashes[i] > rhs.m_varcoeff_hashes[i])
{
return false;
}
}
if(lhs.m_nodalPointsType < rhs.m_nodalPointsType)
{
return true;
}
if(lhs.m_nodalPointsType > rhs.m_nodalPointsType)
{
return false;
}
return false;
}
std::ostream & Nektar::StdRegions::operator<< ( std::ostream &  os,
const IndexMapKey &  rhs 
)

Definition at line 199 of file IndexMapKey.cpp.

References Nektar::StdRegions::IndexMapKey::GetIndexMapType(), and IndexMapTypeMap.

{
os << "IndexMapType: " << IndexMapTypeMap[rhs.GetIndexMapType()]
<< std::endl;
return os;
}
std::ostream & Nektar::StdRegions::operator<< ( std::ostream &  os,
const StdMatrixKey &  rhs 
)

Definition at line 276 of file StdMatrixKey.cpp.

References ConstFactorTypeMap, Nektar::StdRegions::StdMatrixKey::GetBase(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNcoeffs(), Nektar::StdRegions::StdMatrixKey::GetShapeType(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffHashes(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), MatrixTypeMap, Nektar::LibUtilities::ShapeTypeDimMap, Nektar::LibUtilities::ShapeTypeMap, and VarCoeffTypeMap.

{
os << "MatrixType: " << MatrixTypeMap[rhs.GetMatrixType()] << ", ShapeType: "
<< LibUtilities::ShapeTypeMap[rhs.GetShapeType()] << ", Ncoeffs: " << rhs.GetNcoeffs()
<< std::endl;
if(rhs.GetConstFactors().size())
{
os << "Constants: " << endl;
ConstFactorMap::const_iterator x;
for(x = rhs.GetConstFactors().begin(); x != rhs.GetConstFactors().end(); ++x)
{
os << "\t value " << ConstFactorTypeMap[x->first] <<" : " << x->second << endl;
}
}
if(rhs.GetVarCoeffs().size())
{
os << "Variable coefficients: " << endl;
VarCoeffMap::const_iterator x;
unsigned int i = 0;
for (x = rhs.GetVarCoeffs().begin(); x != rhs.GetVarCoeffs().end(); ++x)
{
os << "\t Coeff defined: " << VarCoeffTypeMap[x->first] << endl;
os << "\t Hash: " << rhs.GetVarCoeffHashes()[i++] << endl;
}
}
for(unsigned int i = 0; i < LibUtilities::ShapeTypeDimMap[rhs.GetShapeType()]; ++i)
{
os << rhs.GetBase()[i]->GetBasisKey();
}
return os;
}
bool Nektar::StdRegions::operator== ( const IndexMapKey &  lhs,
const IndexMapKey &  rhs 
)

Definition at line 159 of file IndexMapKey.cpp.

References Nektar::StdRegions::IndexMapKey::m_entityID, Nektar::StdRegions::IndexMapKey::m_indexMapType, Nektar::StdRegions::IndexMapKey::m_orientation, Nektar::StdRegions::IndexMapKey::m_p, Nektar::StdRegions::IndexMapKey::m_q, Nektar::StdRegions::IndexMapKey::m_r, and Nektar::StdRegions::IndexMapKey::m_shapeType.

{
if(lhs.m_indexMapType != rhs.m_indexMapType)
{
return false;
}
if(lhs.m_shapeType != rhs.m_shapeType)
{
return false;
}
if(lhs.m_p != rhs.m_p)
{
return false;
}
if(lhs.m_q != rhs.m_q)
{
return false;
}
if(lhs.m_r != rhs.m_r)
{
return false;
}
if(lhs.m_entityID != rhs.m_entityID)
{
return false;
}
if(lhs.m_orientation != rhs.m_orientation)
{
return false;
}
return true;
}
bool Nektar::StdRegions::operator== ( const StdMatrixKey &  lhs,
const StdMatrixKey &  rhs 
)

Definition at line 195 of file StdMatrixKey.cpp.

References Nektar::StdRegions::StdMatrixKey::m_base, Nektar::StdRegions::StdMatrixKey::m_factors, Nektar::StdRegions::StdMatrixKey::m_matrixType, Nektar::StdRegions::StdMatrixKey::m_ncoeffs, Nektar::StdRegions::StdMatrixKey::m_nodalPointsType, Nektar::StdRegions::StdMatrixKey::m_shapeType, Nektar::StdRegions::StdMatrixKey::m_varcoeff_hashes, Nektar::StdRegions::StdMatrixKey::m_varcoeffs, and Nektar::LibUtilities::ShapeTypeDimMap.

{
if(lhs.m_matrixType != rhs.m_matrixType)
{
return false;
}
if(lhs.m_ncoeffs != rhs.m_ncoeffs)
{
return false;
}
for(unsigned int i = 0; i < LibUtilities::ShapeTypeDimMap[lhs.m_shapeType]; ++i)
{
if(lhs.m_base[i].get() != rhs.m_base[i].get())
{
return false;
}
}
if(lhs.m_factors.size() != rhs.m_factors.size())
{
return false;
}
else
{
ConstFactorMap::const_iterator x, y;
for(x = lhs.m_factors.begin(), y = rhs.m_factors.begin();
x != lhs.m_factors.end(); ++x, ++y)
{
if (x->second != y->second)
{
return false;
}
}
}
if(lhs.m_nodalPointsType != rhs.m_nodalPointsType)
{
return false;
}
if(lhs.m_varcoeffs.size() != rhs.m_varcoeffs.size())
{
return false;
}
for (unsigned int i = 0; i < lhs.m_varcoeff_hashes.size(); ++i)
{
if(lhs.m_varcoeff_hashes[i] != rhs.m_varcoeff_hashes[i])
{
return false;
}
}
VarCoeffMap::const_iterator x;
for (x = lhs.m_varcoeffs.begin(); x != lhs.m_varcoeffs.end(); ++x)
{
VarCoeffMap::const_iterator y;
// Check var coeff is found
if ((y = rhs.m_varcoeffs.find(x->first)) == rhs.m_varcoeffs.end())
{
return false;
}
if (x->second != y->second)
{
return false;
}
}
for (unsigned int i = 0; i < lhs.m_varcoeffs.size(); ++i)
{
if(lhs.m_varcoeff_hashes[i] != rhs.m_varcoeff_hashes[i])
{
return false;
}
}
return true;
}
NekMatrix< NekDouble > & Nektar::StdRegions::SetColumn ( NekMatrix< NekDouble > &  matA,
int  n,
const NekVector< NekDouble > &  x 
)

Definition at line 62 of file StdExpUtil.cpp.

Referenced by Invert().

{
for( int i=0; i<int(matA.GetRows()); ++i )
{
matA(i,n) = x[i];
}
return matA;
}
Array< OneD, NekDouble > Nektar::StdRegions::ToArray ( const NekVector< NekDouble > &  x)

Definition at line 126 of file StdExpUtil.cpp.

References Nektar::NekVector< DataType >::GetPtr(), and GetSize().

{
return Array<OneD, NekDouble>( GetSize(x), x.GetPtr() );
}
NekVector<NekDouble> Nektar::StdRegions::ToVector ( const Array< OneD, const NekDouble > &  x)

Definition at line 121 of file NodalUtil.cpp.

References Nektar::LibUtilities::GetSize().

{
return NekVector<NekDouble>( GetSize(x), x.data() );
}
NekVector<NekDouble> Nektar::StdRegions::ToVector ( const ConstArray< OneD, NekDouble > &  x)

Definition at line 121 of file StdExpUtil.cpp.

References GetSize().

{
return NekVector<NekDouble>( GetSize(x), x.data() );
}
NekVector< NekDouble > Nektar::StdRegions::VectorPower ( const NekVector< NekDouble > &  x,
NekDouble  p 
)

Definition at line 143 of file StdExpUtil.cpp.

References GetSize().

{
int size = GetSize(x);
NekVector<NekDouble> z(size);
for( int i=0; i<size; ++i )
{
z[i] = pow( x[i], p );
}
return z;
}
std::string Nektar::StdRegions::VectorToString ( const NekVector< NekDouble > &  v,
int  precision,
NekDouble  expSigFigs 
)

Definition at line 183 of file StdExpUtil.cpp.

References Nektar::NekVector< DataType >::GetRows(), and MakeRound().

{
stringstream s;
s << setprecision(precision) << "[ ";
int N = int(v.GetRows());
for(int j=0; j<N; ++j )
{
NekDouble x = MakeRound(expSigFigs * v(j)) / expSigFigs;
s << setw(7) << right << x;
if( j < N-1 )
{
s << ", ";
}
}
s << " ]";
return s.str();
}

Variable Documentation

const char* const Nektar::StdRegions::ConstFactorTypeMap[]
Initial value:
{
"FactorLambda",
"FactorTau",
"FactorTime",
"FactorSVVCutoffRatio",
"FactorSVVDiffCoeff",
"FactorGaussVertex",
"FactorGaussEdge"
}

Definition at line 237 of file StdRegions.hpp.

Referenced by Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::MultiRegions::operator<<(), and operator<<().

const char* const Nektar::StdRegions::ElementTypeMap[]
Initial value:
{
"StdSegExp",
"SegExp",
"StdQuadExp",
"StdTriExp",
"StdNodalTriExp",
"QuadExp",
"TriExp",
"NodalTriExp",
"StdHexExp",
"StdPrismExp",
"StdPyrExp",
"StdTetExp",
"StdNodalTetExp",
"HexExp",
"PrismExp",
"PyrExp",
"TetExp",
"NodalTetExp",
}

Definition at line 76 of file StdRegions.hpp.

Referenced by Nektar::Utilities::ProcessJac::Process().

const int Nektar::StdRegions::g_shapenedges[LibUtilities::SIZE_ShapeType] = {0,1,3,4,6,8,9,12}

define list of number of edges corresponding to each ShapeType

Definition at line 49 of file StdExpansion.cpp.

const int Nektar::StdRegions::g_shapenfaces[LibUtilities::SIZE_ShapeType] = {0,0,0,0,4,5,5,6}

define list of number of faces corresponding to each ShapeType

Definition at line 52 of file StdExpansion.cpp.

const int Nektar::StdRegions::g_shapenverts[LibUtilities::SIZE_ShapeType] = {0,2,3,4,4,5,6,8}

define list of number of vertices corresponding to each ShapeType

Definition at line 46 of file StdExpansion.cpp.

const char* const Nektar::StdRegions::IndexMapTypeMap[]
Initial value:
{
"EdgeToElement",
"FaceToElement",
"EdgeInterior",
"FaceInterior",
"Boundary",
"Vertex"
}

Definition at line 259 of file StdRegions.hpp.

Referenced by operator<<().

const char* const Nektar::StdRegions::MatrixTypeMap[]

Definition at line 150 of file StdRegions.hpp.

Referenced by Nektar::MultiRegions::operator<<(), and operator<<().

ConstFactorMap Nektar::StdRegions::NullConstFactorMap
static
VarCoeffMap Nektar::StdRegions::NullVarCoeffMap
static
const char* const Nektar::StdRegions::OrientationMap[]
Initial value:
{
"NoOrientation",
"Fwd",
"Bwd",
"Forwards",
"Backwards",
"Dir1FwdDir1_Dir2FwdDir2",
"Dir1FwdDir1_Dir2BwdDir2",
"Dir1BwdDir1_Dir2FwdDir2",
"Dir1BwdDir1_Dir2BwdDir2",
"Dir1FwdDir2_Dir2FwdDir1",
"Dir1FwdDir2_Dir2BwdDir1",
"Dir1BwdDir2_Dir2FwdDir1",
"Dir1BwdDir2_Dir2BwdDir1"
}

Definition at line 287 of file StdRegions.hpp.

const char* const Nektar::StdRegions::VarCoeffTypeMap[]
Initial value:
{
"VarCoeffMass",
"VarCoeffLaplacian",
"VarCoeffWeakDeriv",
"VarCoeffD00",
"VarCoeffD11",
"VarCoeffD22",
"VarCoeffD01",
"VarCoeffD02",
"VarCoeffD12",
"VarCoeffVelX",
"VarCoeffVelY"
}

Definition at line 210 of file StdRegions.hpp.

Referenced by Nektar::StdRegions::StdMatrixKey::GetVarCoeff(), and operator<<().