Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::PulseWavePropagation:
Collaboration graph
[legend]

Public Member Functions

virtual ~PulseWavePropagation ()
- Public Member Functions inherited from Nektar::PulseWaveSystem
virtual ~PulseWaveSystem ()
 Destructor.
int GetNdomains ()
Array< OneD,
MultiRegions::ExpListSharedPtr
UpdateVessels (void)
void CalcCharacteristicVariables (int omega)
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object.
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
template<class T >
boost::shared_ptr< T > as ()
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name.
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available.
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow.
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
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.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $.
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection.
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection.
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
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.
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
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.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
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.
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.
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.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point.
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time.
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 GetNumElmVelocity ()
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 SetStepsToOne ()
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
SOLVER_UTILS_EXPORT void FwdTransFields ()
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison.

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 Creates an instance of this class.

Static Public Attributes

static std::string className = GetEquationSystemFactory().RegisterCreatorFunction("PulseWavePropagation", PulseWavePropagation::create, "Pulse Wave Propagation equation.")
 Name of class.

Protected Member Functions

 PulseWavePropagation (const LibUtilities::SessionReaderSharedPtr &pSession)
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 ()
virtual void v_GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
virtual void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 DG Pulse Wave Propagation routines: Numerical Flux at interelemental boundaries.
void RiemannSolverUpwind (NekDouble AL, NekDouble uL, NekDouble AR, NekDouble uR, NekDouble &Aflux, NekDouble &uflux, NekDouble A_0, NekDouble beta, NekDouble n)
 Upwinding Riemann solver for interelemental boundaries.
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
- Protected Member Functions inherited from Nektar::PulseWaveSystem
 PulseWaveSystem (const LibUtilities::SessionReaderSharedPtr &m_session)
 Initialises PulseWaveSystem class members.
virtual void v_DoInitialise ()
 Sets up initial conditions.
virtual void v_DoSolve ()
 Solves an unsteady problem.
void LinkSubdomains (Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &fields)
 Links the subdomains.
void BifurcationRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0)
 Riemann Problem for Bifurcation.
void MergingRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0)
 Riemann Problem for Merging Flow.
void JunctionRiemann (Array< OneD, NekDouble > &Au, Array< OneD, NekDouble > &uu, Array< OneD, NekDouble > &beta, Array< OneD, NekDouble > &A_0)
 Riemann Problem for Junction.
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.
NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Compute the L_inf error between fields and a given exact solution.
void WriteVessels (const std::string &outname)
 Write input fields to the given filename.
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)
 Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
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.
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members.
int nocase_cmp (const string &s1, const string &s2)
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time.
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space.
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)
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
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

Array< OneD,
PulseWaveBoundarySharedPtr
m_Boundary
PulseWavePressureAreaSharedPtr m_pressureArea
- Protected Attributes inherited from Nektar::PulseWaveSystem
Array< OneD,
MultiRegions::ExpListSharedPtr
m_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.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme.
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
NekDouble m_epsilon
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used.
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used.
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used.
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state.
std::vector< int > m_intVariables
std::vector< FilterSharedPtrm_filters
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
std::string m_sessionName
 Name of the session.
NekDouble m_time
 Current time of simulation.
NekDouble m_fintime
 Finish time of the simulation.
NekDouble m_timestep
 Time step size.
NekDouble m_lambda
 Lambda constant in real system if one required.
NekDouble m_checktime
 Time between checkpoints.
int m_steps
 Number of steps to take.
int m_checksteps
 Number of steps between checkpoints.
int m_spacedim
 Spatial dimension (>= expansion dim).
int m_expdim
 Expansion dimension.
bool m_SingleMode
 Flag to determine if single homogeneous mode is used.
bool m_HalfMode
 Flag to determine if half homogeneous mode is used.
bool m_MultipleModes
 Flag to determine if use multiple homogenenous modes are used.
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
enum HomogeneousType m_HomogeneousType
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
int m_npointsX
 number of points in X direction (if homogeneous)
