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...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- 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 (bool DeclareField=true)
 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, const Array< OneD, const NekDouble > &input)
 
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 int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool 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...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. 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...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

 Bidomain (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
virtual 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...
 
virtual void v_SetInitialConditions (NekDouble initialtime, bool dumpInitialConditions, const int domain) override
 Sets a custom initial condition. More...
 
virtual 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 NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise () override
 Sets up initial conditions. 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)
 
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...
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
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 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...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth 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_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::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
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...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
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...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity 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_timestepMax = -1.0
 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_pararealIter
 Number of parareal 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_useInitialCondition
 Flag to determine if IC are used. 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. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 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 []
 

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

Constructor & Destructor Documentation

◆ ~Bidomain()

Nektar::Bidomain::~Bidomain ( )
virtual

Desctructor.

Definition at line 166 of file Bidomain.cpp.

167 {
168 }

◆ Bidomain()

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

Constructor.

Definition at line 68 of file Bidomain.cpp.

70  : UnsteadySystem(pSession, pGraph)
71 {
72 }
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 53 of file Bidomain.h.

56  {
59  p->InitObject();
60  return p;
61  }
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::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 176 of file Bidomain.cpp.

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

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, Vmath::Smul(), Vmath::Vadd(), and Vmath::Vcopy().

Referenced by v_InitObject().

◆ 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 420 of file Bidomain.cpp.

423 {
424  int nq = m_fields[0]->GetNpoints();
425  m_cell->TimeIntegrate(inarray, outarray, time);
426  if (m_stimDuration > 0 && time < m_stimDuration)
427  {
428  Array<OneD, NekDouble> x0(nq);
429  Array<OneD, NekDouble> x1(nq);
430  Array<OneD, NekDouble> x2(nq);
431  Array<OneD, NekDouble> result(nq);
432 
433  // get the coordinates
434  m_fields[0]->GetCoords(x0, x1, x2);
435 
437  m_session->GetFunction("Stimulus", "u");
438  ifunc->Evaluate(x0, x1, x2, time, result);
439 
440  Vmath::Vadd(nq, outarray[0], 1, result, 1, outarray[0], 1);
441  }
442  Vmath::Smul(nq, 1.0 / m_capMembrane, outarray[0], 1, outarray[0], 1);
443 }
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:97
NekDouble m_stimDuration
Stimulus current.
Definition: Bidomain.h:110
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129

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

◆ v_GenerateSummary()

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

Prints a summary of the model parameters.

@TODO Update summary

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 457 of file Bidomain.cpp.

