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

A model for cardiac conduction. More...

#include <Bidomain.h>

Inheritance diagram for Nektar::Bidomain:
[legend]

Public Member Functions

virtual ~Bidomain ()
 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 >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

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

Protected Member Functions

 Bidomain (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 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, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_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 NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT 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)
 
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 Attributes

CellModelSharedPtr m_cell
 Cell model. More...
 
NekDouble m_chi
 
NekDouble m_capMembrane
 
NekDouble m_sigmaix
 
NekDouble m_sigmaiy
 
NekDouble m_sigmaiz
 
NekDouble m_sigmaex
 
NekDouble m_sigmaey
 
NekDouble m_sigmaez
 
StdRegions::VarCoeffMap m_vardiffi
 
StdRegions::VarCoeffMap m_vardiffie
 
Array< OneD, Array< OneD, NekDouble > > tmp1
 
Array< OneD, Array< OneD, NekDouble > > tmp2
 
Array< OneD, Array< OneD, NekDouble > > tmp3
 
NekDouble m_stimDuration
 Stimulus current. More...
 

Friends

class MemoryManager< Bidomain >
 

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...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

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

Constructor & Destructor Documentation

◆ ~Bidomain()

Nektar::Bidomain::~Bidomain ( )
virtual

Desctructor.

Definition at line 176 of file Bidomain.cpp.

177  {
178 
179  }

◆ Bidomain()

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

Constructor.

Definition at line 71 of file Bidomain.cpp.

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

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 54 of file Bidomain.h.

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

