35 #ifndef NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
36 #define NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
38 #include <boost/core/ignore_unused.hpp>
57 namespace MultiRegions
64 class GlobalLinSysKey;
103 class ExpList:
public std::enable_shared_from_this<ExpList>
121 const std::vector<unsigned int> &eIDs,
122 const bool DeclareCoeffPhysArrays =
true);
127 const bool DeclareCoeffPhysArrays =
true);
298 const bool PhysSpaceForcing =
true);
353 bool UnShuff =
true);
360 bool UnShuff =
true);
399 std::string var =
"")
405 std::ostream &outfile,
433 std::ostream &outfile,
437 std::ostream &outfile,
439 std::string var =
"v")
487 bool useComm =
true);
552 return v_L2(inarray, soln);
661 inline const std::shared_ptr<LocalRegions::ExpansionVector>
679 bool returnNearestElmt =
false);
690 bool returnNearestElmt =
false);
773 inline std::shared_ptr<ExpList> &
GetTrace();
775 inline std::shared_ptr<AssemblyMapDG> &
GetTraceMap(
void);
812 inline const Array<
OneD,
const SpatialDomains::
820 const std::string varName =
"",
845 std::shared_ptr<ExpList> &result,
846 const bool DeclareCoeffPhysArrays =
true);
864 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
868 std::vector<NekDouble> &HomoLen =
870 bool homoStrips =
false,
871 std::vector<unsigned int> &HomoSIDs =
873 std::vector<unsigned int> &HomoZIDs =
875 std::vector<unsigned int> &HomoYIDs =
897 std::vector<LibUtilities::FieldDefinitionsSharedPtr>
915 std::vector<NekDouble> &fielddata)
925 std::vector<NekDouble> &fielddata,
938 std::vector<NekDouble> &fielddata,
949 const std::shared_ptr<ExpList> &fromExpList,
957 std::vector<NekDouble> &fielddata,
970 return shared_from_this();
974 std::shared_ptr<LibUtilities::SessionReader>
GetSession()
const
1017 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1090 std::shared_ptr<LocalRegions::ExpansionVector>
m_exp;
1134 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1138 const std::shared_ptr<DNekMat> &Gmat,
1149 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1166 return (*m_exp).size();
1186 virtual std::shared_ptr<ExpList> &
v_GetTrace();
1240 const bool PhysSpaceForcing);
1361 bool UnShuff =
true);
1368 bool UnShuff =
true);
1403 std::shared_ptr<ExpList> &result,
1404 const bool DeclareCoeffPhysArrays);
1423 virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1427 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1431 std::vector<NekDouble> &fielddata);
1435 std::vector<NekDouble> &fielddata,
1440 std::vector<NekDouble> &fielddata, std::string &field,
1446 std::string var =
"");
1454 std::ostream &outfile,
1459 std::ostream &outfile,
1493 const std::string &varName,
1494 const std::shared_ptr<ExpList> locExpList);
1502 unsigned int index,
const std::string& variable);
1513 const std::string varName =
"",
1529 "This method is not defined or valid for this class type");
1536 boost::ignore_unused(visc);
1538 "This method is not defined or valid for this class type");
1542 virtual std::shared_ptr<ExpList> &
v_GetPlane(
int n);
1565 return (*
m_exp)[eid]->GetNcoeffs();
1577 for(i= 0; i < (*m_exp).size(); ++i)
1579 returnval = (std::max)(returnval,
1580 (*
m_exp)[i]->EvalBasisNumModesMax());
1595 for(i= 0; i < (*m_exp).size(); ++i)
1598 = (std::max)(returnval[i],(*
m_exp)[i]->EvalBasisNumModesMax());
1615 return (*
m_exp)[eid]->GetTotPoints();
1623 int nbase = (*m_exp)[0]->GetNumBases();
1625 for(
int i = 0; i < (*m_exp).size(); ++i)
1628 for(
int j = 0; j < nbase; ++j)
1630 cnt *= (int)(scale*((*
m_exp)[i]->GetNumPoints(j)));
1681 "Input array does not have correct number of elements.");
1816 const bool PhysSpaceForcing)
1819 v_HelmSolve(inarray, outarray, flags, factors, varcoeff,
1820 varfactors, dirForcing, PhysSpaceForcing);
1836 lambda, coeffstate,dirForcing);
1848 lambda, coeffstate,dirForcing);
2016 "eid is larger than number of elements");
2017 return (*
m_exp)[eid]->GetCoordim();
2026 return (*
m_exp)[0]->GetShapeDimension();
2172 return (*m_exp).size();
2190 inline const std::shared_ptr<LocalRegions::ExpansionVector>
2358 const std::string varName,
2416 std::shared_ptr<ExpList> &result,
2417 const bool DeclareCoeffPhysArrays)
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
#define MULTI_REGIONS_EXPORT
Base class for all multi-elemental spectral/hp expansions.
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
void IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to all local expansion modes .
std::shared_ptr< ExpList > & GetTrace()
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void SetCoeffs(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
void ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract the data in fielddata into the coeffs.
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
std::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
virtual void v_ReadGlobalOptimizationParameters()
void Upwind(const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
void LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve Advection Diffusion Reaction.
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
void GlobalEigenSystem(const std::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
void DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
void Reset()
Reset geometry information and reset matrices.
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
BlockMatrixMapShPtr m_blockMat
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
virtual void v_ExtractPhysToBnd(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
std::shared_ptr< ExpList > & GetPlane(int n)
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
void GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
virtual std::shared_ptr< ExpList > & v_GetPlane(int n)
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n.
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
int GetExpSize(void)
This function returns the number of elements in the expansion.
virtual void v_ClearGlobalLinSysManager(void)
void FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the forward transformation of a function onto the global spectra...
void GeneralGetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, bool homoStrips=false, std::vector< unsigned int > &HomoSIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
void IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to the derivative (in directio...
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associaed with the homogeneous expansion.
virtual void v_LocalToGlobal(bool UseComm)
int GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
void SetHomoLen(const NekDouble lhom)
This function sets the Width of homogeneous direction associaed with the homogeneous expansion.
virtual int v_GetNumElmts(void)
virtual void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
SpatialDomains::MeshGraphSharedPtr GetGraph()
virtual ~ExpList()
The default destructor.
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
const Array< OneD, const int > & GetTraceBndMap(void)
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
virtual void v_FillBndCondFromField()
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
virtual void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
void ClearGlobalLinSysManager(void)
virtual const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
int GetExpIndex(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false)
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
std::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
virtual void v_GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
void HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const MultiRegions::VarFactorsMap &varfactors=MultiRegions::NullVarFactorsMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const bool PhysSpaceForcing=true)
Solve helmholtz problem.
virtual const Array< OneD, const int > & v_GetTraceBndMap()
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
virtual NekDouble v_GetHomoLen(void)
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
NekDouble PhysIntegral(void)
This function integrates a function over the domain consisting of all the elements of the expansion.
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
void ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expans...
void ReadGlobalOptimizationParameters()
void GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
This function calculates the coordinates of all the elemental quadrature points .
const Array< OneD, int > EvalBasisNumModesMaxPerExp(void) const
Returns the vector of the number of modes in the elemental basis order over all elements.
std::shared_ptr< GlobalLinSys > GenGlobalBndLinSys(const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocTo...
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
bool m_physState
The state of the array m_phys.
ExpList()
The default constructor.
virtual void v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
void FillBndCondFromField(void)
Fill Bnd Condition expansion from the values stored in expansion.
std::shared_ptr< LibUtilities::Comm > GetComm()
Returns the comm object.
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
virtual void v_BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetShapeDimension()
This function returns the dimension of the shape of the element eid.
void ExtractElmtDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather tha...
Array< OneD, const unsigned int > GetYIDs(void)
This function returns a vector containing the wave numbers in y-direction associated with the 3D homo...
NekDouble H1(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Calculates the error of the global spectral/hp element approximation.
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
LibUtilities::CommSharedPtr m_comm
Communicator.
std::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
void WriteVtkHeader(std::ostream &outfile)
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
void GenerateElementVector(const int ElementID, const NekDouble scalar1, const NekDouble scalar2, Array< OneD, NekDouble > &outarray)
Generate vector v such that v[i] = scalar1 if i is in the element < ElementID. Otherwise,...
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
void SetPhysState(const bool physState)
This function manually sets whether the array of physical values (implemented as m_phys) is filled o...
void PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function Galerkin projects the physical space points in inarray to outarray where inarray is ass...
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
virtual void v_ExtractElmtToBndPhys(const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
virtual void v_ExtractPhysToBndElmt(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
void CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
void LinearAdvectionDiffusionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve Advection Diffusion Reaction.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
virtual void v_Upwind(const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_SetUpPhysNormals()
virtual void v_GlobalToLocal(void)
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
void ExtractElmtToBndPhys(int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
std::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
void WriteVtkFooter(std::ostream &outfile)
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray, CoeffState coeffstate=eLocal)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual std::shared_ptr< ExpList > & v_GetTrace()
int Get1DScaledTotPoints(const NekDouble scale) const
Returns the total number of qudature points scaled by the factor scale on each 1D direction.
void PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function interpolates the physical space points in inarray to outarray using the same points def...
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
NekDouble L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error with respect to soln of the global spectral/hp element approximat...
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
void IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Directional derivative along a given direction.
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Array< OneD, const NekDouble > HomogeneousEnergy(void)
This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expan...
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Collections::CollectionVector m_collections
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void ApplyGeomInfo()
Apply geometry information to each expansion.
Array< OneD, const unsigned int > GetZIDs(void)
This function returns a vector containing the wave numbers in z-direction associated with the 3D homo...
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
virtual void v_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
LibUtilities::SessionReaderSharedPtr m_session
Session.
virtual void v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata.
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
virtual void v_SetHomoLen(const NekDouble lhom)
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
virtual void v_IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const NekOptimize::GlobalOptParamSharedPtr & GetGlobalOptParam(void)
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
void MultiplyByElmtInvMass(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass ...
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associaed with the homogeneous expansion.
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
int GetPhys_Offset(int n) const
Get the start offset position for a global list of m_phys correspoinding to element n.
void SetModifiedBasis(const bool modbasis)
Set Modified Basis for the stability analysis.
void DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray, CoeffState coeffstate=eLocal)
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract data from raw field data into expansion list.
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
const std::vector< bool > & GetLeftAdjacentFaces(void) const
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
ExpansionType GetExpType(void)
Returns the type of the expansion.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
NekDouble Linf(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error of the global spectral/hp element approximation.
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
void BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the backward transformation of the global spectral/hp element exp...
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
void GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
This function calculates the result of the multiplication of a matrix of type specified by mkey with ...
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Describe a linear system.
Describes a matrix with ordering defined by a local to global map.
std::vector< Collection > CollectionVector
static BasisSharedPtr NullBasisSharedPtr
static std::vector< unsigned int > NullUnsignedIntVector
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
static std::vector< NekDouble > NullNekDoubleVector
std::shared_ptr< Transposition > TranspositionSharedPtr
static Array< OneD, BasisSharedPtr > NullBasisSharedPtr1DArray
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
std::shared_ptr< Expansion > ExpansionSharedPtr
static ExpList NullExpList
An empty ExpList object.
static PeriodicMap NullPeriodicMap
MultiRegions::Direction const DirCartesianMap[]
std::shared_ptr< BlockMatrixMap > BlockMatrixMapShPtr
A shared pointer to a BlockMatrixMap.
static VarFactorsMap NullVarFactorsMap
@ eLocal
Local coefficients.
static ExpListSharedPtr NullExpListSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const Array< OneD, ExpListSharedPtr > NullExpListSharedPtrArray
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
std::map< GlobalMatrixKey, DNekScalBlkMatSharedPtr > BlockMatrixMap
A map between global matrix keys and their associated block matrices.
static const NekDouble kNekUnsetDouble
std::shared_ptr< GlobalOptParam > GlobalOptParamSharedPtr
Pointer to a GlobalOptParam object.
GeomMMF
Principle direction for MMF.
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
std::map< ConstFactorType, NekDouble > ConstFactorMap
static VarCoeffMap NullVarCoeffMap
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
static Array< OneD, NekDouble > NullNekDouble1DArray
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)