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

#include <PulseWavePropagation.h>

Inheritance diagram for Nektar::PulseWavePropagation:
[legend]

Public Member Functions

virtual ~PulseWavePropagation ()
 
Array< OneD, NekDouble > & GetA0 ()
 
Array< OneD, NekDouble > & GetBeta ()
 
Array< OneD, NekDouble > & GetAlpha ()
 
Array< OneD, NekDouble > & GetN ()
 
NekDouble GetRho ()
 
NekDouble GetDomains ()
 
- Public Member Functions inherited from Nektar::PulseWaveSystem
virtual ~PulseWaveSystem ()
 Destructor. More...
 
int GetNdomains ()
 
Array< OneD, MultiRegions::ExpListSharedPtrUpdateVessels (void)
 
- 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

 PulseWavePropagation (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 
void SetPulseWaveBoundaryConditions (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 
virtual void v_InitObject (bool DeclareField=false) override
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
 DG Pulse Wave Propagation routines: More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s) override
 
- Protected Member Functions inherited from Nektar::PulseWaveSystem
 PulseWaveSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises PulseWaveSystem class members. More...
 
virtual void v_DoInitialise () override
 Sets up initial conditions. More...
 
virtual void v_DoSolve () override
 Solves an unsteady problem. More...
 
void LinkSubdomains (Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &fields)
 Links the subdomains. More...
 
void BifurcationRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0, Array< OneD, NekDouble > &alpha)
 Riemann Problem for Bifurcation. More...
 
void MergingRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0, Array< OneD, NekDouble > &alpha)
 Riemann Problem for Merging Flow. More...
 
void InterfaceRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0, Array< OneD, NekDouble > &alpha)
 Riemann Problem for Interface/Junction. More...
 
virtual void v_Output (void) override
 
void CheckPoint_Output (const int n)
 
virtual NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false) override
 Compute the L2 error between fields and a given exact solution. More...
 
virtual NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray) override
 Compute the L_inf error between fields and a given exact solution. More...
 
void WriteVessels (const std::string &outname)
 Write input fields to the given filename. More...
 
void EnforceInterfaceConditions (const Array< OneD, const Array< OneD, NekDouble >> &fields)
 
- 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 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_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
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 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)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
SolverUtils::AdvectionSharedPtr m_advObject
 
Array< OneD, PulseWaveBoundarySharedPtrm_Boundary
 
- Protected Attributes inherited from Nektar::PulseWaveSystem
Array< OneD, MultiRegions::ExpListSharedPtrm_vessels
 
size_t m_nDomains
 
size_t m_currentDomain
 
size_t m_nVariables
 
UpwindTypePulse m_upwindTypePulse
 
Array< OneD, int > m_fieldPhysOffset
 
NekDouble m_rho
 
NekDouble m_pext
 
NekDouble m_C
 
NekDouble m_RT
 
NekDouble m_pout
 
Array< OneD, Array< OneD, NekDouble > > m_A_0
 
Array< OneD, Array< OneD, NekDouble > > m_A_0_trace
 
Array< OneD, Array< OneD, NekDouble > > m_beta
 
Array< OneD, Array< OneD, NekDouble > > m_beta_trace
 
Array< OneD, Array< OneD, NekDouble > > m_gamma
 
Array< OneD, Array< OneD, NekDouble > > m_alpha
 
Array< OneD, Array< OneD, NekDouble > > m_alpha_trace
 
Array< OneD, Array< OneD, NekDouble > > m_trace_fwd_normal
 
std::map< int, SpatialDomains::CompositeMapm_domain
 
std::vector< int > m_domOrder
 
Array< OneD, Array< OneD, NekDouble > > m_pressure
 
PulseWavePressureAreaSharedPtr m_pressureArea
 
bool extraFields = false
 
Array< OneD, Array< OneD, NekDouble > > m_PWV
 
Array< OneD, Array< OneD, NekDouble > > m_W1
 
Array< OneD, Array< OneD, NekDouble > > m_W2
 
std::vector< std::vector< InterfacePointShPtr > > m_vesselIntfcs
 