57  {
59  ::AllocateSharedPtr(pSession, pGraph);
60  p->InitObject();
61  return p;
62  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

◆ DoImplicitSolve()

void Nektar::Bidomain::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 188 of file Bidomain.cpp.

References Nektar::StdRegions::eFactorLambda, m_capMembrane, m_chi, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_vardiffi, m_vardiffie, Nektar::NullFlagList, Vmath::Smul(), Vmath::Vadd(), and Vmath::Vcopy().

Referenced by v_InitObject().

193  {
194  int nvariables = inarray.num_elements();
195  int nq = m_fields[0]->GetNpoints();
196 
197  Array<OneD, NekDouble> grad0(nq), grad1(nq), grad2(nq), grad(nq);
198  Array<OneD, NekDouble> ggrad0(nq), ggrad1(nq), ggrad2(nq), ggrad(nq), temp(nq);
199 
200  // We solve ( \sigma\nabla^2 - HHlambda ) Y[i] = rhs [i]
201  // inarray = input: \hat{rhs} -> output: \hat{Y}
202  // outarray = output: nabla^2 \hat{Y}
203  // where \hat = modal coeffs
204  for (int i = 0; i < nvariables; ++i)
205  {
206  // Only apply diffusion to first variable.
207  if (i > 1) {
208  Vmath::Vcopy(nq, &inarray[i][0], 1, &outarray[i][0], 1);
209  continue;
210  }
211  if (i == 0) {
213  factors[StdRegions::eFactorLambda] = (1.0/lambda)*(m_capMembrane*m_chi);
214  if (m_spacedim==1) {
215  // Take first partial derivative
216  m_fields[i]->PhysDeriv(inarray[1],ggrad0);
217  // Take second partial derivative
218  m_fields[i]->PhysDeriv(0,ggrad0,ggrad0);
219  // Multiply by Intracellular-Conductivity
220  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
221  {
222  Vmath::Smul(nq, m_session->GetParameter("sigmaix"), ggrad0, 1, ggrad0, 1);
223  }
224  // Add partial derivatives together
225  Vmath::Vcopy(nq, ggrad0, 1, ggrad, 1);
226  Vmath::Smul(nq, -1.0, ggrad, 1, ggrad, 1);
227  // Multiply 1.0/timestep/lambda
228  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1, temp, 1);
229  Vmath::Vadd(nq, ggrad, 1, temp, 1, m_fields[i]->UpdatePhys(), 1);
230  // Solve a system of equations with Helmholtz solver and transform
231  // back into physical space.
232  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(), m_fields[i]->UpdateCoeffs(),NullFlagList,factors);
233  m_fields[i]->BwdTrans( m_fields[i]->GetCoeffs(), m_fields[i]->UpdatePhys());
234  m_fields[i]->SetPhysState(true);
235  // Copy the solution vector (required as m_fields must be set).
236  outarray[i] = m_fields[i]->GetPhys();
237  }
238 
239  if (m_spacedim==2) {
240  // Take first partial derivative
241  m_fields[i]->PhysDeriv(inarray[1],ggrad0,ggrad1);
242  // Take second partial derivative
243  m_fields[i]->PhysDeriv(0,ggrad0,ggrad0);
244  m_fields[i]->PhysDeriv(1,ggrad1,ggrad1);
245  // Multiply by Intracellular-Conductivity
246  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
247  {
248  Vmath::Smul(nq, m_session->GetParameter("sigmaix"), ggrad0, 1, ggrad0, 1);
249  Vmath::Smul(nq, m_session->GetParameter("sigmaiy"), ggrad1, 1, ggrad1, 1);
250  }
251  // Add partial derivatives together
252  Vmath::Vadd(nq, ggrad0, 1, ggrad1, 1, ggrad, 1);
253  Vmath::Smul(nq, -1.0, ggrad, 1, ggrad, 1);
254  // Multiply 1.0/timestep/lambda
255  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1, temp, 1);
256  Vmath::Vadd(nq, ggrad, 1, temp, 1, m_fields[i]->UpdatePhys(), 1);
257  // Solve a system of equations with Helmholtz solver and transform
258  // back into physical space.
259  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(), m_fields[i]->UpdateCoeffs(),NullFlagList,factors,m_vardiffi);
260  m_fields[i]->BwdTrans( m_fields[i]->GetCoeffs(), m_fields[i]->UpdatePhys());
261  m_fields[i]->SetPhysState(true);
262  // Copy the solution vector (required as m_fields must be set).
263  outarray[i] = m_fields[i]->GetPhys();
264  }
265 
266  if (m_spacedim==3) {
267  // Take first partial derivative
268  m_fields[i]->PhysDeriv(inarray[1],ggrad0,ggrad1,ggrad2);
269  // Take second partial derivative
270  m_fields[i]->PhysDeriv(0,ggrad0,ggrad0);
271  m_fields[i]->PhysDeriv(1,ggrad1,ggrad1);
272  m_fields[i]->PhysDeriv(2,ggrad2,ggrad2);
273  // Multiply by Intracellular-Conductivity
274  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
275  {
276  Vmath::Smul(nq, m_session->GetParameter("sigmaix"), ggrad0, 1, ggrad0, 1);
277  Vmath::Smul(nq, m_session->GetParameter("sigmaiy"), ggrad1, 1, ggrad1, 1);
278  Vmath::Smul(nq, m_session->GetParameter("sigmaiz"), ggrad2, 1, ggrad2, 1);
279  }
280  // Add partial derivatives together
281  Vmath::Vadd(nq, ggrad0, 1, ggrad1, 1, ggrad, 1);
282  Vmath::Vadd(nq, ggrad2, 1, ggrad, 1, ggrad, 1);
283  Vmath::Smul(nq, -1.0, ggrad, 1, ggrad, 1);
284  // Multiply 1.0/timestep/lambda
285  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1, temp, 1);
286  Vmath::Vadd(nq, ggrad, 1, temp, 1, m_fields[i]->UpdatePhys(), 1);
287  // Solve a system of equations with Helmholtz solver and transform
288  // back into physical space.
289  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(), m_fields[i]->UpdateCoeffs(),NullFlagList,factors,m_vardiffi);
290  m_fields[i]->BwdTrans( m_fields[i]->GetCoeffs(), m_fields[i]->UpdatePhys());
291  m_fields[i]->SetPhysState(true);
292  // Copy the solution vector (required as m_fields must be set).
293  outarray[i] = m_fields[i]->GetPhys();
294  }
295 
296  }
297  if (i == 1) {
299  factors[StdRegions::eFactorLambda] = 0.0;
300  if (m_spacedim==1) {
301  // Take first partial derivative
302  m_fields[i]->PhysDeriv(m_fields[0]->UpdatePhys(),grad0);
303  // Take second derivative
304  m_fields[i]->PhysDeriv(0,grad0,grad0);
305  // Multiply by Intracellular-Conductivity
306  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
307  {
308  Vmath::Smul(nq, m_session->GetParameter("sigmaix"), grad0, 1, grad0, 1);
309  }
310  // and sum terms
311  Vmath::Vcopy(nq, grad0, 1, grad, 1);
312  Vmath::Smul(nq, (-1.0*m_session->GetParameter("sigmaix"))/(m_session->GetParameter("sigmaix")+m_session->GetParameter("sigmaix")), grad, 1, grad, 1);
313  // Now solve Poisson problem for \phi_e
314  m_fields[i]->SetPhys(grad);
315  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(), m_fields[i]->UpdateCoeffs(), NullFlagList, factors);
316  m_fields[i]->BwdTrans( m_fields[i]->GetCoeffs(), m_fields[i]->UpdatePhys());
317  m_fields[i]->SetPhysState(true);
318  // Copy the solution vector (required as m_fields must be set).
319  outarray[i] = m_fields[i]->GetPhys();
320  }
321 
322  if (m_spacedim==2) {
323  // Take first partial derivative
324  m_fields[i]->PhysDeriv(m_fields[0]->UpdatePhys(),grad0,grad1);
325  // Take second derivative
326  m_fields[i]->PhysDeriv(0,grad0,grad0);
327  m_fields[i]->PhysDeriv(1,grad1,grad1);
328  // Multiply by Intracellular-Conductivity
329  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
330  {
331  Vmath::Smul(nq, m_session->GetParameter("sigmaix"), grad0, 1, grad0, 1);
332  Vmath::Smul(nq, m_session->GetParameter("sigmaiy"), grad1, 1, grad1, 1);
333  }
334  // and sum terms
335  Vmath::Vadd(nq, grad0, 1, grad1, 1, grad, 1);
336  Vmath::Smul(nq, -1.0, grad, 1, grad, 1);
337  // Now solve Poisson problem for \phi_e
338  m_fields[i]->SetPhys(grad);
339  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(), m_fields[i]->UpdateCoeffs(), NullFlagList, factors, m_vardiffie);
340  m_fields[i]->BwdTrans( m_fields[i]->GetCoeffs(), m_fields[i]->UpdatePhys());
341  m_fields[i]->SetPhysState(true);
342  // Copy the solution vector (required as m_fields must be set).
343  outarray[i] = m_fields[i]->GetPhys();
344  }
345 
346  if (m_spacedim==3) {
347  // Take first partial derivative
348  m_fields[i]->PhysDeriv(m_fields[0]->UpdatePhys(),grad0,grad1,grad2);
349  // Take second derivative
350  m_fields[i]->PhysDeriv(0,grad0,grad0);
351  m_fields[i]->PhysDeriv(1,grad1,grad1);
352  m_fields[i]->PhysDeriv(2,grad2,grad2);
353  // Multiply by Intracellular-Conductivity
354  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
355  {
356  Vmath::Smul(nq, m_session->GetParameter("sigmaix"), grad0, 1, grad0, 1);
357  Vmath::Smul(nq, m_session->GetParameter("sigmaiy"), grad1, 1, grad1, 1);
358  Vmath::Smul(nq, m_session->GetParameter("sigmaiz"), grad2, 1, grad2, 1);
359  }
360  // and sum terms
361  Vmath::Vadd(nq, grad0, 1, grad1, 1, grad, 1);
362  Vmath::Vadd(nq, grad2, 1, grad, 1, grad, 1);
363  Vmath::Smul(nq, -1.0, grad, 1, grad, 1);
364  // Now solve Poisson problem for \phi_e
365  m_fields[i]->SetPhys(grad);
366  m_fields[i]->HelmSolve(m_fields[i]->GetPhys(), m_fields[i]->UpdateCoeffs(), NullFlagList, factors, m_vardiffie);
367  m_fields[i]->BwdTrans( m_fields[i]->GetCoeffs(), m_fields[i]->UpdatePhys());
368  m_fields[i]->SetPhysState(true);
369  // Copy the solution vector (required as m_fields must be set).
370  outarray[i] = m_fields[i]->GetPhys();
371  }
372 
373  }
374  }
375  }
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
NekDouble m_chi
Definition: Bidomain.h:103
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:216
StdRegions::VarCoeffMap m_vardiffie
Definition: Bidomain.h:106
int m_spacedim
Spatial dimension (>= expansion dim).
StdRegions::VarCoeffMap m_vardiffi
Definition: Bidomain.h:105
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
NekDouble m_capMembrane
Definition: Bidomain.h:103
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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.cpp:302
static FlagList NullFlagList
An empty flag list.

