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

A model for cardiac conduction. More...

#include <Monodomain.h>

Inheritance diagram for Nektar::Monodomain:
[legend]

Public Member Functions

 ~Monodomain () override
 Desctructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT ~UnsteadySystem () override
 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 NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 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 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 NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. 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 GetTime ()
 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 Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
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 void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int 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...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::ALEHelper
virtual ~ALEHelper ()=default
 
virtual SOLVER_UTILS_EXPORT void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
SOLVER_UTILS_EXPORT void InitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity (const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_ALEPreMultiplyMass (Array< OneD, Array< OneD, NekDouble > > &fields)
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMass (Array< OneD, Array< OneD, NekDouble > > &traceNormals, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble time)
 Update m_fields with u^n by multiplying by inverse mass matrix. That's then used in e.g. checkpoint output and L^2 error calculation. More...
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void MoveMesh (const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
 
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace ()
 
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

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

 Monodomain (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
 Solve for the diffusion term. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Computes the reaction terms \(f(u,v)\) and \(g(u,v)\). More...
 
void v_SetInitialConditions (NekDouble initialtime, bool dumpInitialConditions, const int domain) override
 Sets a custom initial condition. More...
 
void v_GenerateSummary (SummaryList &s) override
 Prints a summary of the model parameters. More...
 
- 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 void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_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)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
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...
 
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_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. 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 void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Private Member Functions

void LoadStimuli ()
 

Private Attributes

CellModelSharedPtr m_cell
 Cell model. More...
 
std::vector< StimulusSharedPtrm_stimulus
 
StdRegions::VarCoeffMap m_vardiff
 Variable diffusivity. More...
 
NekDouble m_chi
 
NekDouble m_capMembrane
 
NekDouble m_stimDuration
 Stimulus current. More...
 

Friends

class MemoryManager< Monodomain >
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
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...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. 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_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_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time 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_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_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z. More...
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData 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...
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_ImplicitALESolver = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

A model for cardiac conduction.

Base model of cardiac electrophysiology of the form

\begin{align*} \frac{\partial u}{\partial t} = \nabla^2 u + J_{ion}, \end{align*}

where the reaction term, \(J_{ion}\) is defined by a specific cell model.

This implementation, at present, treats the reaction terms explicitly and the diffusive element implicitly.

Definition at line 48 of file Monodomain.h.

Constructor & Destructor Documentation

◆ ~Monodomain()

Nektar::Monodomain::~Monodomain ( )
override

Desctructor.

Definition at line 302 of file Monodomain.cpp.

303{
304}

◆ Monodomain()

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

Constructor.

Definition at line 70 of file Monodomain.cpp.

72 : UnsteadySystem(pSession, pGraph)
73{
74}
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 54 of file Monodomain.h.

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

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

◆ DoImplicitSolve()

void Nektar::Monodomain::DoImplicitSolve ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
NekDouble  time,
NekDouble  lambda 
)
protected

Solve for the diffusion term.

Parameters
inarrayInput array.
outarrayOutput array.
timeCurrent simulation time.
lambdaTimestep.

Definition at line 312 of file Monodomain.cpp.

316{
317 int nvariables = inarray.size();
318 int nq = m_fields[0]->GetNpoints();
320 // lambda = \Delta t
322
323 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
324 // inarray = input: \hat{rhs} -> output: \hat{Y}
325 // outarray = output: nabla^2 \hat{Y}
326 // where \hat = modal coeffs
327 for (int i = 0; i < nvariables; ++i)
328 {
329 // Multiply 1.0/timestep
330 Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
331 m_fields[i]->UpdatePhys(), 1);
332
333 // Solve a system of equations with Helmholtz solver and transform
334 // back into physical space.
335 m_fields[i]->HelmSolve(m_fields[i]->GetPhys(),
336 m_fields[i]->UpdateCoeffs(), factors, m_vardiff);
337
338 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
339 m_fields[i]->UpdatePhys());
340 m_fields[i]->SetPhysState(true);
341
342 // Copy the solution vector (required as m_fields must be set).
343 outarray[i] = m_fields[i]->GetPhys();
344 }
345}
StdRegions::VarCoeffMap m_vardiff
Variable diffusivity.
Definition: Monodomain.h:103
NekDouble m_capMembrane
Definition: Monodomain.h:106
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References Nektar::StdRegions::eFactorLambda, Nektar::VarcoeffHashingTest::factors, m_capMembrane, m_chi, Nektar::SolverUtils::EquationSystem::m_fields, m_vardiff, and Vmath::Smul().

Referenced by v_InitObject().

