Nektar++
Classes | Typedefs | Enumerations | Functions | Variables
Nektar::StdRegions Namespace Reference

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

Classes

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 std::shared_ptr< StdExpansionStdExpansionSharedPtr
 
typedef std::vector< StdExpansionSharedPtrStdExpansionVector
 
typedef std::shared_ptr< StdExpansion0DStdExpansion0DSharedPtr
 
typedef std::shared_ptr< StdExpansion1DStdExpansion1DSharedPtr
 
typedef std::shared_ptr< StdExpansion2DStdExpansion2DSharedPtr
 
typedef std::shared_ptr< StdExpansion3DStdExpansion3DSharedPtr
 
typedef std::shared_ptr< StdHexExpStdHexExpSharedPtr
 
typedef std::shared_ptr< StdMatrixKeyStdMatrixKeySharedPtr
 
typedef std::shared_ptr< StdNodalPrismExpStdNodalPrismExpSharedPtr
 
typedef std::shared_ptr< StdNodalTetExpStdNodalTetExpSharedPtr
 
typedef std::shared_ptr< StdNodalTriExpStdNodalTriExpSharedPtr
 
typedef std::shared_ptr< StdPointExpStdPointExpSharedPtr
 
typedef std::shared_ptr< StdPrismExpStdPrismExpSharedPtr
 
typedef std::tuple< unsigned int, unsigned int, unsigned int, unsigned int > Mode
 
typedef std::shared_ptr< StdPyrExpStdPyrExpSharedPtr
 
typedef std::shared_ptr< StdQuadExpStdQuadExpSharedPtr
 
typedef std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
 
typedef std::map< ConstFactorType, NekDoubleConstFactorMap
 
typedef ConstFactorMap FactorMap
 
typedef std::shared_ptr< StdSegExpStdSegExpSharedPtr
 
typedef std::shared_ptr< StdTetExpStdTetExpSharedPtr
 
typedef std::shared_ptr< StdTriExpStdTriExpSharedPtr
 

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 , eMassGJP , eInvMass , eLaplacian ,
  eLaplacian00 , eLaplacian01 , eLaplacian02 , eLaplacian10 ,
  eLaplacian11 , eLaplacian12 , eLaplacian20 , eLaplacian21 ,
  eLaplacian22 , eInvLaplacianWithUnityMean , eWeakDeriv0 , eWeakDeriv1 ,
  eWeakDeriv2 , eWeakDirectionalDeriv , eMassLevelCurvature , eLinearAdvectionReaction ,
  eLinearAdvectionDiffusionReaction , eNBasisTrans , eInvNBasisTrans , eBwdTrans ,
  eBwdMat , eIProductWRTBase , eIProductWRTDerivBase0 , eIProductWRTDerivBase1 ,
  eIProductWRTDerivBase2 , eDerivBase0 , eDerivBase1 , eDerivBase2 ,
  eHelmholtz , eHelmholtzGJP , eHybridDGHelmholtz , eInvHybridDGHelmholtz ,
  eHybridDGHelmBndLam , eHybridDGLamToQ0 , eHybridDGLamToQ1 , eHybridDGLamToQ2 ,
  eHybridDGLamToU , eFwdTrans , ePreconR , ePreconRMass ,
  ePreconLinearSpace , ePreconLinearSpaceMass , eInterpGauss , eGaussDG ,
  ePhysInterpToEquiSpaced , eEquiSpacedToCoeffs , eNormDerivOnTrace , SIZE_MatrixType
}
 