◆ DoOdeRhs()

void Nektar::Bidomain::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 378 of file Bidomain.cpp.

References m_capMembrane, m_cell, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, m_stimDuration, Vmath::Smul(), and Vmath::Vadd().

Referenced by v_InitObject().

382  {
383  int nq = m_fields[0]->GetNpoints();
384  m_cell->TimeIntegrate(inarray, outarray, time);
385  if (m_stimDuration > 0 && time < m_stimDuration)
386  {
387  Array<OneD,NekDouble> x0(nq);
388  Array<OneD,NekDouble> x1(nq);
389  Array<OneD,NekDouble> x2(nq);
390  Array<OneD,NekDouble> result(nq);
391 
392  // get the coordinates
393  m_fields[0]->GetCoords(x0,x1,x2);
394 
396  = m_session->GetFunction("Stimulus", "u");
397  ifunc->Evaluate(x0,x1,x2,time, result);
398 
399  Vmath::Vadd(nq, outarray[0], 1, result, 1, outarray[0], 1);
400  }
401  Vmath::Smul(nq, 1.0/m_capMembrane, outarray[0], 1, outarray[0], 1);
402  }
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:216
NekDouble m_stimDuration
Stimulus current.
Definition: Bidomain.h:113
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:131
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:101
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
NekDouble m_capMembrane
Definition: Bidomain.h:103
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.cpp:302

