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

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  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::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 UseGLL)
 
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]
 

Typedef Documentation

◆ ConstFactorMap

Definition at line 430 of file StdRegions.hpp.

◆ FactorMap

Definition at line 434 of file StdRegions.hpp.

◆ StdExpansion0DSharedPtr

Definition at line 73 of file StdExpansion0D.h.

◆ StdExpansion1DSharedPtr

Definition at line 105 of file StdExpansion1D.h.

◆ StdExpansion2DSharedPtr

Definition at line 220 of file StdExpansion2D.h.

◆ StdExpansion3DSharedPtr

Definition at line 45 of file StdExpansion3D.h.

◆ StdExpansionSharedPtr

Definition at line 1794 of file StdExpansion.h.

◆ StdExpansionVector

Definition at line 1795 of file StdExpansion.h.

◆ StdHexExpSharedPtr

Definition at line 228 of file StdHexExp.h.

◆ StdMatrixKeySharedPtr

Definition at line 199 of file StdMatrixKey.h.

◆ StdNodalPrismExpSharedPtr

Definition at line 135 of file StdNodalPrismExp.h.

◆ StdNodalTetExpSharedPtr

Definition at line 135 of file StdNodalTetExp.h.

◆ StdNodalTriExpSharedPtr

Definition at line 167 of file StdNodalTriExp.h.

◆ StdPointExpSharedPtr

Definition at line 142 of file StdPointExp.h.

◆ StdPrismExpSharedPtr

Definition at line 218 of file StdPrismExp.h.

◆ StdPyrExpSharedPtr

Definition at line 214 of file StdPyrExp.h.

◆ StdQuadExpSharedPtr

Definition at line 221 of file StdQuadExp.h.

◆ StdSegExpSharedPtr

Definition at line 217 of file StdSegExp.h.

◆ StdTetExpSharedPtr

Definition at line 227 of file StdTetExp.h.

◆ StdTriExpSharedPtr

Definition at line 218 of file StdTriExp.h.

◆ VarCoeffMap

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

392{
411};

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

46{
47 // eStdPointExp,
49 eSegExp,
54 eTriExp,
61 eHexExp,
63 ePyrExp,
64 eTetExp,
67};

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

81{
83 eMass,
108 eBwdTrans,
109 eBwdMat,
126 eFwdTrans,
127 ePreconR,
132 eGaussDG,
137};
@ eLinearAdvectionDiffusionReactionGJP
Definition: StdRegions.hpp:105

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

438{
440 eFwd,
441 eBwd,
442 eForwards,
444 eDir1FwdDir1_Dir2FwdDir2, // These flags are interpreted as
445 eDir1FwdDir1_Dir2BwdDir2, // taking the second direction to the
446 eDir1BwdDir1_Dir2FwdDir2, // first direction. So Dir1FwdDir2 takes
447 eDir1BwdDir1_Dir2BwdDir2, // direction 2 and makes it backward
448 eDir1FwdDir2_Dir2FwdDir1, // to direction 1 in the mapped face.
449 eDir1FwdDir2_Dir2BwdDir1, // Note be careful not to flip this
450 eDir1BwdDir2_Dir2FwdDir1, // convention especially when using
451 eDir1BwdDir2_Dir2BwdDir1, // transposed mappings.
453};

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

196{
232};

Function Documentation

◆ EvaluateQuadFaceBasisKey()

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

Definition at line 431 of file StdExpansion3D.cpp.