int m_npointsY
 number of points in Y direction (if homogeneous)
int m_npointsZ
 number of points in Z direction (if homogeneous)
int m_HomoDirec
 number of homogenous directions
int m_NumMode
 Mode to use in case of single mode analysis.

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).
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...

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 47 of file PulseWavePropagation.h.

Constructor & Destructor Documentation

Nektar::PulseWavePropagation::~PulseWavePropagation ( )
virtual

Definition at line 78 of file PulseWavePropagation.cpp.

{
}
Nektar::PulseWavePropagation::PulseWavePropagation ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 55 of file PulseWavePropagation.cpp.

: PulseWaveSystem(pSession)
{
}

Member Function Documentation

static EquationSystemSharedPtr Nektar::PulseWavePropagation::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file PulseWavePropagation.h.

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

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

Referenced by v_InitObject().

{
// Just copy over array
for(int i = 0; i < m_nVariables; ++i)
{
Vmath::Vcopy(inarray[i].num_elements(),inarray[i],1,outarray[i],1);
}
}
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, the routine WeakDGAdvection will be called which then calls v_GetFluxVector and v_NumericalFlux implemented in the PulseWavePropagation class.

Definition at line 92 of file PulseWavePropagation.cpp.

References Nektar::PulseWaveSystem::EnforceInterfaceConditions(), Nektar::PulseWaveSystem::m_currentDomain, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::PulseWaveSystem::m_nDomains, Nektar::PulseWaveSystem::m_nVariables, Nektar::PulseWaveSystem::m_vessels, Vmath::Neg(), SetPulseWaveBoundaryConditions(), and Nektar::SolverUtils::EquationSystem::WeakDGAdvection().

Referenced by v_InitObject().

{
int i;
Array<OneD, Array<OneD, NekDouble> > physarray(m_nVariables);
Array<OneD, Array<OneD, NekDouble> > modarray (m_nVariables);
Array<OneD, NekDouble> tmpArray;
int cnt = 0;
// Set up Inflow and Outflow boundary conditions.
SetPulseWaveBoundaryConditions(inarray, outarray, time);
// Set up any interface conditions and write into boundary condition
// do advection evauation in all domains
for(int omega=0; omega < m_nDomains; ++omega)
{
m_currentDomain = omega;
int nq = m_vessels[omega*m_nVariables]->GetTotPoints();
int ncoeffs = m_vessels[omega*m_nVariables]->GetNcoeffs();
for (i = 0; i < m_nVariables; ++i)
{
physarray[i] = inarray[i]+cnt;
modarray[i] = Array<OneD, NekDouble>(ncoeffs);
}
for(i = 0; i < m_nVariables; ++i)
{
m_fields[i] = m_vessels[omega*m_nVariables+ i];
}
WeakDGAdvection(physarray, modarray, true, true);
for(i = 0; i < m_nVariables; ++i)
{
Vmath::Neg(ncoeffs,modarray[i],1);
}
for(i = 0; i < m_nVariables; ++i)
{
m_vessels[omega*m_nVariables+i]->MultiplyByElmtInvMass(modarray[i],modarray[i]);
m_vessels[omega*m_nVariables+i]->BwdTrans(modarray[i],tmpArray = outarray[i]+cnt);
}
cnt += nq;
}
}
void Nektar::PulseWavePropagation::RiemannSolverUpwind ( NekDouble  AL,
NekDouble  uL,
NekDouble  AR,
NekDouble  uR,
NekDouble Aflux,
NekDouble uflux,
NekDouble  A_0,
NekDouble  beta,
NekDouble  n 
)
protected

Upwinding Riemann solver for interelemental boundaries.