◆ v_GenerateSummary()

void Nektar::Bidomain::v_GenerateSummary ( SummaryList s)
protectedvirtual

Prints a summary of the model parameters.

Update summary

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 416 of file Bidomain.cpp.

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

417  {
419 
420  /// @TODO Update summary
421  ASSERTL0(false, "Update the generate summary");
422 //
423 // out << "\tChi : " << m_chi << endl;
424 // out << "\tCm : " << m_capMembrane << endl;
425 // if (m_session->DefinesFunction("IntracellularConductivity", "AnisotropicConductivityX") &&
426 // m_session->GetFunctionType("IntracellularConductivity", "AnisotropicConductivityX") == LibUtilities::eFunctionTypeExpression)
427 // {
428 // out << "\tIntra-Diffusivity-x : "
429 // << m_session->GetFunction("IntracellularConductivity", "AnisotropicConductivityX")->GetExpression()
430 // << endl;
431 // }
432 // if (m_session->DefinesFunction("IntracellularConductivity", "AnisotropicConductivityY") &&
433 // m_session->GetFunctionType("IntracellularConductivity", "AnisotropicConductivityY") == LibUtilities::eFunctionTypeExpression)
434 // {
435 // out << "\tIntra-Diffusivity-y : "
436 // << m_session->GetFunction("IntracellularConductivity", "AnisotropicConductivityY")->GetExpression()
437 // << endl;
438 // }
439 // if (m_session->DefinesFunction("IntracellularConductivity", "AnisotropicConductivityZ") &&
440 // m_session->GetFunctionType("IntracellularConductivity", "AnisotropicConductivityZ") == LibUtilities::eFunctionTypeExpression)
441 // {
442 // out << "\tIntra-Diffusivity-z : "
443 // << m_session->GetFunction("IntracellularConductivity", "AnisotropicConductivityZ")->GetExpression()
444 // << endl;
445 // }
446 // if (m_session->DefinesFunction("ExtracellularConductivity", "AnisotropicConductivityX") &&
447 // m_session->GetFunctionType("ExtracellularConductivity", "AnisotropicConductivityX") == LibUtilities::eFunctionTypeExpression)
448 // {
449 // out << "\tExtra-Diffusivity-x : "
450 // << m_session->GetFunction("ExtracellularConductivity", "AnisotropicConductivityX")->GetExpression()
451 // << endl;
452 // }
453 // if (m_session->DefinesFunction("ExtracellularConductivity", "AnisotropicConductivityY") &&
454 // m_session->GetFunctionType("ExtracellularConductivity", "AnisotropicConductivityY") == LibUtilities::eFunctionTypeExpression)
455 // {
456 // out << "\tExtra-Diffusivity-y : "
457 // << m_session->GetFunction("ExtracellularConductivity", "AnisotropicConductivityY")->GetExpression()
458 // << endl;
459 // }
460 // if (m_session->DefinesFunction("ExtracellularConductivity", "AnisotropicConductivityZ") &&
461 // m_session->GetFunctionType("ExtracellularConductivity", "AnisotropicConductivityZ") == LibUtilities::eFunctionTypeExpression)
462 // {
463 // out << "\tExtra-Diffusivity-z : "
464 // << m_session->GetFunction("ExtracellularConductivity", "AnisotropicConductivityZ")->GetExpression()
465 // << endl;
466 // }
467  m_cell->GenerateSummary(s);
468  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:101

◆ v_InitObject()

void Nektar::Bidomain::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 78 of file Bidomain.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), DoImplicitSolve(), DoOdeRhs(), Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::GetCellModelFactory(), 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_stimDuration, m_vardiffi, m_vardiffie, Nektar::FilterCheckpointCellModel::SetCellModel(), tmp1, tmp2, tmp3, Nektar::SolverUtils::UnsteadySystem::v_InitObject(), and Vmath::Vadd().

