72 m_substeps = pSession->GetParameter(
"Substeps");
78 std::set<enum LibUtilities::ShapeType> s;
79 for (
size_t i = 0; i <
m_field->GetNumElmts(); ++i)
81 s.insert(
m_field->GetExp(i)->DetShapeType());
96 m_nq = pField->GetNcoeffs();
97 int order =
m_field->GetExp(0)->GetBasis(0)->GetNumModes();
107 LibUtilities::eGaussRadauMAlpha1Beta0));
111 LibUtilities::eGaussRadauMAlpha2Beta0));
122 m_nq = pField->GetTotPoints();
131 ASSERTL1(
m_nvar > 0,
"Cell model must have at least 1 variable.");
135 for (
size_t i = 0; i <
m_nvar; ++i)
141 for (
size_t i = 0; i <
m_gates.size(); ++i)
146 if (
m_session->DefinesFunction(
"CellModelInitialConditions"))
170 size_t nvar = inarray.size();
180 for (
size_t k = 0; k < nvar; ++k)
190 for (
size_t k = 0; k < nvar; ++k)
192 for (
size_t i = 0; i <
m_field->GetNumElmts(); ++i)
194 phys_offset =
m_field->GetPhys_Offset(i);
195 coef_offset =
m_field->GetCoeff_Offset(i);
196 if (
m_field->GetExp(0)->DetShapeType() ==
199 m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset,
206 m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset,
240 for (
size_t j = 0; j <
m_gates.size(); ++j)
261 for (
size_t k = 0; k < nvar; ++k)
263 for (
size_t i = 0; i <
m_field->GetNumElmts(); ++i)
265 int phys_offset =
m_field->GetPhys_Offset(i);
266 int coef_offset =
m_field->GetCoeff_Offset(i);
267 if (
m_field->GetExp(0)->DetShapeType() ==
271 m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp = outarray[k] +
277 m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp = outarray[k] +
297 for (
size_t j = 0; j <
m_gates.size(); ++j)
319 for (
size_t i = 0; i <
m_field->GetNumElmts(); ++i)
321 int coef_offset =
m_field->GetCoeff_Offset(i);
325 tmp = outarray + coef_offset);
330 tmp = outarray + coef_offset);
349 const bool root = (
m_session->GetComm()->GetRank() == 0);
350 const std::string fncName =
"CellModelInitialConditions";
361 cout <<
"Cell model initial conditions: " << endl;
365 std::set<std::string> filelist;
366 for (j = 1; j < nvar; ++j)
371 if (
m_session->GetFunctionType(fncName, varName) ==
374 filelist.insert(
m_session->GetFunctionFilename(fncName, varName));
379 typedef std::vector<LibUtilities::FieldDefinitionsSharedPtr> FDef;
380 typedef std::vector<std::vector<NekDouble>> FData;
381 std::map<std::string, FDef> FieldDef;
382 std::map<std::string, FData> FieldData;
385 for (
auto &setIt : filelist)
389 cout <<
" - Reading file: " << setIt << endl;
391 FieldDef[setIt] = FDef(0);
392 FieldData[setIt] = FData(0);
395 fld->Import(setIt, FieldDef[setIt], FieldData[setIt], fieldMetaDataMap);
399 auto iter = fieldMetaDataMap.find(
"Time");
400 if (iter != fieldMetaDataMap.end())
414 if (
m_session->GetFunctionType(fncName, varName) ==
417 const std::string file =
418 m_session->GetFunctionFilename(fncName, varName);
422 cout <<
" - Field " << varName <<
": from file " << file
427 for (
size_t i = 0; i < FieldDef[file].size(); ++i)
430 FieldDef[file][i], FieldData[file][i], varName, coeffs);
437 for (
size_t i = 0; i <
m_field->GetNumElmts(); ++i)
439 int coef_offset =
m_field->GetCoeff_Offset(i);
440 if (
m_field->GetExp(0)->DetShapeType() ==
443 m_nodalTri->ModalToNodal(coeffs + coef_offset,
449 m_nodalTet->ModalToNodal(coeffs + coef_offset,
460 else if (
m_session->GetFunctionType(fncName, varName) ==
464 m_session->GetFunction(fncName, varName);
468 cout <<
" - Field " << varName <<
": " << equ->GetExpression()
472 const size_t nphys =
m_field->GetNpoints();
476 m_field->GetCoords(x0, x1, x2);
484 equ->Evaluate(x0, x1, x2, phys);
485 for (
size_t i = 0; i <
m_field->GetNumElmts(); ++i)
487 int phys_offset =
m_field->GetPhys_Offset(i);
488 int coef_offset =
m_field->GetCoeff_Offset(i);
489 if (
m_field->GetExp(0)->DetShapeType() ==
492 m_field->GetExp(0)->FwdTrans(phys + phys_offset,
499 m_field->GetExp(0)->FwdTrans(phys + phys_offset,
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Array< OneD, Array< OneD, NekDouble > > m_cellSol
Cell model solution variables.
void Initialise()
Initialise the cell model storage and set initial conditions.
void Update(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the derivatives of cell model variables.
Array< OneD, Array< OneD, NekDouble > > m_wsp
Cell model integration workspace.
bool m_useNodal
Flag indicating whether nodal projection in use.
Array< OneD, NekDouble > GetCellSolutionCoeffs(size_t idx)
NekDouble m_lastTime
Timestep for pde model.
size_t m_substeps
Number of substeps to take.
StdRegions::StdNodalTetExpSharedPtr m_nodalTet
Array< OneD, Array< OneD, NekDouble > > m_nodalTmp
Temporary array for nodal projection.
MultiRegions::ExpListSharedPtr m_field
Transmembrane potential field from PDE system.
std::vector< int > m_concentrations
Indices of cell model variables which are concentrations.
std::vector< int > m_gates
Indices of cell model variables which are gates.
Array< OneD, NekDouble > GetCellSolution(size_t idx)
virtual void v_SetInitialConditions()=0
StdRegions::StdNodalTriExpSharedPtr m_nodalTri
StdNodalTri for cell model calculations.
LibUtilities::SessionReaderSharedPtr m_session
Session.
size_t m_nq
Number of physical points.
size_t m_nvar
Number of variables in cell model (inc. transmembrane voltage)
CellModel(const LibUtilities::SessionReaderSharedPtr &pSession, const MultiRegions::ExpListSharedPtr &pField)
void TimeIntegrate(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Time integrate the cell model by one PDE timestep.
Array< OneD, Array< OneD, NekDouble > > m_gates_tau
Storage for gate tau values.
std::string GetCellVarName(size_t idx)
Describes the specification for a Basis.
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Provides a generic Factory class.
Defines a specification for a set of points.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< FieldIO > FieldIOSharedPtr
std::map< std::string, std::string > FieldMetaDataMap
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
@ eNodalTriEvenlySpaced
2D Evenly-spaced points on a Triangle
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eNodalTetEvenlySpaced
3D Evenly-spaced points on a Tetrahedron
@ eFunctionTypeExpression
@ eModified_B
Principle Modified Functions .
@ eModified_C
Principle Modified Functions .
@ eModified_A
Principle Modified Functions .
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
CellModelFactory & GetCellModelFactory()
void Vexp(int n, const T *x, const int incx, T *y, const int incy)
exp y = exp(x)
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.