Nektar++
|
Base class for all multi-elemental spectral/hp expansions. More...
#include <ExpList.h>
Public Member Functions | |
ExpList () | |
The default constructor. | |
ExpList (const LibUtilities::SessionReaderSharedPtr &pSession) | |
The default constructor. | |
ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph) | |
The default constructor. | |
ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true) | |
The copy constructor. | |
virtual | ~ExpList () |
The default destructor. | |
int | GetNcoeffs (void) const |
Returns the total number of local degrees of freedom . | |
int | GetNcoeffs (const int eid) const |
Returns the total number of local degrees of freedom for element eid. | |
ExpansionType | GetExpType (void) |
Returns the type of the expansion. | |
void | SetExpType (ExpansionType Type) |
Returns the type of the expansion. | |
int | EvalBasisNumModesMax (void) const |
Evaulates the maximum number of modes in the elemental basis order over all elements. | |
const Array< OneD, int > | EvalBasisNumModesMaxPerExp (void) const |
Returns the vector of the number of modes in the elemental basis order over all elements. | |
int | GetTotPoints (void) const |
Returns the total number of quadrature points m_npoints . | |
int | GetTotPoints (const int eid) const |
Returns the total number of quadrature points for eid's element . | |
int | GetNpoints (void) const |
Returns the total number of quadrature points m_npoints . | |
int | Get1DScaledTotPoints (const NekDouble scale) const |
Returns the total number of qudature points scaled by the factor scale on each 1D direction. | |
void | SetWaveSpace (const bool wavespace) |
Sets the wave space to the one of the possible configuration true or false. | |
void | SetModifiedBasis (const bool modbasis) |
Set Modified Basis for the stability analysis. | |
void | SetPhys (int i, NekDouble val) |
Set the i th value of m_phys to value val. | |
bool | GetWaveSpace (void) const |
This function returns the third direction expansion condition, which can be in wave space (coefficient) or not It is stored in the variable m_WaveSpace. | |
void | SetPhys (const Array< OneD, const NekDouble > &inarray) |
Fills the array m_phys. | |
void | SetPhysArray (Array< OneD, NekDouble > &inarray) |
Sets the array m_phys. | |
void | SetPhysState (const bool physState) |
This function manually sets whether the array of physical values (implemented as m_phys) is filled or not. | |
bool | GetPhysState (void) const |
This function indicates whether the array of physical values (implemented as m_phys) is filled or not. | |
NekDouble | PhysIntegral (void) |
This function integrates a function over the domain consisting of all the elements of the expansion. | |
NekDouble | PhysIntegral (const Array< OneD, const NekDouble > &inarray) |
This function integrates a function over the domain consisting of all the elements of the expansion. | |
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 . | |
void | IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal) |
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 direction. | |
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 spectral/hp expansion. | |
void | FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal) |
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 matrix. | |
void | MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal) |
void | SmoothField (Array< OneD, NekDouble > &field) |
Smooth a field across elements. | |
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 Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray) |
Solve helmholtz problem. | |
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. | |
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. | |
void | FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
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 expansion. | |
void | BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal) |
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 . | |
void | HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true) |
void | HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true) |
void | DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal) |
void | GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID) |
void | NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID) |
void | ApplyGeomInfo () |
Apply geometry information to each expansion. | |
void | WriteTecplotHeader (std::ofstream &outfile, std::string var="") |
void | WriteTecplotZone (std::ofstream &outfile, int expansion=-1) |
void | WriteTecplotField (std::ofstream &outfile, int expansion=-1) |
void | WriteTecplotConnectivity (std::ofstream &outfile, int expansion=-1) |
void | WriteVtkHeader (std::ofstream &outfile) |
void | WriteVtkFooter (std::ofstream &outfile) |
void | WriteVtkPieceHeader (std::ofstream &outfile, int expansion) |
void | WriteVtkPieceFooter (std::ofstream &outfile, int expansion) |
void | WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var="v") |
int | GetCoordim (int eid) |
This function returns the dimension of the coordinates of the element eid. | |
void | SetCoeff (int i, NekDouble val) |
Set the i th coefficiient in m_coeffs to value val. | |
void | SetCoeffs (int i, NekDouble val) |
Set the i th coefficiient in m_coeffs to value val. | |
void | SetCoeffsArray (Array< OneD, NekDouble > &inarray) |
Set the m_coeffs array to inarray. | |
const Array< OneD, const NekDouble > & | GetCoeffs () const |
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expansion coefficients. | |
void | ImposeDirichletConditions (Array< OneD, NekDouble > &outarray) |
Impose Dirichlet Boundary Conditions onto Array. | |
void | FillBndCondFromField (void) |
Fill Bnd Condition expansion from the values stored in expansion. | |
void | LocalToGlobal (void) |
Put the coefficients into global ordering using m_coeffs. | |
void | GlobalToLocal (void) |
Put the coefficients into local ordering and place in m_coeffs. | |
NekDouble | GetCoeff (int i) |
Get the i th value (coefficient) of m_coeffs. | |
NekDouble | GetCoeffs (int i) |
Get the i th value (coefficient) of m_coeffs. | |
const Array< OneD, const NekDouble > & | GetPhys () const |
This function returns (a reference to) the array (implemented as m_phys) containing the function evaluated at the quadrature points. | |
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. | |
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 approximation. | |
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. | |
NekDouble | Integral (const Array< OneD, const NekDouble > &inarray) |
Array< OneD, const NekDouble > | HomogeneousEnergy (void) |
This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expansion. | |
void | SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc) |
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion. | |
Array< OneD, const unsigned int > | GetZIDs (void) |
This function returns a vector containing the wave numbers in z-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction. | |
LibUtilities::TranspositionSharedPtr | GetTransposition (void) |
This function returns the transposition class associaed with the homogeneous expansion. | |
NekDouble | GetHomoLen (void) |
This function returns the Width of homogeneous direction associaed with the homogeneous expansion. | |
Array< OneD, const unsigned int > | GetYIDs (void) |
This function returns a vector containing the wave numbers in y-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier 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 defined in the expansion but where the number of points are rescaled by 1DScale. | |
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 assumed to be defined in the expansion but where the number of points are rescaled by 1DScale. | |
int | GetExpSize (void) |
This function returns the number of elements in the expansion. | |
int | GetNumElmts (void) |
This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp. | |
const boost::shared_ptr < LocalRegions::ExpansionVector > | GetExp () const |
This function returns the vector of elements in the expansion. | |
LocalRegions::ExpansionSharedPtr & | GetExp (int n) const |
This function returns (a shared pointer to) the local elemental expansion of the element. | |
LocalRegions::ExpansionSharedPtr & | GetExp (const Array< OneD, const NekDouble > &gloCoord) |
This function returns (a shared pointer to) the local elemental expansion containing the arbitrary point given by gloCoord. | |
int | GetExpIndex (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false) |
int | GetExpIndex (const Array< OneD, const NekDouble > &gloCoords, Array< OneD, NekDouble > &locCoords, NekDouble tol=0.0, bool returnNearestElmt=false) |
int | GetCoeff_Offset (int n) const |
Get the start offset position for a global list of m_coeffs correspoinding to element n. | |
int | GetPhys_Offset (int n) const |
Get the start offset position for a global list of m_phys correspoinding to element n. | |
int | GetOffset_Elmt_Id (int n) const |
Get the element id associated with the n th consecutive block of data in m_phys and m_coeffs. | |
Array< OneD, NekDouble > & | UpdateCoeffs () |
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expansion coefficients. | |
Array< OneD, NekDouble > & | UpdatePhys () |
This function returns (a reference to) the array (implemented as m_phys) containing the function evaluated at the quadrature points. | |
void | PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d) |
void | PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) |
This function discretely evaluates the derivative of a function on the domain consisting of all elements of the expansion. | |
void | PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d) |
const Array< OneD, const boost::shared_ptr< ExpList > > & | GetBndCondExpansions () |
boost::shared_ptr< ExpList > & | UpdateBndCondExpansion (int i) |
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 | Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind) |
boost::shared_ptr< ExpList > & | GetTrace () |
boost::shared_ptr < AssemblyMapDG > & | GetTraceMap (void) |
const Array< OneD, const int > & | GetTraceBndMap (void) |
void | GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals) |
void | AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray) |
void | AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray) |
void | AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray) |
void | GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd) |
void | GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd) |
void | ExtractTracePhys (Array< OneD, NekDouble > &outarray) |
void | ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & | GetBndConditions () |
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & | UpdateBndConditions () |
void | EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble) |
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 a vector given by inarray. | |
void | GeneralMatrixOp_IterPerExp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
void | SetUpPhysNormals () |
void | GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID) |
void | GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector) |
const NekOptimize::GlobalOptParamSharedPtr & | GetGlobalOptParam (void) |
map< int, RobinBCInfoSharedPtr > | GetRobinBCInfo () |
void | GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap) |
std::vector < LibUtilities::FieldDefinitionsSharedPtr > | GetFieldDefinitions () |
void | GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef) |
void | AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata) |
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata. | |
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 | 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 than planes in homogeneous case. | |
void | ExtractCoeffsToCoeffs (const boost::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 expansions rather than planes in homogeneous case. | |
void | ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs) |
Extract the data in fielddata into the coeffs. | |
boost::shared_ptr< ExpList > | GetSharedThisPtr () |
Returns a shared pointer to the current object. | |
boost::shared_ptr < LibUtilities::SessionReader > | GetSession () |
Returns the session object. | |
boost::shared_ptr < LibUtilities::Comm > | GetComm () |
Returns the comm object. | |
SpatialDomains::MeshGraphSharedPtr | GetGraph () |
LibUtilities::BasisSharedPtr | GetHomogeneousBasis (void) |
boost::shared_ptr< ExpList > & | GetPlane (int n) |
Public Attributes | |
ExpansionType | m_expType |
Protected Member Functions | |
boost::shared_ptr< DNekMat > | GenGlobalMatrixFull (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap) |
const DNekScalBlkMatSharedPtr | GenBlockMatrix (const GlobalMatrixKey &gkey) |
This function assembles the block diagonal matrix of local matrices of the type mtype. | |
const DNekScalBlkMatSharedPtr & | GetBlockMatrix (const GlobalMatrixKey &gkey) |
void | MultiplyByBlockMatrix (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
boost::shared_ptr< GlobalMatrix > | GenGlobalMatrix (const GlobalMatrixKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap) |
Generates a global matrix from the given key and map. | |
void | GlobalEigenSystem (const boost::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray) |
boost::shared_ptr< GlobalLinSys > | GenGlobalLinSys (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap) |
This operation constructs the global linear system of type mkey. | |
boost::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 LocToGloBaseMap. | |
void | ReadGlobalOptimizationParameters () |
virtual int | v_GetNumElmts (void) |
virtual const Array< OneD, const boost::shared_ptr < ExpList > > & | v_GetBndCondExpansions (void) |
virtual boost::shared_ptr < ExpList > & | v_UpdateBndCondExpansion (int i) |
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) |
virtual void | v_Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind) |
virtual boost::shared_ptr < ExpList > & | v_GetTrace () |
virtual boost::shared_ptr < AssemblyMapDG > & | v_GetTraceMap () |
virtual const Array< OneD, const int > & | v_GetTraceBndMap () |
virtual void | v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals) |
virtual void | v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray) |
virtual void | v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray) |
virtual void | v_AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray) |
virtual void | v_GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd) |
virtual void | v_GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd) |
virtual void | v_ExtractTracePhys (Array< OneD, NekDouble > &outarray) |
virtual void | v_ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
virtual void | v_MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate) |
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 Array< OneD, const NekDouble > &dirForcing) |
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) |
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) |
virtual void | v_ImposeDirichletConditions (Array< OneD, NekDouble > &outarray) |
virtual void | v_FillBndCondFromField () |
virtual void | v_LocalToGlobal (void) |
virtual void | v_GlobalToLocal (void) |
virtual void | v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate) |
virtual void | v_BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
virtual void | v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate) |
virtual void | v_FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
virtual void | v_SmoothField (Array< OneD, NekDouble > &field) |
virtual void | v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate) |
virtual void | v_IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
virtual void | v_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate) |
virtual void | v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray) |
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) |
virtual void | v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d) |
virtual void | v_PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d) |
virtual void | v_HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true) |
virtual void | v_HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true) |
virtual void | v_DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal) |
virtual void | v_GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID) |
virtual void | v_NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID) |
virtual void | v_SetUpPhysNormals () |
virtual void | v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID) |
virtual void | v_ReadGlobalOptimizationParameters () |
virtual std::vector < LibUtilities::FieldDefinitionsSharedPtr > | v_GetFieldDefinitions (void) |
virtual void | v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef) |
virtual void | v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata) |
virtual void | v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs) |
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. | |
virtual void | v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs) |
virtual void | v_WriteTecplotHeader (std::ofstream &outfile, std::string var="") |
virtual void | v_WriteTecplotZone (std::ofstream &outfile, int expansion) |
virtual void | v_WriteTecplotField (std::ofstream &outfile, int expansion) |
virtual void | v_WriteTecplotConnectivity (std::ofstream &outfile, int expansion) |
virtual void | v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion) |
virtual void | v_WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var) |
virtual NekDouble | v_L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray) |
virtual NekDouble | v_Integral (const Array< OneD, const NekDouble > &inarray) |
virtual Array< OneD, const NekDouble > | v_HomogeneousEnergy (void) |
virtual LibUtilities::TranspositionSharedPtr | v_GetTransposition (void) |
virtual NekDouble | v_GetHomoLen (void) |
virtual Array< OneD, const unsigned int > | v_GetZIDs (void) |
virtual Array< OneD, const unsigned int > | v_GetYIDs (void) |
virtual void | v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
virtual void | v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray) |
void | ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList) |
Static Protected Member Functions | |
static SpatialDomains::BoundaryConditionShPtr | GetBoundaryCondition (const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable) |
Protected Attributes | |
LibUtilities::CommSharedPtr | m_comm |
Communicator. | |
LibUtilities::SessionReaderSharedPtr | m_session |
Session. | |
SpatialDomains::MeshGraphSharedPtr | m_graph |
Mesh associated with this expansion list. | |
int | m_ncoeffs |
The total number of local degrees of freedom. m_ncoeffs . | |
int | m_npoints |
Array< OneD, NekDouble > | m_coeffs |
Concatenation of all local expansion coefficients. | |
Array< OneD, NekDouble > | m_phys |
The global expansion evaluated at the quadrature points. | |
bool | m_physState |
The state of the array m_phys. | |
boost::shared_ptr < LocalRegions::ExpansionVector > | m_exp |
The list of local expansions. | |
Array< OneD, int > | m_coeff_offset |
Offset of elemental data into the array m_coeffs. | |
Array< OneD, int > | m_phys_offset |
Offset of elemental data into the array m_phys. | |
Array< OneD, int > | m_offset_elmt_id |
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coeffs and m_phys is associated, i.e. for an array of constant expansion size and single shape elements m_phys[n*m_npoints] is the data related to m_exp[m_offset_elmt_id[n]];. | |
NekOptimize::GlobalOptParamSharedPtr | m_globalOptParam |
BlockMatrixMapShPtr | m_blockMat |
bool | m_WaveSpace |
Private Member Functions | |
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & | v_GetBndConditions () |
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & | v_UpdateBndConditions () |
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) |
virtual map< int, RobinBCInfoSharedPtr > | v_GetRobinBCInfo (void) |
virtual void | v_GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces) |
virtual LibUtilities::BasisSharedPtr | v_GetHomogeneousBasis (void) |
virtual void | v_SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc) |
virtual boost::shared_ptr < ExpList > & | v_GetPlane (int n) |
Base class for all multi-elemental spectral/hp expansions.
All multi-elemental expansions can be considered as the assembly of the various elemental contributions. On a discrete level, this yields,
where is the number of elements and is the local elemental number of expansion modes. As it is the lowest level class, it contains the definition of the common data and common routines to all multi-elemental expansions.
The class stores a vector of expansions, m_exp, (each derived from StdRegions::StdExpansion) which define the constituent components of the domain. The coefficients from these expansions are concatenated in m_coeffs, while the expansion evaluated at the quadrature points is stored in m_phys.
Nektar::MultiRegions::ExpList::ExpList | ( | ) |
The default constructor.
Creates an empty expansion list. The expansion list will typically be populated by a derived class (namely one of MultiRegions::ExpList1D, MultiRegions::ExpList2D or MultiRegions::ExpList3D).
Definition at line 82 of file ExpList.cpp.
References Nektar::MultiRegions::eNoType, and SetExpType().
Nektar::MultiRegions::ExpList::ExpList | ( | const LibUtilities::SessionReaderSharedPtr & | pSession | ) |
The default constructor.
Creates an empty expansion list. The expansion list will typically be populated by a derived class (namely one of MultiRegions::ExpList1D, MultiRegions::ExpList2D or MultiRegions::ExpList3D).
Definition at line 108 of file ExpList.cpp.
References Nektar::MultiRegions::eNoType, and SetExpType().
Nektar::MultiRegions::ExpList::ExpList | ( | const LibUtilities::SessionReaderSharedPtr & | pSession, |
const SpatialDomains::MeshGraphSharedPtr & | pGraph | ||
) |
The default constructor.
Creates an empty expansion list. The expansion list will typically be populated by a derived class (namely one of MultiRegions::ExpList1D, MultiRegions::ExpList2D or MultiRegions::ExpList3D).
Definition at line 135 of file ExpList.cpp.
References Nektar::MultiRegions::eNoType, and SetExpType().
Nektar::MultiRegions::ExpList::ExpList | ( | const ExpList & | in, |
const bool | DeclareCoeffPhysArrays = true |
||
) |
The copy constructor.
Copies an existing expansion list.
in | Source expansion list. |
Definition at line 162 of file ExpList.cpp.
References Nektar::MultiRegions::eNoType, m_coeffs, m_ncoeffs, m_npoints, m_phys, and SetExpType().
|
virtual |
|
inline |
Definition at line 1943 of file ExpList.h.
References v_AddFwdBwdTraceIntegral().
|
inline |
Definition at line 1928 of file ExpList.h.
References v_AddTraceIntegral().
|
inline |
Definition at line 1936 of file ExpList.h.
References v_AddTraceIntegral().
|
inline |
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Definition at line 766 of file ExpList.h.
References v_AppendFieldData().
|
inline |
Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata.
Definition at line 776 of file ExpList.h.
References v_AppendFieldData().
void Nektar::MultiRegions::ExpList::ApplyGeomInfo | ( | ) |
Apply geometry information to each expansion.
Configures geometric info, such as tangent direction, on each expansion.
graph2D | Mesh |
Definition at line 1349 of file ExpList.cpp.
Referenced by Nektar::MultiRegions::DisContField1D::DisContField1D(), and Nektar::MultiRegions::DisContField3D::DisContField3D().
|
inline |
Reimplemented in Nektar::MultiRegions::ContField2D, and Nektar::MultiRegions::ContField1D.
Definition at line 1522 of file ExpList.h.
References v_BwdTrans().
Referenced by Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), and v_WriteTecplotField().
|
inline |
This function elementally evaluates the backward transformation of the global spectral/hp element expansion.
Definition at line 1533 of file ExpList.h.
References v_BwdTrans_IterPerExp().
Referenced by Nektar::MultiRegions::ContField1D::BwdTrans(), Nektar::MultiRegions::ContField2D::BwdTrans(), and Nektar::MultiRegions::ContField3D::v_BwdTrans().
|
inline |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1664 of file ExpList.h.
References v_DealiasedProd().
|
inline |
Evaulates the maximum number of modes in the elemental basis order over all elements.
Evaulates the maximum number of modes in the elemental basis order over all elements
Definition at line 1330 of file ExpList.h.
References m_exp.
Referenced by EvalBasisNumModesMaxPerExp().
|
inline |
Returns the vector of the number of modes in the elemental basis order over all elements.
Definition at line 1347 of file ExpList.h.
References EvalBasisNumModesMax(), and m_exp.
|
inline |
Definition at line 1992 of file ExpList.h.
References v_EvaluateBoundaryConditions().
Referenced by Nektar::MultiRegions::DisContField1D::DisContField1D(), Nektar::MultiRegions::DisContField2D::DisContField2D(), and Nektar::MultiRegions::DisContField3D::DisContField3D().
void Nektar::MultiRegions::ExpList::ExtractCoeffsToCoeffs | ( | const boost::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 expansions rather than planes in homogeneous case.
Definition at line 2022 of file ExpList.cpp.
References v_ExtractCoeffsToCoeffs().
void Nektar::MultiRegions::ExpList::ExtractDataToCoeffs | ( | LibUtilities::FieldDefinitionsSharedPtr & | fielddef, |
std::vector< NekDouble > & | fielddata, | ||
std::string & | field, | ||
Array< OneD, NekDouble > & | coeffs | ||
) |
Extract the data in fielddata into the coeffs.
Definition at line 2013 of file ExpList.cpp.
References v_ExtractDataToCoeffs().
void Nektar::MultiRegions::ExpList::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 than planes in homogeneous case.
|
protected |
Definition at line 1798 of file ExpList.cpp.
References ASSERTL0, Nektar::LibUtilities::FieldIO::Import(), and m_session.
Referenced by Nektar::MultiRegions::DisContField3D::v_EvaluateBoundaryConditions(), and Nektar::MultiRegions::DisContField2D::v_EvaluateBoundaryConditions().
Definition at line 1966 of file ExpList.h.
References v_ExtractTracePhys().
|
inline |
Definition at line 1972 of file ExpList.h.
References v_ExtractTracePhys().
Fill Bnd Condition expansion from the values stored in expansion.
Definition at line 1751 of file ExpList.h.
References v_FillBndCondFromField().
|
inline |
Reimplemented in Nektar::MultiRegions::ContField2D, and Nektar::MultiRegions::ContField1D.
Definition at line 1493 of file ExpList.h.
References v_FwdTrans().
void Nektar::MultiRegions::ExpList::FwdTrans_BndConstrained | ( | const Array< OneD, const NekDouble > & | inarray, |
Array< OneD, NekDouble > & | outarray | ||
) |
Definition at line 531 of file ExpList.cpp.
References m_coeff_offset, and m_phys_offset.
|
inline |
This function elementally evaluates the forward transformation of a function onto the global spectral/hp expansion.
Definition at line 1504 of file ExpList.h.
References v_FwdTrans_IterPerExp().
|
protected |
This function assembles the block diagonal matrix of local matrices of the type mtype.
This function assembles the block diagonal matrix , which is the concatenation of the local matrices of the type mtype, that is
mtype | the type of matrix to be assembled |
scalar | an optional parameter |
constant | an optional parameter |
Definition at line 595 of file ExpList.cpp.
References Nektar::StdRegions::eBwdTrans, Nektar::eDIAGONAL, ErrorUtil::efatal, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eMass, Nektar::LibUtilities::eNoShapeType, Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetShapeType(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_exp, m_offset_elmt_id, m_phys_offset, and NEKERROR.
Referenced by GetBlockMatrix().
void Nektar::MultiRegions::ExpList::GeneralGetFieldDefinitions | ( | std::vector< LibUtilities::FieldDefinitionsSharedPtr > & | fielddef, |
int | NumHomoDir = 0 , |
||
Array< OneD, LibUtilities::BasisSharedPtr > & | HomoBasis = LibUtilities::NullBasisSharedPtr1DArray , |
||
std::vector< NekDouble > & | HomoLen = LibUtilities::NullNekDoubleVector , |
||
std::vector< unsigned int > & | HomoZIDs = LibUtilities::NullUnsignedIntVector , |
||
std::vector< unsigned int > & | HomoYIDs = LibUtilities::NullUnsignedIntVector |
||
) |
Definition at line 1873 of file ExpList.cpp.
References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, and m_exp.
Referenced by v_GetFieldDefinitions().
|
inline |
This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray.
This operation is equivalent to the evaluation of , that is,
where are the local matrices of type specified by the key mkey. The decoupling of the local matrices allows for a local evaluation of the operation. However, rather than a local matrix-vector multiplication, the local operations are evaluated as implemented in the function StdRegions::StdExpansion::GeneralMatrixOp.
mkey | This key uniquely defines the type matrix required for the operation. |
inarray | The vector of size . |
outarray | The resulting vector of size . |
Reimplemented in Nektar::MultiRegions::ContField1D.
Definition at line 2030 of file ExpList.h.
References v_GeneralMatrixOp().
Referenced by Nektar::MultiRegions::ContField3D::GenerateDirBndCondForcing().
void Nektar::MultiRegions::ExpList::GeneralMatrixOp_IterPerExp | ( | const GlobalMatrixKey & | gkey, |
const Array< OneD, const NekDouble > & | inarray, | ||
Array< OneD, NekDouble > & | outarray | ||
) |
Definition at line 746 of file ExpList.cpp.
References Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_coeff_offset, m_globalOptParam, m_offset_elmt_id, m_phys_offset, and MultiplyByBlockMatrix().
Referenced by Nektar::MultiRegions::ContField3D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField1D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField2D::v_GeneralMatrixOp(), and v_GeneralMatrixOp().
|
protected |
Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocToGloBaseMap.
Definition at line 1111 of file ExpList.cpp.
References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::MultiRegions::eSIZE_GlobalSysSolnType, Nektar::MultiRegions::GetGlobalLinSysFactory(), Nektar::MultiRegions::GlobalLinSysKey::GetGlobalSysSolnType(), GetRobinBCInfo(), GetSharedThisPtr(), and Nektar::MultiRegions::GlobalSysSolnTypeMap.
Referenced by Nektar::MultiRegions::DisContField3D::GetGlobalBndLinSys(), Nektar::MultiRegions::DisContField2D::GetGlobalBndLinSys(), and Nektar::MultiRegions::DisContField1D::GetGlobalBndLinSys().
|
protected |
This operation constructs the global linear system of type mkey.
Consider a linear system to be solved. Dependent on the solution method, this function constructs
mkey | A key which uniquely defines the global matrix to be constructed. |
locToGloMap | Contains the mapping array and required information for the transformation from local to global degrees of freedom. |
Definition at line 1092 of file ExpList.cpp.
References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::MultiRegions::eSIZE_GlobalSysSolnType, Nektar::MultiRegions::GetGlobalLinSysFactory(), Nektar::MultiRegions::GlobalLinSysKey::GetGlobalSysSolnType(), GetSharedThisPtr(), and Nektar::MultiRegions::GlobalSysSolnTypeMap.
|
protected |
Generates a global matrix from the given key and map.
Retrieves local matrices from each expansion in the expansion list and combines them together to generate a global matrix system.
mkey | Matrix key for the matrix to be generated. |
locToGloMap | Local to global mapping. |
Definition at line 814 of file ExpList.cpp.
References Nektar::StdRegions::eBwdTrans, ErrorUtil::efatal, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eMass, Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_npoints, m_offset_elmt_id, m_phys_offset, m_session, and NEKERROR.
Referenced by Nektar::MultiRegions::ContField3D::GetGlobalMatrix(), and Nektar::MultiRegions::ContField2D::GetGlobalMatrix().
|
protected |
Definition at line 951 of file ExpList.cpp.
References Nektar::eFULL, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eLaplacian, Nektar::ePOSITIVE_DEFINITE_SYMMETRIC, Nektar::ePOSITIVE_DEFINITE_SYMMETRIC_BANDED, Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), GetRobinBCInfo(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_coeff_offset, m_offset_elmt_id, and m_phys_offset.
Referenced by Nektar::MultiRegions::ContField2D::LinearAdvectionEigs().
|
inline |
Returns the total number of qudature points scaled by the factor scale on each 1D direction.
Definition at line 1377 of file ExpList.h.
References m_exp.
|
inline |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D.
Definition at line 1676 of file ExpList.h.
References v_GetBCValues().
|
protected |
Definition at line 731 of file ExpList.cpp.
References GenBlockMatrix(), Nektar::iterator, and m_blockMat.
Referenced by MultiplyByBlockMatrix(), MultiplyByElmtInvMass(), Nektar::MultiRegions::DisContField3D::v_GeneralMatrixOp(), Nektar::MultiRegions::DisContField2D::v_GeneralMatrixOp(), Nektar::MultiRegions::DisContField3D::v_HelmSolve(), Nektar::MultiRegions::DisContField2D::v_HelmSolve(), and Nektar::MultiRegions::DisContField1D::v_HelmSolve().
|
inline |
Reimplemented in Nektar::MultiRegions::ContField2D, Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::ContField1D, Nektar::MultiRegions::DisContField3DHomogeneous2D, and Nektar::MultiRegions::ContField3D.
Definition at line 1879 of file ExpList.h.
References v_GetBndCondExpansions().
|
inline |
Reimplemented in Nektar::MultiRegions::ContField2D, Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::ContField1D, and Nektar::MultiRegions::DisContField3DHomogeneous2D.
Definition at line 1980 of file ExpList.h.
References v_GetBndConditions().
|
staticprotected |
Definition at line 2544 of file ExpList.cpp.
References ASSERTL1.
Referenced by Nektar::MultiRegions::DisContField2D::FindPeriodicEdges(), Nektar::MultiRegions::DisContField3D::FindPeriodicFaces(), Nektar::MultiRegions::DisContField1D::FindPeriodicVertices(), Nektar::MultiRegions::DisContField3D::GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::DisContField1D::GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::DisContField2D::GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::DisContField1D::SetBoundaryConditionExpansion(), and Nektar::MultiRegions::DisContField3DHomogeneous1D::SetupBoundaryConditions().
|
inline |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, and Nektar::MultiRegions::DisContField3DHomogeneous2D.
Definition at line 2045 of file ExpList.h.
References v_GetBoundaryToElmtMap().
Referenced by Nektar::MultiRegions::DisContField2D::DisContField2D(), Nektar::MultiRegions::DisContField3D::DisContField3D(), Nektar::MultiRegions::DisContField1D::v_GetRobinBCInfo(), Nektar::MultiRegions::DisContField3D::v_GetRobinBCInfo(), and Nektar::MultiRegions::DisContField2D::v_GetRobinBCInfo().
|
inline |
|
inline |
Get the start offset position for a global list of m_coeffs correspoinding to element n.
Definition at line 1832 of file ExpList.h.
References m_coeff_offset.
Referenced by Nektar::MultiRegions::DisContField3DHomogeneous1D::NormVectorIProductWRTBase(), Nektar::MultiRegions::AssemblyMapCG1D::SetUp1DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG::SetUpUniversalC0ContMap(), Nektar::MultiRegions::DisContField3D::v_AddFwdBwdTraceIntegral(), Nektar::MultiRegions::DisContField2D::v_AddFwdBwdTraceIntegral(), Nektar::MultiRegions::DisContField3D::v_AddTraceIntegral(), Nektar::MultiRegions::DisContField1D::v_AddTraceIntegral(), Nektar::MultiRegions::DisContField2D::v_AddTraceIntegral(), and Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy().
|
inline |
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expansion coefficients.
As the function returns a constant reference to a const Array, it is not possible to modify the underlying data of the array m_coeffs. In order to do so, use the function UpdateCoeffs instead.
Definition at line 1740 of file ExpList.h.
References m_coeffs.
Referenced by Nektar::MultiRegions::DisContField3DHomogeneous2D::EvaluateBoundaryConditions(), Nektar::MultiRegions::DisContField3DHomogeneous1D::EvaluateBoundaryConditions(), Nektar::MultiRegions::ExpList1D::PostProcess(), and Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy().
|
inline |
|
inline |
|
inline |
This function returns the dimension of the coordinates of the element eid.
eid | The index of the element to be checked. |
Definition at line 1700 of file ExpList.h.
References ASSERTL2, and m_exp.
Referenced by Nektar::MultiRegions::ExpList2DHomogeneous1D::GetCoords(), GetExp(), GetExpIndex(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_GetCoords(), v_GetCoords(), Nektar::MultiRegions::ExpList2D::v_GetNormals(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_GetNormals(), Nektar::MultiRegions::ExpList1D::v_GetNormals(), Nektar::MultiRegions::ExpList1D::v_Upwind(), Nektar::MultiRegions::ExpList1DHomogeneous2D::v_WriteTecplotZone(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_WriteTecplotZone(), and v_WriteTecplotZone().
|
inline |
This function calculates the coordinates of all the elemental quadrature points .
Reimplemented in Nektar::MultiRegions::ExpList3DHomogeneous1D, Nektar::MultiRegions::ExpList2DHomogeneous1D, Nektar::MultiRegions::ExpList3DHomogeneous2D, and Nektar::MultiRegions::ExpList1DHomogeneous2D.
Definition at line 1597 of file ExpList.h.
References v_GetCoords().
Referenced by v_WriteTecplotZone().
|
inline |
This function returns the vector of elements in the expansion.
Definition at line 1823 of file ExpList.h.
References m_exp.
Referenced by Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::ContField3DHomogeneous2D::ContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::DisContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D(), Nektar::MultiRegions::ExpList2DHomogeneous1D::ExpList2DHomogeneous1D(), Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::GetBCValues(), GetExpIndex(), Nektar::MultiRegions::DisContField3DHomogeneous1D::NormVectorIProductWRTBase(), Nektar::MultiRegions::ExpList1D::PostProcess(), Nektar::MultiRegions::AssemblyMapCG1D::SetUp1DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DGraphC0ContMap(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG::SetUpUniversalC0ContMap(), Nektar::MultiRegions::AssemblyMapDG::SetUpUniversalDGMap(), Nektar::MultiRegions::AssemblyMapDG::SetUpUniversalTraceMap(), v_GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy(), Nektar::MultiRegions::AssemblyMapCG::v_LinearSpaceMap(), and v_WriteTecplotHeader().
|
inline |
LocalRegions::ExpansionSharedPtr & Nektar::MultiRegions::ExpList::GetExp | ( | const Array< OneD, const NekDouble > & | gloCoord | ) |
This function returns (a shared pointer to) the local elemental expansion containing the arbitrary point given by gloCoord.
Definition at line 1188 of file ExpList.cpp.
References ASSERTL0, GetCoordim(), and m_exp.
int Nektar::MultiRegions::ExpList::GetExpIndex | ( | const Array< OneD, const NekDouble > & | gloCoord, |
NekDouble | tol = 0.0 , |
||
bool | returnNearestElmt = false |
||
) |
This function returns the index of the local elemental expansion containing the arbitrary point given by gloCoord.
Definition at line 1209 of file ExpList.cpp.
int Nektar::MultiRegions::ExpList::GetExpIndex | ( | const Array< OneD, const NekDouble > & | gloCoords, |
Array< OneD, NekDouble > & | locCoords, | ||
NekDouble | tol = 0.0 , |
||
bool | returnNearestElmt = false |
||
) |
This function returns the index and the Local Cartesian Coordinates locCoords of the local elemental expansion containing the arbitrary point given by gloCoords.
Definition at line 1220 of file ExpList.cpp.
References GetCoordim(), GetExp(), GetNumElmts(), m_exp, m_graph, Nektar::NekPoint< data_type >::SetX(), Nektar::NekPoint< data_type >::SetY(), Nektar::NekPoint< data_type >::SetZ(), Vmath::Vcopy(), and WARNINGL0.
|
inline |
This function returns the number of elements in the expansion.
Definition at line 1802 of file ExpList.h.
Referenced by Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::ContField3DHomogeneous2D::ContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::DisContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3D::EvaluateHDGPostProcessing(), Nektar::MultiRegions::DisContField2D::EvaluateHDGPostProcessing(), Nektar::MultiRegions::ExpList0D::ExpList0D(), Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), Nektar::MultiRegions::ExpList2DHomogeneous1D::ExpList2DHomogeneous1D(), Nektar::MultiRegions::ExpList3D::ExpList3D(), Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D(), Nektar::MultiRegions::DisContField1D::GetNegatedFluxNormal(), Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), Nektar::MultiRegions::ExpList1D::PeriodicEval(), Nektar::MultiRegions::DisContField3D::SameTypeOfBoundaryConditions(), Nektar::MultiRegions::DisContField2D::SameTypeOfBoundaryConditions(), Nektar::MultiRegions::DisContField3D::v_AddFwdBwdTraceIntegral(), Nektar::MultiRegions::DisContField2D::v_AddFwdBwdTraceIntegral(), Nektar::MultiRegions::DisContField3D::v_AddTraceIntegral(), Nektar::MultiRegions::DisContField1D::v_AddTraceIntegral(), Nektar::MultiRegions::DisContField2D::v_AddTraceIntegral(), Nektar::MultiRegions::DisContField3D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField2D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField1D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField3D::v_GetBoundaryToElmtMap(), Nektar::MultiRegions::DisContField2D::v_GetBoundaryToElmtMap(), Nektar::MultiRegions::DisContField1D::v_GetBoundaryToElmtMap(), Nektar::MultiRegions::DisContField3D::v_GetFwdBwdTracePhys(), Nektar::MultiRegions::DisContField2D::v_GetFwdBwdTracePhys(), Nektar::MultiRegions::DisContField1D::v_GetFwdBwdTracePhys(), Nektar::MultiRegions::DisContField3D::v_HelmSolve(), Nektar::MultiRegions::DisContField2D::v_HelmSolve(), Nektar::MultiRegions::DisContField1D::v_HelmSolve(), Nektar::MultiRegions::ExpList3D::v_PhysGalerkinProjection1DScaled(), Nektar::MultiRegions::ExpList2D::v_PhysGalerkinProjection1DScaled(), Nektar::MultiRegions::ExpList3D::v_PhysInterp1DScaled(), Nektar::MultiRegions::ExpList2D::v_PhysInterp1DScaled(), Nektar::MultiRegions::ExpList3D::v_ReadGlobalOptimizationParameters(), and Nektar::MultiRegions::ExpList2D::v_ReadGlobalOptimizationParameters().
ExpansionType Nektar::MultiRegions::ExpList::GetExpType | ( | void | ) |
Returns the type of the expansion.
Definition at line 189 of file ExpList.cpp.
References m_expType.
|
inline |
Definition at line 751 of file ExpList.h.
References v_GetFieldDefinitions().
|
inline |
Definition at line 757 of file ExpList.h.
References v_GetFieldDefinitions().
|
inline |
Definition at line 1951 of file ExpList.h.
References v_GetFwdBwdTracePhys().
|
inline |
Definition at line 1958 of file ExpList.h.
References v_GetFwdBwdTracePhys().
|
inline |
|
inline |
|
inline |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 839 of file ExpList.h.
References v_GetHomogeneousBasis().
This function returns the Width of homogeneous direction associaed with the homogeneous expansion.
Definition at line 503 of file ExpList.h.
References v_GetHomoLen().
|
inline |
Returns the total number of local degrees of freedom .
Definition at line 1316 of file ExpList.h.
References m_ncoeffs.
Referenced by Nektar::MultiRegions::ContField3D::GenerateDirBndCondForcing(), Nektar::MultiRegions::ContField2D::LaplaceSolve(), Nektar::MultiRegions::ContField3D::v_FillBndCondFromField(), Nektar::MultiRegions::ContField2D::v_FillBndCondFromField(), Nektar::MultiRegions::DisContField3D::v_HelmSolve(), Nektar::MultiRegions::ContField3D::v_HelmSolve(), Nektar::MultiRegions::DisContField2D::v_HelmSolve(), Nektar::MultiRegions::ContField2D::v_HelmSolve(), Nektar::MultiRegions::ContField3D::v_ImposeDirichletConditions(), Nektar::MultiRegions::ContField2D::v_ImposeDirichletConditions(), and Nektar::MultiRegions::ContField2D::v_LinearAdvectionDiffusionReactionSolve().
|
inline |
|
inline |
|
inline |
This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp.
Definition at line 548 of file ExpList.h.
References v_GetNumElmts().
Referenced by GetExpIndex().
|
inline |
Get the element id associated with the n th consecutive block of data in m_phys and m_coeffs.
Definition at line 1848 of file ExpList.h.
References m_offset_elmt_id.
Referenced by Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), Nektar::MultiRegions::AssemblyMapCG1D::SetUp1DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DGraphC0ContMap(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and Nektar::MultiRegions::AssemblyMapDG::SetUpUniversalDGMap().
|
inline |
Definition at line 742 of file ExpList.h.
References v_GetPeriodicEntities().
This function returns (a reference to) the array (implemented as m_phys) containing the function evaluated at the quadrature points.
As the function returns a constant reference to a const Array it is not possible to modify the underlying data of the array m_phys. In order to do so, use the function UpdatePhys instead.
Definition at line 1793 of file ExpList.h.
References m_phys.
Referenced by Nektar::MultiRegions::DisContField3D::v_GetFwdBwdTracePhys(), Nektar::MultiRegions::DisContField2D::v_GetFwdBwdTracePhys(), and Nektar::MultiRegions::DisContField1D::v_GetFwdBwdTracePhys().
|
inline |
Get the start offset position for a global list of m_phys correspoinding to element n.
Definition at line 1840 of file ExpList.h.
References m_phys_offset.
Referenced by Nektar::MultiRegions::DisContField3DHomogeneous1D::GetBCValues(), Nektar::MultiRegions::DisContField3DHomogeneous1D::NormVectorIProductWRTBase(), Nektar::MultiRegions::DisContField3D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField2D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField1D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField3D::v_GetFwdBwdTracePhys(), Nektar::MultiRegions::DisContField2D::v_GetFwdBwdTracePhys(), and Nektar::MultiRegions::DisContField1D::v_GetFwdBwdTracePhys().
|
inline |
This function indicates whether the array of physical values (implemented as m_phys) is filled or not.
Definition at line 1464 of file ExpList.h.
References m_physState.
|
inline |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 844 of file ExpList.h.
References v_GetPlane().
|
inline |
Definition at line 737 of file ExpList.h.
References v_GetRobinBCInfo().
Referenced by GenGlobalBndLinSys(), and GenGlobalMatrixFull().
|
inline |
|
inline |
Returns a shared pointer to the current object.
Definition at line 816 of file ExpList.h.
Referenced by GenGlobalBndLinSys(), and GenGlobalLinSys().
|
inline |
Returns the total number of quadrature points m_npoints .
Definition at line 1366 of file ExpList.h.
References m_npoints.
Referenced by Nektar::MultiRegions::ExpList3DHomogeneous2D::GetCoords(), Nektar::MultiRegions::ExpList2DHomogeneous1D::GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy(), Nektar::MultiRegions::ExpList3DHomogeneous2D::v_L2(), v_WriteTecplotField(), and v_WriteTecplotZone().
|
inline |
|
inline |
Return a reference to the trace space associated with this expansion list.
Definition at line 1907 of file ExpList.h.
References v_GetTrace().
Referenced by Nektar::MultiRegions::DisContField2D::v_AddFwdBwdTraceIntegral(), Nektar::MultiRegions::DisContField1D::v_AddTraceIntegral(), and Nektar::MultiRegions::DisContField2D::v_AddTraceIntegral().
Definition at line 1917 of file ExpList.h.
References v_GetTraceBndMap().
|
inline |
Definition at line 1912 of file ExpList.h.
References v_GetTraceMap().
Referenced by v_GetTraceBndMap().
|
inline |
This function returns the transposition class associaed with the homogeneous expansion.
Definition at line 496 of file ExpList.h.
References v_GetTransposition().
|
inline |
This function returns the third direction expansion condition, which can be in wave space (coefficient) or not It is stored in the variable m_WaveSpace.
Definition at line 1415 of file ExpList.h.
References m_WaveSpace.
This function returns a vector containing the wave numbers in y-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction.
Definition at line 512 of file ExpList.h.
References v_GetYIDs().
This function returns a vector containing the wave numbers in z-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction.
Definition at line 489 of file ExpList.h.
References v_GetZIDs().
Referenced by v_WriteTecplotZone().
|
protected |
Put the coefficients into local ordering and place in m_coeffs.
Definition at line 1761 of file ExpList.h.
References v_GlobalToLocal().
Referenced by Nektar::MultiRegions::ContField1D::BwdTrans(), Nektar::MultiRegions::ContField2D::BwdTrans(), Nektar::MultiRegions::ContField1D::FwdTrans(), Nektar::MultiRegions::ContField2D::FwdTrans(), Nektar::MultiRegions::ContField2D::LaplaceSolve(), Nektar::MultiRegions::ContField1D::MultiplyByInvMassMatrix(), Nektar::MultiRegions::ContField2D::MultiplyByInvMassMatrix(), Nektar::MultiRegions::ContField3D::v_BwdTrans(), Nektar::MultiRegions::ContField3D::v_FwdTrans(), Nektar::MultiRegions::ContField3D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField1D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField2D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField3D::v_HelmSolve(), Nektar::MultiRegions::ContField1D::v_HelmSolve(), Nektar::MultiRegions::ContField2D::v_HelmSolve(), Nektar::MultiRegions::ContField2D::v_LinearAdvectionDiffusionReactionSolve(), Nektar::MultiRegions::ContField2D::v_LinearAdvectionReactionSolve(), and Nektar::MultiRegions::ContField3D::v_MultiplyByInvMassMatrix().
NekDouble Nektar::MultiRegions::ExpList::H1 | ( | const Array< OneD, const NekDouble > & | inarray, |
const Array< OneD, const NekDouble > & | soln = NullNekDouble1DArray |
||
) |
Calculates the error of the global spectral/hp element approximation.
Given a spectral/hp approximation evaluated at the quadrature points (which should be contained in m_phys), this function calculates the error of this approximation with respect to an exact solution. The local distribution of the quadrature points allows an elemental evaluation of this operation through the functions StdRegions::StdExpansion::H1.
The exact solution, also evaluated at the quadrature points, should be contained in the variable m_phys of the ExpList object Sol.
soln | An 1D array, containing the discrete evaluation of the exact solution at the quadrature points. |
Definition at line 1854 of file ExpList.cpp.
References m_comm, m_phys_offset, and Nektar::LibUtilities::ReduceSum.
|
inline |
|
inline |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1651 of file ExpList.h.
References v_HomogeneousBwdTrans().
This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expansion.
Definition at line 473 of file ExpList.h.
References v_HomogeneousEnergy().
|
inline |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1638 of file ExpList.h.
References v_HomogeneousFwdTrans().
|
inline |
Impose Dirichlet Boundary Conditions onto Array.
Definition at line 1745 of file ExpList.h.
References v_ImposeDirichletConditions().
|
inline |
Reimplemented in Nektar::MultiRegions::ContField1D, and Nektar::MultiRegions::ContField2D.
Definition at line 1472 of file ExpList.h.
References v_IProductWRTBase().
Referenced by Nektar::MultiRegions::ContField3D::v_FwdTrans(), Nektar::MultiRegions::DisContField3D::v_HelmSolve(), Nektar::MultiRegions::ContField3D::v_HelmSolve(), Nektar::MultiRegions::DisContField2D::v_HelmSolve(), and Nektar::MultiRegions::DisContField1D::v_HelmSolve().
|
inline |
This function calculates the inner product of a function with respect to all {local} expansion modes .
Definition at line 1483 of file ExpList.h.
References v_IProductWRTBase_IterPerExp().
Referenced by Nektar::MultiRegions::ContField2D::IProductWRTBase(), Nektar::MultiRegions::ContField1D::IProductWRTBase(), v_FwdTrans_IterPerExp(), and Nektar::MultiRegions::ContField3D::v_IProductWRTBase().
void Nektar::MultiRegions::ExpList::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 direction.
dir) | of all {local} expansion modes . |
The operation is evaluated locally for every element by the function StdRegions::StdExpansion::IProductWRTDerivBase.
dir | {0,1} is the direction in which the derivative of the basis should be taken |
inarray | An array of size containing the values of the function at the quadrature points . |
outarray | An array of size used to store the result. |
Definition at line 354 of file ExpList.cpp.
References m_coeff_offset, and m_phys_offset.
|
inline |
This function calculates the error with respect to soln of the global spectral/hp element approximation.
Definition at line 453 of file ExpList.h.
References v_L2().
Referenced by Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), and Nektar::MultiRegions::ExpList3DHomogeneous2D::v_L2().
|
inline |
Solve Advection Diffusion Reaction.
Definition at line 1570 of file ExpList.h.
References v_LinearAdvectionDiffusionReactionSolve().
|
inline |
Solve Advection Diffusion Reaction.
Definition at line 1582 of file ExpList.h.
References v_LinearAdvectionReactionSolve().
NekDouble Nektar::MultiRegions::ExpList::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.
Given a spectral/hp approximation evaluated at the quadrature points (which should be contained in m_phys), this function calculates the error of this approximation with respect to an exact solution. The local distribution of the quadrature points allows an elemental evaluation of this operation through the functions StdRegions::StdExpansion::Linf.
The exact solution, also evaluated at the quadrature points, should be contained in the variable m_phys of the ExpList object Sol.
soln | A 1D array, containing the discrete evaluation of the exact solution at the quadrature points in its array m_phys. |
Definition at line 1660 of file ExpList.cpp.
References m_comm, m_npoints, Nektar::NullNekDouble1DArray, Nektar::LibUtilities::ReduceMax, and Vmath::Vmax().
Put the coefficients into global ordering using m_coeffs.
Reimplemented in Nektar::MultiRegions::ContField1D.
Definition at line 1756 of file ExpList.h.
References v_LocalToGlobal().
Referenced by Nektar::MultiRegions::ContField3D::v_FillBndCondFromField(), Nektar::MultiRegions::ContField2D::v_FillBndCondFromField(), Nektar::MultiRegions::ContField3D::v_HelmSolve(), and Nektar::MultiRegions::ContField2D::v_HelmSolve().
|
protected |
Retrieves the block matrix specified by bkey, and computes .
gkey | GlobalMatrixKey specifying the block matrix to use in the matrix-vector multiply. |
inarray | Input vector . |
outarray | Output vector . |
Definition at line 270 of file ExpList.cpp.
References Nektar::eWrapper, and GetBlockMatrix().
Referenced by GeneralMatrixOp_IterPerExp(), v_BwdTrans_IterPerExp(), and v_IProductWRTBase_IterPerExp().
void Nektar::MultiRegions::ExpList::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 matrix.
The coefficients of the function to be acted upon should be contained in the
inarray. | The resulting coefficients are stored in |
outarray | |
inarray | An array of size containing the inner product. |
Definition at line 482 of file ExpList.cpp.
References Nektar::StdRegions::eInvMass, Nektar::eWrapper, GetBlockMatrix(), and m_ncoeffs.
Referenced by v_FwdTrans_IterPerExp().
|
inline |
Reimplemented in Nektar::MultiRegions::ContField2D, and Nektar::MultiRegions::ContField1D.
Definition at line 1544 of file ExpList.h.
References v_MultiplyByInvMassMatrix().
|
inline |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D.
Definition at line 1687 of file ExpList.h.
References v_NormVectorIProductWRTBase().
|
inline |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1627 of file ExpList.h.
References v_PhysDeriv().
|
inline |
This function discretely evaluates the derivative of a function on the domain consisting of all elements of the expansion.
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1608 of file ExpList.h.
References v_PhysDeriv().
|
inline |
Definition at line 1619 of file ExpList.h.
References v_PhysDeriv().
|
inline |
This function Galerkin projects the physical space points in inarray to outarray where inarray is assumed to be defined in the expansion but where the number of points are rescaled by 1DScale.
Definition at line 533 of file ExpList.h.
References v_PhysGalerkinProjection1DScaled().
This function integrates a function over the domain consisting of all the elements of the expansion.
The integration is evaluated locally, that is
where the integration over the separate elements is done by the function StdRegions::StdExpansion::Integral, which discretely evaluates the integral using Gaussian quadrature.
Note that the array m_phys should be filled with the values of the function at the quadrature points .
Definition at line 222 of file ExpList.cpp.
References ASSERTL1, m_phys, and m_physState.
NekDouble Nektar::MultiRegions::ExpList::PhysIntegral | ( | const Array< OneD, const NekDouble > & | inarray | ) |
This function integrates a function over the domain consisting of all the elements of the expansion.
The integration is evaluated locally, that is
where the integration over the separate elements is done by the function StdRegions::StdExpansion::Integral, which discretely evaluates the integral using Gaussian quadrature.
inarray | An array of size containing the values of the function at the quadrature points . |
Definition at line 247 of file ExpList.cpp.
References m_phys_offset.
|
inline |
This function interpolates the physical space points in inarray to outarray using the same points defined in the expansion but where the number of points are rescaled by 1DScale.
Definition at line 521 of file ExpList.h.
References v_PhysInterp1DScaled().
|
inlineprotected |
Definition at line 992 of file ExpList.h.
References v_ReadGlobalOptimizationParameters().
Referenced by Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), and Nektar::MultiRegions::ExpList3D::ExpList3D().
void Nektar::MultiRegions::ExpList::SetExpType | ( | ExpansionType | Type | ) |
Returns the type of the expansion.
Definition at line 197 of file ExpList.cpp.
References m_expType.
Referenced by ExpList(), Nektar::MultiRegions::ExpList0D::ExpList0D(), Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), Nektar::MultiRegions::ExpList3D::ExpList3D(), Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D(), and Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D().
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition at line 480 of file ExpList.h.
References v_SetHomo1DSpecVanVisc().
|
inline |
Set Modified Basis for the stability analysis.
|
inline |
Fills the array m_phys.
This function fills the array , the evaluation of the expansion at the quadrature points (implemented as m_phys), with the values of the array inarray.
inarray | The array containing the values where m_phys should be filled with. |
Definition at line 1435 of file ExpList.h.
References ASSERTL0, m_npoints, m_phys, m_physState, and Vmath::Vcopy().
|
inline |
This function manually sets whether the array of physical values (implemented as m_phys) is filled or not.
physState | true (=filled) or false (=not filled). |
Definition at line 1455 of file ExpList.h.
References m_physState.
|
inline |
Definition at line 2040 of file ExpList.h.
References v_SetUpPhysNormals().
Referenced by Nektar::MultiRegions::DisContField2D::DisContField2D(), Nektar::MultiRegions::DisContField3D::DisContField3D(), Nektar::MultiRegions::DisContField3D::GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::DisContField2D::GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::DisContField3D::SetUpDG(), Nektar::MultiRegions::DisContField2D::SetUpDG(), and Nektar::MultiRegions::DisContField1D::SetUpDG().
|
inline |
Sets the wave space to the one of the possible configuration true or false.
Definition at line 1407 of file ExpList.h.
References m_WaveSpace.
Smooth a field across elements.
Definition at line 1514 of file ExpList.h.
References v_SmoothField().
|
inline |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, and Nektar::MultiRegions::DisContField3DHomogeneous2D.
Definition at line 1884 of file ExpList.h.
References v_UpdateBndCondExpansion().
|
inline |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, and Nektar::MultiRegions::DisContField3DHomogeneous2D.
Definition at line 1987 of file ExpList.h.
References v_UpdateBndConditions().
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expansion coefficients.
If one wants to get hold of the underlying data without modifying them, rather use the function GetCoeffs instead.
Definition at line 1859 of file ExpList.h.
References m_coeffs.
Referenced by Nektar::MultiRegions::DisContField3DHomogeneous2D::EvaluateBoundaryConditions(), and Nektar::MultiRegions::DisContField3DHomogeneous1D::EvaluateBoundaryConditions().
This function returns (a reference to) the array (implemented as m_phys) containing the function evaluated at the quadrature points.
If one wants to get hold of the underlying data without modifying them, rather use the function GetPhys instead.
Definition at line 1870 of file ExpList.h.
References m_phys, and m_physState.
|
inline |
Definition at line 1889 of file ExpList.h.
References v_Upwind().
|
inline |
Definition at line 1898 of file ExpList.h.
References v_Upwind().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2222 of file ExpList.cpp.
References ASSERTL0.
Referenced by AddFwdBwdTraceIntegral().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D.
Definition at line 2205 of file ExpList.cpp.
References ASSERTL0.
Referenced by AddTraceIntegral().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2214 of file ExpList.cpp.
References ASSERTL0.
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 1985 of file ExpList.cpp.
References m_coeffs.
Referenced by AppendFieldData().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 1990 of file ExpList.cpp.
References m_coeff_offset.
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ContField3D.
Definition at line 2376 of file ExpList.cpp.
References v_BwdTrans_IterPerExp().
Referenced by BwdTrans().
|
protectedvirtual |
Given the elemental coefficients of an expansion, this function evaluates the spectral/hp expansion at the quadrature points . The operation is evaluated locally by the elemental function StdRegions::StdExpansion::BwdTrans.
inarray | An array of size containing the local coefficients . |
outarray | The resulting physical values at the quadrature points will be stored in this array of size . |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 1148 of file ExpList.cpp.
References Nektar::StdRegions::eBwdTrans, m_coeff_offset, m_globalOptParam, m_offset_elmt_id, m_phys_offset, and MultiplyByBlockMatrix().
Referenced by BwdTrans_IterPerExp(), and v_BwdTrans().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 2325 of file ExpList.cpp.
References ASSERTL0.
Referenced by DealiasedProd().
|
privatevirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2513 of file ExpList.cpp.
References ASSERTL0.
Referenced by EvaluateBoundaryConditions().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 2110 of file ExpList.cpp.
References m_coeff_offset, m_ncoeffs, m_offset_elmt_id, and Vmath::Vcopy().
Referenced by ExtractCoeffsToCoeffs().
|
protectedvirtual |
Extract data from raw field data into expansion list.
fielddef | Field definitions. |
fielddata | Data for associated field. |
field | Field variable name. |
coeffs | Resulting coefficient array. |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 2035 of file ExpList.cpp.
References ASSERTL0, Nektar::LibUtilities::GetNumberOfCoefficients(), m_coeff_offset, m_exp, and Vmath::Vcopy().
Referenced by ExtractDataToCoeffs().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2247 of file ExpList.cpp.
References ASSERTL0.
Referenced by ExtractTracePhys().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField1D, Nektar::MultiRegions::DisContField2D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2253 of file ExpList.cpp.
References ASSERTL0.
Reimplemented in Nektar::MultiRegions::DisContField2D.
Definition at line 2357 of file ExpList.cpp.
References ASSERTL0.
Referenced by FillBndCondFromField().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ContField3D.
Definition at line 2383 of file ExpList.cpp.
References v_FwdTrans_IterPerExp().
Referenced by FwdTrans().
|
protectedvirtual |
Given a function defined at the quadrature points, this function determines the transformed elemental coefficients employing a discrete elemental Galerkin projection from physical space to coefficient space. For each element, the operation is evaluated locally by the function StdRegions::StdExpansion::IproductWRTBase followed by a call to #MultiRegions#MultiplyByElmtInvMass.
inarray | An array of size containing the values of the function at the quadrature points . |
outarray | The resulting coefficients will be stored in this array of size . |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 521 of file ExpList.cpp.
References IProductWRTBase_IterPerExp(), m_ncoeffs, and MultiplyByElmtInvMass().
Referenced by FwdTrans_IterPerExp(), and v_FwdTrans().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2398 of file ExpList.cpp.
References GeneralMatrixOp_IterPerExp().
Referenced by GeneralMatrixOp().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D.
Definition at line 2332 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetBCValues().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2141 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetBndCondExpansions().
|
privatevirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2492 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetBndConditions().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField1D, Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField3D, and Nektar::MultiRegions::DisContField3DHomogeneous2D.
Definition at line 2474 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetBoundaryToElmtMap().
|
protectedvirtual |
The operation is evaluated locally by the elemental function StdRegions::StdExpansion::GetCoords.
coord_0 | After calculation, the coordinate will be stored in this array. |
coord_1 | After calculation, the coordinate will be stored in this array. |
coord_2 | After calculation, the coordinate will be stored in this array. |
Reimplemented in Nektar::MultiRegions::ExpList3DHomogeneous1D, Nektar::MultiRegions::ExpList2DHomogeneous1D, Nektar::MultiRegions::ExpList3DHomogeneous2D, and Nektar::MultiRegions::ExpList1DHomogeneous2D.
Definition at line 2418 of file ExpList.cpp.
References ASSERTL0, GetCoordim(), GetExp(), and m_phys_offset.
Referenced by GetCoords().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 1971 of file ExpList.cpp.
Referenced by GetFieldDefinitions().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 1978 of file ExpList.cpp.
References GeneralGetFieldDefinitions().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField1D, Nektar::MultiRegions::DisContField2D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2231 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetFwdBwdTracePhys().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField1D, Nektar::MultiRegions::DisContField2D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2238 of file ExpList.cpp.
References ASSERTL0.
|
inlineprivatevirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1285 of file ExpList.h.
References ASSERTL0, and Nektar::LibUtilities::NullBasisSharedPtr.
Referenced by GetHomogeneousBasis().
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1760 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetHomoLen().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpList1D, Nektar::MultiRegions::ExpList2DHomogeneous1D, Nektar::MultiRegions::ExpList2D, and Nektar::MultiRegions::ExpList0D.
Definition at line 2198 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetNormals().
|
inlineprotectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 999 of file ExpList.h.
Referenced by GetNumElmts().
|
privatevirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField3D, and Nektar::MultiRegions::ExpList3DHomogeneous1D.
Definition at line 2535 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetPeriodicEntities().
|
privatevirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 2558 of file ExpList.cpp.
References ASSERTL0, and Nektar::MultiRegions::NullExpListSharedPtr.
Referenced by GetPlane().
|
privatevirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField3D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3DHomogeneous2D.
Definition at line 2525 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetRobinBCInfo().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField3D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField2D.
Definition at line 2177 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetTrace().
|
protectedvirtual |
Definition at line 2193 of file ExpList.cpp.
References GetTraceMap().
Referenced by GetTraceBndMap().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D, Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField3D, and Nektar::MultiRegions::DisContField1D.
Definition at line 2185 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetTraceMap().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1751 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetTransposition().
|
protectedvirtual |
Definition at line 1777 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetYIDs().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1768 of file ExpList.cpp.
References ASSERTL0.
Referenced by GetZIDs().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField1D, Nektar::MultiRegions::DisContField2D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2270 of file ExpList.cpp.
References ASSERTL0.
Referenced by HelmSolve().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 2315 of file ExpList.cpp.
References ASSERTL0.
Referenced by HomogeneousBwdTrans().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpList3DHomogeneous1D.
Definition at line 1743 of file ExpList.cpp.
References ASSERTL0.
Referenced by HomogeneousEnergy().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 2305 of file ExpList.cpp.
References ASSERTL0.
Referenced by HomogeneousFwdTrans().
|
protectedvirtual |
Definition at line 2349 of file ExpList.cpp.
References ASSERTL0.
Referenced by ImposeDirichletConditions().
|
protectedvirtual |
Definition at line 1729 of file ExpList.cpp.
References m_comm, m_offset_elmt_id, m_phys_offset, and Nektar::LibUtilities::ReduceSum.
Referenced by Integral().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, Nektar::MultiRegions::ExpListHomogeneous2D, and Nektar::MultiRegions::ContField3D.
Definition at line 2390 of file ExpList.cpp.
References v_IProductWRTBase_IterPerExp().
Referenced by IProductWRTBase().
|
protectedvirtual |
The operation is evaluated locally for every element by the function StdRegions::StdExpansion::IProductWRTBase.
inarray | An array of size containing the values of the function at the quadrature points . |
outarray | An array of size used to store the result. |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 300 of file ExpList.cpp.
References Nektar::StdRegions::eIProductWRTBase, m_coeff_offset, m_globalOptParam, m_offset_elmt_id, m_phys_offset, and MultiplyByBlockMatrix().
Referenced by IProductWRTBase_IterPerExp(), and v_IProductWRTBase().
|
protectedvirtual |
Given a spectral/hp approximation evaluated at the quadrature points (which should be contained in m_phys), this function calculates the error of this approximation with respect to an exact solution. The local distribution of the quadrature points allows an elemental evaluation of this operation through the functions StdRegions::StdExpansion::L2.
The exact solution, also evaluated at the quadrature points, should be contained in the variable m_phys of the ExpList object Sol.
Sol | An ExpList, containing the discrete evaluation of the exact solution at the quadrature points in its array m_phys. |
Reimplemented in Nektar::MultiRegions::ExpList3DHomogeneous1D, and Nektar::MultiRegions::ExpList3DHomogeneous2D.
Definition at line 1699 of file ExpList.cpp.
References m_comm, m_phys_offset, Nektar::NullNekDouble1DArray, and Nektar::LibUtilities::ReduceSum.
Referenced by L2().
|
protectedvirtual |
Definition at line 2281 of file ExpList.cpp.
References ASSERTL0.
Referenced by LinearAdvectionDiffusionReactionSolve().
|
protectedvirtual |
Definition at line 2293 of file ExpList.cpp.
References ASSERTL0.
Referenced by LinearAdvectionReactionSolve().
|
protectedvirtual |
Definition at line 2261 of file ExpList.cpp.
References ASSERTL0.
Referenced by MultiplyByInvMassMatrix().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField3DHomogeneous1D.
Definition at line 2340 of file ExpList.cpp.
References ASSERTL0.
Referenced by NormVectorIProductWRTBase().
|
protectedvirtual |
Given a function evaluated at the quadrature points, this function calculates the derivatives , and of the function at the same quadrature points. The local distribution of the quadrature points allows an elemental evaluation of the derivative. This is done by a call to the function StdRegions::StdExpansion::PhysDeriv.
inarray | An array of size containing the values of the function at the quadrature points . |
out_d0 | The discrete evaluation of the derivative will be stored in this array of size . |
out_d1 | The discrete evaluation of the derivative will be stored in this array of size . Note that if no memory is allocated for out_d1, the derivative will not be calculated. |
out_d2 | The discrete evaluation of the derivative will be stored in this array of size . Note that if no memory is allocated for out_d2, the derivative will not be calculated. |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 403 of file ExpList.cpp.
References m_phys_offset.
Referenced by PhysDeriv(), and v_PhysDeriv().
|
protectedvirtual |
Definition at line 429 of file ExpList.cpp.
References Nektar::MultiRegions::DirCartesianMap, and v_PhysDeriv().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D, and Nektar::MultiRegions::ExpListHomogeneous2D.
Definition at line 437 of file ExpList.cpp.
References Nektar::MultiRegions::eN, Nektar::MultiRegions::eS, and m_phys_offset.
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1793 of file ExpList.cpp.
References ASSERTL0.
Referenced by PhysGalerkinProjection1DScaled().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1787 of file ExpList.cpp.
References ASSERTL0.
Referenced by PhysInterp1DScaled().
|
protectedvirtual |
Definition at line 2483 of file ExpList.cpp.
References ASSERTL0.
Referenced by ReadGlobalOptimizationParameters().
|
inlineprivatevirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1293 of file ExpList.h.
References ASSERTL0.
Referenced by SetHomo1DSpecVanVisc().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpList3D.
Definition at line 2466 of file ExpList.cpp.
References ASSERTL0.
Referenced by SetUpPhysNormals().
|
protectedvirtual |
This function smooth a field after some calculaitons which have been done elementally.
field | An array containing the field in physical space |
Reimplemented in Nektar::MultiRegions::ContField3DHomogeneous1D.
Definition at line 553 of file ExpList.cpp.
References ASSERTL0, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, and m_exp.
Referenced by SmoothField().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2149 of file ExpList.cpp.
References ASSERTL0.
Referenced by UpdateBndCondExpansion().
|
privatevirtual |
Reimplemented in Nektar::MultiRegions::DisContField2D, Nektar::MultiRegions::DisContField1D, and Nektar::MultiRegions::DisContField3D.
Definition at line 2503 of file ExpList.cpp.
References ASSERTL0.
Referenced by UpdateBndConditions().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpList1D.
Definition at line 2157 of file ExpList.cpp.
References ASSERTL0.
Referenced by Upwind().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpList1D, Nektar::MultiRegions::ExpList2D, and Nektar::MultiRegions::ExpList0D.
Definition at line 2167 of file ExpList.cpp.
References ASSERTL0.
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpList3DHomogeneous1D.
Definition at line 1487 of file ExpList.cpp.
References ASSERTL0, and m_exp.
Referenced by WriteTecplotConnectivity().
|
protectedvirtual |
Write Tecplot Files Field
outfile | Output file name. |
expansion | Expansion that is considered |
Definition at line 1564 of file ExpList.cpp.
References BwdTrans(), GetTotPoints(), m_coeffs, m_phys, m_phys_offset, and m_physState.
Referenced by WriteTecplotField().
|
protectedvirtual |
Write Tecplot Files Header
outfile | Output file name. |
var | variables names |
Definition at line 1359 of file ExpList.cpp.
References Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, GetExp(), and m_expType.
Referenced by WriteTecplotHeader().
|
protectedvirtual |
Write Tecplot Files Zone
outfile | Output file name. |
expansion | Expansion that is considered |
Reimplemented in Nektar::MultiRegions::ExpList2DHomogeneous1D, Nektar::MultiRegions::ExpList3DHomogeneous2D, and Nektar::MultiRegions::ExpList1DHomogeneous2D.
Definition at line 1393 of file ExpList.cpp.
References ASSERTL0, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, GetCoordim(), GetCoords(), GetTotPoints(), GetZIDs(), m_exp, and m_expType.
Referenced by WriteTecplotZone().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpListHomogeneous1D.
Definition at line 1623 of file ExpList.cpp.
References Nektar::NekConstants::kNekZeroTol, m_phys, and m_phys_offset.
Referenced by WriteVtkPieceData().
|
protectedvirtual |
Reimplemented in Nektar::MultiRegions::ExpList2DHomogeneous1D, Nektar::MultiRegions::ExpList3DHomogeneous1D, Nektar::MultiRegions::ExpList3DHomogeneous2D, and Nektar::MultiRegions::ExpList1DHomogeneous2D.
Definition at line 1612 of file ExpList.cpp.
References ASSERTL0.
Referenced by WriteVtkPieceHeader().
|
inline |
Definition at line 372 of file ExpList.h.
References v_WriteTecplotConnectivity().
|
inline |
Definition at line 366 of file ExpList.h.
References v_WriteTecplotField().
|
inline |
Definition at line 353 of file ExpList.h.
References v_WriteTecplotHeader().
|
inline |
Definition at line 359 of file ExpList.h.
References v_WriteTecplotZone().
void Nektar::MultiRegions::ExpList::WriteVtkFooter | ( | std::ofstream & | outfile | ) |
Definition at line 1606 of file ExpList.cpp.
void Nektar::MultiRegions::ExpList::WriteVtkHeader | ( | std::ofstream & | outfile | ) |
Definition at line 1598 of file ExpList.cpp.
|
inline |
Definition at line 390 of file ExpList.h.
References v_WriteVtkPieceData().
void Nektar::MultiRegions::ExpList::WriteVtkPieceFooter | ( | std::ofstream & | outfile, |
int | expansion | ||
) |
Definition at line 1617 of file ExpList.cpp.
|
inline |
Definition at line 381 of file ExpList.h.
References v_WriteVtkPieceHeader().
|
protected |
Definition at line 946 of file ExpList.h.
Referenced by GetBlockMatrix().
|
protected |
Offset of elemental data into the array m_coeffs.
Definition at line 931 of file ExpList.h.
Referenced by Nektar::MultiRegions::DisContField3D::EvaluateHDGPostProcessing(), Nektar::MultiRegions::DisContField2D::EvaluateHDGPostProcessing(), FwdTrans_BndConstrained(), GeneralMatrixOp_IterPerExp(), GenGlobalMatrixFull(), GetCoeff_Offset(), IProductWRTDerivBase(), Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList0D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList2D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList1D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpListHomogeneous2D::v_AppendFieldData(), Nektar::MultiRegions::ExpListHomogeneous1D::v_AppendFieldData(), v_AppendFieldData(), v_BwdTrans_IterPerExp(), v_ExtractCoeffsToCoeffs(), Nektar::MultiRegions::ExpListHomogeneous2D::v_ExtractDataToCoeffs(), Nektar::MultiRegions::ExpListHomogeneous1D::v_ExtractDataToCoeffs(), v_ExtractDataToCoeffs(), Nektar::MultiRegions::DisContField1D::v_HelmSolve(), and v_IProductWRTBase_IterPerExp().
Concatenation of all local expansion coefficients.
The array of length m_ncoeffs which is the concatenation of the local expansion coefficients over all elements
Definition at line 890 of file ExpList.h.
Referenced by Nektar::MultiRegions::ContField3D::Assemble(), Nektar::MultiRegions::ContField2D::Assemble(), Nektar::MultiRegions::ContField1D::Assemble(), Nektar::MultiRegions::DisContField3D::EvaluateHDGPostProcessing(), Nektar::MultiRegions::DisContField2D::EvaluateHDGPostProcessing(), ExpList(), Nektar::MultiRegions::ExpList0D::ExpList0D(), Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), GetCoeff(), GetCoeffs(), Nektar::MultiRegions::ContField2D::GlobalToLocal(), Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), SetCoeff(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), SetCoeffs(), SetCoeffsArray(), UpdateCoeffs(), Nektar::MultiRegions::ExpListHomogeneous2D::v_AppendFieldData(), Nektar::MultiRegions::ExpListHomogeneous1D::v_AppendFieldData(), v_AppendFieldData(), Nektar::MultiRegions::ContField3D::v_FillBndCondFromField(), Nektar::MultiRegions::ContField2D::v_FillBndCondFromField(), Nektar::MultiRegions::ContField3D::v_GlobalToLocal(), Nektar::MultiRegions::ContField1D::v_GlobalToLocal(), Nektar::MultiRegions::ContField2D::v_GlobalToLocal(), Nektar::MultiRegions::ContField3D::v_LocalToGlobal(), Nektar::MultiRegions::ContField1D::v_LocalToGlobal(), Nektar::MultiRegions::ContField2D::v_LocalToGlobal(), and v_WriteTecplotField().
|
protected |
Communicator.
Definition at line 858 of file ExpList.h.
Referenced by Nektar::MultiRegions::ExpListHomogeneous1D::ExpListHomogeneous1D(), Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D(), Nektar::MultiRegions::ExpListHomogeneous1D::GenHomogeneous1DBlockMatrix(), GetComm(), H1(), Nektar::MultiRegions::ExpListHomogeneous1D::Homogeneous1DTrans(), Linf(), Nektar::MultiRegions::DisContField3D::SameTypeOfBoundaryConditions(), Nektar::MultiRegions::DisContField2D::SameTypeOfBoundaryConditions(), Nektar::MultiRegions::ExpListHomogeneous1D::v_DealiasedProd(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy(), v_Integral(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_L2(), v_L2(), and Nektar::MultiRegions::ExpListHomogeneous1D::v_PhysDeriv().
|
protected |
The list of local expansions.
The (shared pointer to the) vector containing (shared pointers to) all local expansions. The fact that the local expansions are all stored as a (pointer to a) #StdExpansion, the abstract base class for all local expansions, allows a general implementation where most of the routines for the derived classes are defined in the ExpList base class.
Definition at line 928 of file ExpList.h.
Referenced by Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::ContField3DHomogeneous2D::ContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::DisContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D(), EvalBasisNumModesMax(), EvalBasisNumModesMaxPerExp(), Nektar::MultiRegions::ExpList2DHomogeneous1D::ExpList2DHomogeneous1D(), Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D(), GenBlockMatrix(), GeneralGetFieldDefinitions(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D(), Get1DScaledTotPoints(), GetCoordim(), GetExp(), GetExpIndex(), GetNcoeffs(), GetTotPoints(), Nektar::MultiRegions::DisContField2D::IsLeftAdjacentEdge(), Nektar::MultiRegions::DisContField3D::IsLeftAdjacentFace(), Nektar::MultiRegions::DisContField1D::IsLeftAdjacentVertex(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList0D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList2D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList1D::SetCoeffPhysOffsets(), Nektar::MultiRegions::DisContField3D::SetUpDG(), Nektar::MultiRegions::DisContField2D::SetUpDG(), Nektar::MultiRegions::DisContField1D::SetUpDG(), Nektar::MultiRegions::ExpListHomogeneous1D::v_ExtractDataToCoeffs(), v_ExtractDataToCoeffs(), Nektar::MultiRegions::ExpList0D::v_GetNormals(), Nektar::MultiRegions::ExpList2D::v_GetNormals(), Nektar::MultiRegions::ExpList1D::v_GetNormals(), Nektar::MultiRegions::ExpList3D::v_PhysGalerkinProjection1DScaled(), Nektar::MultiRegions::ExpList2D::v_PhysGalerkinProjection1DScaled(), Nektar::MultiRegions::ExpList3D::v_PhysInterp1DScaled(), Nektar::MultiRegions::ExpList2D::v_PhysInterp1DScaled(), Nektar::MultiRegions::ExpList3D::v_ReadGlobalOptimizationParameters(), Nektar::MultiRegions::ExpList2D::v_ReadGlobalOptimizationParameters(), Nektar::MultiRegions::ExpList3D::v_SetUpPhysNormals(), Nektar::MultiRegions::ExpList2D::v_SetUpPhysNormals(), Nektar::MultiRegions::ExpList1D::v_SetUpPhysNormals(), v_SmoothField(), Nektar::MultiRegions::ExpList2D::v_Upwind(), Nektar::MultiRegions::ExpList1D::v_Upwind(), v_WriteTecplotConnectivity(), and v_WriteTecplotZone().
ExpansionType Nektar::MultiRegions::ExpList::m_expType |
Definition at line 850 of file ExpList.h.
Referenced by GetExpType(), SetExpType(), v_WriteTecplotHeader(), and v_WriteTecplotZone().
|
protected |
Definition at line 944 of file ExpList.h.
Referenced by Nektar::MultiRegions::ContField2D::BwdTrans(), Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::ContField3DHomogeneous2D::ContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::DisContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D(), Nektar::MultiRegions::ExpList0D::ExpList0D(), Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList1DHomogeneous2D::ExpList1DHomogeneous2D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), Nektar::MultiRegions::ExpList2DHomogeneous1D::ExpList2DHomogeneous1D(), Nektar::MultiRegions::ExpList3D::ExpList3D(), Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D(), GeneralMatrixOp_IterPerExp(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D(), GetGlobalOptParam(), Nektar::MultiRegions::ContField2D::IProductWRTBase(), Nektar::MultiRegions::ContField3D::v_BwdTrans(), v_BwdTrans_IterPerExp(), Nektar::MultiRegions::ContField3D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField2D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField3D::v_IProductWRTBase(), v_IProductWRTBase_IterPerExp(), Nektar::MultiRegions::ExpList3D::v_ReadGlobalOptimizationParameters(), and Nektar::MultiRegions::ExpList2D::v_ReadGlobalOptimizationParameters().
|
protected |
Mesh associated with this expansion list.
Definition at line 864 of file ExpList.h.
Referenced by Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField2D::FindPeriodicEdges(), Nektar::MultiRegions::DisContField3D::FindPeriodicFaces(), Nektar::MultiRegions::DisContField1D::FindPeriodicVertices(), GetExpIndex(), GetGraph(), Nektar::MultiRegions::DisContField3D::SetUpDG(), Nektar::MultiRegions::DisContField2D::SetUpDG(), Nektar::MultiRegions::DisContField1D::SetUpDG(), Nektar::MultiRegions::DisContField3D::v_GetBoundaryToElmtMap(), and Nektar::MultiRegions::DisContField2D::v_GetBoundaryToElmtMap().
|
protected |
The total number of local degrees of freedom. m_ncoeffs .
Definition at line 868 of file ExpList.h.
Referenced by Nektar::MultiRegions::ContField2D::BwdTrans(), Nektar::MultiRegions::ContField1D::ContField1D(), Nektar::MultiRegions::ContField2D::ContField2D(), Nektar::MultiRegions::ContField3D::ContField3D(), ExpList(), Nektar::MultiRegions::ExpList0D::ExpList0D(), Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), Nektar::MultiRegions::ExpList3D::ExpList3D(), Nektar::MultiRegions::ContField1D::FwdTrans(), GetNcoeffs(), Nektar::MultiRegions::ContField2D::IProductWRTBase(), Nektar::MultiRegions::ContField1D::IProductWRTBase(), Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), Nektar::MultiRegions::ContField2D::LaplaceSolve(), MultiplyByElmtInvMass(), Nektar::MultiRegions::ContField1D::MultiplyByInvMassMatrix(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList0D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList2D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList1D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ContField3D::v_BwdTrans(), Nektar::MultiRegions::ExpListHomogeneous1D::v_ExtractCoeffsToCoeffs(), v_ExtractCoeffsToCoeffs(), v_FwdTrans_IterPerExp(), Nektar::MultiRegions::ContField3D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField1D::v_GeneralMatrixOp(), Nektar::MultiRegions::ContField2D::v_GeneralMatrixOp(), Nektar::MultiRegions::DisContField3D::v_HelmSolve(), Nektar::MultiRegions::DisContField2D::v_HelmSolve(), Nektar::MultiRegions::DisContField1D::v_HelmSolve(), and Nektar::MultiRegions::ContField3D::v_IProductWRTBase().
|
protected |
The total number of quadrature points. m_npoints
Definition at line 873 of file ExpList.h.
Referenced by ExpList(), Nektar::MultiRegions::ExpList0D::ExpList0D(), Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), Nektar::MultiRegions::ExpList3D::ExpList3D(), GenGlobalMatrix(), GetNpoints(), GetTotPoints(), Nektar::MultiRegions::ExpListHomogeneous1D::Homogeneous1DTrans(), Nektar::MultiRegions::ExpListHomogeneous2D::Homogeneous2DTrans(), Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), Nektar::MultiRegions::ContField2D::LinearAdvectionEigs(), Linf(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList0D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList2D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList1D::SetCoeffPhysOffsets(), and SetPhys().
|
protected |
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coeffs and m_phys is associated, i.e. for an array of constant expansion size and single shape elements m_phys[n*m_npoints] is the data related to m_exp[m_offset_elmt_id[n]];.
Definition at line 942 of file ExpList.h.
Referenced by Nektar::MultiRegions::DisContField3D::EvaluateHDGPostProcessing(), Nektar::MultiRegions::DisContField2D::EvaluateHDGPostProcessing(), GenBlockMatrix(), GeneralMatrixOp_IterPerExp(), GenGlobalMatrix(), GenGlobalMatrixFull(), GetOffset_Elmt_Id(), Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList0D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList2D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList1D::SetCoeffPhysOffsets(), v_BwdTrans_IterPerExp(), v_ExtractCoeffsToCoeffs(), Nektar::MultiRegions::DisContField3D::v_HelmSolve(), Nektar::MultiRegions::DisContField2D::v_HelmSolve(), v_Integral(), and v_IProductWRTBase_IterPerExp().
The global expansion evaluated at the quadrature points.
The array of length m_npoints containing the evaluation of at the quadrature points .
Definition at line 907 of file ExpList.h.
Referenced by ExpList(), Nektar::MultiRegions::ExpList0D::ExpList0D(), Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), GetPhys(), Nektar::MultiRegions::DisContField2D::L2_DGDeriv(), PhysIntegral(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), SetPhys(), SetPhysArray(), UpdatePhys(), Nektar::MultiRegions::DisContField3D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField1D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField2D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField3D::v_GetFwdBwdTracePhys(), Nektar::MultiRegions::DisContField2D::v_GetFwdBwdTracePhys(), Nektar::MultiRegions::DisContField1D::v_GetFwdBwdTracePhys(), v_WriteTecplotField(), Nektar::MultiRegions::ExpListHomogeneous1D::v_WriteVtkPieceData(), and v_WriteVtkPieceData().
|
protected |
Offset of elemental data into the array m_phys.
Definition at line 934 of file ExpList.h.
Referenced by FwdTrans_BndConstrained(), GenBlockMatrix(), GeneralMatrixOp_IterPerExp(), GenGlobalMatrix(), GenGlobalMatrixFull(), GetPhys_Offset(), H1(), IProductWRTDerivBase(), PhysIntegral(), Nektar::MultiRegions::ExpList3D::SetCoeffPhys(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList0D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList2D::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList1D::SetCoeffPhysOffsets(), v_BwdTrans_IterPerExp(), v_GetCoords(), Nektar::MultiRegions::ExpList0D::v_GetNormals(), Nektar::MultiRegions::ExpList2D::v_GetNormals(), Nektar::MultiRegions::ExpList1D::v_GetNormals(), v_Integral(), v_IProductWRTBase_IterPerExp(), v_L2(), v_PhysDeriv(), Nektar::MultiRegions::ExpList2D::v_Upwind(), Nektar::MultiRegions::ExpList1D::v_Upwind(), v_WriteTecplotField(), Nektar::MultiRegions::ExpListHomogeneous1D::v_WriteVtkPieceData(), and v_WriteVtkPieceData().
|
protected |
The state of the array m_phys.
Indicates whether the array m_phys, created to contain the evaluation of at the quadrature points , is filled with these values.
Definition at line 916 of file ExpList.h.
Referenced by Nektar::MultiRegions::ExpList2D::ExpList2D(), GetPhysState(), PhysIntegral(), SetPhys(), SetPhysState(), UpdatePhys(), Nektar::MultiRegions::DisContField3D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField1D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField2D::v_ExtractTracePhys(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_ExtractTracePhys(), and v_WriteTecplotField().
|
protected |
Session.
Definition at line 861 of file ExpList.h.
Referenced by Nektar::MultiRegions::ContField1D::ContField1D(), Nektar::MultiRegions::ContField2D::ContField2D(), Nektar::MultiRegions::ContField3D::ContField3D(), Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField1D::DisContField1D(), Nektar::MultiRegions::DisContField2D::DisContField2D(), Nektar::MultiRegions::DisContField3D::DisContField3D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::DisContField3DHomogeneous1D(), Nektar::MultiRegions::ExpList2DHomogeneous1D::ExpList2DHomogeneous1D(), Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D(), ExtractFileBCs(), Nektar::MultiRegions::DisContField2D::FindPeriodicEdges(), Nektar::MultiRegions::DisContField3D::FindPeriodicFaces(), Nektar::MultiRegions::DisContField1D::FindPeriodicVertices(), Nektar::MultiRegions::DisContField3D::GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D(), GenGlobalMatrix(), Nektar::MultiRegions::DisContField1D::GetDomainBCs(), GetSession(), Nektar::MultiRegions::DisContField3DHomogeneous2D::SetupBoundaryConditions(), Nektar::MultiRegions::DisContField3DHomogeneous1D::SetupBoundaryConditions(), Nektar::MultiRegions::DisContField3D::SetUpDG(), Nektar::MultiRegions::DisContField2D::SetUpDG(), Nektar::MultiRegions::DisContField1D::SetUpDG(), Nektar::MultiRegions::ExpList3D::v_ReadGlobalOptimizationParameters(), and Nektar::MultiRegions::ExpList2D::v_ReadGlobalOptimizationParameters().
|
protected |
Definition at line 951 of file ExpList.h.
Referenced by GetWaveSpace(), SetWaveSpace(), Nektar::MultiRegions::ExpListHomogeneous2D::v_BwdTrans(), Nektar::MultiRegions::ExpListHomogeneous1D::v_BwdTrans(), Nektar::MultiRegions::ExpListHomogeneous2D::v_BwdTrans_IterPerExp(), Nektar::MultiRegions::ExpListHomogeneous1D::v_BwdTrans_IterPerExp(), Nektar::MultiRegions::ExpListHomogeneous2D::v_DealiasedProd(), Nektar::MultiRegions::ExpListHomogeneous2D::v_FwdTrans(), Nektar::MultiRegions::ExpListHomogeneous1D::v_FwdTrans(), Nektar::MultiRegions::ExpListHomogeneous2D::v_FwdTrans_IterPerExp(), Nektar::MultiRegions::ExpListHomogeneous1D::v_FwdTrans_IterPerExp(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous2D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::ExpListHomogeneous2D::v_PhysDeriv(), and Nektar::MultiRegions::ExpListHomogeneous1D::v_PhysDeriv().