79  {
81  m_session->LoadParameter("Chi", m_chi);
82  m_session->LoadParameter("Cm", m_capMembrane);
83 
84  std::string vCellModel;
85  m_session->LoadSolverInfo("CELLMODEL", vCellModel, "");
86 
87  ASSERTL0(vCellModel != "", "Cell Model not specified.");
88 
90  m_intVariables.push_back(0);
91  m_intVariables.push_back(1);
92 
93  // Load variable coefficients
94  StdRegions::VarCoeffType varCoeffEnum[3] = {
98  };
99  std::string varName[3] = {
100  "AnisotropicConductivityX",
101  "AnisotropicConductivityY",
102  "AnisotropicConductivityZ"
103  };
104 
105 
106  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
107  {
108  for (int i = 0; i < m_spacedim; ++i)
109  {
110  int nq = m_fields[0]->GetNpoints();
111  Array<OneD,NekDouble> x0(nq);
112  Array<OneD,NekDouble> x1(nq);
113  Array<OneD,NekDouble> x2(nq);
114 
115  // get the coordinates
116  m_fields[0]->GetCoords(x0,x1,x2);
117  tmp1 = Array<OneD, const Array<OneD, NekDouble> >(nq);
118  tmp2 = Array<OneD, const Array<OneD, NekDouble> >(nq);
119  tmp3 = Array<OneD, const Array<OneD, NekDouble> >(nq);
120  tmp1[i] = Array<OneD, NekDouble>(nq);
121  tmp2[i] = Array<OneD, NekDouble>(nq);
122  tmp3[i] = Array<OneD, NekDouble>(nq);
123 
125  = m_session->GetFunction("IntracellularConductivity", varName[i]);
127  = m_session->GetFunction("ExtracellularConductivity", varName[i]);
128  for(int j = 0; j < nq; j++)
129  {
130  tmp1[i][j] = ifunc1->Evaluate(x0[j],x1[j],x2[j],0.0);
131  tmp2[i][j] = ifunc2->Evaluate(x0[j],x1[j],x2[j],0.0);
132  }
133  Vmath::Vadd(nq, tmp1[i], 1, tmp2[i], 1, tmp3[i], 1);
134  m_vardiffi[varCoeffEnum[i]] = tmp1[i];
135  m_vardiffie[varCoeffEnum[i]] = tmp3[i];
136  }
137  }
138 
139 
140  if (m_session->DefinesParameter("StimulusDuration"))
141  {
142  ASSERTL0(m_session->DefinesFunction("Stimulus", "u"),
143  "Stimulus function not defined.");
144  m_session->LoadParameter("StimulusDuration", m_stimDuration);
145  }
146  else
147  {
148  m_stimDuration = 0;
149  }
150 
151 
152  // Search through the loaded filters and pass the cell model to any
153  // CheckpointCellModel filters loaded.
154  for (auto &x : m_filters)
155  {
156  if (x.first == "CheckpointCellModel")
157  {
158  std::shared_ptr<FilterCheckpointCellModel> c
159  = std::dynamic_pointer_cast<FilterCheckpointCellModel>(
160  x.second);
161  c->SetCellModel(m_cell);
162  }
163  }
164 
165  if (!m_explicitDiffusion)
166  {
168  }
170  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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.
std::vector< std::pair< std::string, 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: Bidomain.cpp:378
NekDouble m_chi
Definition: Bidomain.h:103
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
StdRegions::VarCoeffMap m_vardiffie
Definition: Bidomain.h:106
NekDouble m_stimDuration
Stimulus current.
Definition: Bidomain.h:113
int m_spacedim
Spatial dimension (>= expansion dim).
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
StdRegions::VarCoeffMap m_vardiffi
Definition: Bidomain.h:105
CellModelFactory & GetCellModelFactory()
Definition: CellModel.cpp:46
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:131
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:101
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
NekDouble m_capMembrane
Definition: Bidomain.h:103
Array< OneD, Array< OneD, NekDouble > > tmp1
Definition: Bidomain.h:108
Array< OneD, Array< OneD, NekDouble > > tmp2
Definition: Bidomain.h:109
Array< OneD, Array< OneD, NekDouble > > tmp3
Definition: Bidomain.h:110
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: Bidomain.cpp:188
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.cpp:302

◆ v_SetInitialConditions()

void Nektar::Bidomain::v_SetInitialConditions ( NekDouble  initialtime,
bool  dumpInitialConditions,
const int  domain 
)
protectedvirtual

Sets a custom initial condition.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 405 of file Bidomain.cpp.

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

408  {
409  EquationSystem::v_SetInitialConditions(initialtime, dumpInitialConditions, domain);
410  m_cell->Initialise();
411  }
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:101

Friends And Related Function Documentation

◆ MemoryManager< Bidomain >

friend class MemoryManager< Bidomain >
friend

Definition at line 51 of file Bidomain.h.

Member Data Documentation

◆ className

string Nektar::Bidomain::className
static
Initial value:
"Bidomain",
"Bidomain model of cardiac electrophysiology with 3D diffusion.")

