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

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

Classes

struct  cmpop
 
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...
 
class  StdPyrExp
 
class  StdQuadExp
 
class  StdSegExp
 Class representing a segment element in reference space All interface of this class sits in StdExpansion class. More...
 
class  StdTetExp
 
class  StdTriExp
 
struct  VarCoeffEntry
 Representation of a variable coefficient. More...
 

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, VarCoeffEntryVarCoeffMap
 
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 {
  eNoMatrixType , eMass , eMassGJP , eInvMass ,
  eLaplacian , eLaplacian00 , eLaplacian01 , eLaplacian02 ,
  eLaplacian10 , eLaplacian11 , eLaplacian12 , eLaplacian20 ,
  eLaplacian21 , eLaplacian22 , eInvLaplacianWithUnityMean , eWeakDeriv0 ,
  eWeakDeriv1 , eWeakDeriv2 , eWeakDirectionalDeriv , eMassLevelCurvature ,
  eLinearAdvection , eLinearAdvectionReaction , eLinearAdvectionDiffusionReaction , eLinearAdvectionDiffusionReactionGJP ,
  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 ,
  eVarCoeffD01 , eVarCoeffD02 , eVarCoeffD10 , eVarCoeffD11 ,
  eVarCoeffD12 , eVarCoeffD20 , eVarCoeffD21 , eVarCoeffD22 ,
  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)
 
VarCoeffMap RestrictCoeffMap (const VarCoeffMap &m, size_t offset, size_t cnt)
 
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 408 of file StdRegions.hpp.

◆ FactorMap

Definition at line 412 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 83 of file StdExpansion0D.h.

◆ StdExpansion1DSharedPtr

Definition at line 118 of file StdExpansion1D.h.

◆ StdExpansion2DSharedPtr

Definition at line 224 of file StdExpansion2D.h.

◆ StdExpansion3DSharedPtr

Definition at line 49 of file StdExpansion3D.h.

◆ StdExpansionSharedPtr

Definition at line 1775 of file StdExpansion.h.

◆ StdExpansionVector

Definition at line 1776 of file StdExpansion.h.

◆ StdHexExpSharedPtr

Definition at line 255 of file StdHexExp.h.

◆ StdMatrixKeySharedPtr

Definition at line 203 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 174 of file StdNodalTriExp.h.

◆ StdPointExpSharedPtr

Definition at line 152 of file StdPointExp.h.

◆ StdPrismExpSharedPtr

Definition at line 239 of file StdPrismExp.h.

◆ StdPyrExpSharedPtr

Definition at line 258 of file StdPyrExp.h.

◆ StdQuadExpSharedPtr

Definition at line 243 of file StdQuadExp.h.

◆ StdSegExpSharedPtr

Definition at line 267 of file StdSegExp.h.

◆ StdTetExpSharedPtr

Definition at line 255 of file StdTetExp.h.

◆ StdTriExpSharedPtr

Definition at line 232 of file StdTriExp.h.

◆ VarCoeffMap

Definition at line 352 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 369 of file StdRegions.hpp.

370{
389};

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

52{
53 // eStdPointExp,
55 eSegExp,
60 eTriExp,
67 eHexExp,
69 ePyrExp,
70 eTetExp,
73};

◆ MatrixType

Todo:
we need to tidy up matrix construction approach probably using a factory type approach
Enumerator
eNoMatrixType 
eMass 
eMassGJP 
eInvMass 
eLaplacian 
eLaplacian00 
eLaplacian01 
eLaplacian02 
eLaplacian10 
eLaplacian11 
eLaplacian12 
eLaplacian20 
eLaplacian21 
eLaplacian22 
eInvLaplacianWithUnityMean 
eWeakDeriv0 
eWeakDeriv1 
eWeakDeriv2 
eWeakDirectionalDeriv 
eMassLevelCurvature 
eLinearAdvection 
eLinearAdvectionReaction 
eLinearAdvectionDiffusionReaction 
eLinearAdvectionDiffusionReactionGJP 
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 86 of file StdRegions.hpp.

87{
89 eMass,
114 eBwdTrans,
115 eBwdMat,
132 eFwdTrans,
133 ePreconR,
138 eGaussDG,
143};
@ eLinearAdvectionDiffusionReactionGJP
Definition: StdRegions.hpp:111

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