enum  VarCoeffType {
  eVarCoeffMass , eVarCoeffLaplacian , eVarCoeffWeakDeriv , eVarCoeffD00 ,
  eVarCoeffD11 , eVarCoeffD22 , eVarCoeffD01 , eVarCoeffD02 ,
  eVarCoeffD12 , eVarCoeffVelX , eVarCoeffVelY , eVarCoeffVelZ ,
  eVarCoeffMF1x , eVarCoeffMF1y , eVarCoeffMF1z , eVarCoeffMF1Div ,
  eVarCoeffMF1Mag , eVarCoeffMF2x , eVarCoeffMF2y , eVarCoeffMF2z ,
  eVarCoeffMF2Div , eVarCoeffMF2Mag , eVarCoeffMF3x , eVarCoeffMF3y ,
  eVarCoeffMF3z , eVarCoeffMF3Div , eVarCoeffMF3Mag , eVarCoeffMF ,
  eVarCoeffMFDiv , eVarCoeffGmat , eVarCoeffGJPNormVel , SIZE_VarCoeffType
}
 
enum  ConstFactorType {
  eFactorLambda , eFactorCoeffD00 , eFactorCoeffD11 , eFactorCoeffD22 ,
  eFactorCoeffD01 , eFactorCoeffD02 , eFactorCoeffD12 , eFactorTau ,
  eFactorTime , eFactorSVVCutoffRatio , eFactorSVVDiffCoeff , eFactorSVVPowerKerDiffCoeff ,
  eFactorSVVDGKerDiffCoeff , eFactorGaussVertex , eFactorGaussEdge , eFactorGJP ,
  eFactorConst , SIZE_ConstFactorType
}
 
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

