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::BidomainRoth Class Reference

A model for cardiac conduction. More...

#include <BidomainRoth.h>

Inheritance diagram for Nektar::BidomainRoth:
[legend]

Public Member Functions

 ~BidomainRoth () 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...
 

Static Public Member Functions

static SolverUtils::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

 BidomainRoth (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_vardiffi
 
StdRegions::VarCoeffMap m_vardiffe
 
StdRegions::VarCoeffMap m_vardiffie
 
NekDouble m_chi
 
NekDouble m_capMembrane
 
NekDouble m_stimDuration
 Stimulus current. More...
 

Friends

class MemoryManager< BidomainRoth >
 

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_movingFrameVelsxyz
 Moving frame of reference velocities (u, v, w, omega_x, omega_y, omega_z, a_x, a_y, a_z, domega_x, domega_y, domega_z) More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

A model for cardiac conduction.

Definition at line 46 of file BidomainRoth.h.

Constructor & Destructor Documentation

◆ ~BidomainRoth()

Nektar::BidomainRoth::~BidomainRoth ( )
override

Desctructor.

Definition at line 312 of file BidomainRoth.cpp.

313{
314}

◆ BidomainRoth()

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

Constructor.

Definition at line 56 of file BidomainRoth.cpp.

58 : UnsteadySystem(pSession, pGraph)
59{
60}
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 52 of file BidomainRoth.h.

55 {
58 p->InitObject();
59 return p;
60 }
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::BidomainRoth::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 322 of file BidomainRoth.cpp.

326{
327 int nq = m_fields[0]->GetNpoints();
328
329 StdRegions::ConstFactorMap factorsHelmholtz;
330 // lambda = \Delta t
331 factorsHelmholtz[StdRegions::eFactorLambda] =
332 1.0 / lambda * m_chi * m_capMembrane;
333
334 // ------------------------------
335 // Solve Helmholtz problem for Vm
336 // ------------------------------
337 // Multiply 1.0/timestep
338 // Vmath::Vadd(nq, inarray[0], 1, ggrad, 1, m_fields[0]->UpdatePhys(), 1);
339 Vmath::Smul(nq, -factorsHelmholtz[StdRegions::eFactorLambda], inarray[0], 1,
340 m_fields[0]->UpdatePhys(), 1);
341
342 // Solve a system of equations with Helmholtz solver and transform
343 // back into physical space.
344 m_fields[0]->HelmSolve(m_fields[0]->GetPhys(), m_fields[0]->UpdateCoeffs(),
345 factorsHelmholtz, m_vardiffe);
346
347 m_fields[0]->BwdTrans(m_fields[0]->GetCoeffs(), m_fields[0]->UpdatePhys());
348 m_fields[0]->SetPhysState(true);
349
350 // Copy the solution vector (required as m_fields must be set).
351 outarray[0] = m_fields[0]->GetPhys();
352}
StdRegions::VarCoeffMap m_vardiffe
Definition: BidomainRoth.h:101
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
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, m_capMembrane, m_chi, Nektar::SolverUtils::EquationSystem::m_fields, m_vardiffe, and Vmath::Smul().

Referenced by v_InitObject().

◆ DoOdeRhs()

void Nektar::BidomainRoth::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 357 of file BidomainRoth.cpp.

360{
361 int nq = m_fields[0]->GetNpoints();
362
363 // Compute I_ion
364 m_cell->TimeIntegrate(inarray, outarray, time);
365
366 // Compute I_stim
367 for (unsigned int i = 0; i < m_stimulus.size(); ++i)
368 {
369 m_stimulus[i]->Update(outarray, time);
370 }
371
372 Array<OneD, NekDouble> ggrad0(nq), ggrad1(nq), ggrad2(nq), ggrad(nq);
373 StdRegions::ConstFactorMap factorsPoisson;
374 factorsPoisson[StdRegions::eFactorLambda] = 0.0;
375
376 // ----------------------------
377 // Compute \nabla g_i \nabla Vm
378 // ----------------------------
379 m_fields[0]->PhysDeriv(inarray[0], ggrad0, ggrad1, ggrad2);
380 m_fields[0]->PhysDeriv(0, ggrad0, ggrad0);
381 m_fields[0]->PhysDeriv(1, ggrad1, ggrad1);
382 m_fields[0]->PhysDeriv(2, ggrad2, ggrad2);
383 if (m_session->DefinesFunction("IntracellularAnisotropicConductivity") &&
384 m_session->DefinesFunction("ExtracellularAnisotropicConductivity"))
385 {
386 Vmath::Vmul(nq, &m_vardiffi[StdRegions::eVarCoeffD00][0], 1, &ggrad0[0],
387 1, &ggrad0[0], 1);
388 Vmath::Vmul(nq, &m_vardiffi[StdRegions::eVarCoeffD11][0], 1, &ggrad1[0],
389 1, &ggrad1[0], 1);
390 Vmath::Vmul(nq, &m_vardiffi[StdRegions::eVarCoeffD22][0], 1, &ggrad2[0],
391 1, &ggrad2[0], 1);
392 }
393 // Add partial derivatives together
394 Vmath::Vadd(nq, ggrad0, 1, ggrad1, 1, ggrad, 1);
395 Vmath::Vadd(nq, ggrad2, 1, ggrad, 1, ggrad, 1);
396
397 Vmath::Smul(nq, -1.0, ggrad, 1, m_fields[1]->UpdatePhys(), 1);
398
399 // ----------------------------
400 // Solve Poisson problem for Ve
401 // ----------------------------
402 m_fields[1]->HelmSolve(m_fields[1]->GetPhys(), m_fields[1]->UpdateCoeffs(),
403 factorsPoisson, m_vardiffie);
404 m_fields[1]->BwdTrans(m_fields[1]->GetCoeffs(), m_fields[1]->UpdatePhys());
405 m_fields[1]->SetPhysState(true);
406
407 // ------------------------------
408 // Compute Laplacian of Ve (forcing term)
409 // ------------------------------
410 m_fields[1]->PhysDeriv(m_fields[1]->GetPhys(), ggrad0, ggrad1, ggrad2);
411 m_fields[1]->PhysDeriv(0, ggrad0, ggrad0);
412 m_fields[1]->PhysDeriv(1, ggrad1, ggrad1);
413 m_fields[1]->PhysDeriv(2, ggrad2, ggrad2);
414 if (m_session->DefinesFunction("IntracellularAnisotropicConductivity") &&
415 m_session->DefinesFunction("ExtracellularAnisotropicConductivity"))
416 {
417 Vmath::Vmul(nq, &m_vardiffi[StdRegions::eVarCoeffD00][0], 1, &ggrad0[0],
418 1, &ggrad0[0], 1);
419 Vmath::Vmul(nq, &m_vardiffi[StdRegions::eVarCoeffD11][0], 1, &ggrad1[0],
420 1, &ggrad1[0], 1);
421 Vmath::Vmul(nq, &m_vardiffi[StdRegions::eVarCoeffD22][0], 1, &ggrad2[0],
422 1, &ggrad2[0], 1);
423 }
424 // Add partial derivatives together
425 Vmath::Vadd(nq, ggrad0, 1, ggrad1, 1, ggrad, 1);
426 Vmath::Vadd(nq, ggrad2, 1, ggrad, 1, ggrad, 1);
427
428 Vmath::Vadd(nq, ggrad, 1, outarray[0], 1, outarray[0], 1);
429}
StdRegions::VarCoeffMap m_vardiffi
Definition: BidomainRoth.h:100
std::vector< StimulusSharedPtr > m_stimulus
Definition: BidomainRoth.h:98
CellModelSharedPtr m_cell
Cell model.
Definition: BidomainRoth.h:96
StdRegions::VarCoeffMap m_vardiffie
Definition: BidomainRoth.h:102
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
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 Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, m_cell, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, m_stimulus, m_vardiffi, m_vardiffie, Vmath::Smul(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by v_InitObject().

◆ LoadStimuli()

void Nektar::BidomainRoth::LoadStimuli ( )
private

◆ v_GenerateSummary()

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

Prints a summary of the model parameters.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 446 of file BidomainRoth.cpp.

447{
449 m_cell->GenerateSummary(s);
450}
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.

References m_cell, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

◆ v_InitObject()

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

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 65 of file BidomainRoth.cpp.

66{
67 UnsteadySystem::v_InitObject(DeclareField);
68
69 m_session->LoadParameter("Chi", m_chi);
70 m_session->LoadParameter("Cm", m_capMembrane);
71
72 std::string vCellModel;
73 m_session->LoadSolverInfo("CELLMODEL", vCellModel, "");
74
75 ASSERTL0(vCellModel != "", "Cell Model not specified.");
76
78 m_fields[0]);
79
80 m_intVariables.push_back(0);
81
82 // Load variable coefficients
83 StdRegions::VarCoeffType varCoeffEnum[6] = {
87 std::string varCoeffString[6] = {"xx", "xy", "yy", "xz", "yz", "zz"};
88 std::string aniso_var[3] = {"fx", "fy", "fz"};
89
90 const int nq = m_fields[0]->GetNpoints();
91
92 // Allocate storage for variable coeffs and initialize to 1.
93 for (int i = 0, k = 0; i < m_spacedim; ++i)
94 {
95 for (int j = 0; j < i + 1; ++j)
96 {
97 if (i == j)
98 {
99 m_vardiffi[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 1.0);
100 m_vardiffe[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 1.0);
101 m_vardiffie[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 1.0);
102 }
103 else
104 {
105 m_vardiffi[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 0.0);
106 m_vardiffe[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 0.0);
107 m_vardiffie[varCoeffEnum[k]] = Array<OneD, NekDouble>(nq, 0.0);
108 }
109 ++k;
110 }
111 }
112
113 // Apply fibre map f \in [0,1], scale to conductivity range
114 // [o_min,o_max], specified by the session parameters o_min and o_max
115 if (m_session->DefinesFunction("ExtracellularAnisotropicConductivity"))
116 {
117 if (m_session->DefinesCmdLineArgument("verbose"))
118 {
119 cout << "Loading Extracellular Anisotropic Fibre map." << endl;
120 }
121
122 NekDouble o_min = m_session->GetParameter("o_min");
123 NekDouble o_max = m_session->GetParameter("o_max");
124 int k = 0;
125
126 Array<OneD, NekDouble> vTemp_i;
127 Array<OneD, NekDouble> vTemp_j;
128
129 /*
130 * Diffusivity matrix D is upper triangular and defined as
131 * d_00 d_01 d_02
132 * d_11 d_12
133 * d_22
134 *
135 * Given a principle fibre direction _f_ the diffusivity is given
136 * by
137 * d_ij = { D_2 + (D_1 - D_2) f_i f_j if i==j
138 * { (D_1 - D_2) f_i f_j if i!=j
139 *
140 * The vector _f_ is given in terms of the variables fx,fy,fz in the
141 * function AnisotropicConductivity. The values of D_1 and D_2 are
142 * the parameters o_max and o_min, respectively.
143 */
144
145 // Loop through columns of D
146 for (int j = 0; j < m_spacedim; ++j)
147 {
148 ASSERTL0(m_session->DefinesFunction(
149 "ExtracellularAnisotropicConductivity", aniso_var[j]),
150 "Function 'AnisotropicConductivity' not correctly "
151 "defined.");
152
153 GetFunction("ExtracellularAnisotropicConductivity")
154 ->Evaluate(aniso_var[j], vTemp_j);
155
156 // Loop through rows of D
157 for (int i = 0; i < j + 1; ++i)
158 {
159 ASSERTL0(
160 m_session->DefinesFunction(
161 "ExtracellularAnisotropicConductivity", aniso_var[i]),
162 "Function 'ExtracellularAnisotropicConductivity' not "
163 "correctly defined.");
164
165 GetFunction("ExtracellularAnisotropicConductivity")
166 ->Evaluate(aniso_var[i], vTemp_i);
167 Array<OneD, NekDouble> tmp =
168 m_vardiffe[varCoeffEnum[k]].GetValue();
169
170 Vmath::Vmul(nq, vTemp_i, 1, vTemp_j, 1, tmp, 1);
171
172 Vmath::Smul(nq, o_max - o_min, tmp, 1, tmp, 1);
173
174 if (i == j)
175 {
176 Vmath::Sadd(nq, o_min, tmp, 1, tmp, 1);
177 }
178
179 m_vardiffe[varCoeffEnum[k]] = tmp;
180 }
181 }
182 }
183
184 // Apply fibre map f \in [0,1], scale to conductivity range
185 // [o_min,o_max], specified by the session parameters o_min and o_max
186 if (m_session->DefinesFunction("IntracellularAnisotropicConductivity"))
187 {
188 if (m_session->DefinesCmdLineArgument("verbose"))
189 {
190 cout << "Loading Anisotropic Fibre map." << endl;
191 }
192
193 NekDouble o_min = m_session->GetParameter("o_min");
194 NekDouble o_max = m_session->GetParameter("o_max");
195 int k = 0;
196
197 Array<OneD, NekDouble> vTemp_i;
198 Array<OneD, NekDouble> vTemp_j;
199
200 /*
201 * Diffusivity matrix D is upper triangular and defined as
202 * d_00 d_01 d_02
203 * d_11 d_12
204 * d_22
205 *
206 * Given a principle fibre direction _f_ the diffusivity is given
207 * by
208 * d_ij = { D_2 + (D_1 - D_2) f_i f_j if i==j
209 * { (D_1 - D_2) f_i f_j if i!=j
210 *
211 * The vector _f_ is given in terms of the variables fx,fy,fz in the
212 * function AnisotropicConductivity. The values of D_1 and D_2 are
213 * the parameters o_max and o_min, respectively.
214 */
215
216 // Loop through columns of D
217 for (int j = 0; j < m_spacedim; ++j)
218 {
219 ASSERTL0(m_session->DefinesFunction(
220 "IntracellularAnisotropicConductivity", aniso_var[j]),
221 "Function 'IntracellularAnisotropicConductivity' not "
222 "correctly defined.");
223
224 GetFunction("IntracellularAnisotropicConductivity")
225 ->Evaluate(aniso_var[j], vTemp_j);
226
227 // Loop through rows of D
228 for (int i = 0; i < j + 1; ++i)
229 {
230 ASSERTL0(
231 m_session->DefinesFunction(
232 "IntracellularAnisotropicConductivity", aniso_var[i]),
233 "Function 'IntracellularAnisotropicConductivity' not "
234 "correctly defined.");
235 GetFunction("IntracellularAnisotropicConductivity")
236 ->Evaluate(aniso_var[i], vTemp_i);
237
238 Array<OneD, NekDouble> tmp =
239 m_vardiffi[varCoeffEnum[k]].GetValue();
240 Array<OneD, NekDouble> tmp2 =
241 m_vardiffe[varCoeffEnum[k]].GetValue();
242 Vmath::Vmul(nq, vTemp_i, 1, vTemp_j, 1, tmp, 1);
243
244 Vmath::Smul(nq, o_max - o_min, tmp, 1, tmp, 1);
245
246 if (i == j)
247 {
248 Vmath::Sadd(nq, o_min, tmp, 1, tmp, 1);
249 }
250
251 Vmath::Vadd(nq, tmp2, 1, tmp, 1, tmp2, 1);
252
253 m_vardiffi[varCoeffEnum[k]] = tmp;
254 m_vardiffe[varCoeffEnum[k]] = tmp2;
255
256 ++k;
257 }
258 }
259 }
260
261 // Write out conductivity values
262 for (int j = 0, k = 0; j < m_spacedim; ++j)
263 {
264 // Loop through rows of D
265 for (int i = 0; i < j + 1; ++i)
266 {
267 // Transform variable coefficient and write out to file.
268 m_fields[0]->FwdTransLocalElmt(
269 m_vardiffi[varCoeffEnum[k]].GetValue(),
270 m_fields[0]->UpdateCoeffs());
271 std::stringstream filenamei;
272 filenamei << "IConductivity_" << varCoeffString[k] << ".fld";
273 WriteFld(filenamei.str());
274
275 // Transform variable coefficient and write out to file.
276 m_fields[0]->FwdTransLocalElmt(
277 m_vardiffe[varCoeffEnum[k]].GetValue(),
278 m_fields[0]->UpdateCoeffs());
279 std::stringstream filenamee;
280 filenamee << "EConductivity_" << varCoeffString[k] << ".fld";
281 WriteFld(filenamee.str());
282
283 ++k;
284 }
285 }
286
287 // Search through the loaded filters and pass the cell model to any
288 // CheckpointCellModel filters loaded.
289 for (auto &x : m_filters)
290 {
291 if (x.first == "CheckpointCellModel")
292 {
293 std::shared_ptr<FilterCheckpointCellModel> c =
294 std::dynamic_pointer_cast<FilterCheckpointCellModel>(x.second);
295 c->SetCellModel(m_cell);
296 }
297 }
298 // Load stimuli
300
302 {
304 }
307}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Computes the reaction terms and .
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.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
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 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_vardiffe, m_vardiffi, m_vardiffie, Vmath::Sadd(), Vmath::Smul(), Nektar::SolverUtils::UnsteadySystem::v_InitObject(), Vmath::Vadd(), Vmath::Vmul(), and Nektar::SolverUtils::EquationSystem::WriteFld().

◆ v_SetInitialConditions()

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

Sets a custom initial condition.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 434 of file BidomainRoth.cpp.

437{
438 EquationSystem::v_SetInitialConditions(initialtime, dumpInitialConditions,
439 domain);
440 m_cell->Initialise();
441}
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< BidomainRoth >

friend class MemoryManager< BidomainRoth >
friend

Definition at line 1 of file BidomainRoth.h.

Member Data Documentation

◆ className

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

Name of class.

Registers the class with the Factory.

Definition at line 63 of file BidomainRoth.h.

◆ m_capMembrane

NekDouble Nektar::BidomainRoth::m_capMembrane
private

Definition at line 105 of file BidomainRoth.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_cell

CellModelSharedPtr Nektar::BidomainRoth::m_cell
private

Cell model.

Definition at line 96 of file BidomainRoth.h.

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

◆ m_chi

NekDouble Nektar::BidomainRoth::m_chi
private

Definition at line 104 of file BidomainRoth.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_stimDuration

NekDouble Nektar::BidomainRoth::m_stimDuration
private

Stimulus current.

Definition at line 108 of file BidomainRoth.h.

◆ m_stimulus

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

Definition at line 98 of file BidomainRoth.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_vardiffe

StdRegions::VarCoeffMap Nektar::BidomainRoth::m_vardiffe
private

Definition at line 101 of file BidomainRoth.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_vardiffi

StdRegions::VarCoeffMap Nektar::BidomainRoth::m_vardiffi
private

Definition at line 100 of file BidomainRoth.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_vardiffie

StdRegions::VarCoeffMap Nektar::BidomainRoth::m_vardiffie
private

Definition at line 102 of file BidomainRoth.h.

Referenced by DoOdeRhs(), and v_InitObject().