◆ DoOdeRhs()

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

Computes the reaction terms \(f(u,v)\) and \(g(u,v)\).

Definition at line 350 of file Monodomain.cpp.

353{
354 // Compute I_ion
355 m_cell->TimeIntegrate(inarray, outarray, time);
356
357 // Compute I_stim
358 for (unsigned int i = 0; i < m_stimulus.size(); ++i)
359 {
360 m_stimulus[i]->Update(outarray, time);
361 }
362}
CellModelSharedPtr m_cell
Cell model.
Definition: Monodomain.h:98
std::vector< StimulusSharedPtr > m_stimulus
Definition: Monodomain.h:100

References m_cell, and m_stimulus.

Referenced by v_InitObject().

◆ LoadStimuli()

void Nektar::Monodomain::LoadStimuli ( )
private

◆ v_GenerateSummary()

void Nektar::Monodomain::v_GenerateSummary ( SummaryList s)
overrideprotectedvirtual

Prints a summary of the model parameters.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 379 of file Monodomain.cpp.

380{
382 if (m_session->DefinesFunction("d00") &&
383 m_session->GetFunctionType("d00", "intensity") ==
385 {
387 s, "Diffusivity-x",
388 m_session->GetFunction("d00", "intensity")->GetExpression());
389 }
390 if (m_session->DefinesFunction("d11") &&
391 m_session->GetFunctionType("d11", "intensity") ==
393 {
395 s, "Diffusivity-y",
396 m_session->GetFunction("d11", "intensity")->GetExpression());
397 }
398 if (m_session->DefinesFunction("d22") &&
399 m_session->GetFunctionType("d22", "intensity") ==
401 {
403 s, "Diffusivity-z",
404 m_session->GetFunction("d22", "intensity")->GetExpression());
405 }
406 m_cell->GenerateSummary(s);
407}
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47

References Nektar::SolverUtils::AddSummaryItem(), Nektar::LibUtilities::eFunctionTypeExpression, m_cell, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

◆ v_InitObject()

void Nektar::Monodomain::v_InitObject ( bool  DeclareField = true)
overrideprotectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 79 of file Monodomain.cpp.