Name of class.

Registers the class with the Factory.

Definition at line 65 of file Bidomain.h.

◆ m_capMembrane

NekDouble Nektar::Bidomain::m_capMembrane
private

Definition at line 103 of file Bidomain.h.

Referenced by DoImplicitSolve(), DoOdeRhs(), and v_InitObject().

◆ m_cell

CellModelSharedPtr Nektar::Bidomain::m_cell
private

Cell model.

Definition at line 101 of file Bidomain.h.

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

◆ m_chi

NekDouble Nektar::Bidomain::m_chi
private

Definition at line 103 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_sigmaex

NekDouble Nektar::Bidomain::m_sigmaex
private

Definition at line 103 of file Bidomain.h.

◆ m_sigmaey

NekDouble Nektar::Bidomain::m_sigmaey
private

Definition at line 103 of file Bidomain.h.

◆ m_sigmaez

NekDouble Nektar::Bidomain::m_sigmaez
private

Definition at line 103 of file Bidomain.h.

◆ m_sigmaix

NekDouble Nektar::Bidomain::m_sigmaix
private

Definition at line 103 of file Bidomain.h.

◆ m_sigmaiy

NekDouble Nektar::Bidomain::m_sigmaiy
private

Definition at line 103 of file Bidomain.h.

◆ m_sigmaiz

NekDouble Nektar::Bidomain::m_sigmaiz
private

Definition at line 103 of file Bidomain.h.

◆ m_stimDuration

NekDouble Nektar::Bidomain::m_stimDuration
private

Stimulus current.

Definition at line 113 of file Bidomain.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_vardiffi

StdRegions::VarCoeffMap Nektar::Bidomain::m_vardiffi
private

Definition at line 105 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_vardiffie

StdRegions::VarCoeffMap Nektar::Bidomain::m_vardiffie
private

Definition at line 106 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ tmp1

Array<OneD, Array<OneD, NekDouble> > Nektar::Bidomain::tmp1
private

Definition at line 108 of file Bidomain.h.

Referenced by v_InitObject().

◆ tmp2

Array<OneD, Array<OneD, NekDouble> > Nektar::Bidomain::tmp2
private

Definition at line 109 of file Bidomain.h.

Referenced by v_InitObject().

◆ tmp3

Array<OneD, Array<OneD, NekDouble> > Nektar::Bidomain::tmp3
private

Definition at line 110 of file Bidomain.h.

Referenced by v_InitObject().