Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Friends | List of all members
Nektar::LEE Class Reference

#include <LEE.h>

Inheritance diagram for Nektar::LEE:
[legend]

Public Member Functions

virtual ~LEE ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::AcousticSystem
virtual ~AcousticSystem ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble >> &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, const Array< OneD, const NekDouble > &input)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

 LEE (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
virtual void v_AddLinTerm (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray) override
 
virtual void v_InitObject (bool DeclareFields=true) override
 Initialization object for the LEE class. More...
 
virtual void v_GetFluxVector (const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux) override
 Return the flux vector for the LEE equations. More...
 
- Protected Member Functions inherited from Nektar::AcousticSystem
 AcousticSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 Compute the right-hand side. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 Compute the projection and call the method for imposing the boundary conditions in case of discontinuous projection. More...
 
virtual bool v_PreIntegrate (int step) override
 v_PreIntegrate More...
 
virtual void v_Output () override
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 Compute the advection velocity in the standard space for each element of the expansion. More...
 
virtual void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables) override
 
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals ()
 Get the normal vectors. More...
 
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs ()
 Get the locations of the components of the directed fields within the fields array. More...
 
const Array< OneD, const Array< OneD, NekDouble > > & GetBasefieldFwdBwd ()
 Get the baseflow field. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise () override
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 

Private Member Functions

virtual void v_RiemannInvariantBC (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &BfFwd, Array< OneD, Array< OneD, NekDouble >> &physarray) override
 Outflow characteristic boundary conditions for compressible flow problems. More...
 

Friends

class MemoryManager< LEE >
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Protected Attributes inherited from Nektar::AcousticSystem
int m_ip
 indices of the fields More...
 
int m_irho
 
int m_iu
 
bool m_conservative
 we are dealing with a conservative formualtion More...
 
SolverUtils::CouplingSharedPtr m_coupling
 
SolverUtils::AdvectionSharedPtr m_advection
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
Array< OneD, Array< OneD, NekDouble > > m_bfFwdBwd
 
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
 
Array< OneD, Array< OneD, NekDouble > > m_bf
 
std::vector< std::string > m_bfNames
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_timestepMax = -1.0
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_pararealIter
 Number of parareal time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_useInitialCondition
 Flag to determine if IC are used. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 46 of file LEE.h.

Constructor & Destructor Documentation

◆ ~LEE()

Nektar::LEE::~LEE ( )
virtual

Destructor.

Destructor for LEE class.

Definition at line 125 of file LEE.cpp.

126 {
127 }

◆ LEE()

Nektar::LEE::LEE ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
protected

Initialises UnsteadySystem class members.

Definition at line 47 of file LEE.cpp.

49  : UnsteadySystem(pSession, pGraph), AcousticSystem(pSession, pGraph)
50 {
51  m_ip = 0;
52  m_irho = 1;
53  m_iu = 2;
54 
55  m_conservative = true;
56 }
bool m_conservative
we are dealing with a conservative formualtion
AcousticSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
int m_ip
indices of the fields
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

References Nektar::AcousticSystem::m_conservative, Nektar::AcousticSystem::m_ip, Nektar::AcousticSystem::m_irho, and Nektar::AcousticSystem::m_iu.

Member Function Documentation

◆ create()

static EquationSystemSharedPtr Nektar::LEE::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file LEE.h.