LibUtilities::BasisKey EvaluateQuadFaceBasisKey (const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
 
LibUtilities::BasisKey EvaluateTriFaceBasisKey (const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
 
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 char *const ElementTypeMap []
 
const char *const MatrixTypeMap []
 
const char *const VarCoeffTypeMap []
 
static VarCoeffMap NullVarCoeffMap
 
const char *const ConstFactorTypeMap []
 
static ConstFactorMap NullConstFactorMap
 
static FactorMap NullFactorMap
 
const char *const OrientationMap []
 
const int kSVVDGFiltermodesmin = 3
 
const int kSVVDGFiltermodesmax = 11
 
const NekDouble kSVVDGFilter [9][11]
 

Detailed Description

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

Typedef Documentation

◆ ConstFactorMap

Definition at line 282 of file StdRegions.hpp.

◆ FactorMap

Definition at line 286 of file StdRegions.hpp.

◆ Mode

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

Definition at line 47 of file StdPyrExp.h.

◆ StdExpansion0DSharedPtr

Definition at line 86 of file StdExpansion0D.h.

◆ StdExpansion1DSharedPtr

Definition at line 91 of file StdExpansion1D.h.

◆ StdExpansion2DSharedPtr

Definition at line 196 of file StdExpansion2D.h.

◆ StdExpansion3DSharedPtr

Definition at line 49 of file StdExpansion3D.h.

◆ StdExpansionSharedPtr

Definition at line 1632 of file StdExpansion.h.

◆ StdExpansionVector

Definition at line 1633 of file StdExpansion.h.

◆ StdHexExpSharedPtr

Definition at line 254 of file StdHexExp.h.

◆ StdMatrixKeySharedPtr

Definition at line 197 of file StdMatrixKey.h.

◆ StdNodalPrismExpSharedPtr

Definition at line 138 of file StdNodalPrismExp.h.

◆ StdNodalTetExpSharedPtr

Definition at line 138 of file StdNodalTetExp.h.

◆ StdNodalTriExpSharedPtr

Definition at line 172 of file StdNodalTriExp.h.

◆ StdPointExpSharedPtr

Definition at line 130 of file StdPointExp.h.

◆ StdPrismExpSharedPtr

Definition at line 231 of file StdPrismExp.h.

◆ StdPyrExpSharedPtr

Definition at line 248 of file StdPyrExp.h.

◆ StdQuadExpSharedPtr

Definition at line 234 of file StdQuadExp.h.

◆ StdSegExpSharedPtr

Definition at line 47 of file StdSegExp.h.

◆ StdTetExpSharedPtr

Definition at line 247 of file StdTetExp.h.

◆ StdTriExpSharedPtr

Definition at line 235 of file StdTriExp.h.

◆ VarCoeffMap

Definition at line 240 of file StdRegions.hpp.

Enumeration Type Documentation

◆ ConstFactorType

Enumerator
eFactorLambda 
eFactorCoeffD00 
eFactorCoeffD11 
eFactorCoeffD22 
eFactorCoeffD01 
eFactorCoeffD02 
eFactorCoeffD12 
eFactorTau 
eFactorTime 
eFactorSVVCutoffRatio 
eFactorSVVDiffCoeff 
eFactorSVVPowerKerDiffCoeff 
eFactorSVVDGKerDiffCoeff 
eFactorGaussVertex 
eFactorGaussEdge 
eFactorGJP 
eFactorConst 
SIZE_ConstFactorType 

Definition at line 243 of file StdRegions.hpp.

244 {
252  eFactorTau,
253  eFactorTime,
260  eFactorGJP,
261  eFactorConst,
263 };

◆ ElementType

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

Definition at line 50 of file StdRegions.hpp.

51 {
52  // eStdPointExp,
53  eStdSegExp,
54  eSegExp,
56  eStdTriExp,
58  eQuadExp,
59  eTriExp,
61  eStdHexExp,
63  eStdPyrExp,
64  eStdTetExp,
66  eHexExp,
67  ePrismExp,
68  ePyrExp,
69  eTetExp,
72 };

◆ MatrixType

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

Definition at line 84 of file StdRegions.hpp.

85 {
86  eMass,
87  eMassGJP,
88  eInvMass,
89  eLaplacian,
100  eWeakDeriv0,
101  eWeakDeriv1,
102  eWeakDeriv2,
107  eNBasisTrans,
109  eBwdTrans,
110  eBwdMat,
115  eDerivBase0,
116  eDerivBase1,
117  eDerivBase2,
118  eHelmholtz,
127  eFwdTrans,
128  ePreconR,
129  ePreconRMass,
132  eInterpGauss,
133  eGaussDG,
138 };

◆ Orientation

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 289 of file StdRegions.hpp.

290 {
292  eFwd,
293  eBwd,
294  eForwards,
295  eBackwards,
296  eDir1FwdDir1_Dir2FwdDir2, // These flags are interpreted as
297  eDir1FwdDir1_Dir2BwdDir2, // taking the second direction to the
298  eDir1BwdDir1_Dir2FwdDir2, // first direction. So Dir1FwdDir2 takes
299  eDir1BwdDir1_Dir2BwdDir2, // direction 2 and makes it backward
300  eDir1FwdDir2_Dir2FwdDir1, // to direction 1 in the mapped face.
301  eDir1FwdDir2_Dir2BwdDir1, // Note be careful not to flip this
302  eDir1BwdDir2_Dir2FwdDir1, // convention especially when using
303  eDir1BwdDir2_Dir2BwdDir1, // transposed mappings.
305 };

◆ VarCoeffType

Enumerator
eVarCoeffMass 
eVarCoeffLaplacian 
eVarCoeffWeakDeriv 
eVarCoeffD00 
eVarCoeffD11 
eVarCoeffD22 
eVarCoeffD01 
eVarCoeffD02 
eVarCoeffD12 
eVarCoeffVelX 
eVarCoeffVelY 
eVarCoeffVelZ 
eVarCoeffMF1x 
eVarCoeffMF1y 
eVarCoeffMF1z 
eVarCoeffMF1Div 
eVarCoeffMF1Mag 
eVarCoeffMF2x 
eVarCoeffMF2y 
eVarCoeffMF2z 
eVarCoeffMF2Div 
eVarCoeffMF2Mag 
eVarCoeffMF3x 
eVarCoeffMF3y 
eVarCoeffMF3z 
eVarCoeffMF3Div 
eVarCoeffMF3Mag 
eVarCoeffMF 
eVarCoeffMFDiv 
eVarCoeffGmat 
eVarCoeffGJPNormVel 
SIZE_VarCoeffType 

Definition at line 192 of file StdRegions.hpp.

193 {
197  eVarCoeffD00,
198  eVarCoeffD11,
199  eVarCoeffD22,
200  eVarCoeffD01,
201  eVarCoeffD02,
202  eVarCoeffD12,
221  eVarCoeffMF,
226 };

Function Documentation

◆ EvaluateQuadFaceBasisKey()

LibUtilities::BasisKey Nektar::StdRegions::EvaluateQuadFaceBasisKey ( const int  facedir,
const LibUtilities::BasisType  faceDirBasisType,
const int  numpoints,
const int  nummodes 
)

Definition at line 443 of file StdExpansion3D.cpp.

446 {
447  boost::ignore_unused(facedir);
448 
449  switch (faceDirBasisType)
450  {
452  {
453  const LibUtilities::PointsKey pkey(
455  return LibUtilities::BasisKey(LibUtilities::eModified_A, nummodes,
456  pkey);
457  }
460  {
461  const LibUtilities::PointsKey pkey(
462  numpoints + 1, LibUtilities::eGaussLobattoLegendre);
463  return LibUtilities::BasisKey(LibUtilities::eModified_A, nummodes,
464  pkey);
465  }
467  {
468  const LibUtilities::PointsKey pkey(
470  return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange, nummodes,
471  pkey);
472  }
474  {
475  const LibUtilities::PointsKey pkey(
477  return LibUtilities::BasisKey(LibUtilities::eOrtho_A, nummodes,
478  pkey);
479  }
482  {
483  const LibUtilities::PointsKey pkey(
484  numpoints + 1, LibUtilities::eGaussLobattoLegendre);
485  return LibUtilities::BasisKey(LibUtilities::eOrtho_A, nummodes,
486  pkey);
487  }
488  default:
489  {
490  NEKERROR(ErrorUtil::efatal, "expansion type unknown");
491  break;
492  }
493  }
494 
495  // Keep things happy by returning a value.
497 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:52
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:48
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, NEKERROR, and Nektar::LibUtilities::NullBasisKey().

Referenced by Nektar::SpatialDomains::MeshGraph::GetFaceBasisKey(), Nektar::StdRegions::StdHexExp::v_GetTraceBasisKey(), Nektar::StdRegions::StdPrismExp::v_GetTraceBasisKey(), and Nektar::StdRegions::StdPyrExp::v_GetTraceBasisKey().

◆ EvaluateTriFaceBasisKey()

LibUtilities::BasisKey Nektar::StdRegions::EvaluateTriFaceBasisKey ( const int  facedir,
const LibUtilities::BasisType  faceDirBasisType,
const int  numpoints,
const int  nummodes 
)

Definition at line 499 of file StdExpansion3D.cpp.

502 {
503  switch (faceDirBasisType)
504  {
506  {
507  const LibUtilities::PointsKey pkey(
509  return LibUtilities::BasisKey(LibUtilities::eModified_A, nummodes,
510  pkey);
511  }
515  {
516  switch (facedir)
517  {
518  case 0:
519  {
520  const LibUtilities::PointsKey pkey(
521  numpoints + 1, LibUtilities::eGaussLobattoLegendre);
522  return LibUtilities::BasisKey(LibUtilities::eModified_A,
523  nummodes, pkey);
524  }
525  case 1:
526  {
527  // const LibUtilities::PointsKey pkey(
528  // numpoints+1,
529  // LibUtilities::eGaussLobattoLegendre);
530  const LibUtilities::PointsKey pkey(
531  numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
532  return LibUtilities::BasisKey(LibUtilities::eModified_B,
533  nummodes, pkey);
534  }
535  default:
536  {
537 
538  NEKERROR(ErrorUtil::efatal, "invalid value to flag");
539  break;
540  }
541  }
542  break;
543  }
544 
546  {
547  switch (facedir)
548  {
549  case 0:
550  {
551  const LibUtilities::PointsKey pkey(
553  return LibUtilities::BasisKey(LibUtilities::eOrtho_A,
554  nummodes, pkey);
555  }
556  case 1:
557  {
558  const LibUtilities::PointsKey pkey(
559  numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
560  return LibUtilities::BasisKey(LibUtilities::eOrtho_B,
561  nummodes, pkey);
562  }
563  default:
564  {
565  NEKERROR(ErrorUtil::efatal, "invalid value to flag");
566  break;
567  }
568  }
569  break;
570  }
571 
576  {
577  switch (facedir)
578  {
579  case 0:
580  {
581  const LibUtilities::PointsKey pkey(
583  return LibUtilities::BasisKey(LibUtilities::eOrtho_A,
584  nummodes, pkey);
585  }
586  case 1:
587  {
588  const LibUtilities::PointsKey pkey(
589  numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
590  return LibUtilities::BasisKey(LibUtilities::eOrtho_B,
591  nummodes, pkey);
592  }
593  default:
594  {
595  NEKERROR(ErrorUtil::efatal, "invalid value to flag");
596  break;
597  }
598  }
599  break;
600  }
601  default:
602  {
603  NEKERROR(ErrorUtil::efatal, "expansion type unknown");
604  break;
605  }
606  }
607 
608  // Keep things happy by returning a value.
610 }
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:55
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:53

References Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eModifiedPyr_C, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::LibUtilities::eOrthoPyr_C, NEKERROR, and Nektar::LibUtilities::NullBasisKey().

Referenced by Nektar::SpatialDomains::MeshGraph::GetFaceBasisKey(), Nektar::StdRegions::StdPrismExp::v_GetTraceBasisKey(), Nektar::StdRegions::StdPyrExp::v_GetTraceBasisKey(), and Nektar::StdRegions::StdTetExp::v_GetTraceBasisKey().

◆ find()

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

Definition at line 327 of file StdRegions.hpp.

329 {
330  InputIterator val;
331 
332  if (startingpoint == first)
333  {
334  val = find(first, last, value);
335  }
336  else
337  {
338  val = find(startingpoint, last, value);
339  if (val == last)
340  {
341  val = find(first, startingpoint, value);
342  if (val == startingpoint)
343  {
344  val = last;
345  }
346  }
347  }
348  return val;
349 }
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:327

Referenced by Nektar::SolverUtils::ForcingAbsorption::CalcAbsorption(), Nektar::FieldUtils::Iso::Condense(), Nektar::FieldUtils::Field::CreateExp(), Nektar::SolverUtils::FilterThresholdMax::FilterThresholdMax(), Nektar::SolverUtils::FilterThresholdMin::FilterThresholdMin(), Nektar::SolverUtils::Coupling::GenerateVariableMapping(), Nektar::SpatialDomains::MeshGraph::GetCompositeList(), Nektar::LibUtilities::Interpreter::ExpressionEvaluator::GetConstant(), Nektar::FieldUtils::ProcessMapping::GetMapping(), Nektar::MultiRegions::AssemblyCommDG::InitialiseStructure(), Nektar::FieldUtils::ProcessCreateExp::LoadFieldData(), Nektar::FieldUtils::InputFld::Process(), Nektar::FieldUtils::ProcessAddCompositeID::Process(), Nektar::FieldUtils::ProcessAddFld::Process(), Nektar::FieldUtils::ProcessFieldFromString::Process(), Nektar::FieldUtils::ProcessRemoveField::Process(), Nektar::SpatialDomains::Domain::Read(), Nektar::SpatialDomains::BoundaryConditions::ReadBoundaryConditions(), Nektar::SpatialDomains::MeshGraphXml::ReadCurves(), Nektar::SpatialDomains::MeshGraph::ReadExpansionInfo(), Nektar::LibUtilities::SessionReader::ReadVariables(), Nektar::LibUtilities::CsvIO::v_ImportFieldData(), Nektar::SolverUtils::FilterAeroForces::v_Initialise(), Nektar::FilterMovingBody::v_Initialise(), Nektar::Dummy::v_InitObject(), and Nektar::Dummy::v_PostIntegrate().

◆ operator<()

bool Nektar::StdRegions::operator< ( const StdMatrixKey lhs,
const StdMatrixKey rhs 
)

Definition at line 90 of file StdMatrixKey.cpp.

91 {
92  if (lhs.m_matrixType < rhs.m_matrixType)
93  {
94  return true;
95  }
96 
97  if (lhs.m_matrixType > rhs.m_matrixType)
98  {
99  return false;
100  }
101 
102  if (lhs.m_ncoeffs < rhs.m_ncoeffs) // probably do not need this check since
103  // checking the m_base below?
104  {
105  return true;
106  }
107 
108  if (lhs.m_ncoeffs > rhs.m_ncoeffs)
109  {
110  return false;
111  }
112 
113  if (LibUtilities::ShapeTypeDimMap[lhs.m_shapeType] <
114  LibUtilities::ShapeTypeDimMap[rhs.m_shapeType])
115  {
116  return true;
117  }
118 
119  if (LibUtilities::ShapeTypeDimMap[lhs.m_shapeType] >
120  LibUtilities::ShapeTypeDimMap[rhs.m_shapeType])
121  {
122  return false;
123  }
124 
125  for (unsigned int i = 0; i < LibUtilities::ShapeTypeDimMap[lhs.m_shapeType];
126  ++i)
127  {
128  if (lhs.m_base[i].get() < rhs.m_base[i].get())
129  {
130  return true;
131  }
132 
133  if (lhs.m_base[i].get() > rhs.m_base[i].get())
134  {
135  return false;
136  }
137  }
138 
139  if (lhs.m_factors.size() < rhs.m_factors.size())
140  {
141  return true;
142  }
143  else if (lhs.m_factors.size() > rhs.m_factors.size())
144  {
145  return false;
146  }
147  else
148  {
149  for (auto x = lhs.m_factors.begin(), y = rhs.m_factors.begin();
150  x != lhs.m_factors.end(); ++x, ++y)
151  {
152  if (x->second < y->second)
153  {
154  return true;
155  }
156  if (x->second > y->second)
157  {
158  return false;
159  }
160  }
161  }
162 
163  if (lhs.m_varcoeffs.size() < rhs.m_varcoeffs.size())
164  {
165  return true;
166  }
167 
168  if (lhs.m_varcoeffs.size() > rhs.m_varcoeffs.size())
169  {
170  return false;
171  }
172 
173  for (unsigned int i = 0; i < lhs.m_varcoeff_hashes.size(); ++i)
174  {
175  if (lhs.m_varcoeff_hashes[i] < rhs.m_varcoeff_hashes[i])
176  {
177  return true;
178  }
179  if (lhs.m_varcoeff_hashes[i] > rhs.m_varcoeff_hashes[i])
180  {
181  return false;
182  }
183  }
184 
185  if (lhs.m_nodalPointsType < rhs.m_nodalPointsType)
186  {
187  return true;
188  }
189 
190  if (lhs.m_nodalPointsType > rhs.m_nodalPointsType)
191  {
192  return false;
193  }
194 
195  return false;
196 }
constexpr unsigned int ShapeTypeDimMap[SIZE_ShapeType]
Definition: ShapeType.hpp:83

◆ operator<<()

std::ostream & Nektar::StdRegions::operator<< ( std::ostream &  os,
const StdMatrixKey rhs 
)

Definition at line 273 of file StdMatrixKey.cpp.

274 {
275  os << "MatrixType: " << MatrixTypeMap[rhs.GetMatrixType()]
276  << ", ShapeType: " << LibUtilities::ShapeTypeMap[rhs.GetShapeType()]
277  << ", Ncoeffs: " << rhs.GetNcoeffs() << std::endl;
278 
279  if (rhs.GetConstFactors().size())
280  {
281  os << "Constants: " << endl;
282  for (auto &x : rhs.GetConstFactors())
283  {
284  os << "\t value " << ConstFactorTypeMap[x.first] << " : "
285  << x.second << endl;
286  }
287  }
288  if (rhs.GetVarCoeffs().size())
289  {
290  os << "Variable coefficients: " << endl;
291  unsigned int i = 0;
292  for (auto &x : rhs.GetVarCoeffs())
293  {
294  os << "\t Coeff defined: " << VarCoeffTypeMap[x.first] << endl;
295  os << "\t Hash: " << rhs.GetVarCoeffHashes()[i++] << endl;
296  }
297  }
298 
299  for (unsigned int i = 0;
300  i < LibUtilities::ShapeTypeDimMap[rhs.GetShapeType()]; ++i)
301  {
302  os << rhs.GetBase()[i]->GetBasisKey();
303  }
304 
305  return os;
306 }
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:77
const char *const VarCoeffTypeMap[]
Definition: StdRegions.hpp:228
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:265
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:140

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.

◆ operator==()

bool Nektar::StdRegions::operator== ( const StdMatrixKey lhs,
const StdMatrixKey rhs 
)

Definition at line 198 of file StdMatrixKey.cpp.

199 {
200  if (lhs.m_matrixType != rhs.m_matrixType)
201  {
202  return false;
203  }
204 
205  if (lhs.m_ncoeffs != rhs.m_ncoeffs)
206  {
207  return false;
208  }
209 
210  for (unsigned int i = 0; i < LibUtilities::ShapeTypeDimMap[lhs.m_shapeType];
211  ++i)
212  {
213  if (lhs.m_base[i].get() != rhs.m_base[i].get())
214  {
215  return false;
216  }
217  }
218 
219  if (lhs.m_factors.size() != rhs.m_factors.size())
220  {
221  return false;
222  }
223  else
224  {
225  return lhs.m_factors.size() == rhs.m_factors.size() &&
226  std::equal(lhs.m_factors.begin(), lhs.m_factors.end(),
227  rhs.m_factors.begin());
228  }
229 
230  if (lhs.m_nodalPointsType != rhs.m_nodalPointsType)
231  {
232  return false;
233  }
234 
235  if (lhs.m_varcoeffs.size() != rhs.m_varcoeffs.size())
236  {
237  return false;
238  }
239 
240  for (unsigned int i = 0; i < lhs.m_varcoeff_hashes.size(); ++i)
241  {
242  if (lhs.m_varcoeff_hashes[i] != rhs.m_varcoeff_hashes[i])
243  {
244  return false;
245  }
246  }
247 
248  for (auto &x : lhs.m_varcoeffs)
249  {
250  auto y = rhs.m_varcoeffs.find(x.first);
251  // Check var coeff is found
252  if (y == rhs.m_varcoeffs.end())
253  {
254  return false;
255  }
256 
257  if (x.second != y->second)
258  {
259  return false;
260  }
261  }
262  for (unsigned int i = 0; i < lhs.m_varcoeffs.size(); ++i)
263  {
264  if (lhs.m_varcoeff_hashes[i] != rhs.m_varcoeff_hashes[i])
265  {
266  return false;
267  }
268  }
269 
270  return true;
271 }

Variable Documentation

◆ ConstFactorTypeMap

const char* const Nektar::StdRegions::ConstFactorTypeMap[]
Initial value:
= {"FactorLambda",
"FactorCoeffD00",
"FactorCoeffD11",
"FactorCoeffD22",
"FactorCoeffD01",
"FactorCoeffD02",
"FactorCoeffD12",
"FactorTau",
"FactorTime",
"FactorSVVCutoffRatio",
"FactorSVVDiffCoeff",
"FactorSVVPowerKerDiffCoeff",
"FactorSVVDGKerDiffCoeff",
"FactorGaussVertex",
"FactorGaussEdge",
"FactorGJP",
"FactorConstant"}

Definition at line 265 of file StdRegions.hpp.

Referenced by Nektar::Collections::Helmholtz_IterPerExp::CheckFactors(), Nektar::Collections::Helmholtz_MatrixFree::CheckFactors(), export_StdMatrixKey(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::MultiRegions::operator<<(), and operator<<().

◆ ElementTypeMap

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 74 of file StdRegions.hpp.

◆ kSVVDGFilter

const NekDouble Nektar::StdRegions::kSVVDGFilter[9][11]
Initial value:
= {
{0, 0.36212, 1, 0, 0, 0, 0, 0, 0, 0, 0},
{0, 0.70546, 0.078836, 1, 0, 0, 0, 0, 0, 0, 0},
{0, 0, 0.49411, 0.072394, 1, 0, 0, 0, 0, 0, 0},
{0, 0, 0.000073566, 0.40506, 0.094122, 1, 0, 0, 0, 0, 0},
{0, 0, 0, 0.0001422, 0.36863, 0.11815, 1, 0, 0, 0, 0},
{0, 0, 0, 0, 0.00019497, 0.41397, 0.16927, 1, 0, 0, 0},
{0, 0, 0, 0, 0, 0.0009762, 0.12747, 0.13763, 1, 0, 0},
{0, 0, 0, 0, 0, 0, 0.0023592, 0.23683, 0.17196, 1, 0},
{0, 0, 0, 0, 0, 0, 0, 0.0026055, 0.28682, 0.22473, 1}}

Definition at line 355 of file StdRegions.hpp.

Referenced by Nektar::StdRegions::StdHexExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdPrismExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdPyrExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdQuadExp::v_SVVLaplacianFilter(), Nektar::StdRegions::StdTetExp::v_SVVLaplacianFilter(), and Nektar::StdRegions::StdTriExp::v_SVVLaplacianFilter().

◆ kSVVDGFiltermodesmax

const int Nektar::StdRegions::kSVVDGFiltermodesmax = 11

◆ kSVVDGFiltermodesmin

const int Nektar::StdRegions::kSVVDGFiltermodesmin = 3

◆ MatrixTypeMap

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

◆ NullConstFactorMap

ConstFactorMap Nektar::StdRegions::NullConstFactorMap
static

◆ NullFactorMap

FactorMap Nektar::StdRegions::NullFactorMap
static

◆ NullVarCoeffMap

VarCoeffMap Nektar::StdRegions::NullVarCoeffMap
static

◆ OrientationMap

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 307 of file StdRegions.hpp.

◆ VarCoeffTypeMap

const char* const Nektar::StdRegions::VarCoeffTypeMap[]
Initial value:
= {
"VarCoeffMass", "VarCoeffLaplacian", "VarCoeffWeakDeriv",
"VarCoeffD00", "VarCoeffD11", "VarCoeffD22",
"VarCoeffD01", "VarCoeffD02", "VarCoeffD12",
"VarCoeffVelX", "VarCoeffVelY", "VarCoeffVelZ",
"VarCoeffMF1x", "VarCoeffMF1y", "VarCoeffMF1z",
"VarCoeffMF1Div", "VarCoeffMF1Mag", "VarCoeffMF2x",
"VarCoeffMF2y", "VarCoeffMF2z", "VarCoeffMF2Div",
"VarCoeffMF2Mag", "VarCoeffMF3x", "VarCoeffMF3y",
"VarCoeffMF3z", "VarCoeffMF3Div", "VarCoeffMF3Mag",
"VarCoeffMF", "VarCoeffMFDiv", "VarCoeffGmat",
"VarCoeffGJPNormVel"}

Definition at line 228 of file StdRegions.hpp.

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