435{
436
437 switch (faceDirBasisType)
438 {
440 {
441 const LibUtilities::PointsKey pkey(
444 pkey);
445 }
448 {
449 const LibUtilities::PointsKey pkey(
452 pkey);
453 }
455 {
456 const LibUtilities::PointsKey pkey(
458 return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange, nummodes,
459 pkey);
460 }
462 {
463 const LibUtilities::PointsKey pkey(
465 return LibUtilities::BasisKey(LibUtilities::eOrtho_A, nummodes,
466 pkey);
467 }
470 {
471 const LibUtilities::PointsKey pkey(
473 return LibUtilities::BasisKey(LibUtilities::eOrtho_A, nummodes,
474 pkey);
475 }
476 default:
477 {
478 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
479 break;
480 }
481 }
482
483 // Keep things happy by returning a value.
485}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
Describes the specification for a Basis.
Definition: Basis.h:45
Defines a specification for a set of points.
Definition: Points.h:50
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:51
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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,
bool  UseGLL 
)

Definition at line 487 of file StdExpansion3D.cpp.

490{
491 switch (faceDirBasisType)
492 {
494 {
495 const LibUtilities::PointsKey pkey(
498 pkey);
499 }
503 {
504 switch (facedir)
505 {
506 case 0:
507 {
508 const LibUtilities::PointsKey pkey(
511 nummodes, pkey);
512 }
513 case 1:
514 {
515 LibUtilities::PointsKey pkey;
516
517 if (UseGLL)
518 {
519 pkey = LibUtilities::PointsKey(
521 }
522 else
523 {
524 pkey = LibUtilities::PointsKey(
525 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
526 }
527 return LibUtilities::BasisKey(LibUtilities::eModified_B,
528 nummodes, pkey);
529 }
530 default:
531 {
532
533 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
534 break;
535 }
536 }
537 break;
538 }
539
541 {
542 switch (facedir)
543 {
544 case 0:
545 {
546 const LibUtilities::PointsKey pkey(
548 return LibUtilities::BasisKey(LibUtilities::eOrtho_A,
549 nummodes, pkey);
550 }
551 case 1:
552 {
553 const LibUtilities::PointsKey pkey(
554 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
555 return LibUtilities::BasisKey(LibUtilities::eOrtho_B,
556 nummodes, pkey);
557 }
558 default:
559 {
560 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
561 break;
562 }
563 }
564 break;
565 }
566
571 {
572 switch (facedir)
573 {
574 case 0:
575 {
576 const LibUtilities::PointsKey pkey(
578 return LibUtilities::BasisKey(LibUtilities::eOrtho_A,
579 nummodes, pkey);
580 }
581 case 1:
582 {
583 const LibUtilities::PointsKey pkey(
584 numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
585 return LibUtilities::BasisKey(LibUtilities::eOrtho_B,
586 nummodes, pkey);
587 }
588 default:
589 {
590 NEKERROR(ErrorUtil::efatal, "invalid value to flag");
591 break;
592 }
593 }
594 break;
595 }
596 default:
597 {
598 NEKERROR(ErrorUtil::efatal, "expansion type unknown");
599 break;
600 }
601 }
602
603 // Keep things happy by returning a value.
605}
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:53
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:51

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

477{
478 InputIterator val;
479
480 if (startingpoint == first)
481 {
482 val = find(first, last, value);
483 }
484 else
485 {
486 val = find(startingpoint, last, value);
487 if (val == last)
488 {
489 val = find(first, startingpoint, value);
490 if (val == startingpoint)
491 {
492 val = last;
493 }
494 }
495 }
496 return val;
497}
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:475

References find().

Referenced by Nektar::SolverUtils::ForcingAbsorption::CalcAbsorption(), Nektar::FieldUtils::Iso::Condense(), Nektar::SolverUtils::FilterThresholdMax::FilterThresholdMax(), Nektar::SolverUtils::FilterThresholdMin::FilterThresholdMin(), find(), Nektar::SolverUtils::Coupling::GenerateVariableMapping(), Nektar::ParseUtils::GenerateVariableVector(), 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::FieldUtils::Field::ReadFieldDefs(), Nektar::LibUtilities::SessionReader::ReadVariables(), Nektar::LibUtilities::CsvIO::v_ImportPtsFieldData(), 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::MeshGraphIOXml::v_ReadCurves().

◆ operator<()

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

Definition at line 94 of file StdMatrixKey.cpp.

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

◆ operator<<()

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

Definition at line 277 of file StdMatrixKey.cpp.

278{
279 os << "MatrixType: " << MatrixTypeMap[rhs.GetMatrixType()]
280 << ", ShapeType: " << LibUtilities::ShapeTypeMap[rhs.GetShapeType()]
281 << ", Ncoeffs: " << rhs.GetNcoeffs() << std::endl;
282
283 if (rhs.GetConstFactors().size())
284 {
285 os << "Constants: " << endl;
286 for (auto &x : rhs.GetConstFactors())
287 {
288 os << "\t value " << ConstFactorTypeMap[x.first] << " : "
289 << x.second << endl;
290 }
291 }
292 if (rhs.GetVarCoeffs().size())
293 {
294 os << "Variable coefficients: " << endl;
295 unsigned int i = 0;
296 for (auto &x : rhs.GetVarCoeffs())
297 {
298 os << "\t Coeff defined: " << VarCoeffTypeMap[x.first] << endl;
299 os << "\t Hash: " << rhs.GetVarCoeffHashes()[i++] << endl;
300 }
301 }
302
303 for (unsigned int i = 0;
305 {
306 os << rhs.GetBase()[i]->GetBasisKey();
307 }
308
309 return os;
310}
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
Definition: StdMatrixKey.h:103
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:88
const VarCoeffMap & GetVarCoeffs() const
Definition: StdMatrixKey.h:169
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:138
std::vector< std::size_t > GetVarCoeffHashes() const
Definition: StdMatrixKey.h:109
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:75
const char *const VarCoeffTypeMap[]
Definition: StdRegions.hpp:234
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:413
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:139

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

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

◆ RestrictCoeffMap()

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

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

Referenced by export_StdMatrixKey(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::MultiRegions::operator<<(), operator<<(), Nektar::Collections::Helmholtz_IterPerExp::UpdateFactors(), Nektar::Collections::Helmholtz_MatrixFree::UpdateFactors(), Nektar::Collections::LinearAdvectionDiffusionReaction_IterPerExp::UpdateFactors(), Nektar::Collections::LinearAdvectionDiffusionReaction_MatrixFree::UpdateFactors(), Nektar::Collections::PhysInterp1DScaled_Helper::UpdateFactors(), and Nektar::Collections::PhysInterp1DScaled_MatrixFree::UpdateFactors().

◆ 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 69 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 503 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

Definition at line 376 of file StdRegions.hpp.

Referenced by Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::LocalRegions::Expansion2D::CreateMatrix(), Nektar::LocalRegions::Expansion3D::CreateMatrix(), export_ExpList(), export_MatrixKey(), export_StdMatrixKey(), Nektar::Collections::Helmholtz_IterPerExp::Helmholtz_IterPerExp(), Nektar::Collections::Helmholtz_MatrixFree::Helmholtz_MatrixFree(), Nektar::Collections::Helmholtz_NoCollection::Helmholtz_NoCollection(), Nektar::Collections::LinearAdvectionDiffusionReaction_IterPerExp::LinearAdvectionDiffusionReaction_IterPerExp(), Nektar::Collections::LinearAdvectionDiffusionReaction_MatrixFree::LinearAdvectionDiffusionReaction_MatrixFree(), Nektar::Collections::LinearAdvectionDiffusionReaction_NoCollection::LinearAdvectionDiffusionReaction_NoCollection(), Nektar::StdRegions::StdNodalPrismExp::ModalToNodal(), Nektar::StdRegions::StdNodalTetExp::ModalToNodal(), Nektar::StdRegions::StdNodalTriExp::ModalToNodal(), Nektar::StdRegions::StdNodalPrismExp::NodalToModal(), Nektar::StdRegions::StdNodalTetExp::NodalToModal(), Nektar::StdRegions::StdNodalPrismExp::NodalToModalTranspose(), Nektar::StdRegions::StdNodalTetExp::NodalToModalTranspose(), Nektar::StdRegions::StdNodalTriExp::NodalToModalTranspose(), Nektar::Collections::Helmholtz_NoCollection::operator()(), Nektar::Collections::LinearAdvectionDiffusionReaction_NoCollection::operator()(), Nektar::LocalRegions::NodalTriExp::v_FwdTrans(), Nektar::StdRegions::StdNodalPrismExp::v_FwdTrans(), Nektar::StdRegions::StdNodalTetExp::v_FwdTrans(), Nektar::StdRegions::StdNodalTriExp::v_FwdTrans(), Nektar::StdRegions::StdNodalTriExp::v_NodalToModal(), Nektar::VCSImplicit::v_SolvePressure(), Nektar::VCSWeakPressure::v_SolvePressure(), and Nektar::VelocityCorrectionScheme::v_SolveViscous().

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

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