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 | 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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Monodomain:
Collaboration graph
[legend]

Public Member Functions

virtual ~Monodomain ()
 Desctructor. More...
 
- 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 >
boost::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 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 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. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. 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, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
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) $. More...
 
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. More...
 
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. More...
 
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. 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 ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to 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::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 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 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 SetSteps (const int steps)
 
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 std::string &s1, const std::string &s2)
 Perform a case-insensitive string comparison. More...
 
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)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 

Protected Member Functions

 Monodomain (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
virtual void v_InitObject ()
 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...
 
virtual void v_SetInitialConditions (NekDouble initialtime, bool dumpInitialConditions, const int domain)
 Sets a custom initial condition. More...
 
virtual void v_GenerateSummary (SummaryList &s)
 Prints a summary of the model parameters. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. 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 void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
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. 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_SteadyStateCheck (int step)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
 
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)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const std::string &s1, const std::string &s2)
 
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_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 void v_Output (void)
 
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)
 

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

- 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...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information. More...
 
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...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
std::map< std::string, Array
< OneD, Array< OneD, float > > > 
m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
std::map< std::string, Array
< OneD, Array< OneD, unsigned
int > > > 
m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 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...
 
std::set< std::string > m_loadedFields
 
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, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq 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...
 

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 50 of file Monodomain.h.

Constructor & Destructor Documentation

Nektar::Monodomain::~Monodomain ( )
virtual

Desctructor.

Definition at line 317 of file Monodomain.cpp.

318  {
319 
320  }
Nektar::Monodomain::Monodomain ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Constructor.

Definition at line 73 of file Monodomain.cpp.

75  : UnsteadySystem(pSession)
76  {
77  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.

Member Function Documentation

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

Creates an instance of this class.

Definition at line 56 of file Monodomain.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

58  {
60  p->InitObject();
61  return p;
62  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
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 329 of file Monodomain.cpp.

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

Referenced by v_InitObject().

334  {
335  int nvariables = inarray.num_elements();
336  int nq = m_fields[0]->GetNpoints();
338  // lambda = \Delta t
339  factors[StdRegions::eFactorLambda] = 1.0/lambda*m_chi*m_capMembrane;
340 
341  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
342  // inarray = input: \hat{rhs} -> output: \hat{Y}
343  // outarray = output: nabla^2 \hat{Y}
344  // where \hat = modal coeffs
345  for (int i = 0; i < nvariables; ++i)
346  {
347  // Multiply 1.0/timestep
348  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
349  m_fields[i]->UpdatePhys(), 1);
350 
351  // Solve a system of equations with Helmholtz solver and transform
352  // back into physical space.
353  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(),
354  m_fields[i]->UpdateCoeffs(), NullFlagList,
355  factors, m_vardiff);
356 
357  m_fields[i]->BwdTrans( m_fields[i]->GetCoeffs(),
358  m_fields[i]->UpdatePhys());
359  m_fields[i]->SetPhysState(true);
360 
361  // Copy the solution vector (required as m_fields must be set).
362  outarray[i] = m_fields[i]->GetPhys();
363  }
364  }
StdRegions::VarCoeffMap m_vardiff
Variable diffusivity.
Definition: Monodomain.h:105
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
NekDouble m_capMembrane
Definition: Monodomain.h:108
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
static FlagList NullFlagList
An empty flag list.
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 370 of file Monodomain.cpp.

References m_cell, and m_stimulus.

Referenced by v_InitObject().

374  {
375  // Compute I_ion
376  m_cell->TimeIntegrate(inarray, outarray, time);
377 
378  // Compute I_stim
379  for (unsigned int i = 0; i < m_stimulus.size(); ++i)
380  {
381  m_stimulus[i]->Update(outarray, time);
382  }
383  }
CellModelSharedPtr m_cell
Cell model.
Definition: Monodomain.h:100
std::vector< StimulusSharedPtr > m_stimulus
Definition: Monodomain.h:102
void Nektar::Monodomain::LoadStimuli ( )
private
void Nektar::Monodomain::v_GenerateSummary ( SummaryList s)
protectedvirtual

Prints a summary of the model parameters.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 403 of file Monodomain.cpp.

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

404  {
406  if (m_session->DefinesFunction("d00") &&
407  m_session->GetFunctionType("d00", "intensity")
409  {
410  AddSummaryItem(s, "Diffusivity-x",
411  m_session->GetFunction("d00", "intensity")->GetExpression());
412  }
413  if (m_session->DefinesFunction("d11") &&
414  m_session->GetFunctionType("d11", "intensity")
416  {
417  AddSummaryItem(s, "Diffusivity-y",
418  m_session->GetFunction("d11", "intensity")->GetExpression());
419  }
420  if (m_session->DefinesFunction("d22") &&
421  m_session->GetFunctionType("d22", "intensity")
423  {
424  AddSummaryItem(s, "Diffusivity-z",
425  m_session->GetFunction("d22", "intensity")->GetExpression());
426  }
427  m_cell->GenerateSummary(s);
428  }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
CellModelSharedPtr m_cell
Cell model.
Definition: Monodomain.h:100
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void Nektar::Monodomain::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 83 of file Monodomain.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), DoImplicitSolve(), DoOdeRhs(), Nektar::SolverUtils::EquationSystem::EvaluateFunction(), Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD01, Nektar::StdRegions::eVarCoeffD02, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD12, Nektar::StdRegions::eVarCoeffD22, Nektar::GetCellModelFactory(), 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(), Nektar::FilterCellHistoryPoints::SetCellModel(), Nektar::FilterCheckpointCellModel::SetCellModel(), Vmath::Smul(), Nektar::SolverUtils::UnsteadySystem::v_InitObject(), Vmath::Vmul(), and Nektar::SolverUtils::EquationSystem::WriteFld().

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

Sets a custom initial condition.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 389 of file Monodomain.cpp.

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

392  {
394  dumpInitialConditions,
395  domain);
396  m_cell->Initialise();
397  }
CellModelSharedPtr m_cell
Cell model.
Definition: Monodomain.h:100
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)

Friends And Related Function Documentation

friend class MemoryManager< Monodomain >
friend

Definition at line 53 of file Monodomain.h.

Member Data Documentation

string Nektar::Monodomain::className
static
Initial value:
"Monodomain",
"Monodomain model of cardiac electrophysiology.")

Name of class.

Registers the class with the Factory.

Definition at line 65 of file Monodomain.h.

NekDouble Nektar::Monodomain::m_capMembrane
private

Definition at line 108 of file Monodomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

CellModelSharedPtr Nektar::Monodomain::m_cell
private

Cell model.

Definition at line 100 of file Monodomain.h.

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

NekDouble Nektar::Monodomain::m_chi
private

Definition at line 107 of file Monodomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

NekDouble Nektar::Monodomain::m_stimDuration
private

Stimulus current.

Definition at line 111 of file Monodomain.h.

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

Definition at line 102 of file Monodomain.h.

Referenced by DoOdeRhs(), and v_InitObject().

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

Variable diffusivity.

Definition at line 105 of file Monodomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().