Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
Nektar::Bidomain Class Reference

A model for cardiac conduction. More...

#include <Bidomain.h>

Inheritance diagram for Nektar::Bidomain:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Bidomain:
Collaboration graph
[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 >
boost::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $. More...
 
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
 
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const
Array< OneD, int > 
GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const
std::string 
GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetNumElmVelocity ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison. More...
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

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

Static Public Attributes

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

Protected Member Functions

 Bidomain (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
virtual void v_InitObject ()
 Init object for UnsteadySystem class. More...
 
void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
 Solve for the diffusion term. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Computes the reaction terms $f(u,v)$ and $g(u,v)$. More...
 
virtual void v_SetInitialConditions (NekDouble initialtime, bool dumpInitialConditions, const int domain)
 Sets a custom initial condition. More...
 
virtual void v_GenerateSummary (SummaryList &s)
 Prints a summary of the model parameters. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file. More...
 
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_SteadyStateCheck (int step)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const string &s1, const string &s2)
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
 
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Private 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...
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
map< std::string, Array< OneD,
Array< OneD, float > > > 
m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD,
Array< OneD, unsigned int > > > 
m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
std::set< std::string > m_loadedFields
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
int m_NumMode
 Mode to use in case of single mode analysis. More...
 

Detailed Description

A model for cardiac conduction.

Base model of cardiac electrophysiology of the form

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

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

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

Definition at line 49 of file Bidomain.h.

Constructor & Destructor Documentation

Nektar::Bidomain::~Bidomain ( )
virtual

Desctructor.

Definition at line 177 of file Bidomain.cpp.

178  {
179 
180  }
Nektar::Bidomain::Bidomain ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Constructor.

Definition at line 70 of file Bidomain.cpp.

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

Member Function Documentation

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

Creates an instance of this class.

Definition at line 55 of file Bidomain.h.

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

57  {
59  p->InitObject();
60  return p;
61  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::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 189 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().

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

383  {
384  int nq = m_fields[0]->GetNpoints();
385  m_cell->TimeIntegrate(inarray, outarray, time);
386  if (m_stimDuration > 0 && time < m_stimDuration)
387  {
388  Array<OneD,NekDouble> x0(nq);
389  Array<OneD,NekDouble> x1(nq);
390  Array<OneD,NekDouble> x2(nq);
391  Array<OneD,NekDouble> result(nq);
392 
393  // get the coordinates
394  m_fields[0]->GetCoords(x0,x1,x2);
395 
397  = m_session->GetFunction("Stimulus", "u");
398  ifunc->Evaluate(x0,x1,x2,time, result);
399 
400  Vmath::Vadd(nq, outarray[0], 1, result, 1, outarray[0], 1);
401  }
402  Vmath::Smul(nq, 1.0/m_capMembrane, outarray[0], 1, outarray[0], 1);
403  }
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
NekDouble m_stimDuration
Stimulus current.
Definition: Bidomain.h:111
boost::shared_ptr< Equation > EquationSharedPtr
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:99
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
NekDouble m_capMembrane
Definition: Bidomain.h:101
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:285
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 417 of file Bidomain.cpp.

References ASSERTL0, and m_cell.

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

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 76 of file Bidomain.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::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, and Vmath::Vadd().

77  {
78  UnsteadySystem::v_InitObject();
79  m_session->LoadParameter("Chi", m_chi);
80  m_session->LoadParameter("Cm", m_capMembrane);
81 
82  std::string vCellModel;
83  m_session->LoadSolverInfo("CELLMODEL", vCellModel, "");
84 
85  ASSERTL0(vCellModel != "", "Cell Model not specified.");
86 
88  m_intVariables.push_back(0);
89  m_intVariables.push_back(1);
90 
91  // Load variable coefficients
92  StdRegions::VarCoeffType varCoeffEnum[3] = {
96  };
97  std::string varName[3] = {
98  "AnisotropicConductivityX",
99  "AnisotropicConductivityY",
100  "AnisotropicConductivityZ"
101  };
102 
103 
104  if (m_session->DefinesFunction("IntracellularConductivity") && m_session->DefinesFunction("ExtracellularConductivity"))
105  {
106  for (int i = 0; i < m_spacedim; ++i)
107  {
108  int nq = m_fields[0]->GetNpoints();
109  Array<OneD,NekDouble> x0(nq);
110  Array<OneD,NekDouble> x1(nq);
111  Array<OneD,NekDouble> x2(nq);
112 
113  // get the coordinates
114  m_fields[0]->GetCoords(x0,x1,x2);
118  tmp1[i] = Array<OneD, NekDouble>(nq);
119  tmp2[i] = Array<OneD, NekDouble>(nq);
120  tmp3[i] = Array<OneD, NekDouble>(nq);
121 
123  = m_session->GetFunction("IntracellularConductivity", varName[i]);
125  = m_session->GetFunction("ExtracellularConductivity", varName[i]);
126  for(int j = 0; j < nq; j++)
127  {
128  tmp1[i][j] = ifunc1->Evaluate(x0[j],x1[j],x2[j],0.0);
129  tmp2[i][j] = ifunc2->Evaluate(x0[j],x1[j],x2[j],0.0);
130  }
131  Vmath::Vadd(nq, tmp1[i], 1, tmp2[i], 1, tmp3[i], 1);
132  m_vardiffi[varCoeffEnum[i]] = tmp1[i];
133  m_vardiffie[varCoeffEnum[i]] = tmp3[i];
134  }
135  }
136 
137 
138  if (m_session->DefinesParameter("StimulusDuration"))
139  {
140  ASSERTL0(m_session->DefinesFunction("Stimulus", "u"),
141  "Stimulus function not defined.");
142  m_session->LoadParameter("StimulusDuration", m_stimDuration);
143  }
144  else
145  {
146  m_stimDuration = 0;
147  }
148 
149 
150  // Search through the loaded filters and pass the cell model to any
151  // CheckpointCellModel filters loaded.
152  int k = 0;
153  const LibUtilities::FilterMap& f = m_session->GetFilters();
154  LibUtilities::FilterMap::const_iterator x;
155  for (x = f.begin(); x != f.end(); ++x, ++k)
156  {
157  if (x->first == "CheckpointCellModel")
158  {
159  boost::shared_ptr<FilterCheckpointCellModel> c
160  = boost::dynamic_pointer_cast<FilterCheckpointCellModel>(
161  m_filters[k]);
162  c->SetCellModel(m_cell);
163  }
164  }
165 
166  if (!m_explicitDiffusion)
167  {
169  }
171  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
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:379
void SetCellModel(CellModelSharedPtr &pCellModel)
NekDouble m_chi
Definition: Bidomain.h:101
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
StdRegions::VarCoeffMap m_vardiffie
Definition: Bidomain.h:104
NekDouble m_stimDuration
Stimulus current.
Definition: Bidomain.h:111
std::vector< std::pair< std::string, FilterParams > > FilterMap
Definition: SessionReader.h:66
int m_spacedim
Spatial dimension (>= expansion dim).
boost::shared_ptr< Equation > EquationSharedPtr
StdRegions::VarCoeffMap m_vardiffi
Definition: Bidomain.h:103
CellModelFactory & GetCellModelFactory()
Definition: CellModel.cpp:45
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:99
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
NekDouble m_capMembrane
Definition: Bidomain.h:101
std::vector< FilterSharedPtr > m_filters
Array< OneD, Array< OneD, NekDouble > > tmp1
Definition: Bidomain.h:106
Array< OneD, Array< OneD, NekDouble > > tmp2
Definition: Bidomain.h:107
Array< OneD, Array< OneD, NekDouble > > tmp3
Definition: Bidomain.h:108
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:189
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:285
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 406 of file Bidomain.cpp.

References m_cell.

409  {
410  EquationSystem::v_SetInitialConditions(initialtime, dumpInitialConditions, domain);
411  m_cell->Initialise();
412  }
CellModelSharedPtr m_cell
Cell model.
Definition: Bidomain.h:99

Friends And Related Function Documentation

friend class MemoryManager< Bidomain >
friend

Definition at line 52 of file Bidomain.h.

Member Data Documentation

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

NekDouble Nektar::Bidomain::m_capMembrane
private

Definition at line 101 of file Bidomain.h.

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

CellModelSharedPtr Nektar::Bidomain::m_cell
private

Cell model.

Definition at line 99 of file Bidomain.h.

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

NekDouble Nektar::Bidomain::m_chi
private

Definition at line 101 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

NekDouble Nektar::Bidomain::m_sigmaex
private

Definition at line 101 of file Bidomain.h.

NekDouble Nektar::Bidomain::m_sigmaey
private

Definition at line 101 of file Bidomain.h.

NekDouble Nektar::Bidomain::m_sigmaez
private

Definition at line 101 of file Bidomain.h.

NekDouble Nektar::Bidomain::m_sigmaix
private

Definition at line 101 of file Bidomain.h.

NekDouble Nektar::Bidomain::m_sigmaiy
private

Definition at line 101 of file Bidomain.h.

NekDouble Nektar::Bidomain::m_sigmaiz
private

Definition at line 101 of file Bidomain.h.

NekDouble Nektar::Bidomain::m_stimDuration
private

Stimulus current.

Definition at line 111 of file Bidomain.h.

Referenced by DoOdeRhs(), and v_InitObject().

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

Definition at line 103 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

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

Definition at line 104 of file Bidomain.h.

Referenced by DoImplicitSolve(), and v_InitObject().

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

Definition at line 106 of file Bidomain.h.

Referenced by v_InitObject().

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

Definition at line 107 of file Bidomain.h.

Referenced by v_InitObject().

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

Definition at line 108 of file Bidomain.h.

Referenced by v_InitObject().