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  // This is disabled for now as it causes problems at high order.
90  // m_useNodal = true;
91  }
92 
93  // ---------------------------
94  // Move to nodal points
95  if (m_useNodal)
96  {
97  m_nq = pField->GetNcoeffs();
98  int order = m_field->GetExp(0)->GetBasis(0)->GetNumModes();
99 
100  // Set up a nodal tri
110 
115  }
116  else
117  {
118  m_nq = pField->GetTotPoints();
119  }
120  }
121 
122 
123  /**
124  * Initialise the cell model. Allocate workspace and variable storage.
125  */
127  {
128  ASSERTL1(m_nvar > 0, "Cell model must have at least 1 variable.");
129 
130  m_cellSol = Array<OneD, Array<OneD, NekDouble> >(m_nvar);
131  m_wsp = Array<OneD, Array<OneD, NekDouble> >(m_nvar);
132  for (unsigned int i = 0; i < m_nvar; ++i)
133  {
134  m_cellSol[i] = Array<OneD, NekDouble>(m_nq);
135  m_wsp[i] = Array<OneD, NekDouble>(m_nq);
136  }
137  m_gates_tau = Array<OneD, Array<OneD, NekDouble> >(m_gates.size());
138  for (unsigned int i = 0; i < m_gates.size(); ++i)
139  {
140  m_gates_tau[i] = Array<OneD, NekDouble>(m_nq);
141  }
142 
143  if (m_session->DefinesFunction("CellModelInitialConditions"))
144  {
145  LoadCellModel();
146  }
147  else
148  {
150  }
151  }
152 
153  /**
154  * Integrates the cell model for one PDE time-step. Cell model is
155  * sub-stepped.
156  *
157  * Ion concentrations and membrane potential are integrated using forward
158  * Euler, while gating variables are integrated using the Rush-Larsen
159  * scheme.
160  */
162  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
163  Array<OneD, Array<OneD, NekDouble> > &outarray,
164  const NekDouble time)
165  {
166  int phys_offset = 0;
167  int coef_offset = 0;
168  int nvar = inarray.num_elements();
169  Array<OneD, NekDouble> tmp;
170 
171  // ---------------------------
172  // Check nodal temp array set up
173  if (m_useNodal)
174  {
175  if (!m_nodalTmp.num_elements())
176  {
177  m_nodalTmp = Array<OneD, Array<OneD, NekDouble> >(nvar);
178  for (unsigned int k = 0; k < nvar; ++k)
179  {
180  m_nodalTmp[k] = Array<OneD, NekDouble>(m_nq);
181  }
182  }
183 
184  // Move to nodal points
185  Array<OneD, NekDouble> tmpCoeffs(max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
186 
187  for (unsigned int k = 0; k < nvar; ++k)
188  {
189  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
190  {
191  phys_offset = m_field->GetPhys_Offset(i);
192  coef_offset = m_field->GetCoeff_Offset(i);
193  if (m_field->GetExp(0)->DetShapeType() == LibUtilities::eTriangle)
194  {
195  m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset, tmpCoeffs);
196  m_nodalTri->ModalToNodal(tmpCoeffs, tmp=m_nodalTmp[k]+coef_offset);
197  }
198  else
199  {
200  m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset, tmpCoeffs);
201  m_nodalTet->ModalToNodal(tmpCoeffs, tmp=m_nodalTmp[k]+coef_offset);
202  }
203  }
204  }
205  // Copy new transmembrane potential into cell model
206  Vmath::Vcopy(m_nq, m_nodalTmp[0], 1, m_cellSol[0], 1);
207  }
208  else
209  {
210  // Copy new transmembrane potential into cell model
211  Vmath::Vcopy(m_nq, inarray[0], 1, m_cellSol[0], 1);
212  }
213  // -------------------------
214 
215  NekDouble delta_t = (time - m_lastTime)/m_substeps;
216 
217 
218  // Perform substepping
219  for (unsigned int i = 0; i < m_substeps - 1; ++i)
220  {
221  Update(m_cellSol, m_wsp, time);
222  // Voltage
223  Vmath::Svtvp(m_nq, delta_t, m_wsp[0], 1, m_cellSol[0], 1, m_cellSol[0], 1);
224  // Ion concentrations
225  for (unsigned int j = 0; j < m_concentrations.size(); ++j)
226  {
228  }
229  // Gating variables: Rush-Larsen scheme
230  for (unsigned int j = 0; j < m_gates.size(); ++j)
231  {
232  Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
233  Vmath::Vexp(m_nq, m_gates_tau[j], 1, m_gates_tau[j], 1);
234  Vmath::Vsub(m_nq, m_cellSol[m_gates[j]], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
235  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);
236  }
237  }
238 
239  // Perform final cell model step
240  Update(m_cellSol, m_wsp, time);
241 
242  // Output dV/dt from last step but integrate remaining cell model vars
243  // Transform cell model I_total from nodal to modal space
244  if (m_useNodal)
245  {
246  Array<OneD, NekDouble> tmpCoeffs(max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
247 
248  for (unsigned int k = 0; k < nvar; ++k)
249  {
250  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
251  {
252  int phys_offset = m_field->GetPhys_Offset(i);
253  int coef_offset = m_field->GetCoeff_Offset(i);
254  if (m_field->GetExp(0)->DetShapeType() == LibUtilities::eTriangle)
255  {
256  m_nodalTri->NodalToModal(m_wsp[k]+coef_offset, tmpCoeffs);
257  m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp=outarray[k] + phys_offset);
258  }
259  else
260  {
261  m_nodalTet->NodalToModal(m_wsp[k]+coef_offset, tmpCoeffs);
262  m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp=outarray[k] + phys_offset);
263  }
264  }
265  }
266  }
267  else
268  {
269  Vmath::Vcopy(m_nq, m_wsp[0], 1, outarray[0], 1);
270  }
271 
272  // Ion concentrations
273  for (unsigned int j = 0; j < m_concentrations.size(); ++j)
274  {
276  }
277 
278  // Gating variables: Rush-Larsen scheme
279  for (unsigned int j = 0; j < m_gates.size(); ++j)
280  {
281  Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
282  Vmath::Vexp(m_nq, m_gates_tau[j], 1, m_gates_tau[j], 1);
283  Vmath::Vsub(m_nq, m_cellSol[m_gates[j]], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
284  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);
285  }
286 
287  m_lastTime = time;
288  }
289 
290  Array<OneD, NekDouble> CellModel::GetCellSolutionCoeffs(unsigned int idx)
291  {
292  ASSERTL0(idx < m_nvar, "Index out of range for cell model.");
293 
294  Array<OneD, NekDouble> outarray(m_field->GetNcoeffs());
295  Array<OneD, NekDouble> tmp;
296 
297  if (m_useNodal)
298  {
299  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
300  {
301  int coef_offset = m_field->GetCoeff_Offset(i);
302  if (m_field->GetExp(0)->DetShapeType() == LibUtilities::eTriangle)
303  {
304  m_nodalTri->NodalToModal(m_cellSol[idx]+coef_offset, tmp=outarray+coef_offset);
305  }
306  else
307  {
308  m_nodalTet->NodalToModal(m_cellSol[idx]+coef_offset, tmp=outarray+coef_offset);
309  }
310  }
311  }
312  else
313  {
314  m_field->FwdTrans_IterPerExp(m_cellSol[idx], outarray);
315  }
316 
317  return outarray;
318  }
319 
320  Array<OneD, NekDouble> CellModel::GetCellSolution(unsigned int idx)
321  {
322  return m_cellSol[idx];
323  }
324 
326  {
327  const bool root = (m_session->GetComm()->GetRank() == 0);
328  const std::string fncName = "CellModelInitialConditions";
329  std::string varName;
330  Array<OneD, NekDouble> coeffs(m_field->GetNcoeffs());
331  Array<OneD, NekDouble> tmp;
332 
333  SpatialDomains::MeshGraphSharedPtr vGraph = m_field->GetGraph();
334 
335  if (root)
336  {
337  cout << "Cell model initial conditions: " << endl;
338  }
339 
340  // Load each cell model variable
341  // j=0 and j=1 are for transmembrane or intra/extra-cellular volt.
342  Vmath::Zero(m_nq, m_cellSol[0], 1);
343  for(int j = 1; j < m_cellSol.num_elements(); ++j)
344  {
345  // Get the name of the jth variable
346  varName = GetCellVarName(j);
347 
348  // Check if this variable is defined in a file or analytically
349  if (m_session->GetFunctionType(fncName, varName) ==
351  {
352  const std::string file =
353  m_session->GetFunctionFilename(fncName, varName);
354 
355  if (root)
356  {
357  cout << " - Field " << varName << ": from file "
358  << file << endl;
359  }
360 
361  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
362  std::vector<std::vector<NekDouble> > FieldData;
363 
364  // Read the restart file containing this variable
367  fld->Import(file, FieldDef, FieldData);
368 
369  LibUtilities::FieldMetaDataMap fieldMetaDataMap;
371  fld->ImportFieldMetaData(file,fieldMetaDataMap);
372  iter = fieldMetaDataMap.find("Time");
373  if(iter != fieldMetaDataMap.end())
374  {
375  m_lastTime = boost::lexical_cast<NekDouble>(iter->second);
376  }
377 
378  // Extract the data into the modal coefficients
379  for(int i = 0; i < FieldDef.size(); ++i)
380  {
381  m_field->ExtractDataToCoeffs(FieldDef[i],
382  FieldData[i],
383  varName,
384  coeffs);
385  }
386 
387  // If using nodal cell model then we do a modal->nodal transform
388  // otherwise we do a backward transform onto physical points.
389  if (m_useNodal)
390  {
391  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
392  {
393  int coef_offset = m_field->GetCoeff_Offset(i);
394  if (m_field->GetExp(0)->DetShapeType() ==
396  {
397  m_nodalTri->ModalToNodal(coeffs+coef_offset,
398  tmp=m_cellSol[j]+coef_offset);
399  }
400  else
401  {
402  m_nodalTet->ModalToNodal(coeffs+coef_offset,
403  tmp=m_cellSol[j]+coef_offset);
404  }
405  }
406  }
407  else
408  {
409  m_field->BwdTrans(coeffs, m_cellSol[j]);
410  }
411  }
412  else if (m_session->GetFunctionType(fncName, varName) ==
414  {
416  m_session->GetFunction(fncName, varName);
417 
418  if (root)
419  {
420  cout << " - Field " << varName << ": "
421  << equ->GetExpression() << endl;
422  }
423 
424  const unsigned int nphys = m_field->GetNpoints();
425  Array<OneD, NekDouble> x0(nphys);
426  Array<OneD, NekDouble> x1(nphys);
427  Array<OneD, NekDouble> x2(nphys);
428  m_field->GetCoords(x0,x1,x2);
429 
430  if (m_useNodal)
431  {
432  Array<OneD, NekDouble> phys(nphys);
433  Array<OneD, NekDouble> tmpCoeffs(max(m_nodalTri->GetNcoeffs(), m_nodalTet->GetNcoeffs()));
434 
435  equ->Evaluate(x0, x1, x2, phys);
436  for (unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
437  {
438  int phys_offset = m_field->GetPhys_Offset(i);
439  int coef_offset = m_field->GetCoeff_Offset(i);
440  if (m_field->GetExp(0)->DetShapeType() ==
442  {
443  m_field->GetExp(0)->FwdTrans(
444  phys + phys_offset, tmpCoeffs);
445  m_nodalTri->ModalToNodal(
446  tmpCoeffs,
447  tmp = m_cellSol[j] + coef_offset);
448  }
449  else
450  {
451  m_field->GetExp(0)->FwdTrans(
452  phys + phys_offset, tmpCoeffs);
453  m_nodalTet->ModalToNodal(
454  tmpCoeffs,
455  tmp = m_cellSol[j] + coef_offset);
456  }
457  }
458  }
459  else
460  {
461  equ->Evaluate(x0, x1, x2, m_cellSol[j]);
462  }
463  }
464  }
465  }
466 }