80{
81 UnsteadySystem::v_InitObject(DeclareField);
82
83 m_session->LoadParameter("Chi", m_chi);
84 m_session->LoadParameter("Cm", m_capMembrane);
85
86 std::string vCellModel;
87 m_session->LoadSolverInfo("CELLMODEL", vCellModel, "");
88
89 ASSERTL0(vCellModel != "", "Cell Model not specified.");
90
92 m_fields[0]);
93
94 m_intVariables.push_back(0);
95
96 // Load variable coefficients
97 StdRegions::VarCoeffType varCoeffEnum[6] = {
101 std::string varCoeffString[6] = {"xx", "xy", "yy", "xz", "yz", "zz"};
102 std::string aniso_var[3] = {"fx", "fy", "fz"};
103
104 const int nq = m_fields[0]->GetNpoints();
105 const int nVarDiffCmpts = m_spacedim * (m_spacedim + 1) / 2;
106
107 // Allocate storage for variable coeffs and initialize to 1.
108 for (int i = 0, k = 0; i < m_spacedim; ++i)
109 {
110 for (int j = 0; j < i + 1; ++j)
111 {
112 if (i == j)
113 {
114 m_vardiff[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 1.0);
115 }
116 else
117 {
118 m_vardiff[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 0.0);
119 }
120 ++k;
121 }
122 }
123
124 // Apply fibre map f \in [0,1], scale to conductivity range
125 // [o_min,o_max], specified by the session parameters o_min and o_max
126 if (m_session->DefinesFunction("AnisotropicConductivity"))
127 {
128 if (m_session->DefinesCmdLineArgument("verbose"))
129 {
130 cout << "Loading Anisotropic Fibre map." << endl;
131 }
132
133 NekDouble o_min = m_session->GetParameter("o_min");
134 NekDouble o_max = m_session->GetParameter("o_max");
135 int k = 0;
136
137 Array<OneD, NekDouble> vTemp_i;
138 Array<OneD, NekDouble> vTemp_j;
139
140 /*
141 * Diffusivity matrix D is upper triangular and defined as
142 * d_00 d_01 d_02
143 * d_11 d_12
144 * d_22
145 *
146 * Given a principle fibre direction _f_ the diffusivity is given
147 * by
148 * d_ij = { D_2 + (D_1 - D_2) f_i f_j if i==j
149 * { (D_1 - D_2) f_i f_j if i!=j
150 *
151 * The vector _f_ is given in terms of the variables fx,fy,fz in the
152 * function AnisotropicConductivity. The values of D_1 and D_2 are
153 * the parameters o_max and o_min, respectively.
154 */
155
156 // Loop through columns of D
157 for (int j = 0; j < m_spacedim; ++j)
158 {
159 ASSERTL0(m_session->DefinesFunction("AnisotropicConductivity",
160 aniso_var[j]),
161 "Function 'AnisotropicConductivity' not correctly "
162 "defined.");
163 GetFunction("AnisotropicConductivity")
164 ->Evaluate(aniso_var[j], vTemp_j);
165
166 // Loop through rows of D
167 for (int i = 0; i < j + 1; ++i)
168 {
169 ASSERTL0(m_session->DefinesFunction("AnisotropicConductivity",
170 aniso_var[i]),
171 "Function 'AnisotropicConductivity' not correctly "
172 "defined.");
173 GetFunction("AnisotropicConductivity")
174 ->Evaluate(aniso_var[i], vTemp_i);
175
176 Array<OneD, NekDouble> tmp(vTemp_i.size());
177
178 Vmath::Vmul(nq, vTemp_i, 1, vTemp_j, 1, tmp, 1);
179 Vmath::Smul(nq, o_max - o_min, tmp, 1, tmp, 1);
180
181 if (i == j)
182 {
183 Vmath::Sadd(nq, o_min, tmp, 1, tmp, 1);
184 }
185
186 m_vardiff[varCoeffEnum[k]] = tmp;
187 ++k;
188 }
189 }
190 }
191 else
192 {
193 // Otherwise apply isotropic conductivity value (o_max) to
194 // diagonal components of tensor
195 NekDouble o_max = m_session->GetParameter("o_max");
196 for (int i = 0; i < nVarDiffCmpts; ++i)
197 {
198 Array<OneD, NekDouble> tmp(m_vardiff[varCoeffEnum[i]].GetValue());
199 Vmath::Smul(nq, o_max, tmp, 1, tmp, 1);
200 m_vardiff[varCoeffEnum[i]] = tmp;
201 }
202 }
203
204 // Scale by scar map (range 0->1) derived from intensity map
205 // (range d_min -> d_max)
206 if (m_session->DefinesFunction("IsotropicConductivity"))
207 {
208 if (m_session->DefinesCmdLineArgument("verbose"))
209 {
210 cout << "Loading Isotropic Conductivity map." << endl;
211 }
212
213 const std::string varName = "intensity";
214 Array<OneD, NekDouble> vTemp;
215 GetFunction("IsotropicConductivity")->Evaluate(varName, vTemp);
216
217 // If the d_min and d_max parameters are defined, then we need to
218 // rescale the isotropic conductivity to convert from the source
219 // domain (e.g. late-gad intensity) to conductivity
220 if (m_session->DefinesParameter("d_min") ||
221 m_session->DefinesParameter("d_max"))
222 {
223 const NekDouble f_min = m_session->GetParameter("d_min");
224 const NekDouble f_max = m_session->GetParameter("d_max");
225 const NekDouble scar_min = 0.0;
226 const NekDouble scar_max = 1.0;
227
228 // Threshold based on d_min, d_max
229 for (int j = 0; j < nq; ++j)
230 {
231 vTemp[j] = (vTemp[j] < f_min ? f_min : vTemp[j]);
232 vTemp[j] = (vTemp[j] > f_max ? f_max : vTemp[j]);
233 }
234
235 // Rescale to s \in [0,1] (0 maps to d_max, 1 maps to d_min)
236 Vmath::Sadd(nq, -f_min, vTemp, 1, vTemp, 1);
237 Vmath::Smul(nq, -1.0 / (f_max - f_min), vTemp, 1, vTemp, 1);
238 Vmath::Sadd(nq, 1.0, vTemp, 1, vTemp, 1);
239 Vmath::Smul(nq, scar_max - scar_min, vTemp, 1, vTemp, 1);
240 Vmath::Sadd(nq, scar_min, vTemp, 1, vTemp, 1);
241 }
242
243 // Scale anisotropic conductivity values
244 for (int i = 0; i < nVarDiffCmpts; ++i)
245 {
246 Array<OneD, NekDouble> tmp = m_vardiff[varCoeffEnum[i]].GetValue();
247 Vmath::Vmul(nq, vTemp, 1, tmp, 1, tmp, 1);
248 m_vardiff[varCoeffEnum[i]] = tmp;
249 }
250 }
251
252 // Write out conductivity values
253 for (int j = 0, k = 0; j < m_spacedim; ++j)
254 {
255 // Loop through rows of D
256 for (int i = 0; i < j + 1; ++i)
257 {
258 // Transform variable coefficient and write out to file.
259 m_fields[0]->FwdTransLocalElmt(
260 m_vardiff[varCoeffEnum[k]].GetValue(),
261 m_fields[0]->UpdateCoeffs());
262 std::stringstream filename;
263 filename << "Conductivity_" << varCoeffString[k] << ".fld";
264 WriteFld(filename.str());
265
266 ++k;
267 }
268 }
269
270 // Search through the loaded filters and pass the cell model to any
271 // CheckpointCellModel filters loaded.
272 for (auto &x : m_filters)
273 {
274 if (x.first == "CheckpointCellModel")
275 {
276 std::shared_ptr<FilterCheckpointCellModel> c =
277 std::dynamic_pointer_cast<FilterCheckpointCellModel>(x.second);
278 c->SetCellModel(m_cell);
279 }
280 if (x.first == "CellHistoryPoints")
281 {
282 std::shared_ptr<FilterCellHistoryPoints> c =
283 std::dynamic_pointer_cast<FilterCellHistoryPoints>(x.second);
284 c->SetCellModel(m_cell);
285 }
286 }
287
288 // Load stimuli
290
292 {
294 }
297}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Solve for the diffusion term.
Definition: Monodomain.cpp:312
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Computes the reaction terms and .
Definition: Monodomain.cpp:350
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
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.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
static std::vector< StimulusSharedPtr > LoadStimuli(const LibUtilities::SessionReaderSharedPtr &pSession, const MultiRegions::ExpListSharedPtr &pField)
Definition: Stimulus.cpp:89
CellModelFactory & GetCellModelFactory()
Definition: CellModel.cpp:46
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), Nektar::SolverUtils::UnsteadySystem::DoDummyProjection(), DoImplicitSolve(), DoOdeRhs(), Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD01, Nektar::StdRegions::eVarCoeffD02, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD12, Nektar::StdRegions::eVarCoeffD22, Nektar::GetCellModelFactory(), Nektar::SolverUtils::EquationSystem::GetFunction(), Nektar::Stimulus::LoadStimuli(), m_capMembrane, m_cell, m_chi, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_filters, Nektar::SolverUtils::UnsteadySystem::m_intVariables, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_stimulus, m_vardiff, Vmath::Sadd(), Vmath::Smul(), Nektar::SolverUtils::UnsteadySystem::v_InitObject(), Vmath::Vmul(), and Nektar::SolverUtils::EquationSystem::WriteFld().