55  {
57  MemoryManager<LEE>::AllocateSharedPtr(pSession, pGraph);
58  p->InitObject();
59  return p;
60  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ v_AddLinTerm()

void Nektar::LEE::v_AddLinTerm ( const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::AcousticSystem.

Definition at line 193 of file LEE.cpp.

195 {
196  int nq = GetTotPoints();
197 
198  Array<OneD, Array<OneD, NekDouble>> linTerm(m_spacedim + 2);
199  for (int i = 0; i < m_spacedim + 2; ++i)
200  {
201  if (i == 1)
202  {
203  // skip rho
204  continue;
205  }
206 
207  linTerm[i] = Array<OneD, NekDouble>(nq);
208  }
209 
210  Array<OneD, const NekDouble> c0sq = m_bf[0];
211  Array<OneD, const NekDouble> rho0 = m_bf[1];
212  Array<OneD, const NekDouble> gamma = m_bf[m_iu + m_spacedim];
213 
214  Array<OneD, NekDouble> gammaMinOne(nq);
215  Vmath::Sadd(nq, -1.0, gamma, 1, gammaMinOne, 1);
216 
217  Array<OneD, NekDouble> p0(nq);
218  Vmath::Vmul(nq, c0sq, 1, rho0, 1, p0, 1);
219  Vmath::Vdiv(nq, p0, 1, gamma, 1, p0, 1);
220 
221  Array<OneD, Array<OneD, const NekDouble>> u0(m_spacedim);
222  for (int i = 0; i < m_spacedim; ++i)
223  {
224  u0[i] = m_bf[2 + i];
225  }
226 
227  Array<OneD, const NekDouble> p = inarray[0];
228  Array<OneD, const NekDouble> rho = inarray[1];
229 
230  Array<OneD, Array<OneD, const NekDouble>> ru(m_spacedim);
231  for (int i = 0; i < m_spacedim; ++i)
232  {
233  ru[i] = inarray[2 + i];
234  }
235 
236  Array<OneD, NekDouble> grad(nq);
237  Array<OneD, NekDouble> tmp1(nq);
238 
239  // p
240  {
241  Vmath::Zero(nq, linTerm[m_ip], 1);
242  // (1-gamma) ( ru_j / rho0 * dp0/dx_j - p * du0_j/dx_j )
243  for (int j = 0; j < m_spacedim; ++j)
244  {
245  // ru_j / rho0 * dp0/dx_j
246  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], p0, grad);
247  Vmath::Vmul(nq, grad, 1, ru[j], 1, tmp1, 1);
248  Vmath::Vdiv(nq, tmp1, 1, rho0, 1, tmp1, 1);
249  // p * du0_j/dx_j - ru_j / rho0 * dp0/dx_j
250  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], u0[j],
251  grad);
252  Vmath::Vvtvm(nq, grad, 1, p, 1, tmp1, 1, tmp1, 1);
253  // (gamma-1) (p * du0_j/dx_j - ru_j / rho0 * dp0/dx_j)
254  Vmath::Vvtvp(nq, gammaMinOne, 1, tmp1, 1, linTerm[m_ip], 1,
255  linTerm[m_ip], 1);
256  }
257  }
258 
259  // rho has no linTerm
260 
261  // ru_i
262  for (int i = 0; i < m_spacedim; ++i)
263  {
264  Vmath::Zero(nq, linTerm[m_iu + i], 1);
265  // du0_i/dx_j * (u0_j * rho + ru_j)
266  for (int j = 0; j < m_spacedim; ++j)
267  {
268  // d u0_i / d x_j
269  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], u0[i],
270  grad);
271  // u0_j * rho + ru_j
272  Vmath::Vvtvp(nq, u0[j], 1, rho, 1, ru[j], 1, tmp1, 1);
273  // du0_i/dx_j * (u0_j * rho + ru_j)
274  Vmath::Vvtvp(nq, grad, 1, tmp1, 1, linTerm[m_iu + i], 1,
275  linTerm[m_iu + i], 1);
276  }
277  }
278 
279  Array<OneD, NekDouble> tmpC(GetNcoeffs());
280  for (int i = 0; i < m_spacedim + 2; ++i)
281  {
282  if (i == 1)
283  {
284  // skip rho
285  continue;
286  }
287 
288  m_fields[0]->FwdTrans(linTerm[i], tmpC);
289  m_fields[0]->BwdTrans(tmpC, linTerm[i]);
290 
291  Vmath::Vadd(nq, outarray[i], 1, linTerm[i], 1, outarray[i], 1);
292  }
293 }
Array< OneD, Array< OneD, NekDouble > > m_bf
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
SOLVER_UTILS_EXPORT int GetTotPoints()
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:91
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359
void Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.cpp:598
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:284
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384

References Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::AcousticSystem::m_bf, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::AcousticSystem::m_ip, Nektar::AcousticSystem::m_iu, Nektar::SolverUtils::EquationSystem::m_spacedim, CellMLToNektar.cellml_metadata::p, Vmath::Sadd(), Vmath::Vadd(), Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vvtvm(), Vmath::Vvtvp(), and Vmath::Zero().

◆ v_GetFluxVector()

void Nektar::LEE::v_GetFluxVector ( const Array< OneD, Array< OneD, NekDouble >> &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &  flux 
)
overrideprotectedvirtual

Return the flux vector for the LEE equations.

Parameters
physfieldFields.
fluxResulting flux. flux[eq][dir][pt]

Implements Nektar::AcousticSystem.

Definition at line 135 of file LEE.cpp.

138 {
139  int nq = physfield[0].size();
140 
141  ASSERTL1(flux[0].size() == m_spacedim,
142  "Dimension of flux array and velocity array do not match");
143 
144  Array<OneD, const NekDouble> c0sq = m_bf[0];
145  Array<OneD, Array<OneD, const NekDouble>> u0(m_spacedim);
146  for (int i = 0; i < m_spacedim; ++i)
147  {
148  u0[i] = m_bf[2 + i];
149  }
150 
151  Array<OneD, const NekDouble> p = physfield[m_ip];
152  Array<OneD, const NekDouble> rho = physfield[m_irho];
153  Array<OneD, Array<OneD, const NekDouble>> ru(m_spacedim);
154  for (int i = 0; i < m_spacedim; ++i)
155  {
156  ru[i] = physfield[m_iu + i];
157  }
158 
159  // F_{adv,p',j} = c0^2 * ru_j + u0_j * p
160  for (int j = 0; j < m_spacedim; ++j)
161  {
162  int i = 0;
163  Vmath::Vvtvvtp(nq, c0sq, 1, ru[j], 1, u0[j], 1, p, 1, flux[i][j], 1);
164  }
165 
166  // F_{adv,rho',j} = u0_j * rho' + ru_j
167  for (int j = 0; j < m_spacedim; ++j)
168  {
169  int i = 1;
170  // u0_j * rho' + ru_j
171  Vmath::Vvtvp(nq, u0[j], 1, rho, 1, ru[j], 1, flux[i][j], 1);
172  }
173 
174  for (int i = 0; i < m_spacedim; ++i)
175  {
176  // F_{adv,u'_i,j} = ru_i * u0_j + delta_ij * p
177  for (int j = 0; j < m_spacedim; ++j)
178  {
179  // ru_i * u0_j
180  Vmath::Vmul(nq, ru[i], 1, u0[j], 1, flux[m_iu + i][j], 1);
181 
182  // kronecker delta
183  if (i == j)
184  {
185  // delta_ij + p
186  Vmath::Vadd(nq, p, 1, flux[m_iu + i][j], 1, flux[m_iu + i][j],
187  1);
188  }
189  }
190  }
191 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:692

References ASSERTL1, Nektar::AcousticSystem::m_bf, Nektar::AcousticSystem::m_ip, Nektar::AcousticSystem::m_irho, Nektar::AcousticSystem::m_iu, Nektar::SolverUtils::EquationSystem::m_spacedim, CellMLToNektar.cellml_metadata::p, Vmath::Vadd(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

Referenced by v_InitObject().

◆ v_InitObject()

void Nektar::LEE::v_InitObject ( bool  DeclareFields = true)
overrideprotectedvirtual

Initialization object for the LEE class.

Reimplemented from Nektar::AcousticSystem.

Definition at line 61 of file LEE.cpp.

62 {
63  AcousticSystem::v_InitObject(DeclareFields);
64 
65  m_bfNames.push_back("gamma");
66 
67  // Initialize basefield again
68  m_bf = Array<OneD, Array<OneD, NekDouble>>(m_bfNames.size());
69  for (int i = 0; i < m_bf.size(); ++i)
70  {
71  m_bf[i] = Array<OneD, NekDouble>(GetTotPoints());
72  }
73  GetFunction("Baseflow", m_fields[0], true)
74  ->Evaluate(m_bfNames, m_bf, m_time);
75 
76  // Define the normal velocity fields
77  m_bfFwdBwd = Array<OneD, Array<OneD, NekDouble>>(2 * m_bfNames.size());
78  for (int i = 0; i < m_bfFwdBwd.size(); i++)
79  {
80  m_bfFwdBwd[i] = Array<OneD, NekDouble>(GetTraceNpoints(), 0.0);
81  }
82 
83  string riemName;
84  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
85  if (boost::to_lower_copy(riemName) == "characteristics" ||
86  boost::to_lower_copy(riemName) == "leeupwind" ||
87  boost::to_lower_copy(riemName) == "upwind")
88  {
89  riemName = "LEEUpwind";
90  }
91  if (boost::to_lower_copy(riemName) == "laxfriedrichs")
92  {
93  riemName = "LEELaxFriedrichs";
94  }
96  riemName, m_session);
97  m_riemannSolver->SetVector("N", &LEE::GetNormals, this);
98  m_riemannSolver->SetVector("basefieldFwdBwd", &LEE::GetBasefieldFwdBwd,
99  this);
100  m_riemannSolver->SetAuxVec("vecLocs", &LEE::GetVecLocs, this);
101 
102  // Set up advection operator
103  string advName;
104  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
105  m_advection =
107  m_advection->SetFluxVector(&LEE::v_GetFluxVector, this);
108  m_advection->SetRiemannSolver(m_riemannSolver);
109  m_advection->InitObject(m_session, m_fields);
110 
112  {
115  }
116  else
117  {
118  ASSERTL0(false, "Implicit LEE not set up.");
119  }
120 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
Get the locations of the components of the directed fields within the fields array.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
Get the normal vectors.
std::vector< std::string > m_bfNames
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the projection and call the method for imposing the boundary conditions in case of discontinu...
SolverUtils::AdvectionSharedPtr m_advection
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the right-hand side.
Array< OneD, Array< OneD, NekDouble > > m_bfFwdBwd
virtual void v_InitObject(bool DeclareFields=true) override
Initialization object for the AcousticSystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetBasefieldFwdBwd()
Get the baseflow field.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
virtual void v_GetFluxVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux) override
Return the flux vector for the LEE equations.
Definition: LEE.cpp:135
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
SOLVER_UTILS_EXPORT int GetTraceNpoints()
NekDouble m_time
Current time of simulation.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), Nektar::AcousticSystem::DoOdeProjection(), Nektar::AcousticSystem::DoOdeRhs(), Nektar::SolverUtils::GetAdvectionFactory(), Nektar::AcousticSystem::GetBasefieldFwdBwd(), Nektar::SolverUtils::EquationSystem::GetFunction(), Nektar::AcousticSystem::GetNormals(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::AcousticSystem::GetVecLocs(), Nektar::AcousticSystem::m_advection, Nektar::AcousticSystem::m_bf, Nektar::AcousticSystem::m_bfFwdBwd, Nektar::AcousticSystem::m_bfNames, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::AcousticSystem::m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_time, v_GetFluxVector(), and Nektar::AcousticSystem::v_InitObject().

◆ v_RiemannInvariantBC()

void Nektar::LEE::v_RiemannInvariantBC ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble >> &  Fwd,
Array< OneD, Array< OneD, NekDouble >> &  BfFwd,
Array< OneD, Array< OneD, NekDouble >> &  physarray 
)
overrideprivatevirtual

