Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Nektar::CellModel Class Referenceabstract

Cell model base class. More...

#include <CellModel.h>

Inheritance diagram for Nektar::CellModel:
[legend]

Public Member Functions

 CellModel (const LibUtilities::SessionReaderSharedPtr &pSession, const MultiRegions::ExpListSharedPtr &pField)
 
virtual ~CellModel ()
 
void Initialise ()
 Initialise the cell model storage and set initial conditions. More...
 
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. More...
 
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. More...
 
void GenerateSummary (SummaryList &s)
 Print a summary of the cell model. More...
 
size_t GetNumCellVariables ()
 
std::string GetCellVarName (size_t idx)
 
Array< OneD, NekDoubleGetCellSolutionCoeffs (size_t idx)
 
Array< OneD, NekDoubleGetCellSolution (size_t idx)
 

Protected Member Functions

virtual void v_Update (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)=0
 
virtual void v_GenerateSummary (SummaryList &s)=0
 
virtual std::string v_GetCellVarName (size_t idx)
 
virtual void v_SetInitialConditions ()=0
 
void LoadCellModel ()
 

Protected Attributes

LibUtilities::SessionReaderSharedPtr m_session
 Session. More...
 
MultiRegions::ExpListSharedPtr m_field
 Transmembrane potential field from PDE system. More...
 
size_t m_nq
 Number of physical points. More...
 
size_t m_nvar
 Number of variables in cell model (inc. transmembrane voltage) More...
 
NekDouble m_lastTime
 Timestep for pde model. More...
 
size_t m_substeps
 Number of substeps to take. More...
 
Array< OneD, Array< OneD, NekDouble > > m_cellSol
 Cell model solution variables. More...
 
Array< OneD, Array< OneD, NekDouble > > m_wsp
 Cell model integration workspace. More...
 
bool m_useNodal
 Flag indicating whether nodal projection in use. More...
 
StdRegions::StdNodalTriExpSharedPtr m_nodalTri
 StdNodalTri for cell model calculations. More...
 
StdRegions::StdNodalTetExpSharedPtr m_nodalTet
 
Array< OneD, Array< OneD, NekDouble > > m_nodalTmp
 Temporary array for nodal projection. More...
 
std::vector< int > m_concentrations
 Indices of cell model variables which are concentrations. More...
 
std::vector< int > m_gates
 Indices of cell model variables which are gates. More...
 
Array< OneD, Array< OneD, NekDouble > > m_gates_tau
 Storage for gate tau values. More...
 

Detailed Description

Cell model base class.

The CellModel class and derived classes implement a range of cell model ODE systems. A cell model comprises a system of ion concentration variables and zero or more gating variables. Gating variables are time-integrated using the Rush-Larsen method and for each variable y, the corresponding y_inf and tau_y value is computed by Update(). The tau values are stored in separate storage to inarray/outarray, m_gates_tau.

Definition at line 65 of file CellModel.h.

Constructor & Destructor Documentation

◆ CellModel()

Nektar::CellModel::CellModel ( const LibUtilities::SessionReaderSharedPtr pSession,
const MultiRegions::ExpListSharedPtr pField 
)

Cell model base class constructor.

Definition at line 66 of file CellModel.cpp.