416{
418 eFwd,
419 eBwd,
420 eForwards,
422 eDir1FwdDir1_Dir2FwdDir2, // These flags are interpreted as
423 eDir1FwdDir1_Dir2BwdDir2, // taking the second direction to the
424 eDir1BwdDir1_Dir2FwdDir2, // first direction. So Dir1FwdDir2 takes
425 eDir1BwdDir1_Dir2BwdDir2, // direction 2 and makes it backward
426 eDir1FwdDir2_Dir2FwdDir1, // to direction 1 in the mapped face.
427 eDir1FwdDir2_Dir2BwdDir1, // Note be careful not to flip this
428 eDir1BwdDir2_Dir2FwdDir1, // convention especially when using
429 eDir1BwdDir2_Dir2BwdDir1, // transposed mappings.
431};

◆ VarCoeffType

Enumerator
eVarCoeffMass 
eVarCoeffLaplacian 
eVarCoeffWeakDeriv 
eVarCoeffD00 
eVarCoeffD01 
eVarCoeffD02 
eVarCoeffD10 
eVarCoeffD11 
eVarCoeffD12 
eVarCoeffD20 
eVarCoeffD21 
eVarCoeffD22 
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 201 of file StdRegions.hpp.

202{
238};

Function Documentation

◆ EvaluateQuadFaceBasisKey()

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

Definition at line 453 of file StdExpansion3D.cpp.