Riemann solver for upwinding at an interface between two elements. Uses the characteristic variables for calculating the upwinded state $(A_u,u_u)$ from the left $(A_L,u_L)$ and right state $(A_R,u_R)$. Returns the upwinded flux ${F}^u$ needed for the weak formulation (1). Details can be found in "Pulse wave propagation in the human vascular system", section 3.3

Definition at line 322 of file PulseWavePropagation.cpp.

References ASSERTL1, Nektar::PulseWaveSystem::m_pext, and Nektar::PulseWaveSystem::m_rho.

Referenced by v_NumericalFlux().

{
Array<OneD, NekDouble> W(2);
Array<OneD, NekDouble> upwindedphysfield(2);
NekDouble cL = 0.0;
NekDouble cR = 0.0;
NekDouble rho = m_rho;
NekDouble pext = m_pext;
NekDouble p = 0.0;
NekDouble p_t = 0.0;
// Compute the wave speeds. The use of the normal here allows
// for the definition of the characteristics to be inverted
// (and hence the left and right state) if n is in the -ve
// x-direction. This means we end up with the positive
// defintion of the flux which has to therefore be multiplied
// by the normal at the end of the methods This is a bit of a
// mind twister but is efficient from a coding perspective.
cL = sqrt(beta*sqrt(AL)/(2*rho))*n;
cR = sqrt(beta*sqrt(AR)/(2*rho))*n;
ASSERTL1(fabs(cL+cR) > fabs(uL+uR),"Conditions are not sub-sonic");
// If upwinding from left and right for subsonic domain
// then know characteristics immediately
W[0] = uL + 4*cL;
W[1] = uR - 4*cR;
// Calculate conservative variables from characteristics
NekDouble w0mw1 = 0.25*(W[0]-W[1]);
NekDouble fac = rho/(2*beta);
w0mw1 *= w0mw1; // squared
w0mw1 *= w0mw1; // fourth power
fac *= fac; // squared
upwindedphysfield[0]= w0mw1*fac;
upwindedphysfield[1]= 0.5*(W[0] + W[1]);
// Compute the fluxes multipled by the normal.
Aflux = upwindedphysfield[0] * upwindedphysfield[1]*n;
p = pext + beta*(sqrt(upwindedphysfield[0]) - sqrt(A_0));
p_t = 0.5*(upwindedphysfield[1]*upwindedphysfield[1]) + p/rho;
uflux = p_t*n;
}
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 162 of file PulseWavePropagation.cpp.

References Nektar::LibUtilities::NekFactory< tKey, tBase, >::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().

{
int omega;
Array<OneD, MultiRegions::ExpListSharedPtr> vessel(2);
int offset=0;
//This will be moved to the RCR boundary condition once factory is setup
if (time == 0)
{
m_Boundary = Array<OneD,PulseWaveBoundarySharedPtr>(2*m_nDomains);
for(omega = 0; omega < m_nDomains; ++omega)
{
vessel[0] = m_vessels[2*omega];
for(int j = 0; j < 2; ++j)
{
std::string BCType = vessel[0]->GetBndConditions()[j]->
GetBndTypeAsString(vessel[0]->GetBndConditions()[j]->GetUserDefined());
}
}
}
// Loop over all vessesls and set boundary conditions
for(omega = 0; omega < m_nDomains; ++omega)
{
for(int n = 0; n < 2; ++n)
{
m_Boundary[2*omega+n]->DoBoundary(inarray,m_A_0,m_beta,time,omega,offset,n);
}
offset += m_vessels[2*omega]->GetTotPoints();
}
}
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 376 of file PulseWavePropagation.cpp.