458 {
460 
461  /// @TODO Update summary
462  ASSERTL0(false, "Update the generate summary");
463  //
464  // out << "\tChi : " << m_chi << endl;
465  // out << "\tCm : " << m_capMembrane << endl;
466  // if (m_session->DefinesFunction("IntracellularConductivity",
467  // "AnisotropicConductivityX") &&
468  // m_session->GetFunctionType("IntracellularConductivity",
469  // "AnisotropicConductivityX") ==
470  // LibUtilities::eFunctionTypeExpression)
471  // {
472  // out << "\tIntra-Diffusivity-x : "
473  // << m_session->GetFunction("IntracellularConductivity",
474  // "AnisotropicConductivityX")->GetExpression()
475  // << endl;
476  // }
477  // if (m_session->DefinesFunction("IntracellularConductivity",
478  // "AnisotropicConductivityY") &&
479  // m_session->GetFunctionType("IntracellularConductivity",
480  // "AnisotropicConductivityY") ==
481  // LibUtilities::eFunctionTypeExpression)
482  // {
483  // out << "\tIntra-Diffusivity-y : "
484  // << m_session->GetFunction("IntracellularConductivity",
485  // "AnisotropicConductivityY")->GetExpression()
486  // << endl;
487  // }
488  // if (m_session->DefinesFunction("IntracellularConductivity",
489  // "AnisotropicConductivityZ") &&
490  // m_session->GetFunctionType("IntracellularConductivity",
491  // "AnisotropicConductivityZ") ==
492  // LibUtilities::eFunctionTypeExpression)
493  // {
494  // out << "\tIntra-Diffusivity-z : "
495  // << m_session->GetFunction("IntracellularConductivity",
496  // "AnisotropicConductivityZ")->GetExpression()
497  // << endl;
498  // }
499  // if (m_session->DefinesFunction("ExtracellularConductivity",
500  // "AnisotropicConductivityX") &&
501  // m_session->GetFunctionType("ExtracellularConductivity",
502  // "AnisotropicConductivityX") ==
503  // LibUtilities::eFunctionTypeExpression)
504  // {
505  // out << "\tExtra-Diffusivity-x : "
506  // << m_session->GetFunction("ExtracellularConductivity",
507  // "AnisotropicConductivityX")->GetExpression()
508  // << endl;
509  // }
510  // if (m_session->DefinesFunction("ExtracellularConductivity",
511  // "AnisotropicConductivityY") &&
512  // m_session->GetFunctionType("ExtracellularConductivity",
513  // "AnisotropicConductivityY") ==
514  // LibUtilities::eFunctionTypeExpression)
515  // {
516  // out << "\tExtra-Diffusivity-y : "
517  // << m_session->GetFunction("ExtracellularConductivity",
518  // "AnisotropicConductivityY")->GetExpression()
519  // << endl;
520  // }
521  // if (m_session->DefinesFunction("ExtracellularConductivity",
522  // "AnisotropicConductivityZ") &&
523  // m_session->GetFunctionType("ExtracellularConductivity",
524  // "AnisotropicConductivityZ") ==
525  // LibUtilities::eFunctionTypeExpression)
526  // {
527  // out << "\tExtra-Diffusivity-z : "
528  // << m_session->GetFunction("ExtracellularConductivity",
529  // "AnisotropicConductivityZ")->GetExpression()
530  // << endl;
531  // }
532  m_cell->GenerateSummary(s);
533 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.

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

◆ v_InitObject()

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

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 74 of file Bidomain.cpp.

75 {
76  UnsteadySystem::v_InitObject(DeclareField);
77  m_session->LoadParameter("Chi", m_chi);
78  m_session->LoadParameter("Cm", m_capMembrane);
79 
80  std::string vCellModel;
81  m_session->LoadSolverInfo("CELLMODEL", vCellModel, "");
82 
83  ASSERTL0(vCellModel != "", "Cell Model not specified.");
84 
86  m_fields[0]);
87  m_intVariables.push_back(0);
88  m_intVariables.push_back(1);
89 
90  // Load variable coefficients
94  std::string varName[3] = {"AnisotropicConductivityX",
95  "AnisotropicConductivityY",
96  "AnisotropicConductivityZ"};
97 
98  if (m_session->DefinesFunction("IntracellularConductivity") &&
99  m_session->DefinesFunction("ExtracellularConductivity"))
100  {
101  for (int i = 0; i < m_spacedim; ++i)
102  {
103  int nq = m_fields[0]->GetNpoints();
104  Array<OneD, NekDouble> x0(nq);
105  Array<OneD, NekDouble> x1(nq);
106  Array<OneD, NekDouble> x2(nq);
107 
108  // get the coordinates
109  m_fields[0]->GetCoords(x0, x1, x2);
110  tmp1 = Array<OneD, const Array<OneD, NekDouble>>(nq);
111  tmp2 = Array<OneD, const Array<OneD, NekDouble>>(nq);
112  tmp3 = Array<OneD, const Array<OneD, NekDouble>>(nq);
113  tmp1[i] = Array<OneD, NekDouble>(nq);
114  tmp2[i] = Array<OneD, NekDouble>(nq);
115  tmp3[i] = Array<OneD, NekDouble>(nq);
116 
118  m_session->GetFunction("IntracellularConductivity", varName[i]);
120  m_session->GetFunction("ExtracellularConductivity", varName[i]);
121  for (int j = 0; j < nq; j++)
122  {
123  tmp1[i][j] = ifunc1->Evaluate(x0[j], x1[j], x2[j], 0.0);
124  tmp2[i][j] = ifunc2->Evaluate(x0[j], x1[j], x2[j], 0.0);
125  }
126  Vmath::Vadd(nq, tmp1[i], 1, tmp2[i], 1, tmp3[i], 1);
127  m_vardiffi[varCoeffEnum[i]] = tmp1[i];
128  m_vardiffie[varCoeffEnum[i]] = tmp3[i];
129  }
130  }
131 
132  if (m_session->DefinesParameter("StimulusDuration"))
133  {
134  ASSERTL0(m_session->DefinesFunction("Stimulus", "u"),
135  "Stimulus function not defined.");
136  m_session->LoadParameter("StimulusDuration", m_stimDuration);
137  }
138  else
139  {
140  m_stimDuration = 0;
141  }
142 
143  // Search through the loaded filters and pass the cell model to any
144  // CheckpointCellModel filters loaded.
145  for (auto &x : m_filters)
146  {
147  if (x.first == "CheckpointCellModel")
148  {
149  std::shared_ptr<FilterCheckpointCellModel> c =
150  std::dynamic_pointer_cast<FilterCheckpointCellModel>(x.second);
151  c->SetCellModel(m_cell);
152  }
153  }
154 
155  if (!m_explicitDiffusion)
156  {
158  }
161 }
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:176
Array< OneD, Array< OneD, NekDouble > > tmp3
Definition: Bidomain.h:107
Array< OneD, Array< OneD, NekDouble > > tmp2
Definition: Bidomain.h:106
Array< OneD, Array< OneD, NekDouble > > tmp1
Definition: Bidomain.h:105
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:420
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
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
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
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.
CellModelFactory & GetCellModelFactory()
Definition: CellModel.cpp:46

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::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, tmp1, tmp2, tmp3, Nektar::SolverUtils::UnsteadySystem::v_InitObject(), and Vmath::Vadd().