456{
457 boost::ignore_unused(facedir);
458
459 switch (faceDirBasisType)
460 {
462 {
463 const LibUtilities::PointsKey pkey(
466 pkey);
467 }
470 {
471 const LibUtilities::PointsKey pkey(
474 pkey);
475 }
477 {
478 const LibUtilities::PointsKey pkey(
480 return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange, nummodes,
481 pkey);
482 }
484 {
485 const LibUtilities::PointsKey pkey(
487 return LibUtilities::BasisKey(LibUtilities::eOrtho_A, nummodes,
488 pkey);
489 }
492 {
493 const LibUtilities::PointsKey pkey(
495 return LibUtilities::BasisKey(LibUtilities::eOrtho_A, nummodes,
496 pkey);
497 }
498 default:
499 {
500 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
501 break;
502 }
503 }
504
505 // Keep things happy by returning a value.
507}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
Describes the specification for a Basis.
Definition: Basis.h:47
Defines a specification for a set of points.
Definition: Points.h:55
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::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 509 of file StdExpansion3D.cpp.

512{
513 switch (faceDirBasisType)
514 {
516 {
517 const LibUtilities::PointsKey pkey(
520 pkey);
521 }
525 {
526 switch (facedir)
527 {
528 case 0:
529 {
530 const LibUtilities::PointsKey pkey(
533 nummodes, pkey);
534 }
535 case 1:
536 {
537 // const LibUtilities::PointsKey pkey(
538 // numpoints+1,
539 // LibUtilities::eGaussLobattoLegendre);
540 const LibUtilities::PointsKey pkey(
541 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
542 return LibUtilities::BasisKey(LibUtilities::eModified_B,
543 nummodes, pkey);
544 }
545 default:
546 {
547
548 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
549 break;
550 }
551 }
552 break;
553 }
554
556 {
557 switch (facedir)
558 {
559 case 0:
560 {
561 const LibUtilities::PointsKey pkey(
563 return LibUtilities::BasisKey(LibUtilities::eOrtho_A,
564 nummodes, pkey);
565 }
566 case 1:
567 {
568 const LibUtilities::PointsKey pkey(
569 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
570 return LibUtilities::BasisKey(LibUtilities::eOrtho_B,
571 nummodes, pkey);
572 }
573 default:
574 {
575 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
576 break;
577 }
578 }
579 break;
580 }
581
586 {
587 switch (facedir)
588 {
589 case 0:
590 {
591 const LibUtilities::PointsKey pkey(
593 return LibUtilities::BasisKey(LibUtilities::eOrtho_A,
594 nummodes, pkey);
595 }
596 case 1:
597 {
598 const LibUtilities::PointsKey pkey(
599 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
600 return LibUtilities::BasisKey(LibUtilities::eOrtho_B,
601 nummodes, pkey);
602 }
603 default:
604 {
605 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
606 break;
607 }
608 }
609 break;
610 }
611 default:
612 {
613 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
614 break;
615 }
616 }
617
618 // Keep things happy by returning a value.
620}
@ 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::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 453 of file StdRegions.hpp.

455{
456 InputIterator val;
457
458 if (startingpoint == first)
459 {
460 val = find(first, last, value);
461 }
462 else
463 {
464 val = find(startingpoint, last, value);
465 if (val == last)
466 {
467 val = find(first, startingpoint, value);
468 if (val == startingpoint)
469 {
470 val = last;
471 }
472 }
473 }
474 return val;
475}
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:453

References find().

Referenced by Nektar::SolverUtils::ForcingAbsorption::CalcAbsorption(), Nektar::FieldUtils::Iso::Condense(), Nektar::FieldUtils::Field::CreateExp(), Nektar::SolverUtils::FilterThresholdMax::FilterThresholdMax(), Nektar::SolverUtils::FilterThresholdMin::FilterThresholdMin(), find(), Nektar::SolverUtils::Coupling::GenerateVariableMapping(), Nektar::ParseUtils::GenerateVariableSet(), Nektar::SpatialDomains::MeshGraph::GetCompositeList(), Nektar::LibUtilities::Interpreter::ExpressionEvaluator::GetConstant(), Nektar::FieldUtils::ProcessMapping::GetMapping(), Nektar::MultiRegions::AssemblyCommDG::InitialiseStructure(), Nektar::FieldUtils::ProcessCreateExp::LoadFieldData(), Nektar::SpatialDomains::BoundaryConditions::ReadBoundaryConditions(), 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(), Nektar::Dummy::v_PostIntegrate(), Nektar::FieldUtils::InputFld::v_Process(), Nektar::FieldUtils::ProcessAddCompositeID::v_Process(), Nektar::FieldUtils::ProcessAddFld::v_Process(), Nektar::FieldUtils::ProcessFieldFromString::v_Process(), Nektar::FieldUtils::ProcessRemoveField::v_Process(), and Nektar::SpatialDomains::MeshGraphXml::v_ReadCurves().

◆ operator<()

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

Definition at line 96 of file StdMatrixKey.cpp.

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

◆ operator<<()

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

Definition at line 279 of file StdMatrixKey.cpp.

280{
281 os << "MatrixType: " << MatrixTypeMap[rhs.GetMatrixType()]
282 << ", ShapeType: " << LibUtilities::ShapeTypeMap[rhs.GetShapeType()]
283 << ", Ncoeffs: " << rhs.GetNcoeffs() << std::endl;
284
285 if (rhs.GetConstFactors().size())
286 {
287 os << "Constants: " << endl;
288 for (auto &x : rhs.GetConstFactors())
289 {
290 os << "\t value " << ConstFactorTypeMap[x.first] << " : "
291 << x.second << endl;
292 }
293 }
294 if (rhs.GetVarCoeffs().size())
295 {
296 os << "Variable coefficients: " << endl;
297 unsigned int i = 0;
298 for (auto &x : rhs.GetVarCoeffs())
299 {
300 os << "\t Coeff defined: " << VarCoeffTypeMap[x.first] << endl;
301 os << "\t Hash: " << rhs.GetVarCoeffHashes()[i++] << endl;
302 }
303 }
304
305 for (unsigned int i = 0;
307 {
308 os << rhs.GetBase()[i]->GetBasisKey();
309 }
310
311 return os;
312}
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
Definition: StdMatrixKey.h:107
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:92
const VarCoeffMap & GetVarCoeffs() const
Definition: StdMatrixKey.h:173
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:87
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:142
std::vector< std::size_t > GetVarCoeffHashes() const
Definition: StdMatrixKey.h:113
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:79
const char *const VarCoeffTypeMap[]
Definition: StdRegions.hpp:240
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:391
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:145

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 204 of file StdMatrixKey.cpp.

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

◆ RestrictCoeffMap()

VarCoeffMap Nektar::StdRegions::RestrictCoeffMap ( const VarCoeffMap m,
size_t  offset,
size_t  cnt 
)
inline

Definition at line 355 of file StdRegions.hpp.

357{
358 VarCoeffMap ret;
359
360 for (auto &x : m)
361 {
362 ret[x.first] =
363 Array<OneD, NekDouble>(cnt, x.second.GetValue() + offset);
364 }
365
366 return ret;
367}
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:352

Referenced by Nektar::MultiRegions::ExpList::GenBlockMatrix(), Nektar::MultiRegions::ExpList::GeneralMatrixOp(), Nektar::MultiRegions::ExpList::GenGlobalMatrix(), Nektar::MultiRegions::ExpList::GenGlobalMatrixFull(), and Nektar::MultiRegions::GlobalLinSys::GetBlockMatrixKey().

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 391 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 75 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 481 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 433 of file StdRegions.hpp.

◆ VarCoeffTypeMap

const char* const Nektar::StdRegions::VarCoeffTypeMap[]
Initial value:
= {
"VarCoeffMass", "VarCoeffLaplacian", "VarCoeffWeakDeriv",
"VarCoeffD00", "VarCoeffD01", "VarCoeffD02",
"VarCoeffD10", "VarCoeffD11", "VarCoeffD12",
"VarCoeffD20", "VarCoeffD21", "VarCoeffD22",
"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 240 of file StdRegions.hpp.

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