◆ v_SetInitialConditions()

void Nektar::Monodomain::v_SetInitialConditions ( NekDouble  initialtime,
bool  dumpInitialConditions,
const int  domain 
)
overrideprotectedvirtual

Sets a custom initial condition.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 367 of file Monodomain.cpp.

370{
371 EquationSystem::v_SetInitialConditions(initialtime, dumpInitialConditions,
372 domain);
373 m_cell->Initialise();
374}
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)

References m_cell, and Nektar::SolverUtils::EquationSystem::v_SetInitialConditions().

Friends And Related Function Documentation

◆ MemoryManager< Monodomain >

friend class MemoryManager< Monodomain >
friend

Definition at line 1 of file Monodomain.h.

Member Data Documentation

◆ className

string Nektar::Monodomain::className
static
Initial value:
=
"Monodomain", Monodomain::create,
"Monodomain model of cardiac electrophysiology.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: Monodomain.h:54
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Registers the class with the Factory.

Definition at line 65 of file Monodomain.h.

◆ m_capMembrane

NekDouble Nektar::Monodomain::m_capMembrane
private

Definition at line 106 of file Monodomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_cell

CellModelSharedPtr Nektar::Monodomain::m_cell
private

Cell model.

Definition at line 98 of file Monodomain.h.

Referenced by DoOdeRhs(), v_GenerateSummary(), v_InitObject(), and v_SetInitialConditions().

◆ m_chi

NekDouble Nektar::Monodomain::m_chi
private

Definition at line 105 of file Monodomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_stimDuration

NekDouble Nektar::Monodomain::m_stimDuration
private

Stimulus current.

Definition at line 109 of file Monodomain.h.

◆ m_stimulus

std::vector<StimulusSharedPtr> Nektar::Monodomain::m_stimulus
private

Definition at line 100 of file Monodomain.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_vardiff

StdRegions::VarCoeffMap Nektar::Monodomain::m_vardiff
private

Variable diffusivity.

Definition at line 103 of file Monodomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().