◆ v_SetInitialConditions()

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

Sets a custom initial condition.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 445 of file Bidomain.cpp.

448 {
449  EquationSystem::v_SetInitialConditions(initialtime, dumpInitialConditions,
450  domain);
451  m_cell->Initialise();
452 }
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< Bidomain >

friend class MemoryManager< Bidomain >
friend

Definition at line 1 of file Bidomain.h.

Member Data Documentation

◆ className

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

Name of class.

Registers the class with the Factory.

Definition at line 64 of file Bidomain.h.

◆ m_capMembrane

NekDouble Nektar::Bidomain::m_capMembrane
private

Definition at line 99 of file Bidomain.h.

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

◆ m_cell

CellModelSharedPtr Nektar::Bidomain::m_cell
private

Cell model.

Definition at line 97 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 99 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_sigmaex

NekDouble Nektar::Bidomain::m_sigmaex
private

Definition at line 99 of file Bidomain.h.

◆ m_sigmaey

NekDouble Nektar::Bidomain::m_sigmaey
private

Definition at line 100 of file Bidomain.h.

◆ m_sigmaez

NekDouble Nektar::Bidomain::m_sigmaez
private

Definition at line 100 of file Bidomain.h.

◆ m_sigmaix

NekDouble Nektar::Bidomain::m_sigmaix
private

Definition at line 99 of file Bidomain.h.

◆ m_sigmaiy

NekDouble Nektar::Bidomain::m_sigmaiy
private

Definition at line 99 of file Bidomain.h.

◆ m_sigmaiz

NekDouble Nektar::Bidomain::m_sigmaiz
private

Definition at line 99 of file Bidomain.h.

◆ m_stimDuration

NekDouble Nektar::Bidomain::m_stimDuration
private

Stimulus current.

Definition at line 110 of file Bidomain.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_vardiffi

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

Definition at line 102 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_vardiffie

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

Definition at line 103 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ tmp1

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

Definition at line 105 of file Bidomain.h.

Referenced by v_InitObject().

◆ tmp2

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

Definition at line 106 of file Bidomain.h.

Referenced by v_InitObject().

◆ tmp3

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

Definition at line 107 of file Bidomain.h.

Referenced by v_InitObject().