std::vector< std::vector< InterfacePointShPtr > > m_bifurcations
 
std::vector< std::vector< InterfacePointShPtr > > m_mergingJcts
 
- 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...
 

Friends

class MemoryManager< PulseWavePropagation >
 

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...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Set up the routines based on the weak formulation from "Computational Modelling of 1D blood flow with variable mechanical properties" by S. J. Sherwin et al. The weak formulation (1) reads: \( \sum_{e=1}^{N_{el}} \left[ \left( \frac{\partial \mathbf{U}^{\delta} }{\partial t} , \mathbf{\psi}^{\delta} \right)_{\Omega_e} - \left( \frac{\partial \mathbf{F(\mathbf{U})}^{\delta} } {\partial x}, \mathbf{\psi}^{\delta} \right)_{\Omega_e} + \left[ \mathbf{\psi}^{\delta} \cdot \{ \mathbf{F}^u - \mathbf{F}(\mathbf{U}^{\delta}) \} \right]_{x_e^l}^{x_eû} \right] = 0 \)

Definition at line 50 of file PulseWavePropagation.h.

Constructor & Destructor Documentation

◆ ~PulseWavePropagation()

Nektar::PulseWavePropagation::~PulseWavePropagation ( )
virtual

Definition at line 138 of file PulseWavePropagation.cpp.

139 {
140 }

◆ PulseWavePropagation()

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

Definition at line 65 of file PulseWavePropagation.cpp.

68  : PulseWaveSystem(pSession, pGraph)
69 {
70 }
PulseWaveSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises PulseWaveSystem class members.

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 56 of file PulseWavePropagation.h.

59  {
62  pGraph);
63  p->InitObject();
64  return p;
65  }
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.

◆ DoOdeProjection()

void Nektar::PulseWavePropagation::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const NekDouble  time 
)
protected

Definition at line 203 of file PulseWavePropagation.cpp.

206 {
207  boost::ignore_unused(time);
208 
209  // Just copy over array
210  for (size_t i = 0; i < m_nVariables; ++i)
211  {
212  Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
213  }
214 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Nektar::PulseWaveSystem::m_nVariables, and Vmath::Vcopy().

Referenced by v_InitObject().

◆ DoOdeRhs()

void Nektar::PulseWavePropagation::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const NekDouble  time 
)
protected

Computes the right hand side of (1). The RHS is everything except the term that contains the time derivative \(\frac{\partial \mathbf{U}}{\partial t}\). In case of a Discontinuous Galerkin projection, m_advObject->Advect will be called

Definition at line 150 of file PulseWavePropagation.cpp.

