Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Private Member Functions | List of all members
Nektar::MultiRegions::ExpList3D Class Reference

Abstraction of a three-dimensional multi-elemental expansion which is merely a collection of local expansions. More...

#include <ExpList3D.h>

Inheritance diagram for Nektar::MultiRegions::ExpList3D:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::ExpList3D:
Collaboration graph
[legend]

Public Member Functions

 ExpList3D ()
 Default constructor.
 ExpList3D (const ExpList3D &In)
 Copy constructor.
 ExpList3D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &TBa, const LibUtilities::BasisKey &TBb, const LibUtilities::BasisKey &TBc, const LibUtilities::BasisKey &HBa, const LibUtilities::BasisKey &HBb, const LibUtilities::BasisKey &HBc, const SpatialDomains::MeshGraphSharedPtr &graph3D, const LibUtilities::PointsType TetNb=LibUtilities::SIZE_PointsType)
 ExpList3D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable="DefaultVar")
 Sets up a list of local expansions based on an input mesh.
 ExpList3D (const SpatialDomains::ExpansionMap &expansions)
 Sets up a list of local expansions based on an expansion vector.
virtual ~ExpList3D ()
 Destructor.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 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 $N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m$.
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 $=Q_{\mathrm{tot}}$.
int GetTotPoints (const int eid) const
 Returns the total number of quadrature points for eid's element $=Q_{\mathrm{tot}}$.
int GetNpoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
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 $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
bool GetPhysState (void) const
 This function indicates whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
