Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CellModel.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File CellModel.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Cell model base class.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
39 
41 //#include <LibUtilities/LinearAlgebra/Blas.hpp>
42 
43 namespace Nektar
44 {
46  {
47  typedef Loki::SingletonHolder<CellModelFactory,
48  Loki::CreateUsingNew,
49  Loki::NoDestroy > Type;
50  return Type::Instance();
51  }
52 
53  /**
54  * @class CellModel
55  *
56  * The CellModel class and derived classes implement a range of cell model
57  * ODE systems. A cell model comprises a system of ion concentration
58  * variables and zero or more gating variables. Gating variables are
59  * time-integrated using the Rush-Larsen method and for each variable y,
60  * the corresponding y_inf and tau_y value is computed by Update(). The tau
61  * values are stored in separate storage to inarray/outarray, #m_gates_tau.
62  */
63 
64  /**
65  * Cell model base class constructor.
66  */
68  const MultiRegions::ExpListSharedPtr& pField)
69  {
70  m_session = pSession;
71  m_field = pField;
72  m_lastTime = 0.0;
73  m_substeps = pSession->GetParameter("Substeps");
74  m_nvar = 0;
75  m_useNodal = false;
76 
77  // Number of points in nodal space is the number of coefficients
78  // in modified basis
79  std::set<enum LibUtilities::ShapeType> s;
80  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
81  {
82  s.insert(m_field->GetExp(i)->DetShapeType());
83  }
84 
85  // Use nodal projection if only triangles
86  if (s.size() == 1 && (s.count(LibUtilities::eTriangle) == 1 ||
87  s.count(LibUtilities::eTetrahedron) == 1))
88  {
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
109 
114  }
115  else
116  {
117  m_nq = pField->GetTotPoints();
118  }
119  }
120 
121 
122  /**
123  * Initialise the cell model. Allocate workspace and variable storage.
124  */
126  {
127  ASSERTL1(m_nvar > 0, "Cell model must have at least 1 variable.");
128 
129  m_cellSol = Array<OneD, Array<OneD, NekDouble> >(m_nvar);
130  m_wsp = Array<OneD, Array<OneD, NekDouble> >(m_nvar);
131  for (unsigned int i = 0; i < m_nvar; ++i)
132  {
133  m_cellSol[i] = Array<OneD, NekDouble>(m_nq);
134  m_wsp[i] = Array<OneD, NekDouble>(m_nq);
135  }
136  m_gates_tau = Array<OneD, Array<OneD, NekDouble> >(m_gates.size());
137  for (unsigned int i = 0; i < m_gates.size(); ++i)
138  {
139  m_gates_tau[i] = Array<OneD, NekDouble>(m_nq);
140  }
141 
142  if (m_session->DefinesFunction("CellModelInitialConditions"))
143  {
144  LoadCellModel();
145  }
146  else
147  {
149  }
150  }
151 
152  /**
153  * Integrates the cell model for one PDE time-step. Cell model is
154  * sub-stepped.
155  *
156  * Ion concentrations and membrane potential are integrated using forward
157  * Euler, while gating variables are integrated using the Rush-Larsen
158  * scheme.
159  */
161  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
162  Array<OneD, Array<OneD, NekDouble> > &outarray,
163  const NekDouble time)
164  {
165  int phys_offset = 0;
166  int coef_offset = 0;
167  int nvar = inarray.num_elements();
168  Array<OneD, NekDouble> tmp;
169 
170  // ---------------------------
171  // Check nodal temp array set up
172  if (m_useNodal)
173  {
174  if (!m_nodalTmp.num_elements())
175  {
176  m_nodalTmp = Array<OneD, Array<OneD, NekDouble> >(nvar);
177  for (unsigned int k = 0; k < nvar; ++k)
178  {
179  m_nodalTmp[k] = Array<OneD, NekDouble>(m_nq);
180  }
181  }
182 
183  // Move to nodal points
184  Array<OneD, NekDouble> tmpCoeffs(max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
185 
186  for (unsigned int k = 0; k < nvar; ++k)
187  {
188  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
189  {
190  phys_offset = m_field->GetPhys_Offset(i);
191  coef_offset = m_field->GetCoeff_Offset(i);
192  if (m_field->GetExp(0)->DetShapeType() == LibUtilities::eTriangle)
193  {
194  m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset, tmpCoeffs);
195  m_nodalTri->ModalToNodal(tmpCoeffs, tmp=m_nodalTmp[k]+coef_offset);
196  }
197  else
198  {
199  m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset, tmpCoeffs);
200  m_nodalTet->ModalToNodal(tmpCoeffs, tmp=m_nodalTmp[k]+coef_offset);
201  }
202  }
203  }
204  // Copy new transmembrane potential into cell model
205  Vmath::Vcopy(m_nq, m_nodalTmp[0], 1, m_cellSol[0], 1);
206  }
207  else
208  {
209  // Copy new transmembrane potential into cell model
210  Vmath::Vcopy(m_nq, inarray[0], 1, m_cellSol[0], 1);
211  }
212  // -------------------------
213 
214  NekDouble delta_t = (time - m_lastTime)/m_substeps;
215 
216 
217  // Perform substepping
218  for (unsigned int i = 0; i < m_substeps - 1; ++i)
219  {
220  Update(m_cellSol, m_wsp, time);
221  // Voltage
222  Vmath::Svtvp(m_nq, delta_t, m_wsp[0], 1, m_cellSol[0], 1, m_cellSol[0], 1);
223  // Ion concentrations
224  for (unsigned int j = 0; j < m_concentrations.size(); ++j)
225  {
227  }
228  // Gating variables: Rush-Larsen scheme
229  for (unsigned int j = 0; j < m_gates.size(); ++j)
230  {
231  Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
232  Vmath::Vexp(m_nq, m_gates_tau[j], 1, m_gates_tau[j], 1);
233  Vmath::Vsub(m_nq, m_cellSol[m_gates[j]], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
234  Vmath::Vvtvp(m_nq, m_cellSol[m_gates[j]], 1, m_gates_tau[j], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
235  }
236  }
237 
238  // Perform final cell model step
239  Update(m_cellSol, m_wsp, time);
240 
241  // Output dV/dt from last step but integrate remaining cell model vars
242  // Transform cell model I_total from nodal to modal space
243  if (m_useNodal)
244  {
245  Array<OneD, NekDouble> tmpCoeffs(max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
246 
247  for (unsigned int k = 0; k < nvar; ++k)
248  {
249  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
250  {
251  int phys_offset = m_field->GetPhys_Offset(i);
252  int coef_offset = m_field->GetCoeff_Offset(i);
253  if (m_field->GetExp(0)->DetShapeType() == LibUtilities::eTriangle)
254  {
255  m_nodalTri->NodalToModal(m_wsp[k]+coef_offset, tmpCoeffs);
256  m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp=outarray[k] + phys_offset);
257  }
258  else
259  {
260  m_nodalTet->NodalToModal(m_wsp[k]+coef_offset, tmpCoeffs);
261  m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp=outarray[k] + phys_offset);
262  }
263  }
264  }
265  }
266  else
267  {
268  Vmath::Vcopy(m_nq, m_wsp[0], 1, outarray[0], 1);
269  }
270 
271  // Ion concentrations
272  for (unsigned int j = 0; j < m_concentrations.size(); ++j)
273  {
275  }
276 
277  // Gating variables: Rush-Larsen scheme
278  for (unsigned int j = 0; j < m_gates.size(); ++j)
279  {
280  Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
281  Vmath::Vexp(m_nq, m_gates_tau[j], 1, m_gates_tau[j], 1);
282  Vmath::Vsub(m_nq, m_cellSol[m_gates[j]], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
283  Vmath::Vvtvp(m_nq, m_cellSol[m_gates[j]], 1, m_gates_tau[j], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
284  }
285 
286  m_lastTime = time;
287  }
288 
289  Array<OneD, NekDouble> CellModel::GetCellSolutionCoeffs(unsigned int idx)
290  {
291  ASSERTL0(idx < m_nvar, "Index out of range for cell model.");
292 
293  Array<OneD, NekDouble> outarray(m_field->GetNcoeffs());
294  Array<OneD, NekDouble> tmp;
295 
296  if (m_useNodal)
297  {
298  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
299  {
300  int coef_offset = m_field->GetCoeff_Offset(i);
301  if (m_field->GetExp(0)->DetShapeType() == LibUtilities::eTriangle)
302  {
303  m_nodalTri->NodalToModal(m_cellSol[idx]+coef_offset, tmp=outarray+coef_offset);
304  }
305  else
306  {
307  m_nodalTet->NodalToModal(m_cellSol[idx]+coef_offset, tmp=outarray+coef_offset);
308  }
309  }
310  }
311  else
312  {
313  m_field->FwdTrans_IterPerExp(m_cellSol[idx], outarray);
314  }
315 
316  return outarray;
317  }
318 
319  Array<OneD, NekDouble> CellModel::GetCellSolution(unsigned int idx)
320  {
321  return m_cellSol[idx];
322  }
323 
325  {
326  const bool root = (m_session->GetComm()->GetRank() == 0);
327  const std::string fncName = "CellModelInitialConditions";
328  std::string varName;
329  Array<OneD, NekDouble> coeffs(m_field->GetNcoeffs());
330  Array<OneD, NekDouble> tmp;
331 
332  SpatialDomains::MeshGraphSharedPtr vGraph = m_field->GetGraph();
333 
334  if (root)
335  {
336  cout << "Cell model initial conditions: " << endl;
337  }
338 
339  // Load each cell model variable
340  // j=0 and j=1 are for transmembrane or intra/extra-cellular volt.
341  Vmath::Zero(m_nq, m_cellSol[0], 1);
342  for(int j = 1; j < m_cellSol.num_elements(); ++j)
343  {
344  // Get the name of the jth variable
345  varName = GetCellVarName(j);
346 
347  // Check if this variable is defined in a file or analytically
348  if (m_session->GetFunctionType(fncName, varName) ==
350  {
351  const std::string file =
352  m_session->GetFunctionFilename(fncName, varName);
353 
354  if (root)
355  {
356  cout << " - Field " << varName << ": from file "
357  << file << endl;
358  }
359 
360  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
361  std::vector<std::vector<NekDouble> > FieldData;
362 
363  // Read the restart file containing this variable
366  fld->Import(file, FieldDef, FieldData);
367 
368  LibUtilities::FieldMetaDataMap fieldMetaDataMap;
370  fld->ImportFieldMetaData(file,fieldMetaDataMap);
371  iter = fieldMetaDataMap.find("Time");
372  if(iter != fieldMetaDataMap.end())
373  {
374  m_lastTime = boost::lexical_cast<NekDouble>(iter->second);
375  }
376 
377  // Extract the data into the modal coefficients
378  for(int i = 0; i < FieldDef.size(); ++i)
379  {
380  m_field->ExtractDataToCoeffs(FieldDef[i],
381  FieldData[i],
382  varName,
383  coeffs);
384  }
385 
386  // If using nodal cell model then we do a modal->nodal transform
387  // otherwise we do a backward transform onto physical points.
388  if (m_useNodal)
389  {
390  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
391  {
392  int coef_offset = m_field->GetCoeff_Offset(i);
393  if (m_field->GetExp(0)->DetShapeType() ==
395  {
396  m_nodalTri->ModalToNodal(coeffs+coef_offset,
397  tmp=m_cellSol[j]+coef_offset);
398  }
399  else
400  {
401  m_nodalTet->ModalToNodal(coeffs+coef_offset,
402  tmp=m_cellSol[j]+coef_offset);
403  }
404  }
405  }
406  else
407  {
408  m_field->BwdTrans(coeffs, m_cellSol[j]);
409  }
410  }
411  else if (m_session->GetFunctionType(fncName, varName) ==
413  {
415  m_session->GetFunction(fncName, varName);
416 
417  if (root)
418  {
419  cout << " - Field " << varName << ": "
420  << equ->GetExpression() << endl;
421  }
422 
423  const unsigned int nphys = m_field->GetNpoints();
424  Array<OneD, NekDouble> x0(nphys);
425  Array<OneD, NekDouble> x1(nphys);
426  Array<OneD, NekDouble> x2(nphys);
427  m_field->GetCoords(x0,x1,x2);
428 
429  if (m_useNodal)
430  {
431  Array<OneD, NekDouble> phys(nphys);
432  Array<OneD, NekDouble> tmpCoeffs(max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
433 
434  equ->Evaluate(x0, x1, x2, phys);
435  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
436  {
437  int phys_offset = m_field->GetPhys_Offset(i);
438  int coef_offset = m_field->GetCoeff_Offset(i);
439  if (m_field->GetExp(0)->DetShapeType() ==
441  {
442  m_field->GetExp(0)->FwdTrans(
443  phys + phys_offset, tmpCoeffs);
444  m_nodalTri->ModalToNodal(
445  tmpCoeffs,
446  tmp = m_cellSol[j] + coef_offset);
447  }
448  else
449  {
450  m_field->GetExp(0)->FwdTrans(
451  phys + phys_offset, tmpCoeffs);
452  m_nodalTet->ModalToNodal(
453  tmpCoeffs,
454  tmp = m_cellSol[j] + coef_offset);
455  }
456  }
457  }
458  else
459  {
460  equ->Evaluate(x0, x1, x2, m_cellSol[j]);
461  }
462  }
463  }
464  }
465 }