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 > & GetN ()
 
NekDouble GetRho ()
 
NekDouble GetPext ()
 
- Public Member Functions inherited from Nektar::PulseWaveSystem
virtual ~PulseWaveSystem ()
 Destructor. More...
 
int GetNdomains ()
 
Array< OneD, MultiRegions::ExpListSharedPtrUpdateVessels (void)
 
void CalcCharacteristicVariables (int omega)
 
- 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...
 
- 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 ()
 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, Array< OneD, NekDouble > &output)
 
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 void SetTime (const NekDouble time)
 
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...
 

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...
 

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 ()
 
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)
 
- 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 ()
 Sets up initial conditions. More...
 
virtual void v_DoSolve ()
 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)
 Riemann Problem for Bifurcation. More...
 
void MergingRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0)
 Riemann Problem for Merging Flow. More...
 
void JunctionRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0)
 Riemann Problem for Junction. More...
 
virtual void v_Output (void)
 
void CheckPoint_Output (const int n)
 
NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 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 void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. 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 ()
 
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...
 
- 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
 
PulseWavePressureAreaSharedPtr m_pressureArea
 
- Protected Attributes inherited from Nektar::PulseWaveSystem
Array< OneD, MultiRegions::ExpListSharedPtrm_vessels
 
int m_nDomains
 
int m_currentDomain
 
int 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_trace_fwd_normal
 
std::vector< std::vector< InterfacePointShPtr > > m_vesselJcts
 
std::vector< std::vector< InterfacePointShPtr > > m_bifurcations
 
std::vector< std::vector< InterfacePointShPtr > > m_mergingJcts
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information. More...
 
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::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
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...
 
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...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
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_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
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_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_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...
 
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...
 
- 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 48 of file PulseWavePropagation.h.

Constructor & Destructor Documentation

◆ ~PulseWavePropagation()

Nektar::PulseWavePropagation::~PulseWavePropagation ( )
virtual

Definition at line 119 of file PulseWavePropagation.cpp.

120 {
121 }

◆ PulseWavePropagation()

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

Definition at line 62 of file PulseWavePropagation.cpp.

65  : PulseWaveSystem(pSession, pGraph)
66 {
67 }
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 54 of file PulseWavePropagation.h.

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

57  {
59  ::AllocateSharedPtr(pSession, pGraph);
60  p->InitObject();
61  return p;
62  }
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.

◆ 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 181 of file PulseWavePropagation.cpp.

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

Referenced by v_InitObject().