NekDouble PhysIntegral (void)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
NekDouble PhysIntegral (const Array< OneD, const NekDouble > &inarray)
 This function integrates a function $f(\boldsymbol{x})$ 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 $f(\boldsymbol{x})$ with respect to all {local} expansion modes $\phi_n^e(\boldsymbol{x})$.
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 $f(\boldsymbol{x})$ 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 $u(\boldsymbol{x})$ 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 $\boldsymbol{x}_i$.
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 $\boldsymbol{\hat{u}}_l$ (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 $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ evaluated at the quadrature points.
NekDouble Linf (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the $L_\infty$ 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 $L_2$ 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 $H^1$ error of the global spectral/hp element approximation.
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
Array< OneD, const NekDoubleHomogeneousEnergy (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::ExpansionSharedPtrGetExp (int n) const
 This function returns (a shared pointer to) the local elemental expansion of the $n^{\mathrm{th}}$ element.
LocalRegions::ExpansionSharedPtrGetExp (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 $\boldsymbol{\hat{u}}_l$ (implemented as m_coeffs) containing all local expansion coefficients.
Array< OneD, NekDouble > & UpdatePhys ()
 This function returns (a reference to) the array $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ 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 $f(\boldsymbol{x})$ 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, RobinBCInfoSharedPtrGetRobinBCInfo ()
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< ExpListGetSharedThisPtr ()
 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)

Protected Member Functions

virtual void v_SetUpPhysNormals ()
 Set up the normals on each expansion.
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList
boost::shared_ptr< DNekMatGenGlobalMatrixFull (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 DNekScalBlkMatSharedPtrGetBlockMatrix (const GlobalMatrixKey &gkey)
void MultiplyByBlockMatrix (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
boost::shared_ptr< GlobalMatrixGenGlobalMatrix (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< GlobalLinSysGenGlobalLinSys (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 This operation constructs the global linear system of type mkey.
boost::shared_ptr< GlobalLinSysGenGlobalBndLinSys (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_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)

Private Member Functions

void SetCoeffPhys (void)
 Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys.
virtual void v_ReadGlobalOptimizationParameters ()
virtual void v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
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)

Additional Inherited Members

- Public Attributes inherited from Nektar::MultiRegions::ExpList
ExpansionType m_expType
- Static Protected Member Functions inherited from Nektar::MultiRegions::ExpList
static
SpatialDomains::BoundaryConditionShPtr 
GetBoundaryCondition (const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
- Protected Attributes inherited from Nektar::MultiRegions::ExpList
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 $=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l$.
int m_npoints
Array< OneD, NekDoublem_coeffs
 Concatenation of all local expansion coefficients.
Array< OneD, NekDoublem_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

Detailed Description

Abstraction of a three-dimensional multi-elemental expansion which is merely a collection of local expansions.

Definition at line 49 of file ExpList3D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::ExpList3D::ExpList3D ( )

Default constructor.

Definition at line 52 of file ExpList3D.cpp.

References Nektar::MultiRegions::e3D, and Nektar::MultiRegions::ExpList::SetExpType().

: ExpList()
{
}
Nektar::MultiRegions::ExpList3D::ExpList3D ( const ExpList3D In)

Copy constructor.

Definition at line 57 of file ExpList3D.cpp.

References Nektar::MultiRegions::e3D, and Nektar::MultiRegions::ExpList::SetExpType().

: ExpList(In)
{
}
Nektar::MultiRegions::ExpList3D::ExpList3D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey TBa,
const LibUtilities::BasisKey TBb,
const LibUtilities::BasisKey TBc,
const LibUtilities::BasisKey HBa,
const LibUtilities::BasisKey HBb,
const LibUtilities::BasisKey HBc,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const LibUtilities::PointsType  TetNb = LibUtilities::SIZE_PointsType 
)

Definition at line 67 of file ExpList3D.cpp.

References ASSERTL0, Nektar::MultiRegions::e3D, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LibUtilities::StdSegData::getNumberOfCoefficients(), Nektar::LibUtilities::BasisKey::GetNumModes(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::ReadGlobalOptimizationParameters(), SetCoeffPhys(), Nektar::MultiRegions::ExpList::SetExpType(), and Nektar::LibUtilities::SIZE_PointsType.

:
ExpList(pSession,graph3D)
{
const SpatialDomains::ExpansionMap &expansions = graph3D->GetExpansions();
SpatialDomains::ExpansionMap::const_iterator expIt;
for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
{
if((TetGeom = boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
{
{
// Ntet = MemoryManager<LocalRegions::NodalTetExp>::AllocateSharedPtr(TetBa,TetBb,TetBc,TetNb,TetGeom);
// (*m_exp).push_back(Ntet);
}
else
{
(*m_exp).push_back(tet);
}
}
/*
else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>(expansions[i]->m_geomShPtr)))
{
prism = MemoryManager<LocalRegions::PrismExp>::AllocateSharedPtr(Ba,Bb,Bc,PrismGeom);
(*m_exp).push_back(prism);
m_ncoeffs += StdRegions::StdPrismData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
}
else if((PyrGeom = boost::dynamic_pointer_cast<SpatialDomains::PyrGeom>(expansions[i]->m_geomShPtr)))
{
pyramid = MemoryManager<LocalRegions::PyrExp>::AllocateSharedPtr(Ba,Bb,Bc,PyrGeom);
(*m_exp).push_back(pyramid);
m_ncoeffs += StdRegions::StdPyrData::getNumberOfCoefficients(Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes());
m_npoints += Ba.GetNumPoints()*Bb.GetNumPoints()*Bc.GetNumPoints();
}
*/
else if((HexGeom = boost::dynamic_pointer_cast<SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
{
(*m_exp).push_back(hex);
m_npoints += HBa.GetNumPoints()*HBb.GetNumPoints()*HBc.GetNumPoints();
}
else
{
ASSERTL0(false,"dynamic cast to a proper Geometry3D failed");
}
}
// Setup Default optimisation information.
int nel = GetExpSize();
}
Nektar::MultiRegions::ExpList3D::ExpList3D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const std::string &  variable = "DefaultVar" 
)

Sets up a list of local expansions based on an input mesh.

Given a mesh graph3D, containing information about the domain and the spectral/hp element expansion, this constructor fills the list of local expansions {m_exp} with the proper expansions, calculates the total number of quadrature points $\boldsymbol{x}_i$ and the local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys.

Parameters
graph3DA mesh, containing information about the domain and the spectral/hp element expansion.

Definition at line 169 of file ExpList3D.cpp.

References ASSERTL0, Nektar::MultiRegions::e3D, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::ReadGlobalOptimizationParameters(), SetCoeffPhys(), and Nektar::MultiRegions::ExpList::SetExpType().

:
ExpList(pSession,graph3D)
{
const SpatialDomains::ExpansionMap &expansions
= graph3D->GetExpansions(variable);
SpatialDomains::ExpansionMap::const_iterator expIt;
for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
{
if((TetGeom = boost::dynamic_pointer_cast<
SpatialDomains::TetGeom>(expIt->second->m_geomShPtr)))
{
= expIt->second->m_basisKeyVector[0];
= expIt->second->m_basisKeyVector[1];
= expIt->second->m_basisKeyVector[2];
{
ASSERTL0(false,"LocalRegions::NodalTetExp is not "
"implemented yet");
}
else
{
::AllocateSharedPtr(TetBa,TetBb,TetBc,
TetGeom);
(*m_exp).push_back(tet);
}
}
else if((PrismGeom = boost::dynamic_pointer_cast<SpatialDomains
::PrismGeom>(expIt->second->m_geomShPtr)))
{
= expIt->second->m_basisKeyVector[0];
= expIt->second->m_basisKeyVector[1];
= expIt->second->m_basisKeyVector[2];
::AllocateSharedPtr(PrismBa,PrismBb,
PrismBc,PrismGeom);
(*m_exp).push_back(prism);
}
else if((PyrGeom = boost::dynamic_pointer_cast<
SpatialDomains::PyrGeom>(expIt->second->m_geomShPtr)))
{
= expIt->second->m_basisKeyVector[0];
= expIt->second->m_basisKeyVector[1];
= expIt->second->m_basisKeyVector[2];
::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
PyrGeom);
(*m_exp).push_back(pyramid);
}
else if((HexGeom = boost::dynamic_pointer_cast<
SpatialDomains::HexGeom>(expIt->second->m_geomShPtr)))
{
= expIt->second->m_basisKeyVector[0];
= expIt->second->m_basisKeyVector[1];
= expIt->second->m_basisKeyVector[2];
::AllocateSharedPtr(HexBa,HexBb,HexBc,
HexGeom);
(*m_exp).push_back(hex);
}
else
{
ASSERTL0(false,"dynamic cast to a proper Geometry3D "
"failed");
}
}
// Setup Default optimisation information.
int nel = GetExpSize();
}
Nektar::MultiRegions::ExpList3D::ExpList3D ( const SpatialDomains::ExpansionMap expansions)

Sets up a list of local expansions based on an expansion vector.

Given an expansion vector expansions, containing information about the domain and the spectral/hp element expansion, this constructor fills the list of local expansions {m_exp} with the proper expansions, calculates the total number of quadrature points $\boldsymbol{x}_i$ and the local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys.

Parameters
expansionsAn expansion vector, containing information about the domain and the spectral/hp element expansion.

Definition at line 292 of file ExpList3D.cpp.

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::e3D, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_globalOptParam, SetCoeffPhys(), and Nektar::MultiRegions::ExpList::SetExpType().

:
{
for(int i = 0; i < expansions.size(); ++i)
{
SpatialDomains::ExpansionMap::const_iterator expmap = expansions.find(i);
ASSERTL1(expmap != expansions.end(), "Unable to find expansion.");
const SpatialDomains::ExpansionShPtr exp = expmap->second;
if((TetGeom = boost::dynamic_pointer_cast<
SpatialDomains::TetGeom>(exp->m_geomShPtr)))
{
= exp->m_basisKeyVector[0];
= exp->m_basisKeyVector[1];
= exp->m_basisKeyVector[2];
{
ASSERTL0(false,"LocalRegions::NodalTetExp is not "
"implemented yet");
}
else
{
::AllocateSharedPtr(TetBa,TetBb,TetBc,
TetGeom);
(*m_exp).push_back(tet);
}
}
else if((PrismGeom = boost::dynamic_pointer_cast<
SpatialDomains::PrismGeom>(exp->m_geomShPtr)))
{
= exp->m_basisKeyVector[0];
= exp->m_basisKeyVector[1];
= exp->m_basisKeyVector[2];
::AllocateSharedPtr(PrismBa,PrismBb,
PrismBc,PrismGeom);
(*m_exp).push_back(prism);
}
else if((PyrGeom = boost::dynamic_pointer_cast<
SpatialDomains::PyrGeom>(exp->m_geomShPtr)))
{
= exp->m_basisKeyVector[0];
= exp->m_basisKeyVector[1];
= exp->m_basisKeyVector[2];
::AllocateSharedPtr(PyrBa,PyrBb,PyrBc,
PyrGeom);
(*m_exp).push_back(pyramid);
}
else if((HexGeom = boost::dynamic_pointer_cast<
SpatialDomains::HexGeom>(exp->m_geomShPtr)))
{
= exp->m_basisKeyVector[0];
= exp->m_basisKeyVector[1];
= exp->m_basisKeyVector[2];
::AllocateSharedPtr(HexBa,HexBb,HexBc,
HexGeom);
(*m_exp).push_back(hex);
}
else
{
ASSERTL0(false,"dynamic cast to a proper Geometry3D "
"failed");
}
}
// Setup Default optimisation information.
int nel = GetExpSize();
}
Nektar::MultiRegions::ExpList3D::~ExpList3D ( )
virtual

Destructor.

Definition at line 62 of file ExpList3D.cpp.

{
}

Member Function Documentation

void Nektar::MultiRegions::ExpList3D::SetCoeffPhys ( void  )
private

Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys.

Set up the storage for the concatenated list of coefficients and physical evaluations at the quadrature points. Each expansion (local element) is processed in turn to determine the number of coefficients and physical data points it contributes to the domain. Three arrays, m_coeff_offset, m_phys_offset and m_offset_elmt_id, are also initialised and updated to store the data offsets of each element in the m_coeffs and m_phys arrays, and the element id that each consecutive block is associated respectively.

Definition at line 411 of file ExpList3D.cpp.

References Nektar::MultiRegions::ExpList::m_coeff_offset, Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_exp, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_offset_elmt_id, Nektar::MultiRegions::ExpList::m_phys, and Nektar::MultiRegions::ExpList::m_phys_offset.

Referenced by ExpList3D().

{
int i;
// Set up offset information and array sizes
m_coeff_offset = Array<OneD,int>(m_exp->size());
m_phys_offset = Array<OneD,int>(m_exp->size());
m_offset_elmt_id = Array<OneD,int>(m_exp->size());
for(i = 0; i < m_exp->size(); ++i)
{
m_phys_offset [i] = m_npoints;
m_offset_elmt_id[i] = i;
m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
m_npoints += (*m_exp)[i]->GetTotPoints();
}
m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
m_phys = Array<OneD, NekDouble>(m_npoints);
}
void Nektar::MultiRegions::ExpList3D::v_PhysGalerkinProjection1DScaled ( const NekDouble  scale,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
privatevirtual

Definition at line 581 of file ExpList3D.cpp.

References Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_exp, and Nektar::LibUtilities::PhysGalerkinProject3D().

{
int cnt,cnt1;
cnt = cnt1 = 0;
for(int i = 0; i < GetExpSize(); ++i)
{
// get new points key
int pt0 = (*m_exp)[i]->GetNumPoints(0);
int pt1 = (*m_exp)[i]->GetNumPoints(1);
int pt2 = (*m_exp)[i]->GetNumPoints(2);
int npt0 = (int) pt0*scale;
int npt1 = (int) pt1*scale;
int npt2 = (int) pt2*scale;
LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
// Project points;
newPointsKey1,
newPointsKey2,
&inarray[cnt],
(*m_exp)[i]->GetBasis(0)->GetPointsKey(),
(*m_exp)[i]->GetBasis(1)->GetPointsKey(),
(*m_exp)[i]->GetBasis(2)->GetPointsKey(),
&outarray[cnt1]);
cnt += npt0*npt1*npt2;
cnt1 += pt0*pt1*pt2;
}
}
void Nektar::MultiRegions::ExpList3D::v_PhysInterp1DScaled ( const NekDouble  scale,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
privatevirtual

Definition at line 547 of file ExpList3D.cpp.

References Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LibUtilities::Interp3D(), and Nektar::MultiRegions::ExpList::m_exp.

{
int cnt,cnt1;
cnt = cnt1 = 0;
for(int i = 0; i < GetExpSize(); ++i)
{
// get new points key
int pt0 = (*m_exp)[i]->GetNumPoints(0);
int pt1 = (*m_exp)[i]->GetNumPoints(1);
int pt2 = (*m_exp)[i]->GetNumPoints(2);
int npt0 = (int) pt0*scale;
int npt1 = (int) pt1*scale;
int npt2 = (int) pt2*scale;
LibUtilities::PointsKey newPointsKey0(npt0,(*m_exp)[i]->GetPointsType(0));
LibUtilities::PointsKey newPointsKey1(npt1,(*m_exp)[i]->GetPointsType(1));
LibUtilities::PointsKey newPointsKey2(npt2,(*m_exp)[i]->GetPointsType(2));
// Interpolate points;
LibUtilities::Interp3D((*m_exp)[i]->GetBasis(0)->GetPointsKey(),
(*m_exp)[i]->GetBasis(1)->GetPointsKey(),
(*m_exp)[i]->GetBasis(2)->GetPointsKey(),
&inarray[cnt], newPointsKey0,
newPointsKey1, newPointsKey2,
&outarray[cnt1]);
cnt += pt0*pt1*pt2;
cnt1 += npt0*npt1*npt2;
}
}
void Nektar::MultiRegions::ExpList3D::v_ReadGlobalOptimizationParameters ( )
privatevirtual

Definition at line 435 of file ExpList3D.cpp.

References ASSERTL0, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eTetrahedron, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_exp, Nektar::MultiRegions::ExpList::m_globalOptParam, and Nektar::MultiRegions::ExpList::m_session.

{
Array<OneD, int> NumShape(4,0);
for(int i = 0; i < GetExpSize(); ++i)
{
switch ((*m_exp)[i]->DetShapeType())
{
case LibUtilities::eTetrahedron: NumShape[0]++; break;
case LibUtilities::ePyramid: NumShape[1]++; break;
case LibUtilities::ePrism: NumShape[2]++; break;
case LibUtilities::eHexahedron: NumShape[3]++; break;
default:
ASSERTL0(false, "Unknown expansion type.");
break;
}
}
int three = 3;
}
void Nektar::MultiRegions::ExpList3D::v_SetUpPhysNormals ( )
protectedvirtual

Set up the normals on each expansion.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 535 of file ExpList3D.cpp.

References Nektar::MultiRegions::ExpList::m_exp.

{
int i, j;
for (i = 0; i < m_exp->size(); ++i)
{
for (j = 0; j < (*m_exp)[i]->GetNfaces(); ++j)
{
(*m_exp)[i]->ComputeFaceNormal(j);
}
}
}
void Nektar::MultiRegions::ExpList3D::v_WriteVtkPieceHeader ( std::ofstream &  outfile,
int  expansion 
)
privatevirtual

Definition at line 458 of file ExpList3D.cpp.

{
int i,j,k;
int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
int nquad2 = (*m_exp)[expansion]->GetNumPoints(2);
int ntot = nquad0*nquad1*nquad2;
int ntotminus = (nquad0-1)*(nquad1-1)*(nquad2-1);
Array<OneD,NekDouble> coords[3];
coords[0] = Array<OneD,NekDouble>(ntot);
coords[1] = Array<OneD,NekDouble>(ntot);
coords[2] = Array<OneD,NekDouble>(ntot);
(*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
outfile << " <Piece NumberOfPoints=\""
<< ntot << "\" NumberOfCells=\""
<< ntotminus << "\">" << endl;
outfile << " <Points>" << endl;
outfile << " <DataArray type=\"Float64\" "
<< "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
outfile << " ";
for (i = 0; i < ntot; ++i)
{
for (j = 0; j < 3; ++j)
{
outfile << setprecision(8) << scientific
<< (float)coords[j][i] << " ";
}
outfile << endl;
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " </Points>" << endl;
outfile << " <Cells>" << endl;
outfile << " <DataArray type=\"Int32\" "
<< "Name=\"connectivity\" format=\"ascii\">" << endl;
for (i = 0; i < nquad0-1; ++i)
{
for (j = 0; j < nquad1-1; ++j)
{
for (k = 0; k < nquad2-1; ++k)
{
outfile << k*nquad0*nquad1 + j*nquad0 + i << " "
<< k*nquad0*nquad1 + j*nquad0 + i + 1 << " "
<< k*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
<< k*nquad0*nquad1 + (j+1)*nquad0 + i << " "
<< (k+1)*nquad0*nquad1 + j*nquad0 + i << " "
<< (k+1)*nquad0*nquad1 + j*nquad0 + i + 1 << " "
<< (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i + 1 << " "
<< (k+1)*nquad0*nquad1 + (j+1)*nquad0 + i << " " << endl;
}
}
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " <DataArray type=\"Int32\" "
<< "Name=\"offsets\" format=\"ascii\">" << endl;
for (i = 0; i < ntotminus; ++i)
{
outfile << i*8+8 << " ";
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " <DataArray type=\"UInt8\" "
<< "Name=\"types\" format=\"ascii\">" << endl;
for (i = 0; i < ntotminus; ++i)
{
outfile << "12 ";
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " </Cells>" << endl;
outfile << " <PointData>" << endl;
}