153 {
154  size_t i;
155 
156  Array<OneD, Array<OneD, NekDouble>> physarray(m_nVariables);
157 
158  // Dummy array for WeakDG advection
159  Array<OneD, Array<OneD, NekDouble>> advVel(m_spacedim);
160 
161  // Output array for advection
162  Array<OneD, Array<OneD, NekDouble>> out(m_nVariables);
163 
164  size_t cnt = 0;
165 
166  // Set up Inflow and Outflow boundary conditions.
167  SetPulseWaveBoundaryConditions(inarray, outarray, time);
168 
169  // Set up any interface conditions and write into boundary condition
171 
172  // do advection evaluation in all domains
173  for (size_t omega = 0; omega < m_nDomains; ++omega)
174  {
175  LibUtilities::Timer timer;
176  m_currentDomain = omega;
177  size_t nq = m_vessels[omega * m_nVariables]->GetTotPoints();
178 
179  timer.Start();
180  for (i = 0; i < m_nVariables; ++i)
181  {
182  physarray[i] = inarray[i] + cnt;
183  out[i] = outarray[i] + cnt;
184  }
185 
186  for (i = 0; i < m_nVariables; ++i)
187  {
188  m_fields[i] = m_vessels[omega * m_nVariables + i];
189  }
190 
191  m_advObject->Advect(m_nVariables, m_fields, advVel, physarray, out,
192  time);
193  for (i = 0; i < m_nVariables; ++i)
194  {
195  Vmath::Neg(nq, out[i], 1);
196  }
197  timer.Stop();
198  timer.AccumulateRegion("PulseWavePropagation:_DoOdeRHS", 1);
199  cnt += nq;
200  }
201 }
SolverUtils::AdvectionSharedPtr m_advObject
void SetPulseWaveBoundaryConditions(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
void EnforceInterfaceConditions(const Array< OneD, const Array< OneD, NekDouble >> &fields)
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::PulseWaveSystem::EnforceInterfaceConditions(), m_advObject, Nektar::PulseWaveSystem::m_currentDomain, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::PulseWaveSystem::m_nDomains, Nektar::PulseWaveSystem::m_nVariables, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::PulseWaveSystem::m_vessels, Vmath::Neg(), SetPulseWaveBoundaryConditions(), Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by v_InitObject().

◆ GetA0()

Array< OneD, NekDouble > & Nektar::PulseWavePropagation::GetA0 ( )

Definition at line 387 of file PulseWavePropagation.cpp.

388 {
390 }
Array< OneD, Array< OneD, NekDouble > > m_A_0_trace

References Nektar::PulseWaveSystem::m_A_0_trace, and Nektar::PulseWaveSystem::m_currentDomain.

Referenced by v_InitObject().

◆ GetAlpha()

Array< OneD, NekDouble > & Nektar::PulseWavePropagation::GetAlpha ( )

Definition at line 397 of file PulseWavePropagation.cpp.

398 {
400 }
Array< OneD, Array< OneD, NekDouble > > m_alpha_trace

References Nektar::PulseWaveSystem::m_alpha_trace, and Nektar::PulseWaveSystem::m_currentDomain.

Referenced by v_InitObject().

◆ GetBeta()

Array< OneD, NekDouble > & Nektar::PulseWavePropagation::GetBeta ( )

Definition at line 392 of file PulseWavePropagation.cpp.

393 {
395 }
Array< OneD, Array< OneD, NekDouble > > m_beta_trace

References Nektar::PulseWaveSystem::m_beta_trace, and Nektar::PulseWaveSystem::m_currentDomain.

Referenced by v_InitObject().

◆ GetDomains()

NekDouble Nektar::PulseWavePropagation::GetDomains ( )

Definition at line 412 of file PulseWavePropagation.cpp.

413 {
414  return m_nDomains;
415 }

References Nektar::PulseWaveSystem::m_nDomains.

Referenced by v_InitObject().

◆ GetFluxVector()

void Nektar::PulseWavePropagation::GetFluxVector ( const Array< OneD, Array< OneD, NekDouble >> &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &  flux 
)
protected

DG Pulse Wave Propagation routines:

Calculates the second term of the weak form (1): \( \left( \frac{\partial \mathbf{F(\mathbf{U})}^{\delta} }{\partial x}, \mathbf{\psi}^{\delta} \right)_{\Omega_e} \) The variables of the system are $\mathbf{U} = [A,u]^T$ physfield[0] = A physfield[1] = u flux[0] = F[0] = A*u flux[1] = F[1] = u^2/2 + p/rho

Definition at line 311 of file PulseWavePropagation.cpp.

314 {
315  size_t nq = m_vessels[m_currentDomain * m_nVariables]->GetTotPoints();
316  NekDouble domain = m_currentDomain;
317  m_pressure[domain] = Array<OneD, NekDouble>(nq);
318  Array<OneD, NekDouble> dAUdx(nq);
319  NekDouble viscoelasticGradient = 0.0;
320 
321  LibUtilities::Timer timer;
322 
323  for (size_t j = 0; j < nq; ++j)
324  {
325  timer.Start();
326  flux[0][0][j] = physfield[0][j] * physfield[1][j];
327  timer.Stop();
328  timer.AccumulateRegion("PulseWavePropagation:GetFluxVector-flux", 3);
329  }
330 
331  // d/dx of AU, for the viscoelastic tube law and extra fields
332  m_fields[0]->PhysDeriv(flux[0][0], dAUdx);
333 
334  for (size_t j = 0; j < nq; ++j)
335  {
336  if ((j == 0) || (j == nq - 1))
337  {
338  viscoelasticGradient = dAUdx[j];
339  }
340  else
341  {
342  viscoelasticGradient = (dAUdx[j] + dAUdx[j + 1]) / 2;
343  }
344 
345  m_pressureArea->GetPressure(m_pressure[domain][j], m_beta[domain][j],
346  physfield[0][j], m_A_0[domain][j],
347  viscoelasticGradient, m_gamma[domain][j],
348  m_alpha[domain][j]);
349 
350  flux[1][0][j] = physfield[1][j] * physfield[1][j] / 2 +
351  m_pressure[domain][j] / m_rho;
352  }
353 
354  m_session->MatchSolverInfo("OutputExtraFields", "True", extraFields, true);
355 
356  if (extraFields)
357  {
358  /*
359  Calculates wave speed and characteristic variables.
360 
361  Ideally this should be moved to PulseWaveSystem, but it's easiest to
362  implement here.
363  */
364  size_t counter = 0;
365 
366  m_PWV[domain] = Array<OneD, NekDouble>(nq);
367  m_W1[domain] = Array<OneD, NekDouble>(nq);
368  m_W2[domain] = Array<OneD, NekDouble>(nq);
369 
370  for (size_t j = 0; j < nq; ++j)
371  {
372  m_pressureArea->GetC(m_PWV[domain][j], m_beta[domain][j],
373  physfield[0][counter + j], m_A_0[domain][j],
374  m_alpha[domain][j]);
375  m_pressureArea->GetW1(m_W1[domain][j], physfield[1][counter + j],
376  m_beta[domain][j], physfield[0][counter + j],
377  m_A_0[domain][j], m_alpha[domain][j]);
378  m_pressureArea->GetW2(m_W2[domain][j], physfield[1][counter + j],
379  m_beta[domain][j], physfield[0][counter + j],
380  m_A_0[domain][j], m_alpha[domain][j]);
381  }
382 
383  counter += nq;
384  }
385 }
Array< OneD, Array< OneD, NekDouble > > m_A_0
Array< OneD, Array< OneD, NekDouble > > m_W2
PulseWavePressureAreaSharedPtr m_pressureArea
Array< OneD, Array< OneD, NekDouble > > m_W1
Array< OneD, Array< OneD, NekDouble > > m_pressure
Array< OneD, Array< OneD, NekDouble > > m_gamma
Array< OneD, Array< OneD, NekDouble > > m_alpha
Array< OneD, Array< OneD, NekDouble > > m_PWV
Array< OneD, Array< OneD, NekDouble > > m_beta
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
double NekDouble

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::PulseWaveSystem::extraFields, Nektar::PulseWaveSystem::m_A_0, Nektar::PulseWaveSystem::m_alpha, Nektar::PulseWaveSystem::m_beta, Nektar::PulseWaveSystem::m_currentDomain, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::PulseWaveSystem::m_gamma, Nektar::PulseWaveSystem::m_nVariables, Nektar::PulseWaveSystem::m_pressure, Nektar::PulseWaveSystem::m_pressureArea, Nektar::PulseWaveSystem::m_PWV, Nektar::PulseWaveSystem::m_rho, Nektar::SolverUtils::EquationSystem::m_session, Nektar::PulseWaveSystem::m_vessels, Nektar::PulseWaveSystem::m_W1, Nektar::PulseWaveSystem::m_W2, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by v_InitObject().

◆ GetN()

Array< OneD, NekDouble > & Nektar::PulseWavePropagation::GetN ( )

Definition at line 402 of file PulseWavePropagation.cpp.

403 {
405 }
Array< OneD, Array< OneD, NekDouble > > m_trace_fwd_normal

References Nektar::PulseWaveSystem::m_currentDomain, and Nektar::PulseWaveSystem::m_trace_fwd_normal.

Referenced by v_InitObject().

◆ GetRho()

NekDouble Nektar::PulseWavePropagation::GetRho ( )

Definition at line 407 of file PulseWavePropagation.cpp.

408 {
409  return m_rho;
410 }

References Nektar::PulseWaveSystem::m_rho.

Referenced by v_InitObject().

◆ SetPulseWaveBoundaryConditions()

void Nektar::PulseWavePropagation::SetPulseWaveBoundaryConditions ( const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const NekDouble  time 
)
protected
Does the projection between ... space and the ... space. Also checks for

Q-inflow boundary conditions at the inflow of the current arterial segment and applies the Q-inflow if specified

Definition at line 221 of file PulseWavePropagation.cpp.

225 {
226  boost::ignore_unused(outarray);
227 
228  size_t omega;
229 
230  Array<OneD, MultiRegions::ExpListSharedPtr> vessel(2);
231 
232  size_t offset = 0;
233 
234  // This will be moved to the RCR boundary condition once factory is setup
235  if (time == 0)
236  {
237  m_Boundary = Array<OneD, PulseWaveBoundarySharedPtr>(2 * m_nDomains);
238 
239  for (omega = 0; omega < m_nDomains; ++omega)
240  {
241  vessel[0] = m_vessels[2 * omega];
242  vessel[1] = m_vessels[2 * omega + 1];
243 
244  for (size_t j = 0; j < 2; ++j)
245  {
246  std::string BCType;
247 
248  if (j < vessel[0]->GetBndConditions().size())
249  {
250  BCType = vessel[0]->GetBndConditions()[j]->GetUserDefined();
251  }
252 
253  // if no condition given define it to be NoUserDefined
254  if (BCType.empty() || BCType == "Interface")
255  {
256  BCType = "NoUserDefined";
257  }
258 
259  m_Boundary[2 * omega + j] = GetBoundaryFactory().CreateInstance(
261 
262  // turn on time dependent BCs
263  if (BCType == "Q-inflow")
264  {
265  vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
266  }
267  else if (BCType == "A-inflow")
268  {
269  vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
270  }
271  else if (BCType == "U-inflow")
272  {
273  vessel[1]->GetBndConditions()[j]->SetIsTimeDependent(true);
274  }
275  else if (BCType == "RCR-terminal")
276  {
277  vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
278  }
279  }
280  }
281  }
282 
283  SetBoundaryConditions(time);
284 
285  // Loop over all vessels and set boundary conditions
286  LibUtilities::Timer timer;
287  for (omega = 0; omega < m_nDomains; ++omega)
288  {
289  timer.Start();
290  for (size_t n = 0; n < 2; ++n)
291  {
292  m_Boundary[2 * omega + n]->DoBoundary(
293  inarray, m_A_0, m_beta, m_alpha, time, omega, offset, n);
294  }
295 
296  offset += m_vessels[2 * omega]->GetTotPoints();
297  timer.Stop();
298  timer.AccumulateRegion("PulseWavePropagation:_SetBCs", 1);
299  }
300 }
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
Array< OneD, PulseWaveBoundarySharedPtr > m_Boundary
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
BoundaryFactory & GetBoundaryFactory()

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::GetBoundaryFactory(), Nektar::PulseWaveSystem::m_A_0, Nektar::PulseWaveSystem::m_alpha, Nektar::PulseWaveSystem::m_beta, m_Boundary, Nektar::PulseWaveSystem::m_nDomains, Nektar::PulseWaveSystem::m_pressureArea, Nektar::SolverUtils::EquationSystem::m_session, Nektar::PulseWaveSystem::m_vessels, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by DoOdeRhs().

◆ v_GenerateSummary()

void Nektar::PulseWavePropagation::v_GenerateSummary ( SolverUtils::SummaryList s)
overrideprotectedvirtual

Print summary routine, calls virtual routine reimplemented in UnsteadySystem

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 421 of file PulseWavePropagation.cpp.

422 {
424 }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.

References Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

◆ v_InitObject()

void Nektar::PulseWavePropagation::v_InitObject ( bool  DeclareField = false)
overrideprotectedvirtual

Initialisation routine for multidomain solver. Sets up the expansions for every arterial segment (m_vessels) and for one complete field m_outfield which is needed to write the postprocessing output. Also determines which upwind strategy is used (currently only upwinding scheme available) and reads blodd flow specific parameters from the inputfile

Reimplemented from Nektar::PulseWaveSystem.

Definition at line 72 of file PulseWavePropagation.cpp.

73 {
74  boost::ignore_unused(DeclareField);
75 
76  // Will set up an array of vessels/fields in PulseWaveSystem::v_InitObject
77  // so set DeclareField to false so that the fields are not set up in
78  // EquationSystem unnecessarily. Note the number of fields in Equation
79  // system is related to the number of variables. The number of vessels is
80  // therefore held in PulwWaveSystem.
82 
83  if (m_session->DefinesSolverInfo("PressureArea"))
84  {
86  m_session->GetSolverInfo("PressureArea"), m_vessels, m_session);
87  }
88  else
89  {
91  "Beta", m_vessels, m_session);
92  }
93 
95  {
98  }
99  else
100  {
101  ASSERTL0(false, "Implicit Pulse Wave Propagation not set up.");
102  }
103 
104  // Create advection object
105  string advName;
106  string riemName;
107  switch (m_upwindTypePulse)
108  {
109  case eUpwindPulse:
110  {
111  advName = "WeakDG";
112  riemName = "UpwindPulse";
113  }
114  break;
115  default:
116  {
117  ASSERTL0(false, "populate switch statement for upwind flux");
118  }
119  break;
120  }
121  m_advObject =
123  m_advObject->SetFluxVector(&PulseWavePropagation::GetFluxVector, this);
125  riemName, m_session);
126  m_riemannSolver->SetScalar("A0", &PulseWavePropagation::GetA0, this);
127  m_riemannSolver->SetScalar("beta", &PulseWavePropagation::GetBeta, this);
128  m_riemannSolver->SetScalar("alpha", &PulseWavePropagation::GetAlpha, this);
129  m_riemannSolver->SetScalar("N", &PulseWavePropagation::GetN, this);
130  m_riemannSolver->SetParam("rho", &PulseWavePropagation::GetRho, this);
132  this);
133 
134  m_advObject->SetRiemannSolver(m_riemannSolver);
135  m_advObject->InitObject(m_session, m_fields);
136 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
Array< OneD, NekDouble > & GetAlpha()
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
DG Pulse Wave Propagation routines:
Array< OneD, NekDouble > & GetN()
Array< OneD, NekDouble > & GetA0()
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Array< OneD, NekDouble > & GetBeta()
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
virtual void v_InitObject(bool DeclareField=false) override
UpwindTypePulse m_upwindTypePulse
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()
PressureAreaFactory & GetPressureAreaFactory()
@ eUpwindPulse
simple upwinding scheme

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::eUpwindPulse, GetA0(), Nektar::SolverUtils::GetAdvectionFactory(), GetAlpha(), GetBeta(), GetDomains(), GetFluxVector(), GetN(), Nektar::GetPressureAreaFactory(), GetRho(), Nektar::SolverUtils::GetRiemannSolverFactory(), m_advObject, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::PulseWaveSystem::m_pressureArea, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::PulseWaveSystem::m_upwindTypePulse, Nektar::PulseWaveSystem::m_vessels, and Nektar::PulseWaveSystem::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< PulseWavePropagation >

friend class MemoryManager< PulseWavePropagation >
friend

Definition at line 1 of file PulseWavePropagation.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 68 of file PulseWavePropagation.h.

◆ m_advObject

SolverUtils::AdvectionSharedPtr Nektar::PulseWavePropagation::m_advObject
protected

Definition at line 103 of file PulseWavePropagation.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_Boundary

Array<OneD, PulseWaveBoundarySharedPtr> Nektar::PulseWavePropagation::m_Boundary
protected

Definition at line 105 of file PulseWavePropagation.h.

Referenced by SetPulseWaveBoundaryConditions().

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::PulseWavePropagation::m_riemannSolver
protected

Definition at line 102 of file PulseWavePropagation.h.

Referenced by v_InitObject().