68{
69 m_session = pSession;
70 m_field = pField;
71 m_lastTime = 0.0;
72 m_substeps = pSession->GetParameter("Substeps");
73 m_nvar = 0;
74 m_useNodal = false;
75
76 // Number of points in nodal space is the number of coefficients
77 // in modified basis
78 std::set<enum LibUtilities::ShapeType> s;
79 for (size_t i = 0; i < m_field->GetNumElmts(); ++i)
80 {
81 s.insert(m_field->GetExp(i)->DetShapeType());
82 }
83
84 // Use nodal projection if only triangles
85 if (s.size() == 1 && (s.count(LibUtilities::eTriangle) == 1 ||
86 s.count(LibUtilities::eTetrahedron) == 1))
87 {
88 // This is disabled for now as it causes problems at high order.
89 // m_useNodal = true;
90 }
91
92 // ---------------------------
93 // Move to nodal points
94 if (m_useNodal)
95 {
96 m_nq = pField->GetNcoeffs();
97 int order = m_field->GetExp(0)->GetBasis(0)->GetNumModes();
98
99 // Set up a nodal tri
100 LibUtilities::BasisKey B0(
102 LibUtilities::PointsKey(order,
104 LibUtilities::BasisKey B1(
106 LibUtilities::PointsKey(order,
107 LibUtilities::eGaussRadauMAlpha1Beta0));
108 LibUtilities::BasisKey B2(
110 LibUtilities::PointsKey(order,
111 LibUtilities::eGaussRadauMAlpha2Beta0));
112
113 m_nodalTri =
116 m_nodalTet =
119 }
120 else
121 {
122 m_nq = pField->GetTotPoints();
123 }
124}
bool m_useNodal
Flag indicating whether nodal projection in use.
Definition: CellModel.h:131
NekDouble m_lastTime
Timestep for pde model.
Definition: CellModel.h:121
size_t m_substeps
Number of substeps to take.
Definition: CellModel.h:123
StdRegions::StdNodalTetExpSharedPtr m_nodalTet
Definition: CellModel.h:134
MultiRegions::ExpListSharedPtr m_field
Transmembrane potential field from PDE system.
Definition: CellModel.h:115
StdRegions::StdNodalTriExpSharedPtr m_nodalTri
StdNodalTri for cell model calculations.
Definition: CellModel.h:133
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: CellModel.h:113
size_t m_nq
Number of physical points.
Definition: CellModel.h:117
size_t m_nvar
Number of variables in cell model (inc. transmembrane voltage)
Definition: CellModel.h:119
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
@ eNodalTriEvenlySpaced
2D Evenly-spaced points on a Triangle
Definition: PointsType.h:83
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eNodalTetEvenlySpaced
3D Evenly-spaced points on a Tetrahedron
Definition: PointsType.h:84
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eNodalTetEvenlySpaced, Nektar::LibUtilities::eNodalTriEvenlySpaced, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, m_field, m_lastTime, m_nodalTet, m_nodalTri, m_nq, m_nvar, m_session, m_substeps, and m_useNodal.

◆ ~CellModel()

virtual Nektar::CellModel::~CellModel ( )
inlinevirtual

Definition at line 71 of file CellModel.h.

72 {
73 }

Member Function Documentation

◆ GenerateSummary()

void Nektar::CellModel::GenerateSummary ( SummaryList s)
inline

Print a summary of the cell model.

Definition at line 92 of file CellModel.h.

93 {
95 }
virtual void v_GenerateSummary(SummaryList &s)=0

References v_GenerateSummary().

◆ GetCellSolution()

Array< OneD, NekDouble > Nektar::CellModel::GetCellSolution ( size_t  idx)

Definition at line 342 of file CellModel.cpp.

343{
344 return m_cellSol[idx];
345}
Array< OneD, Array< OneD, NekDouble > > m_cellSol
Cell model solution variables.
Definition: CellModel.h:126

References m_cellSol.

◆ GetCellSolutionCoeffs()

Array< OneD, NekDouble > Nektar::CellModel::GetCellSolutionCoeffs ( size_t  idx)

Definition at line 310 of file CellModel.cpp.

311{
312 ASSERTL0(idx < m_nvar, "Index out of range for cell model.");
313
314 Array<OneD, NekDouble> outarray(m_field->GetNcoeffs());
315 Array<OneD, NekDouble> tmp;
316
317 if (m_useNodal)
318 {
319 for (size_t i = 0; i < m_field->GetNumElmts(); ++i)
320 {
321 int coef_offset = m_field->GetCoeff_Offset(i);
322 if (m_field->GetExp(0)->DetShapeType() == LibUtilities::eTriangle)
323 {
324 m_nodalTri->NodalToModal(m_cellSol[idx] + coef_offset,
325 tmp = outarray + coef_offset);
326 }
327 else
328 {
329 m_nodalTet->NodalToModal(m_cellSol[idx] + coef_offset,
330 tmp = outarray + coef_offset);
331 }
332 }
333 }
334 else
335 {
336 m_field->FwdTransLocalElmt(m_cellSol[idx], outarray);
337 }
338
339 return outarray;
340}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

References ASSERTL0, Nektar::LibUtilities::eTriangle, m_cellSol, m_field, m_nodalTet, m_nodalTri, m_nvar, and m_useNodal.

◆ GetCellVarName()

std::string Nektar::CellModel::GetCellVarName ( size_t  idx)
inline

Definition at line 102 of file CellModel.h.

103 {
104 return v_GetCellVarName(idx);
105 }
virtual std::string v_GetCellVarName(size_t idx)
Definition: CellModel.h:152

References v_GetCellVarName().

Referenced by LoadCellModel().

◆ GetNumCellVariables()

size_t Nektar::CellModel::GetNumCellVariables ( )
inline

Definition at line 97 of file CellModel.h.

98 {
99 return m_nvar;
100 }

References m_nvar.

◆ Initialise()

void Nektar::CellModel::Initialise ( )

Initialise the cell model storage and set initial conditions.

Initialise the cell model. Allocate workspace and variable storage.

Definition at line 129 of file CellModel.cpp.

130{
131 ASSERTL1(m_nvar > 0, "Cell model must have at least 1 variable.");
132
133 m_cellSol = Array<OneD, Array<OneD, NekDouble>>(m_nvar);
134 m_wsp = Array<OneD, Array<OneD, NekDouble>>(m_nvar);
135 for (size_t i = 0; i < m_nvar; ++i)
136 {
137 m_cellSol[i] = Array<OneD, NekDouble>(m_nq);
138 m_wsp[i] = Array<OneD, NekDouble>(m_nq);
139 }
140 m_gates_tau = Array<OneD, Array<OneD, NekDouble>>(m_gates.size());
141 for (size_t i = 0; i < m_gates.size(); ++i)
142 {
143 m_gates_tau[i] = Array<OneD, NekDouble>(m_nq);
144 }
145
146 if (m_session->DefinesFunction("CellModelInitialConditions"))
147 {
149 }
150 else
151 {
153 }
154}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Array< OneD, Array< OneD, NekDouble > > m_wsp
Cell model integration workspace.
Definition: CellModel.h:128
std::vector< int > m_gates
Indices of cell model variables which are gates.
Definition: CellModel.h:141
virtual void v_SetInitialConditions()=0
Array< OneD, Array< OneD, NekDouble > > m_gates_tau
Storage for gate tau values.
Definition: CellModel.h:143

References ASSERTL1, LoadCellModel(), m_cellSol, m_gates, m_gates_tau, m_nq, m_nvar, m_session, m_wsp, and v_SetInitialConditions().

◆ LoadCellModel()

void Nektar::CellModel::LoadCellModel ( )
protected

Definition at line 347 of file CellModel.cpp.

348{
349 const bool root = (m_session->GetComm()->GetRank() == 0);
350 const std::string fncName = "CellModelInitialConditions";
351 const size_t nvar = m_cellSol[0].size();
352 std::string varName;
353 Array<OneD, NekDouble> coeffs(m_field->GetNcoeffs());
354 Array<OneD, NekDouble> tmp;
355 size_t j = 0;
356
357 SpatialDomains::MeshGraphSharedPtr vGraph = m_field->GetGraph();
358
359 if (root)
360 {
361 cout << "Cell model initial conditions: " << endl;
362 }
363
364 // First determine all the files we need to load
365 std::set<std::string> filelist;
366 for (j = 1; j < nvar; ++j)
367 {
368 // Get the name of the jth variable
369 varName = GetCellVarName(j);
370
371 if (m_session->GetFunctionType(fncName, varName) ==
373 {
374 filelist.insert(m_session->GetFunctionFilename(fncName, varName));
375 }
376 }
377
378 // Read files
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;
383 LibUtilities::FieldMetaDataMap fieldMetaDataMap;
384
385 for (auto &setIt : filelist)
386 {
387 if (root)
388 {
389 cout << " - Reading file: " << setIt << endl;
390 }
391 FieldDef[setIt] = FDef(0);
392 FieldData[setIt] = FData(0);
395 fld->Import(setIt, FieldDef[setIt], FieldData[setIt], fieldMetaDataMap);
396 }
397
398 // Get time of checkpoint from file if available
399 auto iter = fieldMetaDataMap.find("Time");
400 if (iter != fieldMetaDataMap.end())
401 {
402 m_lastTime = std::stod(iter->second);
403 }
404
405 // Load each cell model variable
406 // j=0 and j=1 are for transmembrane or intra/extra-cellular volt.
407 Vmath::Zero(m_nq, m_cellSol[0], 1);
408 for (j = 1; j < m_cellSol.size(); ++j)
409 {
410 // Get the name of the jth variable
411 varName = GetCellVarName(j);
412
413 // Check if this variable is defined in a file or analytically
414 if (m_session->GetFunctionType(fncName, varName) ==
416 {
417 const std::string file =
418 m_session->GetFunctionFilename(fncName, varName);
419
420 if (root)
421 {
422 cout << " - Field " << varName << ": from file " << file
423 << endl;
424 }
425
426 // Extract the data into the modal coefficients
427 for (size_t i = 0; i < FieldDef[file].size(); ++i)
428 {
429 m_field->ExtractDataToCoeffs(
430 FieldDef[file][i], FieldData[file][i], varName, coeffs);
431 }
432
433 // If using nodal cell model then we do a modal->nodal transform
434 // otherwise we do a backward transform onto physical points.
435 if (m_useNodal)
436 {
437 for (size_t i = 0; i < m_field->GetNumElmts(); ++i)
438 {
439 int coef_offset = m_field->GetCoeff_Offset(i);
440 if (m_field->GetExp(0)->DetShapeType() ==
442 {
443 m_nodalTri->ModalToNodal(coeffs + coef_offset,
444 tmp = m_cellSol[j] +
445 coef_offset);
446 }
447 else
448 {
449 m_nodalTet->ModalToNodal(coeffs + coef_offset,
450 tmp = m_cellSol[j] +
451 coef_offset);
452 }
453 }
454 }
455 else
456 {
457 m_field->BwdTrans(coeffs, m_cellSol[j]);
458 }
459 }
460 else if (m_session->GetFunctionType(fncName, varName) ==
462 {
464 m_session->GetFunction(fncName, varName);
465
466 if (root)
467 {
468 cout << " - Field " << varName << ": " << equ->GetExpression()
469 << endl;
470 }
471
472 const size_t nphys = m_field->GetNpoints();
473 Array<OneD, NekDouble> x0(nphys);
474 Array<OneD, NekDouble> x1(nphys);
475 Array<OneD, NekDouble> x2(nphys);
476 m_field->GetCoords(x0, x1, x2);
477
478 if (m_useNodal)
479 {
480 Array<OneD, NekDouble> phys(nphys);
481 Array<OneD, NekDouble> tmpCoeffs(
482 max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
483
484 equ->Evaluate(x0, x1, x2, phys);
485 for (size_t i = 0; i < m_field->GetNumElmts(); ++i)
486 {
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() ==
491 {
492 m_field->GetExp(0)->FwdTrans(phys + phys_offset,
493 tmpCoeffs);
494 m_nodalTri->ModalToNodal(tmpCoeffs, tmp = m_cellSol[j] +
495 coef_offset);
496 }
497 else
498 {
499 m_field->GetExp(0)->FwdTrans(phys + phys_offset,
500 tmpCoeffs);
501 m_nodalTet->ModalToNodal(tmpCoeffs, tmp = m_cellSol[j] +
502 coef_offset);
503 }
504 }
505 }
506 else
507 {
508 equ->Evaluate(x0, x1, x2, m_cellSol[j]);
509 }
510 }
511 }
512}
std::string GetCellVarName(size_t idx)
Definition: CellModel.h:102
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:224
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::LibUtilities::FieldIO::CreateForFile(), Nektar::LibUtilities::eFunctionTypeExpression, Nektar::LibUtilities::eFunctionTypeFile, Nektar::LibUtilities::eTriangle, GetCellVarName(), m_cellSol, m_field, m_lastTime, m_nodalTet, m_nodalTri, m_nq, m_session, m_useNodal, and Vmath::Zero().

Referenced by Initialise().

◆ TimeIntegrate()

void Nektar::CellModel::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.

Integrates the cell model for one PDE time-step. Cell model is sub-stepped.

Ion concentrations and membrane potential are integrated using forward Euler, while gating variables are integrated using the Rush-Larsen scheme.

Definition at line 164 of file CellModel.cpp.

167{
168 int phys_offset = 0;
169 int coef_offset = 0;
170 size_t nvar = inarray.size();
171 Array<OneD, NekDouble> tmp;
172
173 // ---------------------------
174 // Check nodal temp array set up
175 if (m_useNodal)
176 {
177 if (!m_nodalTmp.size())
178 {
179 m_nodalTmp = Array<OneD, Array<OneD, NekDouble>>(nvar);
180 for (size_t k = 0; k < nvar; ++k)
181 {
182 m_nodalTmp[k] = Array<OneD, NekDouble>(m_nq);
183 }
184 }
185
186 // Move to nodal points
187 Array<OneD, NekDouble> tmpCoeffs(
188 max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
189
190 for (size_t k = 0; k < nvar; ++k)
191 {
192 for (size_t i = 0; i < m_field->GetNumElmts(); ++i)
193 {
194 phys_offset = m_field->GetPhys_Offset(i);
195 coef_offset = m_field->GetCoeff_Offset(i);
196 if (m_field->GetExp(0)->DetShapeType() ==
198 {
199 m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset,
200 tmpCoeffs);
201 m_nodalTri->ModalToNodal(tmpCoeffs,
202 tmp = m_nodalTmp[k] + coef_offset);
203 }
204 else
205 {
206 m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset,
207 tmpCoeffs);
208 m_nodalTet->ModalToNodal(tmpCoeffs,
209 tmp = m_nodalTmp[k] + coef_offset);
210 }
211 }
212 }
213 // Copy new transmembrane potential into cell model
214 Vmath::Vcopy(m_nq, m_nodalTmp[0], 1, m_cellSol[0], 1);
215 }
216 else
217 {
218 // Copy new transmembrane potential into cell model
219 Vmath::Vcopy(m_nq, inarray[0], 1, m_cellSol[0], 1);
220 }
221 // -------------------------
222
223 NekDouble delta_t = (time - m_lastTime) / m_substeps;
224
225 // Perform substepping
226 for (size_t i = 0; i < m_substeps - 1; ++i)
227 {
228 Update(m_cellSol, m_wsp, time);
229 // Voltage
230 Vmath::Svtvp(m_nq, delta_t, m_wsp[0], 1, m_cellSol[0], 1, m_cellSol[0],
231 1);
232 // Ion concentrations
233 for (size_t j = 0; j < m_concentrations.size(); ++j)
234 {
235 Vmath::Svtvp(m_nq, delta_t, m_wsp[m_concentrations[j]], 1,
238 }
239 // Gating variables: Rush-Larsen scheme
240 for (size_t j = 0; j < m_gates.size(); ++j)
241 {
242 Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
245 m_cellSol[m_gates[j]], 1);
247 m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
248 }
249 }
250
251 // Perform final cell model step
252 Update(m_cellSol, m_wsp, time);
253
254 // Output dV/dt from last step but integrate remaining cell model vars
255 // Transform cell model I_total from nodal to modal space
256 if (m_useNodal)
257 {
258 Array<OneD, NekDouble> tmpCoeffs(
259 max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
260
261 for (size_t k = 0; k < nvar; ++k)
262 {
263 for (size_t i = 0; i < m_field->GetNumElmts(); ++i)
264 {
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() ==
269 {
270 m_nodalTri->NodalToModal(m_wsp[k] + coef_offset, tmpCoeffs);
271 m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp = outarray[k] +
272 phys_offset);
273 }
274 else
275 {
276 m_nodalTet->NodalToModal(m_wsp[k] + coef_offset, tmpCoeffs);
277 m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp = outarray[k] +
278 phys_offset);
279 }
280 }
281 }
282 }
283 else
284 {
285 Vmath::Vcopy(m_nq, m_wsp[0], 1, outarray[0], 1);
286 }
287
288 // Ion concentrations
289 for (size_t j = 0; j < m_concentrations.size(); ++j)
290 {
291 Vmath::Svtvp(m_nq, delta_t, m_wsp[m_concentrations[j]], 1,
294 }
295
296 // Gating variables: Rush-Larsen scheme
297 for (size_t j = 0; j < m_gates.size(); ++j)
298 {
299 Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
302 m_cellSol[m_gates[j]], 1);
304 m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
305 }
306
307 m_lastTime = time;
308}
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.
Definition: CellModel.h:84
Array< OneD, Array< OneD, NekDouble > > m_nodalTmp
Temporary array for nodal projection.
Definition: CellModel.h:136
std::vector< int > m_concentrations
Indices of cell model variables which are concentrations.
Definition: CellModel.h:139
double NekDouble
void Vexp(int n, const T *x, const int incx, T *y, const int incy)
exp y = exp(x)
Definition: Vmath.hpp:315
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.
Definition: Vmath.hpp:396
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
Definition: Vmath.hpp:366
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.hpp:154
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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.
Definition: Vmath.hpp:220

References Nektar::LibUtilities::eTriangle, m_cellSol, m_concentrations, m_field, m_gates, m_gates_tau, m_lastTime, m_nodalTet, m_nodalTmp, m_nodalTri, m_nq, m_substeps, m_useNodal, m_wsp, Vmath::Sdiv(), Vmath::Svtvp(), Update(), Vmath::Vcopy(), Vmath::Vexp(), Vmath::Vsub(), and Vmath::Vvtvp().

◆ Update()

void Nektar::CellModel::Update ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
inline

Compute the derivatives of cell model variables.

Definition at line 84 of file CellModel.h.

87 {
88 v_Update(inarray, outarray, time);
89 }
virtual void v_Update(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)=0

References v_Update().

Referenced by TimeIntegrate().

◆ v_GenerateSummary()

virtual void Nektar::CellModel::v_GenerateSummary ( SummaryList s)
protectedpure virtual

◆ v_GetCellVarName()

virtual std::string Nektar::CellModel::v_GetCellVarName ( size_t  idx)
inlineprotectedvirtual

Reimplemented in Nektar::CourtemancheRamirezNattel98, and Nektar::FentonKarma.

Definition at line 152 of file CellModel.h.

153 {
154 return "Var" + std::to_string(idx);
155 }

Referenced by GetCellVarName().

◆ v_SetInitialConditions()

virtual void Nektar::CellModel::v_SetInitialConditions ( )
protectedpure virtual

◆ v_Update()

virtual void Nektar::CellModel::v_Update ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protectedpure virtual

Member Data Documentation

◆ m_cellSol

Array<OneD, Array<OneD, NekDouble> > Nektar::CellModel::m_cellSol
protected

◆ m_concentrations

std::vector<int> Nektar::CellModel::m_concentrations
protected

◆ m_field

MultiRegions::ExpListSharedPtr Nektar::CellModel::m_field
protected

Transmembrane potential field from PDE system.

Definition at line 115 of file CellModel.h.

Referenced by CellModel(), GetCellSolutionCoeffs(), LoadCellModel(), and TimeIntegrate().

◆ m_gates

std::vector<int> Nektar::CellModel::m_gates
protected

◆ m_gates_tau

Array<OneD, Array<OneD, NekDouble> > Nektar::CellModel::m_gates_tau
protected

◆ m_lastTime

NekDouble Nektar::CellModel::m_lastTime
protected

Timestep for pde model.

Definition at line 121 of file CellModel.h.

Referenced by CellModel(), LoadCellModel(), and TimeIntegrate().

◆ m_nodalTet

StdRegions::StdNodalTetExpSharedPtr Nektar::CellModel::m_nodalTet
protected

Definition at line 134 of file CellModel.h.

Referenced by CellModel(), GetCellSolutionCoeffs(), LoadCellModel(), and TimeIntegrate().

◆ m_nodalTmp

Array<OneD, Array<OneD, NekDouble> > Nektar::CellModel::m_nodalTmp
protected

Temporary array for nodal projection.

Definition at line 136 of file CellModel.h.

Referenced by TimeIntegrate().

◆ m_nodalTri

StdRegions::StdNodalTriExpSharedPtr Nektar::CellModel::m_nodalTri
protected

StdNodalTri for cell model calculations.

Definition at line 133 of file CellModel.h.

Referenced by CellModel(), GetCellSolutionCoeffs(), LoadCellModel(), and TimeIntegrate().

◆ m_nq

size_t Nektar::CellModel::m_nq
protected

◆ m_nvar

size_t Nektar::CellModel::m_nvar
protected

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::CellModel::m_session
protected

Session.

Definition at line 113 of file CellModel.h.

Referenced by CellModel(), Initialise(), and LoadCellModel().

◆ m_substeps

size_t Nektar::CellModel::m_substeps
protected

Number of substeps to take.

Definition at line 123 of file CellModel.h.

Referenced by CellModel(), and TimeIntegrate().

◆ m_useNodal

bool Nektar::CellModel::m_useNodal
protected

Flag indicating whether nodal projection in use.

Definition at line 131 of file CellModel.h.

Referenced by CellModel(), GetCellSolutionCoeffs(), LoadCellModel(), and TimeIntegrate().

◆ m_wsp

Array<OneD, Array<OneD, NekDouble> > Nektar::CellModel::m_wsp
protected

Cell model integration workspace.

Definition at line 128 of file CellModel.h.

Referenced by Initialise(), and TimeIntegrate().