184 {
185  // Just copy over array
186  for (int i = 0; i < m_nVariables; ++i)
187  {
188  Vmath::Vcopy(inarray[i].num_elements(), inarray[i], 1, outarray[i], 1);
189  }
190 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ 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 131 of file PulseWavePropagation.cpp.

References 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(), and SetPulseWaveBoundaryConditions().

Referenced by v_InitObject().

134 {
135  int i;
136 
137  Array<OneD, Array<OneD, NekDouble>> physarray(m_nVariables);
138 
139  // Dummy array for WeakDG advection
140  Array<OneD, Array<OneD, NekDouble>> advVel(m_spacedim);
141 
142  // Output array for advection
143  Array<OneD, Array<OneD, NekDouble>> out(m_nVariables);
144 
145  int cnt = 0;
146 
147  // Set up Inflow and Outflow boundary conditions.
148  SetPulseWaveBoundaryConditions(inarray, outarray, time);
149 
150  // Set up any interface conditions and write into boundary condition
152 
153  // do advection evauation in all domains
154  for (int omega = 0; omega < m_nDomains; ++omega)
155  {
156  m_currentDomain = omega;
157  int nq = m_vessels[omega * m_nVariables]->GetTotPoints();
158 
159  for (i = 0; i < m_nVariables; ++i)
160  {
161  physarray[i] = inarray[i] + cnt;
162  out[i] = outarray[i] + cnt;
163  }
164 
165  for (i = 0; i < m_nVariables; ++i)
166  {
167  m_fields[i] = m_vessels[omega * m_nVariables + i];
168  }
169 
170  m_advObject->Advect(m_nVariables, m_fields, advVel, physarray, out,
171  time);
172  for (i = 0; i < m_nVariables; ++i)
173  {
174  Vmath::Neg(nq, out[i], 1);
175  }
176 
177  cnt += nq;
178  }
179 }
void SetPulseWaveBoundaryConditions(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
SolverUtils::AdvectionSharedPtr m_advObject
void EnforceInterfaceConditions(const Array< OneD, const Array< OneD, NekDouble > > &fields)
int m_spacedim
Spatial dimension (>= expansion dim).
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels

◆ GetA0()

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

Definition at line 299 of file PulseWavePropagation.cpp.

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

Referenced by v_InitObject().

300 {
302 }
Array< OneD, Array< OneD, NekDouble > > m_A_0_trace

◆ GetBeta()

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

Definition at line 304 of file PulseWavePropagation.cpp.

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

Referenced by v_InitObject().

305 {
307 }
Array< OneD, Array< OneD, NekDouble > > m_beta_trace

◆ 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 ${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 p-A-relationship: p = p_ext + beta*(sqrt(A)-sqrt(A_0))

Definition at line 276 of file PulseWavePropagation.cpp.

References ASSERTL0, Nektar::PulseWaveSystem::m_A_0, Nektar::PulseWaveSystem::m_beta, Nektar::PulseWaveSystem::m_currentDomain, Nektar::PulseWaveSystem::m_nVariables, Nektar::PulseWaveSystem::m_pext, Nektar::PulseWaveSystem::m_rho, Nektar::PulseWaveSystem::m_vessels, and CellMLToNektar.cellml_metadata::p.

Referenced by v_InitObject().

279 {
280  int nq = m_vessels[m_currentDomain * m_nVariables]->GetTotPoints();
281  NekDouble p = 0.0;
282  NekDouble p_t = 0.0;
283 
284  for (int j = 0; j < nq; j++)
285  {
286  flux[0][0][j] = physfield[0][j] * physfield[1][j];
287 
288  ASSERTL0(physfield[0][j] >= 0, "Negative A not allowed.");
289 
290  p = m_pext +
291  m_beta[m_currentDomain][j] *
292  (sqrt(physfield[0][j]) - sqrt(m_A_0[m_currentDomain][j]));
293 
294  p_t = (physfield[1][j] * physfield[1][j]) / 2 + p / m_rho;
295  flux[1][0][j] = p_t;
296  }
297 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< OneD, Array< OneD, NekDouble > > m_A_0
double NekDouble
Array< OneD, Array< OneD, NekDouble > > m_beta
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels

◆ GetN()

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

Definition at line 309 of file PulseWavePropagation.cpp.

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

Referenced by v_InitObject().

310 {
312 }
Array< OneD, Array< OneD, NekDouble > > m_trace_fwd_normal

◆ GetPext()

NekDouble Nektar::PulseWavePropagation::GetPext ( )

Definition at line 319 of file PulseWavePropagation.cpp.

References Nektar::PulseWaveSystem::m_pext.

Referenced by v_InitObject().

320 {
321  return m_pext;
322 }

◆ GetRho()

NekDouble Nektar::PulseWavePropagation::GetRho ( )

Definition at line 314 of file PulseWavePropagation.cpp.

References Nektar::PulseWaveSystem::m_rho.

Referenced by v_InitObject().

315 {
316  return m_rho;
317 }

◆ 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 197 of file PulseWavePropagation.cpp.

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

Referenced by DoOdeRhs().

201 {
202  int omega;
203 
204  Array<OneD, MultiRegions::ExpListSharedPtr> vessel(2);
205 
206  int offset = 0;
207 
208  // This will be moved to the RCR boundary condition once factory is setup
209  if (time == 0)
210  {
211  m_Boundary = Array<OneD, PulseWaveBoundarySharedPtr>(2 * m_nDomains);
212 
213  for (omega = 0; omega < m_nDomains; ++omega)
214  {
215  vessel[0] = m_vessels[2 * omega];
216  vessel[1] = m_vessels[2 * omega + 1];
217 
218  for (int j = 0; j < 2; ++j)
219  {
220  std::string BCType =
221  vessel[0]->GetBndConditions()[j]->GetUserDefined();
222  if (BCType.empty()) // if not condition given define it to be
223  // NoUserDefined
224  {
225  BCType = "NoUserDefined";
226  }
227 
228  m_Boundary[2 * omega + j] = GetBoundaryFactory().CreateInstance(
230 
231  // turn on time depedent BCs
232  if (BCType == "Q-inflow")
233  {
234  vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
235  }
236  if (BCType == "A-inflow")
237  {
238  vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
239  }
240  if (BCType == "U-inflow")
241  {
242  vessel[1]->GetBndConditions()[j]->SetIsTimeDependent(true);
243  }
244  else if (BCType == "RCR-terminal")
245  {
246  vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
247  }
248  }
249  }
250  }
251 
252  SetBoundaryConditions(time);
253 
254  // Loop over all vessesls and set boundary conditions
255  for (omega = 0; omega < m_nDomains; ++omega)
256  {
257  for (int n = 0; n < 2; ++n)
258  {
259  m_Boundary[2 * omega + n]->DoBoundary(inarray, m_A_0, m_beta, time,
260  omega, offset, n);
261  }
262  offset += m_vessels[2 * omega]->GetTotPoints();
263  }
264 }
BoundaryFactory & GetBoundaryFactory()
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
Array< OneD, Array< OneD, NekDouble > > m_A_0
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
Array< OneD, Array< OneD, NekDouble > > m_beta
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, PulseWaveBoundarySharedPtr > m_Boundary
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
PulseWavePressureAreaSharedPtr m_pressureArea

◆ v_GenerateSummary()

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

Print summary routine, calls virtual routine reimplemented in UnsteadySystem

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 328 of file PulseWavePropagation.cpp.

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

329 {
331 }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.

◆ v_InitObject()

void Nektar::PulseWavePropagation::v_InitObject ( )
protectedvirtual

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

Gets the Material Properties of each arterial segment
specified in the inputfile from section MaterialProperties

Also gets the Area at static equilibrium A_0 specified in the inputfile.

Having found these points also extract the values at the trace points and the normal direction consistent with the left adjacent definition of Fwd and Bwd

Reimplemented from Nektar::PulseWaveSystem.

Definition at line 69 of file PulseWavePropagation.cpp.

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(), GetBeta(), GetFluxVector(), GetN(), GetPext(), Nektar::GetPressureAreaFactory(), GetRho(), Nektar::SolverUtils::GetRiemannSolverFactory(), m_advObject, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_ode, m_pressureArea, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::PulseWaveSystem::m_upwindTypePulse, Nektar::PulseWaveSystem::m_vessels, and Nektar::PulseWaveSystem::v_InitObject().

70 {
72 
74  "Lymphatic", m_vessels, m_session);
75  m_pressureArea->DoPressure();
76 
78  {
81  }
82  else
83  {
84  ASSERTL0(false, "Implicit Pulse Wave Propagation not set up.");
85  }
86 
87  // Create advection object
88  string advName;
89  string riemName;
90  switch (m_upwindTypePulse)
91  {
92  case eUpwindPulse:
93  {
94  advName = "WeakDG";
95  riemName = "UpwindPulse";
96  }
97  break;
98  default:
99  {
100  ASSERTL0(false, "populate switch statement for upwind flux");
101  }
102  break;
103  }
104  m_advObject =
106  m_advObject->SetFluxVector(&PulseWavePropagation::GetFluxVector, this);
108  riemName, m_session);
109  m_riemannSolver->SetScalar("A0", &PulseWavePropagation::GetA0, this);
110  m_riemannSolver->SetScalar("beta", &PulseWavePropagation::GetBeta, this);
111  m_riemannSolver->SetScalar("N", &PulseWavePropagation::GetN, this);
112  m_riemannSolver->SetParam("rho", &PulseWavePropagation::GetRho, this);
113  m_riemannSolver->SetParam("pext", &PulseWavePropagation::GetPext, this);
114 
115  m_advObject->SetRiemannSolver(m_riemannSolver);
116  m_advObject->InitObject(m_session, m_fields);
117 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
UpwindTypePulse m_upwindTypePulse
Array< OneD, NekDouble > & GetA0()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
SolverUtils::AdvectionSharedPtr m_advObject
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
DG Pulse Wave Propagation routines:
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Array< OneD, NekDouble > & GetN()
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
RiemannSolverFactory & GetRiemannSolverFactory()
Array< OneD, NekDouble > & GetBeta()
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
PressureAreaFactory & GetPressureAreaFactory()
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
PulseWavePressureAreaSharedPtr m_pressureArea
simple upwinding scheme

Friends And Related Function Documentation

◆ MemoryManager< PulseWavePropagation >

friend class MemoryManager< PulseWavePropagation >
friend

Definition at line 51 of file PulseWavePropagation.h.

Member Data Documentation

◆ className

string Nektar::PulseWavePropagation::className
static
Initial value:
=
"PulseWavePropagation", PulseWavePropagation::create,
"Pulse Wave Propagation equation.")

Name of class.

Definition at line 65 of file PulseWavePropagation.h.

◆ m_advObject

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

Definition at line 99 of file PulseWavePropagation.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_Boundary

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

Definition at line 101 of file PulseWavePropagation.h.

Referenced by SetPulseWaveBoundaryConditions().

◆ m_pressureArea

PulseWavePressureAreaSharedPtr Nektar::PulseWavePropagation::m_pressureArea
protected

Definition at line 103 of file PulseWavePropagation.h.

Referenced by SetPulseWaveBoundaryConditions(), and v_InitObject().

◆ m_riemannSolver

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

Definition at line 98 of file PulseWavePropagation.h.

Referenced by v_InitObject().