Outflow characteristic boundary conditions for compressible flow problems.

Implements Nektar::AcousticSystem.

Definition at line 299 of file LEE.cpp.

303 {
304  int id1, id2, nBCEdgePts;
305  int nVariables = physarray.size();
306 
307  const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
308 
309  int eMax = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
310 
311  for (int e = 0; e < eMax; ++e)
312  {
313  nBCEdgePts = m_fields[0]
314  ->GetBndCondExpansions()[bcRegion]
315  ->GetExp(e)
316  ->GetTotPoints();
317  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
318  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[cnt + e]);
319 
320  // Calculate (v.n)
321  Array<OneD, NekDouble> RVn(nBCEdgePts, 0.0);
322  for (int i = 0; i < m_spacedim; ++i)
323  {
324  Vmath::Vvtvp(nBCEdgePts, &Fwd[m_iu + i][id2], 1,
325  &m_traceNormals[i][id2], 1, &RVn[0], 1, &RVn[0], 1);
326  }
327 
328  // Calculate (v0.n)
329  Array<OneD, NekDouble> RVn0(nBCEdgePts, 0.0);
330  for (int i = 0; i < m_spacedim; ++i)
331  {
332  Vmath::Vvtvp(nBCEdgePts, &BfFwd[2 + i][id2], 1,
333  &m_traceNormals[i][id2], 1, &RVn0[0], 1, &RVn0[0], 1);
334  }
335 
336  for (int i = 0; i < nBCEdgePts; ++i)
337  {
338  NekDouble c = sqrt(BfFwd[0][id2 + i]);
339 
340  NekDouble h1, h4, h5;
341 
342  if (RVn0[i] > 0)
343  {
344  // rho - p / c^2
345  h1 = Fwd[m_irho][id2 + i] - Fwd[m_ip][id2 + i] / (c * c);
346  }
347  else
348  {
349  h1 = 0.0;
350  }
351 
352  if (RVn0[i] - c > 0)
353  {
354  // ru / 2 - p / (2*c)
355  h4 = RVn[i] / 2 - Fwd[m_ip][id2 + i] / (2 * c);
356  }
357  else
358  {
359  h4 = 0.0;
360  }
361 
362  if (RVn0[i] + c > 0)
363  {
364  // ru / 2 + p / (2*c)
365  h5 = RVn[i] / 2 + Fwd[m_ip][id2 + i] / (2 * c);
366  }
367  else
368  {
369  h5 = 0.0;
370  }
371 
372  // compute conservative variables
373  // p = c0*(h5-h4)
374  // rho = h1 + (h5-h4)/c0
375  // ru = h4+h5
376  Fwd[m_ip][id2 + i] = c * (h5 - h4);
377  Fwd[m_irho][id2 + i] = h1 + (h5 - h4) / c;
378  NekDouble RVnNew = h4 + h5;
379 
380  // adjust velocity pert. according to new value
381  // here we just omit the wall parallel velocity components, i.e
382  // setting them to zero. This is equivalent to setting the two
383  // vorticity characteristics h2 and h3 to zero. Mathematically,
384  // this is only legitimate for incoming characteristics. However,
385  // as h2 and h3 are convected by the flow, the value we precribe at
386  // an the boundary for putgoing characteristics does not matter.
387  // This implementation saves a few operations and is more robust
388  // for mixed in/outflow boundaries and at the boundaries edges.
389  for (int j = 0; j < m_spacedim; ++j)
390  {
391  Fwd[m_iu + j][id2 + i] = RVnNew * m_traceNormals[j][id2 + i];
392  }
393  }
394 
395  // Copy boundary adjusted values into the boundary expansion
396  for (int i = 0; i < nVariables; ++i)
397  {
398  Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
399  &(m_fields[i]
400  ->GetBndCondExpansions()[bcRegion]
401  ->UpdatePhys())[id1],
402  1);
403  }
404  }
405 }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::AcousticSystem::m_ip, Nektar::AcousticSystem::m_irho, Nektar::AcousticSystem::m_iu, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_traceNormals, tinysimd::sqrt(), Vmath::Vcopy(), and Vmath::Vvtvp().

Friends And Related Function Documentation

◆ MemoryManager< LEE >

friend class MemoryManager< LEE >
friend

Definition at line 1 of file LEE.h.

Member Data Documentation

◆ className

string Nektar::LEE::className
static
Initial value:
"LEE", LEE::create, "Linearized Euler Equations")
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: LEE.h:52
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 62 of file LEE.h.