void Nektar::PulseWavePropagation::v_GetFluxVector ( const int  i,
Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protectedvirtual

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

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 218 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, and Nektar::PulseWaveSystem::m_vessels.

{
int nq = m_vessels[m_currentDomain*m_nVariables]->GetTotPoints();
NekDouble p = 0.0;
NekDouble p_t = 0.0;
switch (i)
{
case 0: // Flux for A equation
{
for (int j = 0; j < nq; j++)
{
flux[0][j] = physfield[0][j]*physfield[1][j];
}
}
break;
case 1: // Flux for u equation
{
for (int j = 0; j < nq; j++)
{
ASSERTL0(physfield[0][j]>=0,"Negative A not allowed.");
(sqrt(physfield[0][j]) - sqrt(m_A_0[m_currentDomain][j]));
p_t = (physfield[1][j]*physfield[1][j])/2 + p/m_rho;
flux[0][j] = p_t;
}
}
break;
default:
ASSERTL0(false,"GetFluxVector: illegal vector index");
break;
}
}
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 60 of file PulseWavePropagation.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::GetPressureAreaFactory(), Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::UnsteadySystem::m_ode, m_pressureArea, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::PulseWaveSystem::m_vessels.

void Nektar::PulseWavePropagation::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numflux 
)
protectedvirtual

DG Pulse Wave Propagation routines: Numerical Flux at interelemental boundaries.

Calculates the third term of the weak form (1): numerical flux at boundary $ \left[ \mathbf{\psi}^{\delta} \cdot \{ \mathbf{F}^u - \mathbf{F}(\mathbf{U}^{\delta}) \} \right]_{x_e^l}^{x_eû} $

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 262 of file PulseWavePropagation.cpp.

References ASSERTL0, Nektar::eUpwindPulse, Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::PulseWaveSystem::m_A_0_trace, Nektar::PulseWaveSystem::m_beta_trace, Nektar::PulseWaveSystem::m_currentDomain, Nektar::PulseWaveSystem::m_nVariables, Nektar::PulseWaveSystem::m_trace_fwd_normal, Nektar::PulseWaveSystem::m_upwindTypePulse, Nektar::PulseWaveSystem::m_vessels, and RiemannSolverUpwind().

{
int i;
int nTracePts = GetTraceTotPoints();
Array<OneD, Array<OneD, NekDouble> > Fwd(m_nVariables);
Array<OneD, Array<OneD, NekDouble> > Bwd(m_nVariables);
for (i = 0; i < m_nVariables; ++i)
{
Fwd[i] = Array<OneD, NekDouble>(nTracePts);
Bwd[i] = Array<OneD, NekDouble>(nTracePts);
}
// Get the physical values at the trace
for (i = 0; i < m_nVariables; ++i)
{
m_vessels[m_currentDomain*m_nVariables+ i]->
GetFwdBwdTracePhys(physfield[i],Fwd[i],Bwd[i]);
}
// Solve the upwinding Riemann problem within one arterial
// segment by calling the upwinding Riemann solver implemented
// in this file
NekDouble Aflux, uflux;
for (i = 0; i < nTracePts; ++i)
{
{
{
RiemannSolverUpwind(Fwd[0][i],Fwd[1][i],Bwd[0][i],Bwd[1][i],
Aflux, uflux, m_A_0_trace[m_currentDomain][i],
}
break;
default:
{
ASSERTL0(false,"populate switch statement for upwind flux");
}
break;
}
numflux[0][i] = Aflux;
numflux[1][i] = uflux;
}
}

Friends And Related Function Documentation

friend class MemoryManager< PulseWavePropagation >
friend

Definition at line 50 of file PulseWavePropagation.h.

Member Data Documentation

string Nektar::PulseWavePropagation::className = GetEquationSystemFactory().RegisterCreatorFunction("PulseWavePropagation", PulseWavePropagation::create, "Pulse Wave Propagation equation.")
static

Name of class.

Definition at line 61 of file PulseWavePropagation.h.

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

Definition at line 95 of file PulseWavePropagation.h.

Referenced by SetPulseWaveBoundaryConditions().

PulseWavePressureAreaSharedPtr Nektar::PulseWavePropagation::m_pressureArea
protected

Definition at line 97 of file PulseWavePropagation.h.

Referenced by SetPulseWaveBoundaryConditions(), and v_InitObject().