Nektar++
Winslow99.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File Winslow99.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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Winslow 1999 cell model
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iostream>
36 #include <string>
37 //#include <LibUtilities/BasicUtils/Vmath.hpp>
39 
40 namespace Nektar
41 {
42  std::string Winslow99::className
44  "Winslow99",
46  "Winslow 1999 cell model.");
47 
48 
49  /**
50  *
51  */
54  const MultiRegions::ExpListSharedPtr& pField):
55  CellModel(pSession, pField)
56  {
57  m_nvar = 33;
58  m_gates.push_back(1);
59  m_gates.push_back(2);
60  m_gates.push_back(3);
61  m_gates.push_back(4);
62  m_gates.push_back(5);
63  m_gates.push_back(6);
64  m_gates.push_back(7);
65  m_concentrations.push_back(8);
66  m_concentrations.push_back(9);
67  m_concentrations.push_back(10);
68  m_concentrations.push_back(11);
69  m_concentrations.push_back(12);
70  m_concentrations.push_back(13);
71  m_concentrations.push_back(14);
72  m_concentrations.push_back(15);
73  m_concentrations.push_back(16);
74  m_concentrations.push_back(17);
75  m_concentrations.push_back(18);
76  m_concentrations.push_back(19);
77  m_gates.push_back(20);
78  m_concentrations.push_back(21);
79  m_concentrations.push_back(22);
80  m_concentrations.push_back(23);
81  m_concentrations.push_back(24);
82  m_concentrations.push_back(25);
83  m_concentrations.push_back(26);
84  m_concentrations.push_back(27);
85  m_concentrations.push_back(28);
86  m_concentrations.push_back(29);
87  m_concentrations.push_back(30);
88  m_concentrations.push_back(31);
89  m_concentrations.push_back(32);
90  }
91 
92 
93 
94 
95 
97  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
98  Array<OneD, Array<OneD, NekDouble> >&outarray,
99  const NekDouble time)
100  {
101  int nq = m_nq;
102  for (unsigned int i = 0; i < nq; ++i)
103  {
104 
105  // Inputs:
106  // Time units: millisecond
107  NekDouble var_chaste_interface__membrane__V = inarray[0][i];
108  // Units: millivolt; Initial value: -96.1638
109  NekDouble var_chaste_interface__fast_sodium_current_m_gate__m = inarray[1][i];
110  // Units: dimensionless; Initial value: 0.0328302
111  NekDouble var_chaste_interface__fast_sodium_current_h_gate__h = inarray[2][i];
112  // Units: dimensionless; Initial value: 0.988354
113  NekDouble var_chaste_interface__fast_sodium_current_j_gate__j = inarray[3][i];
114  // Units: dimensionless; Initial value: 0.99254
115  NekDouble var_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr = inarray[4][i];
116  // Units: dimensionless; Initial value: 0.51
117  NekDouble var_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks = inarray[5][i];
118  // Units: dimensionless; Initial value: 0.264
119  NekDouble var_chaste_interface__transient_outward_potassium_current_X_to1_gate__X_to1 = inarray[6][i];
120  // Units: dimensionless; Initial value: 2.63
121  NekDouble var_chaste_interface__transient_outward_potassium_current_Y_to1_gate__Y_to1 = inarray[7][i];
122  // Units: dimensionless; Initial value: 0.99
123  NekDouble var_chaste_interface__L_type_Ca_current__O = inarray[8][i];
124  // Units: dimensionless; Initial value: 9.84546e-21
125  NekDouble var_chaste_interface__L_type_Ca_current__O_Ca = inarray[9][i];
126  // Units: dimensionless; Initial value: 0
127  NekDouble var_chaste_interface__L_type_Ca_current__C0 = inarray[10][i];
128  // Units: dimensionless; Initial value: 0.997208
129  NekDouble var_chaste_interface__L_type_Ca_current__C1 = inarray[11][i];
130  // Units: dimensionless; Initial value: 6.38897e-5
131  NekDouble var_chaste_interface__L_type_Ca_current__C2 = inarray[12][i];
132  // Units: dimensionless; Initial value: 1.535e-9
133  NekDouble var_chaste_interface__L_type_Ca_current__C3 = inarray[13][i];
134  // Units: dimensionless; Initial value: 1.63909e-14
135  NekDouble var_chaste_interface__L_type_Ca_current__C4 = inarray[14][i];
136  // Units: dimensionless; Initial value: 6.56337e-20
137  NekDouble var_chaste_interface__L_type_Ca_current__C_Ca0 = inarray[15][i];
138  // Units: dimensionless; Initial value: 0.00272826
139  NekDouble var_chaste_interface__L_type_Ca_current__C_Ca1 = inarray[16][i];
140  // Units: dimensionless; Initial value: 6.99215e-7
141  NekDouble var_chaste_interface__L_type_Ca_current__C_Ca2 = inarray[17][i];
142  // Units: dimensionless; Initial value: 6.71989e-11
143  NekDouble var_chaste_interface__L_type_Ca_current__C_Ca3 = inarray[18][i];
144  // Units: dimensionless; Initial value: 2.87031e-15
145  NekDouble var_chaste_interface__L_type_Ca_current__C_Ca4 = inarray[19][i];
146  // Units: dimensionless; Initial value: 4.59752e-20
147  NekDouble var_chaste_interface__L_type_Ca_current_y_gate__y = inarray[20][i];
148  // Units: dimensionless; Initial value: 0.798
149  NekDouble var_chaste_interface__RyR_channel__P_O1 = inarray[21][i];
150  // Units: dimensionless; Initial value: 0
151  NekDouble var_chaste_interface__RyR_channel__P_O2 = inarray[22][i];
152  // Units: dimensionless; Initial value: 0
153  NekDouble var_chaste_interface__RyR_channel__P_C1 = inarray[23][i];
154  // Units: dimensionless; Initial value: 0.47
155  NekDouble var_chaste_interface__RyR_channel__P_C2 = inarray[24][i];
156  // Units: dimensionless; Initial value: 0.53
157  NekDouble var_chaste_interface__intracellular_Ca_fluxes__HTRPNCa = inarray[25][i];
158  // Units: millimolar; Initial value: 0.98
159  NekDouble var_chaste_interface__intracellular_Ca_fluxes__LTRPNCa = inarray[26][i];
160  // Units: millimolar; Initial value: 0.078
161  NekDouble var_chaste_interface__intracellular_ion_concentrations__Nai = inarray[27][i];
162  // Units: millimolar; Initial value: 10
163  NekDouble var_chaste_interface__intracellular_ion_concentrations__Cai = inarray[28][i];
164  // Units: millimolar; Initial value: 0.00008
165  NekDouble var_chaste_interface__intracellular_ion_concentrations__Ki = inarray[29][i];
166  // Units: millimolar; Initial value: 157.8
167  NekDouble var_chaste_interface__intracellular_ion_concentrations__Ca_ss = inarray[30][i];
168  // Units: millimolar; Initial value: 0.00011
169  NekDouble var_chaste_interface__intracellular_ion_concentrations__Ca_JSR = inarray[31][i];
170  // Units: millimolar; Initial value: 0.257
171  NekDouble var_chaste_interface__intracellular_ion_concentrations__Ca_NSR = inarray[32][i];
172  // Units: millimolar; Initial value: 0.257
173 
174 
175  // Mathematics
176  NekDouble d_dt_chaste_interface__membrane__V;
177  const NekDouble var_membrane__R = 8.314472; // joule_per_mole_kelvin
178  const NekDouble var_membrane__T = 310.0; // kelvin
179  const NekDouble var_membrane__F = 96.4853415; // coulomb_per_millimole
180  const NekDouble var_fast_sodium_current__j = var_chaste_interface__fast_sodium_current_j_gate__j; // dimensionless
181  const NekDouble var_fast_sodium_current__h = var_chaste_interface__fast_sodium_current_h_gate__h; // dimensionless
182  const NekDouble var_fast_sodium_current__g_Na = 12.8; // milliS_per_microF
183  const NekDouble var_fast_sodium_current__m = var_chaste_interface__fast_sodium_current_m_gate__m; // dimensionless
184  const NekDouble var_fast_sodium_current__V = var_chaste_interface__membrane__V; // millivolt
185  const NekDouble var_fast_sodium_current__R = var_membrane__R; // joule_per_mole_kelvin
186  const NekDouble var_fast_sodium_current__F = var_membrane__F; // coulomb_per_millimole
187  const NekDouble var_standard_ionic_concentrations__Nao = 138.0; // millimolar
188  const NekDouble var_fast_sodium_current__Nao = var_standard_ionic_concentrations__Nao; // millimolar
189  const NekDouble var_fast_sodium_current__Nai = var_chaste_interface__intracellular_ion_concentrations__Nai; // millimolar
190  const NekDouble var_fast_sodium_current__T = var_membrane__T; // kelvin
191  const NekDouble var_fast_sodium_current__E_Na = ((var_fast_sodium_current__R * var_fast_sodium_current__T) / var_fast_sodium_current__F) * log(var_fast_sodium_current__Nao / var_fast_sodium_current__Nai); // millivolt
192  const NekDouble var_fast_sodium_current__i_Na = var_fast_sodium_current__g_Na * pow(var_fast_sodium_current__m, 3.0) * var_fast_sodium_current__h * var_fast_sodium_current__j * (var_fast_sodium_current__V - var_fast_sodium_current__E_Na); // microA_per_microF
193  const NekDouble var_L_type_Ca_current__O = var_chaste_interface__L_type_Ca_current__O; // dimensionless
194  const NekDouble var_L_type_Ca_current__F = var_membrane__F; // coulomb_per_millimole
195  const NekDouble var_L_type_Ca_current__P_Ca = 0.0003125; // cm_per_second
196  const NekDouble var_standard_ionic_concentrations__Cao = 2.0; // millimolar
197  const NekDouble var_L_type_Ca_current__Cao = var_standard_ionic_concentrations__Cao; // millimolar
198  const NekDouble var_L_type_Ca_current__V = var_chaste_interface__membrane__V; // millivolt
199  const NekDouble var_L_type_Ca_current__T = var_membrane__T; // kelvin
200  const NekDouble var_L_type_Ca_current__R = var_membrane__R; // joule_per_mole_kelvin
201  const NekDouble var_L_type_Ca_current__i_Ca_max = ((((var_L_type_Ca_current__P_Ca / (1.0 * 1.0)) * 4.0 * var_L_type_Ca_current__V * pow(var_L_type_Ca_current__F, 2.0) * 1000.0) / (var_L_type_Ca_current__R * var_L_type_Ca_current__T)) * ((0.001 * exp((2.0 * var_L_type_Ca_current__V * var_L_type_Ca_current__F) / (var_L_type_Ca_current__R * var_L_type_Ca_current__T))) - (0.341 * var_L_type_Ca_current__Cao))) / (exp((2.0 * var_L_type_Ca_current__V * var_L_type_Ca_current__F) / (var_L_type_Ca_current__R * var_L_type_Ca_current__T)) - 1.0); // microA_per_microF
202  const NekDouble var_L_type_Ca_current__y = var_chaste_interface__L_type_Ca_current_y_gate__y; // dimensionless
203  const NekDouble var_L_type_Ca_current__O_Ca = var_chaste_interface__L_type_Ca_current__O_Ca; // dimensionless
204  const NekDouble var_L_type_Ca_current__i_Ca = var_L_type_Ca_current__i_Ca_max * var_L_type_Ca_current__y * (var_L_type_Ca_current__O + var_L_type_Ca_current__O_Ca); // microA_per_microF
205  const NekDouble var_L_type_Ca_current__P_K = 5.79e-07; // cm_per_second
206  const NekDouble var_L_type_Ca_current__i_Ca_half = -0.265; // microA_per_microF
207  const NekDouble var_L_type_Ca_current__p_prime_k = var_L_type_Ca_current__P_K / (1.0 + (var_L_type_Ca_current__i_Ca_max / var_L_type_Ca_current__i_Ca_half)); // cm_per_second
208  const NekDouble var_standard_ionic_concentrations__Ko = 4.0; // millimolar
209  const NekDouble var_L_type_Ca_current__Ko = var_standard_ionic_concentrations__Ko; // millimolar
210  const NekDouble var_L_type_Ca_current__Ki = var_chaste_interface__intracellular_ion_concentrations__Ki; // millimolar
211  const NekDouble var_L_type_Ca_current__i_Ca_K = ((((var_L_type_Ca_current__p_prime_k / (1.0 * 1.0)) * var_L_type_Ca_current__y * (var_L_type_Ca_current__O + var_L_type_Ca_current__O_Ca) * var_L_type_Ca_current__V * pow(var_L_type_Ca_current__F, 2.0)) / (var_L_type_Ca_current__R * var_L_type_Ca_current__T)) * ((var_L_type_Ca_current__Ki * exp((var_L_type_Ca_current__V * var_L_type_Ca_current__F) / (var_L_type_Ca_current__R * var_L_type_Ca_current__T))) - var_L_type_Ca_current__Ko)) / (exp((var_L_type_Ca_current__V * var_L_type_Ca_current__F) / (var_L_type_Ca_current__R * var_L_type_Ca_current__T)) - 1.0); // microA_per_microF
212  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__Ko = var_standard_ionic_concentrations__Ko; // millimolar
213  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__f_Ko = sqrt(var_rapid_activating_delayed_rectifiyer_K_current__Ko / 4.0); // dimensionless
214  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__Ki = var_chaste_interface__intracellular_ion_concentrations__Ki; // millimolar
215  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__R = var_membrane__R; // joule_per_mole_kelvin
216  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__F = var_membrane__F; // coulomb_per_millimole
217  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__T = var_membrane__T; // kelvin
218  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__E_K = ((var_rapid_activating_delayed_rectifiyer_K_current__R * var_rapid_activating_delayed_rectifiyer_K_current__T) / var_rapid_activating_delayed_rectifiyer_K_current__F) * log(var_rapid_activating_delayed_rectifiyer_K_current__Ko / var_rapid_activating_delayed_rectifiyer_K_current__Ki); // millivolt
219  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__V = var_chaste_interface__membrane__V; // millivolt
220  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__R_V = 1.0 / (1.0 + (1.4945 * exp(0.0446 * var_rapid_activating_delayed_rectifiyer_K_current__V))); // dimensionless
221  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__X_kr = var_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr; // dimensionless
222  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__g_Kr = 0.0034; // milliS_per_microF
223  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__i_Kr = var_rapid_activating_delayed_rectifiyer_K_current__g_Kr * var_rapid_activating_delayed_rectifiyer_K_current__f_Ko * var_rapid_activating_delayed_rectifiyer_K_current__R_V * var_rapid_activating_delayed_rectifiyer_K_current__X_kr * (var_rapid_activating_delayed_rectifiyer_K_current__V - var_rapid_activating_delayed_rectifiyer_K_current__E_K); // microA_per_microF
224  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__g_Ks = 0.0027134; // milliS_per_microF
225  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__Ko = var_standard_ionic_concentrations__Ko; // millimolar
226  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__Nao = var_standard_ionic_concentrations__Nao; // millimolar
227  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__Ki = var_chaste_interface__intracellular_ion_concentrations__Ki; // millimolar
228  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__Nai = var_chaste_interface__intracellular_ion_concentrations__Nai; // millimolar
229  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__R = var_membrane__R; // joule_per_mole_kelvin
230  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__F = var_membrane__F; // coulomb_per_millimole
231  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__T = var_membrane__T; // kelvin
232  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__E_Ks = ((var_slow_activating_delayed_rectifiyer_K_current__R * var_slow_activating_delayed_rectifiyer_K_current__T) / var_slow_activating_delayed_rectifiyer_K_current__F) * log((var_slow_activating_delayed_rectifiyer_K_current__Ko + (0.01833 * var_slow_activating_delayed_rectifiyer_K_current__Nao)) / (var_slow_activating_delayed_rectifiyer_K_current__Ki + (0.01833 * var_slow_activating_delayed_rectifiyer_K_current__Nai))); // millivolt
233  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__V = var_chaste_interface__membrane__V; // millivolt
234  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__X_ks = var_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks; // dimensionless
235  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__i_Ks = var_slow_activating_delayed_rectifiyer_K_current__g_Ks * pow(var_slow_activating_delayed_rectifiyer_K_current__X_ks, 2.0) * (var_slow_activating_delayed_rectifiyer_K_current__V - var_slow_activating_delayed_rectifiyer_K_current__E_Ks); // microA_per_microF
236  const NekDouble var_transient_outward_potassium_current__X_to1 = var_chaste_interface__transient_outward_potassium_current_X_to1_gate__X_to1; // dimensionless
237  const NekDouble var_transient_outward_potassium_current__Y_to1 = var_chaste_interface__transient_outward_potassium_current_Y_to1_gate__Y_to1; // dimensionless
238  const NekDouble var_transient_outward_potassium_current__g_to1 = 0.23815; // milliS_per_microF
239  const NekDouble var_transient_outward_potassium_current__E_K = var_rapid_activating_delayed_rectifiyer_K_current__E_K; // millivolt
240  const NekDouble var_transient_outward_potassium_current__V = var_chaste_interface__membrane__V; // millivolt
241  const NekDouble var_transient_outward_potassium_current__i_to1 = var_transient_outward_potassium_current__g_to1 * var_transient_outward_potassium_current__X_to1 * var_transient_outward_potassium_current__Y_to1 * (var_transient_outward_potassium_current__V - var_transient_outward_potassium_current__E_K); // microA_per_microF
242  const NekDouble var_time_independent_potassium_current__Ko = var_standard_ionic_concentrations__Ko; // millimolar
243  const NekDouble var_time_independent_potassium_current__E_K = var_rapid_activating_delayed_rectifiyer_K_current__E_K; // millivolt
244  const NekDouble var_time_independent_potassium_current__F = var_membrane__F; // coulomb_per_millimole
245  const NekDouble var_time_independent_potassium_current_K1_gate__F = var_time_independent_potassium_current__F; // coulomb_per_millimole
246  const NekDouble var_time_independent_potassium_current__V = var_chaste_interface__membrane__V; // millivolt
247  const NekDouble var_time_independent_potassium_current_K1_gate__V = var_time_independent_potassium_current__V; // millivolt
248  const NekDouble var_time_independent_potassium_current__T = var_membrane__T; // kelvin
249  const NekDouble var_time_independent_potassium_current_K1_gate__T = var_time_independent_potassium_current__T; // kelvin
250  const NekDouble var_time_independent_potassium_current_K1_gate__E_K = var_time_independent_potassium_current__E_K; // millivolt
251  const NekDouble var_time_independent_potassium_current__R = var_membrane__R; // joule_per_mole_kelvin
252  const NekDouble var_time_independent_potassium_current_K1_gate__R = var_time_independent_potassium_current__R; // joule_per_mole_kelvin
253  const NekDouble var_time_independent_potassium_current_K1_gate__K1_infinity_V = 1.0 / (2.0 + exp(((1.5 * var_time_independent_potassium_current_K1_gate__F) / (var_time_independent_potassium_current_K1_gate__R * var_time_independent_potassium_current_K1_gate__T)) * (var_time_independent_potassium_current_K1_gate__V - var_time_independent_potassium_current_K1_gate__E_K))); // dimensionless
254  const NekDouble var_time_independent_potassium_current__K1_infinity_V = var_time_independent_potassium_current_K1_gate__K1_infinity_V; // dimensionless
255  const NekDouble var_time_independent_potassium_current__g_K1 = 2.8; // milliS_per_microF
256  const NekDouble var_time_independent_potassium_current__K_mK1 = 13.0; // millimolar
257  const NekDouble var_time_independent_potassium_current__i_K1 = ((var_time_independent_potassium_current__g_K1 * var_time_independent_potassium_current__K1_infinity_V * var_time_independent_potassium_current__Ko) / (var_time_independent_potassium_current__Ko + var_time_independent_potassium_current__K_mK1)) * (var_time_independent_potassium_current__V - var_time_independent_potassium_current__E_K); // microA_per_microF
258  const NekDouble var_plateau_potassium_current__g_Kp = 0.002216; // milliS_per_microF
259  const NekDouble var_plateau_potassium_current__V = var_chaste_interface__membrane__V; // millivolt
260  const NekDouble var_plateau_potassium_current_Kp_gate__V = var_plateau_potassium_current__V; // millivolt
261  const NekDouble var_plateau_potassium_current_Kp_gate__Kp_V = 1.0 / (1.0 + exp((7.488 - var_plateau_potassium_current_Kp_gate__V) / 5.98)); // dimensionless
262  const NekDouble var_plateau_potassium_current__Kp_V = var_plateau_potassium_current_Kp_gate__Kp_V; // dimensionless
263  const NekDouble var_plateau_potassium_current__E_K = var_rapid_activating_delayed_rectifiyer_K_current__E_K; // millivolt
264  const NekDouble var_plateau_potassium_current__i_Kp = var_plateau_potassium_current__g_Kp * var_plateau_potassium_current__Kp_V * (var_plateau_potassium_current__V - var_plateau_potassium_current__E_K); // microA_per_microF
265  const NekDouble var_Na_Ca_exchanger__Nao = var_standard_ionic_concentrations__Nao; // millimolar
266  const NekDouble var_Na_Ca_exchanger__K_sat = 0.2; // dimensionless
267  const NekDouble var_Na_Ca_exchanger__K_mNa = 87.5; // millimolar
268  const NekDouble var_Na_Ca_exchanger__Nai = var_chaste_interface__intracellular_ion_concentrations__Nai; // millimolar
269  const NekDouble var_Na_Ca_exchanger__K_NaCa = 0.3; // microA_per_microF
270  const NekDouble var_Na_Ca_exchanger__V = var_chaste_interface__membrane__V; // millivolt
271  const NekDouble var_Na_Ca_exchanger__T = var_membrane__T; // kelvin
272  const NekDouble var_Na_Ca_exchanger__Cao = var_standard_ionic_concentrations__Cao; // millimolar
273  const NekDouble var_Na_Ca_exchanger__eta = 0.35; // dimensionless
274  const NekDouble var_Na_Ca_exchanger__K_mCa = 1.38; // millimolar
275  const NekDouble var_Na_Ca_exchanger__F = var_membrane__F; // coulomb_per_millimole
276  const NekDouble var_Na_Ca_exchanger__R = var_membrane__R; // joule_per_mole_kelvin
277  const NekDouble var_Na_Ca_exchanger__Cai = var_chaste_interface__intracellular_ion_concentrations__Cai; // millimolar
278  const NekDouble var_Na_Ca_exchanger__i_NaCa = ((var_Na_Ca_exchanger__K_NaCa * 5000.0) / ((pow(var_Na_Ca_exchanger__K_mNa, 3.0) + pow(var_Na_Ca_exchanger__Nao, 3.0)) * (var_Na_Ca_exchanger__K_mCa + var_Na_Ca_exchanger__Cao) * (1.0 + (var_Na_Ca_exchanger__K_sat * exp(((var_Na_Ca_exchanger__eta - 1.0) * var_Na_Ca_exchanger__V * var_Na_Ca_exchanger__F) / (var_Na_Ca_exchanger__R * var_Na_Ca_exchanger__T)))))) * ((exp((var_Na_Ca_exchanger__eta * var_Na_Ca_exchanger__V * var_Na_Ca_exchanger__F) / (var_Na_Ca_exchanger__R * var_Na_Ca_exchanger__T)) * pow(var_Na_Ca_exchanger__Nai, 3.0) * var_Na_Ca_exchanger__Cao) - (exp(((var_Na_Ca_exchanger__eta - 1.0) * var_Na_Ca_exchanger__V * var_Na_Ca_exchanger__F) / (var_Na_Ca_exchanger__R * var_Na_Ca_exchanger__T)) * pow(var_Na_Ca_exchanger__Nao, 3.0) * var_Na_Ca_exchanger__Cai)); // microA_per_microF
279  const NekDouble var_sodium_potassium_pump__Ko = var_standard_ionic_concentrations__Ko; // millimolar
280  const NekDouble var_sodium_potassium_pump__K_mNai = 10.0; // millimolar
281  const NekDouble var_sodium_potassium_pump__Nai = var_chaste_interface__intracellular_ion_concentrations__Nai; // millimolar
282  const NekDouble var_sodium_potassium_pump__V = var_chaste_interface__membrane__V; // millivolt
283  const NekDouble var_sodium_potassium_pump__F = var_membrane__F; // coulomb_per_millimole
284  const NekDouble var_sodium_potassium_pump__T = var_membrane__T; // kelvin
285  const NekDouble var_sodium_potassium_pump__Nao = var_standard_ionic_concentrations__Nao; // millimolar
286  const NekDouble var_sodium_potassium_pump__sigma = (1.0 / 7.0) * (exp(var_sodium_potassium_pump__Nao / 67.3) - 1.0); // dimensionless
287  const NekDouble var_sodium_potassium_pump__R = var_membrane__R; // joule_per_mole_kelvin
288  const NekDouble var_sodium_potassium_pump__f_NaK = 1.0 / (1.0 + (0.1245 * exp(((-0.1) * var_sodium_potassium_pump__V * var_sodium_potassium_pump__F) / (var_sodium_potassium_pump__R * var_sodium_potassium_pump__T))) + (0.0365 * var_sodium_potassium_pump__sigma * exp(((-var_sodium_potassium_pump__V) * var_sodium_potassium_pump__F) / (var_sodium_potassium_pump__R * var_sodium_potassium_pump__T)))); // dimensionless
289  const NekDouble var_sodium_potassium_pump__I_NaK = 0.693; // microA_per_microF
290  const NekDouble var_sodium_potassium_pump__K_mKo = 1.5; // millimolar
291  const NekDouble var_sodium_potassium_pump__i_NaK = (((var_sodium_potassium_pump__I_NaK * var_sodium_potassium_pump__f_NaK) / (1.0 + pow(var_sodium_potassium_pump__K_mNai / var_sodium_potassium_pump__Nai, 1.5))) * var_sodium_potassium_pump__Ko) / (var_sodium_potassium_pump__Ko + var_sodium_potassium_pump__K_mKo); // microA_per_microF
292  const NekDouble var_sarcolemmal_calcium_pump__I_pCa = 0.05; // microA_per_microF
293  const NekDouble var_sarcolemmal_calcium_pump__Cai = var_chaste_interface__intracellular_ion_concentrations__Cai; // millimolar
294  const NekDouble var_sarcolemmal_calcium_pump__K_mpCa = 5e-05; // millimolar
295  const NekDouble var_sarcolemmal_calcium_pump__i_p_Ca = (var_sarcolemmal_calcium_pump__I_pCa * var_sarcolemmal_calcium_pump__Cai) / (var_sarcolemmal_calcium_pump__K_mpCa + var_sarcolemmal_calcium_pump__Cai); // microA_per_microF
296  const NekDouble var_calcium_background_current__R = var_membrane__R; // joule_per_mole_kelvin
297  const NekDouble var_calcium_background_current__Cai = var_chaste_interface__intracellular_ion_concentrations__Cai; // millimolar
298  const NekDouble var_calcium_background_current__F = var_membrane__F; // coulomb_per_millimole
299  const NekDouble var_calcium_background_current__T = var_membrane__T; // kelvin
300  const NekDouble var_calcium_background_current__Cao = var_standard_ionic_concentrations__Cao; // millimolar
301  const NekDouble var_calcium_background_current__E_Ca = ((var_calcium_background_current__R * var_calcium_background_current__T) / (2.0 * var_calcium_background_current__F)) * log(var_calcium_background_current__Cao / var_calcium_background_current__Cai); // millivolt
302  const NekDouble var_calcium_background_current__g_Cab = 0.0003842; // milliS_per_microF
303  const NekDouble var_calcium_background_current__V = var_chaste_interface__membrane__V; // millivolt
304  const NekDouble var_calcium_background_current__i_Ca_b = var_calcium_background_current__g_Cab * (var_calcium_background_current__V - var_calcium_background_current__E_Ca); // microA_per_microF
305  const NekDouble var_sodium_background_current__g_Nab = 0.0031; // milliS_per_microF
306  const NekDouble var_sodium_background_current__V = var_chaste_interface__membrane__V; // millivolt
307  const NekDouble var_sodium_background_current__E_Na = var_fast_sodium_current__E_Na; // millivolt
308  const NekDouble var_sodium_background_current__i_Na_b = var_sodium_background_current__g_Nab * (var_sodium_background_current__V - var_sodium_background_current__E_Na); // microA_per_microF
309  const NekDouble var_fast_sodium_current_m_gate__V = var_fast_sodium_current__V; // millivolt
310  const NekDouble var_fast_sodium_current_m_gate__beta_m = 80.0 * exp((-var_fast_sodium_current_m_gate__V) / 11.0); // per_second
311  const NekDouble var_fast_sodium_current_m_gate__E0_m = var_fast_sodium_current_m_gate__V + 47.13; // millivolt
312  const NekDouble var_fast_sodium_current_m_gate__alpha_m = (fabs(var_fast_sodium_current_m_gate__E0_m) < 1e-05) ? (1000.0 / (0.1 - (0.005 * var_fast_sodium_current_m_gate__E0_m))) : ((320.0 * var_fast_sodium_current_m_gate__E0_m) / (1.0 - exp((-0.1) * var_fast_sodium_current_m_gate__E0_m))); // per_second
313  const NekDouble var_fast_sodium_current_m_gate__m = var_fast_sodium_current__m; // dimensionless
314  const NekDouble var_fast_sodium_current_m_gate__d_m_d_environment__time = (var_fast_sodium_current_m_gate__V >= (-90.0)) ? ((var_fast_sodium_current_m_gate__alpha_m * (1.0 - var_fast_sodium_current_m_gate__m)) - (var_fast_sodium_current_m_gate__beta_m * var_fast_sodium_current_m_gate__m)) : 0.0; // per_second
315  const NekDouble var_fast_sodium_current__fast_sodium_current_m_gate__d_m_d_environment__time = var_fast_sodium_current_m_gate__d_m_d_environment__time; // per_second
316  const NekDouble var_fast_sodium_current_h_gate__V = var_fast_sodium_current__V; // millivolt
317  const NekDouble var_fast_sodium_current_h_gate__beta_h = (var_fast_sodium_current_h_gate__V < (-40.0)) ? ((3560.0 * exp(0.079 * var_fast_sodium_current_h_gate__V)) + (310000.0 * exp(0.35 * var_fast_sodium_current_h_gate__V))) : (1000.0 / (0.13 * (1.0 + exp((var_fast_sodium_current_h_gate__V + 10.66) / (-11.1))))); // per_second
318  const NekDouble var_fast_sodium_current_h_gate__alpha_h = (var_fast_sodium_current_h_gate__V < (-40.0)) ? (135.0 * exp((80.0 + var_fast_sodium_current_h_gate__V) / (-6.8))) : 0.0; // per_second
319  const NekDouble var_fast_sodium_current_h_gate__h = var_fast_sodium_current__h; // dimensionless
320  const NekDouble var_fast_sodium_current_h_gate__d_h_d_environment__time = (var_fast_sodium_current_h_gate__alpha_h * (1.0 - var_fast_sodium_current_h_gate__h)) - (var_fast_sodium_current_h_gate__beta_h * var_fast_sodium_current_h_gate__h); // per_second
321  const NekDouble var_fast_sodium_current__fast_sodium_current_h_gate__d_h_d_environment__time = var_fast_sodium_current_h_gate__d_h_d_environment__time; // per_second
322  const NekDouble var_fast_sodium_current_j_gate__V = var_fast_sodium_current__V; // millivolt
323  const NekDouble var_fast_sodium_current_j_gate__alpha_j = (var_fast_sodium_current_j_gate__V < (-40.0)) ? ((1000.0 * (-((127140.0 * exp(0.2444 * var_fast_sodium_current_j_gate__V)) + (3.474e-05 * exp((-0.04391) * var_fast_sodium_current_j_gate__V)))) * (var_fast_sodium_current_j_gate__V + 37.78)) / (1.0 + exp(0.311 * (var_fast_sodium_current_j_gate__V + 79.23)))) : 0.0; // per_second
324  const NekDouble var_fast_sodium_current_j_gate__beta_j = (var_fast_sodium_current_j_gate__V < (-40.0)) ? ((121.2 * exp((-0.01052) * var_fast_sodium_current_j_gate__V)) / (1.0 + exp((-0.1378) * (var_fast_sodium_current_j_gate__V + 40.14)))) : ((300.0 * exp((-2.535e-07) * var_fast_sodium_current_j_gate__V)) / (1.0 + exp((-0.1) * (var_fast_sodium_current_j_gate__V + 32.0)))); // per_second
325  const NekDouble var_fast_sodium_current_j_gate__j = var_fast_sodium_current__j; // dimensionless
326  const NekDouble var_fast_sodium_current_j_gate__d_j_d_environment__time = (var_fast_sodium_current_j_gate__alpha_j * (1.0 - var_fast_sodium_current_j_gate__j)) - (var_fast_sodium_current_j_gate__beta_j * var_fast_sodium_current_j_gate__j); // per_second
327  const NekDouble var_fast_sodium_current__fast_sodium_current_j_gate__d_j_d_environment__time = var_fast_sodium_current_j_gate__d_j_d_environment__time; // per_second
328  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__tau_factor = 1.0; // dimensionless
329  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__V = var_rapid_activating_delayed_rectifiyer_K_current__V; // millivolt
330  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__K21 = exp((-7.677) - (0.0128 * var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__V)); // dimensionless
331  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__K12 = exp((-5.495) + (0.1691 * var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__V)); // dimensionless
332  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__tau_X_kr = (0.001 / (var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__K12 + var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__K21)) + (var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__tau_factor * 0.027); // second
333  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr = var_rapid_activating_delayed_rectifiyer_K_current__X_kr; // dimensionless
334  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr_inf = var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__K12 / (var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__K12 + var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__K21); // dimensionless
335  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time = (var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr_inf - var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr) / var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__tau_X_kr; // per_second
336  const NekDouble var_rapid_activating_delayed_rectifiyer_K_current__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time = var_rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time; // per_second
337  const NekDouble var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__V = var_slow_activating_delayed_rectifiyer_K_current__V; // millivolt
338  const NekDouble var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__tau_X_ks = 0.001 / (((7.19e-05 * (var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__V - 10.0)) / (1.0 - exp((-0.148) * (var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__V - 10.0)))) + ((0.000131 * (var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__V - 10.0)) / (exp(0.0687 * (var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__V - 10.0)) - 1.0))); // second
339  const NekDouble var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks_infinity = 1.0 / (1.0 + exp((-(var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__V - 24.7)) / 13.6)); // dimensionless
340  const NekDouble var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks = var_slow_activating_delayed_rectifiyer_K_current__X_ks; // dimensionless
341  const NekDouble var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time = (var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks_infinity - var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks) / var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__tau_X_ks; // per_second
342  const NekDouble var_slow_activating_delayed_rectifiyer_K_current__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time = var_slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time; // per_second
343  const NekDouble var_transient_outward_potassium_current_X_to1_gate__V = var_transient_outward_potassium_current__V; // millivolt
344  const NekDouble var_transient_outward_potassium_current_X_to1_gate__alpha_X_to1 = 45.16 * exp(0.03577 * var_transient_outward_potassium_current_X_to1_gate__V); // per_second
345  const NekDouble var_transient_outward_potassium_current_X_to1_gate__X_to1 = var_transient_outward_potassium_current__X_to1; // dimensionless
346  const NekDouble var_transient_outward_potassium_current_X_to1_gate__beta_X_to1 = 98.9 * exp((-0.06237) * var_transient_outward_potassium_current_X_to1_gate__V); // per_second
347  const NekDouble var_transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time = (var_transient_outward_potassium_current_X_to1_gate__alpha_X_to1 * (1.0 - var_transient_outward_potassium_current_X_to1_gate__X_to1)) - (var_transient_outward_potassium_current_X_to1_gate__beta_X_to1 * var_transient_outward_potassium_current_X_to1_gate__X_to1); // per_second
348  const NekDouble var_transient_outward_potassium_current__transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time = var_transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time; // per_second
349  const NekDouble var_transient_outward_potassium_current_Y_to1_gate__V = var_transient_outward_potassium_current__V; // millivolt
350  const NekDouble var_transient_outward_potassium_current_Y_to1_gate__alpha_Y_to1 = (5.415 * exp((-(var_transient_outward_potassium_current_Y_to1_gate__V + 33.5)) / 5.0)) / (1.0 + (0.051335 * exp((-(var_transient_outward_potassium_current_Y_to1_gate__V + 33.5)) / 5.0))); // per_second
351  const NekDouble var_transient_outward_potassium_current_Y_to1_gate__Y_to1 = var_transient_outward_potassium_current__Y_to1; // dimensionless
352  const NekDouble var_transient_outward_potassium_current_Y_to1_gate__beta_Y_to1 = (5.415 * exp((var_transient_outward_potassium_current_Y_to1_gate__V + 33.5) / 5.0)) / (1.0 + (0.051335 * exp((var_transient_outward_potassium_current_Y_to1_gate__V + 33.5) / 5.0))); // per_second
353  const NekDouble var_transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time = (var_transient_outward_potassium_current_Y_to1_gate__alpha_Y_to1 * (1.0 - var_transient_outward_potassium_current_Y_to1_gate__Y_to1)) - (var_transient_outward_potassium_current_Y_to1_gate__beta_Y_to1 * var_transient_outward_potassium_current_Y_to1_gate__Y_to1); // per_second
354  const NekDouble var_transient_outward_potassium_current__transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time = var_transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time; // per_second
355  const NekDouble var_L_type_Ca_current__alpha = 400.0 * exp((var_L_type_Ca_current__V + 2.0) / 10.0); // per_second
356  const NekDouble var_L_type_Ca_current__beta = 50.0 * exp((-(var_L_type_Ca_current__V + 2.0)) / 13.0); // per_second
357  const NekDouble var_L_type_Ca_current__Ca_ss = var_chaste_interface__intracellular_ion_concentrations__Ca_ss; // millimolar
358  const NekDouble var_L_type_Ca_current__gamma = (103.75 * var_L_type_Ca_current__Ca_ss) / 1.0; // per_second
359  const NekDouble var_L_type_Ca_current__a = 2.0; // dimensionless
360  const NekDouble var_L_type_Ca_current__alpha_a = var_L_type_Ca_current__alpha * var_L_type_Ca_current__a; // per_second
361  const NekDouble var_L_type_Ca_current__b = 2.0; // dimensionless
362  const NekDouble var_L_type_Ca_current__beta_b = var_L_type_Ca_current__beta / var_L_type_Ca_current__b; // per_second
363  const NekDouble var_L_type_Ca_current__g = 2000.0; // per_second
364  const NekDouble var_L_type_Ca_current__f = 300.0; // per_second
365  const NekDouble var_L_type_Ca_current__gprime = 7000.0; // per_second
366  const NekDouble var_L_type_Ca_current__fprime = 7.0; // per_second
367  const NekDouble var_L_type_Ca_current__omega = 10.0; // per_second
368  const NekDouble var_L_type_Ca_current__C0 = var_chaste_interface__L_type_Ca_current__C0; // dimensionless
369  const NekDouble var_L_type_Ca_current__C1 = var_chaste_interface__L_type_Ca_current__C1; // dimensionless
370  const NekDouble var_L_type_Ca_current__C2 = var_chaste_interface__L_type_Ca_current__C2; // dimensionless
371  const NekDouble var_L_type_Ca_current__C3 = var_chaste_interface__L_type_Ca_current__C3; // dimensionless
372  const NekDouble var_L_type_Ca_current__C4 = var_chaste_interface__L_type_Ca_current__C4; // dimensionless
373  const NekDouble var_L_type_Ca_current__C_Ca0 = var_chaste_interface__L_type_Ca_current__C_Ca0; // dimensionless
374  const NekDouble var_L_type_Ca_current__C_Ca1 = var_chaste_interface__L_type_Ca_current__C_Ca1; // dimensionless
375  const NekDouble var_L_type_Ca_current__C_Ca2 = var_chaste_interface__L_type_Ca_current__C_Ca2; // dimensionless
376  const NekDouble var_L_type_Ca_current__C_Ca3 = var_chaste_interface__L_type_Ca_current__C_Ca3; // dimensionless
377  const NekDouble var_L_type_Ca_current__C_Ca4 = var_chaste_interface__L_type_Ca_current__C_Ca4; // dimensionless
378  const NekDouble var_L_type_Ca_current__d_O_d_environment__time = (var_L_type_Ca_current__f * var_L_type_Ca_current__C4) - (var_L_type_Ca_current__g * var_L_type_Ca_current__O); // per_second
379  const NekDouble var_L_type_Ca_current__d_O_Ca_d_environment__time = (var_L_type_Ca_current__fprime * var_L_type_Ca_current__C_Ca4) - (var_L_type_Ca_current__gprime * var_L_type_Ca_current__O_Ca); // per_second
380  const NekDouble var_L_type_Ca_current__d_C0_d_environment__time = ((var_L_type_Ca_current__beta * var_L_type_Ca_current__C1) + (var_L_type_Ca_current__omega * var_L_type_Ca_current__C_Ca0)) - (((4.0 * var_L_type_Ca_current__alpha) + var_L_type_Ca_current__gamma) * var_L_type_Ca_current__C0); // per_second
381  const NekDouble var_L_type_Ca_current__d_C1_d_environment__time = ((4.0 * var_L_type_Ca_current__alpha * var_L_type_Ca_current__C0) + (2.0 * var_L_type_Ca_current__beta * var_L_type_Ca_current__C2) + ((var_L_type_Ca_current__omega / var_L_type_Ca_current__b) * var_L_type_Ca_current__C_Ca1)) - ((var_L_type_Ca_current__beta + (3.0 * var_L_type_Ca_current__alpha) + (var_L_type_Ca_current__gamma * var_L_type_Ca_current__a)) * var_L_type_Ca_current__C1); // per_second
382  const NekDouble var_L_type_Ca_current__d_C2_d_environment__time = ((3.0 * var_L_type_Ca_current__alpha * var_L_type_Ca_current__C1) + (3.0 * var_L_type_Ca_current__beta * var_L_type_Ca_current__C3) + ((var_L_type_Ca_current__omega / pow(var_L_type_Ca_current__b, 2.0)) * var_L_type_Ca_current__C_Ca2)) - (((var_L_type_Ca_current__beta * 2.0) + (2.0 * var_L_type_Ca_current__alpha) + (var_L_type_Ca_current__gamma * pow(var_L_type_Ca_current__a, 2.0))) * var_L_type_Ca_current__C2); // per_second
383  const NekDouble var_L_type_Ca_current__d_C3_d_environment__time = ((2.0 * var_L_type_Ca_current__alpha * var_L_type_Ca_current__C2) + (4.0 * var_L_type_Ca_current__beta * var_L_type_Ca_current__C4) + ((var_L_type_Ca_current__omega / pow(var_L_type_Ca_current__b, 3.0)) * var_L_type_Ca_current__C_Ca3)) - (((var_L_type_Ca_current__beta * 3.0) + var_L_type_Ca_current__alpha + (var_L_type_Ca_current__gamma * pow(var_L_type_Ca_current__a, 3.0))) * var_L_type_Ca_current__C3); // per_second
384  const NekDouble var_L_type_Ca_current__d_C4_d_environment__time = ((var_L_type_Ca_current__alpha * var_L_type_Ca_current__C3) + (var_L_type_Ca_current__g * var_L_type_Ca_current__O) + ((var_L_type_Ca_current__omega / pow(var_L_type_Ca_current__b, 4.0)) * var_L_type_Ca_current__C_Ca4)) - (((var_L_type_Ca_current__beta * 4.0) + var_L_type_Ca_current__f + (var_L_type_Ca_current__gamma * pow(var_L_type_Ca_current__a, 4.0))) * var_L_type_Ca_current__C4); // per_second
385  const NekDouble var_L_type_Ca_current__d_C_Ca0_d_environment__time = ((var_L_type_Ca_current__beta_b * var_L_type_Ca_current__C_Ca1) + (var_L_type_Ca_current__gamma * var_L_type_Ca_current__C0)) - (((4.0 * var_L_type_Ca_current__alpha_a) + var_L_type_Ca_current__omega) * var_L_type_Ca_current__C_Ca0); // per_second
386  const NekDouble var_L_type_Ca_current__d_C_Ca1_d_environment__time = ((4.0 * var_L_type_Ca_current__alpha_a * var_L_type_Ca_current__C_Ca0) + (2.0 * var_L_type_Ca_current__beta_b * var_L_type_Ca_current__C_Ca2) + (var_L_type_Ca_current__gamma * var_L_type_Ca_current__a * var_L_type_Ca_current__C1)) - ((var_L_type_Ca_current__beta_b + (3.0 * var_L_type_Ca_current__alpha_a) + (var_L_type_Ca_current__omega / var_L_type_Ca_current__b)) * var_L_type_Ca_current__C_Ca1); // per_second
387  const NekDouble var_L_type_Ca_current__d_C_Ca2_d_environment__time = ((3.0 * var_L_type_Ca_current__alpha_a * var_L_type_Ca_current__C_Ca1) + (3.0 * var_L_type_Ca_current__beta_b * var_L_type_Ca_current__C_Ca3) + (var_L_type_Ca_current__gamma * pow(var_L_type_Ca_current__a, 2.0) * var_L_type_Ca_current__C2)) - (((var_L_type_Ca_current__beta_b * 2.0) + (2.0 * var_L_type_Ca_current__alpha_a) + (var_L_type_Ca_current__omega / pow(var_L_type_Ca_current__b, 2.0))) * var_L_type_Ca_current__C_Ca2); // per_second
388  const NekDouble var_L_type_Ca_current__d_C_Ca3_d_environment__time = ((2.0 * var_L_type_Ca_current__alpha_a * var_L_type_Ca_current__C_Ca2) + (4.0 * var_L_type_Ca_current__beta_b * var_L_type_Ca_current__C_Ca4) + (var_L_type_Ca_current__gamma * pow(var_L_type_Ca_current__a, 3.0) * var_L_type_Ca_current__C3)) - (((var_L_type_Ca_current__beta_b * 3.0) + var_L_type_Ca_current__alpha_a + (var_L_type_Ca_current__omega / pow(var_L_type_Ca_current__b, 3.0))) * var_L_type_Ca_current__C_Ca3); // per_second
389  const NekDouble var_L_type_Ca_current__d_C_Ca4_d_environment__time = ((var_L_type_Ca_current__alpha_a * var_L_type_Ca_current__C_Ca3) + (var_L_type_Ca_current__gprime * var_L_type_Ca_current__O_Ca) + (var_L_type_Ca_current__gamma * pow(var_L_type_Ca_current__a, 4.0) * var_L_type_Ca_current__C4)) - (((var_L_type_Ca_current__beta_b * 4.0) + var_L_type_Ca_current__fprime + (var_L_type_Ca_current__omega / pow(var_L_type_Ca_current__b, 4.0))) * var_L_type_Ca_current__C_Ca4); // per_second
390  const NekDouble var_L_type_Ca_current_y_gate__y = var_L_type_Ca_current__y; // dimensionless
391  const NekDouble var_L_type_Ca_current_y_gate__V = var_L_type_Ca_current__V; // millivolt
392  const NekDouble var_L_type_Ca_current_y_gate__y_infinity = (0.8 / (1.0 + exp((var_L_type_Ca_current_y_gate__V + 12.5) / 5.0))) + 0.2; // dimensionless
393  const NekDouble var_L_type_Ca_current_y_gate__tau_y = (20.0 + (600.0 / (1.0 + exp((var_L_type_Ca_current_y_gate__V + 20.0) / 9.5)))) / 1000.0; // second
394  const NekDouble var_L_type_Ca_current_y_gate__d_y_d_environment__time = (var_L_type_Ca_current_y_gate__y_infinity - var_L_type_Ca_current_y_gate__y) / var_L_type_Ca_current_y_gate__tau_y; // per_second
395  const NekDouble var_L_type_Ca_current__L_type_Ca_current_y_gate__d_y_d_environment__time = var_L_type_Ca_current_y_gate__d_y_d_environment__time; // per_second
396  const NekDouble var_RyR_channel__P_O2 = var_chaste_interface__RyR_channel__P_O2; // dimensionless
397  const NekDouble var_RyR_channel__Ca_ss = var_chaste_interface__intracellular_ion_concentrations__Ca_ss; // millimolar
398  const NekDouble var_RyR_channel__P_O1 = var_chaste_interface__RyR_channel__P_O1; // dimensionless
399  const NekDouble var_RyR_channel__v1 = 1800.0; // per_second
400  const NekDouble var_RyR_channel__Ca_JSR = var_chaste_interface__intracellular_ion_concentrations__Ca_JSR; // millimolar
401  const NekDouble var_RyR_channel__J_rel = var_RyR_channel__v1 * (var_RyR_channel__P_O1 + var_RyR_channel__P_O2) * (var_RyR_channel__Ca_JSR - var_RyR_channel__Ca_ss); // millimolar_per_second
402  const NekDouble var_RyR_channel__k_a_plus = 1.215e+13; // millimolar4_per_second
403  const NekDouble var_RyR_channel__k_a_minus = 576.0; // per_second
404  const NekDouble var_RyR_channel__k_b_plus = 4050000000.0; // millimolar3_per_second
405  const NekDouble var_RyR_channel__k_b_minus = 1930.0; // per_second
406  const NekDouble var_RyR_channel__k_c_plus = 100.0; // per_second
407  const NekDouble var_RyR_channel__k_c_minus = 0.8; // per_second
408  const NekDouble var_RyR_channel__P_C1 = var_chaste_interface__RyR_channel__P_C1; // dimensionless
409  const NekDouble var_RyR_channel__P_C2 = var_chaste_interface__RyR_channel__P_C2; // dimensionless
410  const NekDouble var_RyR_channel__n = 4.0; // dimensionless
411  const NekDouble var_RyR_channel__m = 3.0; // dimensionless
412  const NekDouble var_RyR_channel__d_P_O1_d_environment__time = ((var_RyR_channel__k_a_plus * pow(var_RyR_channel__Ca_ss, var_RyR_channel__n) * var_RyR_channel__P_C1) - ((var_RyR_channel__k_a_minus * var_RyR_channel__P_O1) + (var_RyR_channel__k_b_plus * pow(var_RyR_channel__Ca_ss, var_RyR_channel__m) * var_RyR_channel__P_O1) + (var_RyR_channel__k_c_plus * var_RyR_channel__P_O1))) + (var_RyR_channel__k_b_minus * var_RyR_channel__P_O2) + (var_RyR_channel__k_c_minus * var_RyR_channel__P_C2); // per_second
413  const NekDouble var_RyR_channel__d_P_O2_d_environment__time = (var_RyR_channel__k_b_plus * pow(var_RyR_channel__Ca_ss, var_RyR_channel__m) * var_RyR_channel__P_O1) - (var_RyR_channel__k_b_minus * var_RyR_channel__P_O2); // per_second
414  const NekDouble var_RyR_channel__d_P_C1_d_environment__time = ((-var_RyR_channel__k_a_plus) * pow(var_RyR_channel__Ca_ss, var_RyR_channel__n) * var_RyR_channel__P_C1) + (var_RyR_channel__k_a_minus * var_RyR_channel__P_O1); // per_second
415  const NekDouble var_RyR_channel__d_P_C2_d_environment__time = (var_RyR_channel__k_c_plus * var_RyR_channel__P_O1) - (var_RyR_channel__k_c_minus * var_RyR_channel__P_C2); // per_second
416  const NekDouble var_SERCA2a_pump__Cai = var_chaste_interface__intracellular_ion_concentrations__Cai; // millimolar
417  const NekDouble var_SERCA2a_pump__N_fb = 1.2; // dimensionless
418  const NekDouble var_SERCA2a_pump__K_fb = 0.000168; // millimolar
419  const NekDouble var_SERCA2a_pump__fb = pow(var_SERCA2a_pump__Cai / var_SERCA2a_pump__K_fb, var_SERCA2a_pump__N_fb); // dimensionless
420  const NekDouble var_SERCA2a_pump__Vmaxf = 0.0813; // millimolar_per_second
421  const NekDouble var_SERCA2a_pump__K_SR = 1.0; // dimensionless
422  const NekDouble var_SERCA2a_pump__Ca_NSR = var_chaste_interface__intracellular_ion_concentrations__Ca_NSR; // millimolar
423  const NekDouble var_SERCA2a_pump__K_rb = 3.29; // millimolar
424  const NekDouble var_SERCA2a_pump__N_rb = 1.0; // dimensionless
425  const NekDouble var_SERCA2a_pump__rb = pow(var_SERCA2a_pump__Ca_NSR / var_SERCA2a_pump__K_rb, var_SERCA2a_pump__N_rb); // dimensionless
426  const NekDouble var_SERCA2a_pump__Vmaxr = 0.318; // millimolar_per_second
427  const NekDouble var_SERCA2a_pump__J_up = (var_SERCA2a_pump__K_SR * ((var_SERCA2a_pump__Vmaxf * var_SERCA2a_pump__fb) - (var_SERCA2a_pump__Vmaxr * var_SERCA2a_pump__rb))) / (1.0 + var_SERCA2a_pump__fb + var_SERCA2a_pump__rb); // millimolar_per_second
428  const NekDouble var_intracellular_Ca_fluxes__Ca_NSR = var_chaste_interface__intracellular_ion_concentrations__Ca_NSR; // millimolar
429  const NekDouble var_intracellular_Ca_fluxes__Ca_JSR = var_chaste_interface__intracellular_ion_concentrations__Ca_JSR; // millimolar
430  const NekDouble var_intracellular_Ca_fluxes__tau_tr = 0.0005747; // second
431  const NekDouble var_intracellular_Ca_fluxes__J_tr = (var_intracellular_Ca_fluxes__Ca_NSR - var_intracellular_Ca_fluxes__Ca_JSR) / var_intracellular_Ca_fluxes__tau_tr; // millimolar_per_second
432  const NekDouble var_intracellular_Ca_fluxes__Cai = var_chaste_interface__intracellular_ion_concentrations__Cai; // millimolar
433  const NekDouble var_intracellular_Ca_fluxes__tau_xfer = 0.0267; // second
434  const NekDouble var_intracellular_Ca_fluxes__Ca_ss = var_chaste_interface__intracellular_ion_concentrations__Ca_ss; // millimolar
435  const NekDouble var_intracellular_Ca_fluxes__J_xfer = (var_intracellular_Ca_fluxes__Ca_ss - var_intracellular_Ca_fluxes__Cai) / var_intracellular_Ca_fluxes__tau_xfer; // millimolar_per_second
436  const NekDouble var_intracellular_Ca_fluxes__k_htrpn_minus = 0.066; // per_second
437  const NekDouble var_intracellular_Ca_fluxes__k_htrpn_plus = 20000.0; // per_millimolar_second
438  const NekDouble var_intracellular_Ca_fluxes__HTRPNCa = var_chaste_interface__intracellular_Ca_fluxes__HTRPNCa; // millimolar
439  const NekDouble var_intracellular_Ca_fluxes__d_HTRPNCa_d_environment__time = (var_intracellular_Ca_fluxes__k_htrpn_plus * var_intracellular_Ca_fluxes__Cai * (1.0 - var_intracellular_Ca_fluxes__HTRPNCa)) - (var_intracellular_Ca_fluxes__k_htrpn_minus * var_intracellular_Ca_fluxes__HTRPNCa); // 'millimole per litre per second'
440  const NekDouble var_intracellular_Ca_fluxes__J_HTRPNCa = var_intracellular_Ca_fluxes__d_HTRPNCa_d_environment__time; // millimolar_per_second
441  const NekDouble var_intracellular_Ca_fluxes__LTRPN_tot = 0.07; // dimensionless
442  const NekDouble var_intracellular_Ca_fluxes__LTRPNCa = var_chaste_interface__intracellular_Ca_fluxes__LTRPNCa; // millimolar
443  const NekDouble var_intracellular_Ca_fluxes__k_ltrpn_minus = 40.0; // per_second
444  const NekDouble var_intracellular_Ca_fluxes__k_ltrpn_plus = 40000.0; // per_millimolar_second
445  const NekDouble var_intracellular_Ca_fluxes__d_LTRPNCa_d_environment__time = (var_intracellular_Ca_fluxes__k_ltrpn_plus * var_intracellular_Ca_fluxes__Cai * (1.0 - var_intracellular_Ca_fluxes__LTRPNCa)) - (var_intracellular_Ca_fluxes__k_ltrpn_minus * var_intracellular_Ca_fluxes__LTRPNCa); // 'millimole per litre per second'
446  const NekDouble var_intracellular_Ca_fluxes__J_LTRPNCa = var_intracellular_Ca_fluxes__d_LTRPNCa_d_environment__time; // millimolar_per_second
447  const NekDouble var_intracellular_Ca_fluxes__HTRPN_tot = 0.14; // dimensionless
448  const NekDouble var_intracellular_Ca_fluxes__J_trpn = (var_intracellular_Ca_fluxes__HTRPN_tot * var_intracellular_Ca_fluxes__J_HTRPNCa) + (var_intracellular_Ca_fluxes__LTRPN_tot * var_intracellular_Ca_fluxes__J_LTRPNCa); // millimolar_per_second
449  const NekDouble var_intracellular_ion_concentrations__Cai = var_chaste_interface__intracellular_ion_concentrations__Cai; // millimolar
450  const NekDouble var_intracellular_ion_concentrations__Ca_ss = var_chaste_interface__intracellular_ion_concentrations__Ca_ss; // millimolar
451  const NekDouble var_intracellular_ion_concentrations__Ca_JSR = var_chaste_interface__intracellular_ion_concentrations__Ca_JSR; // millimolar
452  const NekDouble var_intracellular_ion_concentrations__A_cap = 0.0001534; // cm2
453  const NekDouble var_intracellular_ion_concentrations__V_myo = 2.584e-05; // micro_litre
454  const NekDouble var_intracellular_ion_concentrations__V_JSR = 1.6e-07; // micro_litre
455  const NekDouble var_intracellular_ion_concentrations__V_NSR = 2.1e-06; // micro_litre
456  const NekDouble var_intracellular_ion_concentrations__V_SS = 1.2e-09; // micro_litre
457  const NekDouble var_intracellular_ion_concentrations__K_mCMDN = 0.00238; // millimolar
458  const NekDouble var_intracellular_ion_concentrations__K_mEGTA = 0.00015; // millimolar
459  const NekDouble var_intracellular_ion_concentrations__K_mCSQN = 0.8; // millimolar
460  const NekDouble var_intracellular_ion_concentrations__CMDN_tot = 0.05; // millimolar
461  const NekDouble var_intracellular_ion_concentrations__EGTA_tot = 0.0; // millimolar
462  const NekDouble var_intracellular_ion_concentrations__CSQN_tot = 15.0; // millimolar
463  const NekDouble var_intracellular_ion_concentrations__beta_i = 1.0 / (1.0 + ((var_intracellular_ion_concentrations__CMDN_tot * var_intracellular_ion_concentrations__K_mCMDN) / pow(var_intracellular_ion_concentrations__K_mCMDN + var_intracellular_ion_concentrations__Cai, 2.0)) + ((var_intracellular_ion_concentrations__EGTA_tot * var_intracellular_ion_concentrations__K_mEGTA) / pow(var_intracellular_ion_concentrations__K_mEGTA + var_intracellular_ion_concentrations__Cai, 2.0))); // dimensionless
464  const NekDouble var_intracellular_ion_concentrations__beta_SS = 1.0 / (1.0 + ((var_intracellular_ion_concentrations__CMDN_tot * var_intracellular_ion_concentrations__K_mCMDN) / pow(var_intracellular_ion_concentrations__K_mCMDN + var_intracellular_ion_concentrations__Ca_ss, 2.0)) + ((var_intracellular_ion_concentrations__EGTA_tot * var_intracellular_ion_concentrations__K_mEGTA) / pow(var_intracellular_ion_concentrations__K_mEGTA + var_intracellular_ion_concentrations__Ca_ss, 2.0))); // dimensionless
465  const NekDouble var_intracellular_ion_concentrations__beta_JSR = 1.0 / (1.0 + ((var_intracellular_ion_concentrations__CSQN_tot * var_intracellular_ion_concentrations__K_mCSQN) / pow(var_intracellular_ion_concentrations__K_mCSQN + var_intracellular_ion_concentrations__Ca_JSR, 2.0))); // dimensionless
466  const NekDouble var_intracellular_ion_concentrations__F = var_membrane__F; // coulomb_per_millimole
467  const NekDouble var_intracellular_ion_concentrations__i_Na = var_fast_sodium_current__i_Na; // microA_per_microF
468  const NekDouble var_intracellular_ion_concentrations__i_Ca = var_L_type_Ca_current__i_Ca; // microA_per_microF
469  const NekDouble var_intracellular_ion_concentrations__i_Na_b = var_sodium_background_current__i_Na_b; // microA_per_microF
470  const NekDouble var_intracellular_ion_concentrations__i_NaCa = var_Na_Ca_exchanger__i_NaCa; // microA_per_microF
471  const NekDouble var_intracellular_ion_concentrations__i_NaK = var_sodium_potassium_pump__i_NaK; // microA_per_microF
472  const NekDouble var_intracellular_ion_concentrations__i_Ca_K = var_L_type_Ca_current__i_Ca_K; // microA_per_microF
473  const NekDouble var_intracellular_ion_concentrations__i_Kr = var_rapid_activating_delayed_rectifiyer_K_current__i_Kr; // microA_per_microF
474  const NekDouble var_intracellular_ion_concentrations__i_Ks = var_slow_activating_delayed_rectifiyer_K_current__i_Ks; // microA_per_microF
475  const NekDouble var_intracellular_ion_concentrations__i_K1 = var_time_independent_potassium_current__i_K1; // microA_per_microF
476  const NekDouble var_intracellular_ion_concentrations__i_Kp = var_plateau_potassium_current__i_Kp; // microA_per_microF
477  const NekDouble var_intracellular_ion_concentrations__i_to1 = var_transient_outward_potassium_current__i_to1; // microA_per_microF
478  const NekDouble var_intracellular_ion_concentrations__i_p_Ca = var_sarcolemmal_calcium_pump__i_p_Ca; // microA_per_microF
479  const NekDouble var_intracellular_ion_concentrations__i_Ca_b = var_calcium_background_current__i_Ca_b; // microA_per_microF
480  const NekDouble var_intracellular_ion_concentrations__J_up = var_SERCA2a_pump__J_up; // millimolar_per_second
481  const NekDouble var_intracellular_ion_concentrations__J_rel = var_RyR_channel__J_rel; // millimolar_per_second
482  const NekDouble var_intracellular_ion_concentrations__J_xfer = var_intracellular_Ca_fluxes__J_xfer; // millimolar_per_second
483  const NekDouble var_intracellular_ion_concentrations__J_trpn = var_intracellular_Ca_fluxes__J_trpn; // millimolar_per_second
484  const NekDouble var_intracellular_ion_concentrations__J_tr = var_intracellular_Ca_fluxes__J_tr; // millimolar_per_second
485  const NekDouble var_intracellular_ion_concentrations__d_Nai_d_environment__time = ((-0.0) * (var_intracellular_ion_concentrations__i_Na + var_intracellular_ion_concentrations__i_Na_b + (var_intracellular_ion_concentrations__i_NaCa * 3.0) + (var_intracellular_ion_concentrations__i_NaK * 3.0)) * var_intracellular_ion_concentrations__A_cap * 1.0) / (var_intracellular_ion_concentrations__V_myo * var_intracellular_ion_concentrations__F); // 'millimole per litre per second'
486  const NekDouble var_intracellular_ion_concentrations__d_Cai_d_environment__time = var_intracellular_ion_concentrations__beta_i * ((var_intracellular_ion_concentrations__J_xfer - (var_intracellular_ion_concentrations__J_up + var_intracellular_ion_concentrations__J_trpn)) + ((((2.0 * var_intracellular_ion_concentrations__i_NaCa) - (var_intracellular_ion_concentrations__i_p_Ca + var_intracellular_ion_concentrations__i_Ca_b)) * var_intracellular_ion_concentrations__A_cap * 1.0) / (2.0 * var_intracellular_ion_concentrations__V_myo * var_intracellular_ion_concentrations__F))); // 'millimole per litre per second'
487  const NekDouble var_intracellular_ion_concentrations__d_Ki_d_environment__time = ((-0.0) * (var_intracellular_ion_concentrations__i_Ca_K + var_intracellular_ion_concentrations__i_Kr + var_intracellular_ion_concentrations__i_Ks + var_intracellular_ion_concentrations__i_K1 + var_intracellular_ion_concentrations__i_Kp + var_intracellular_ion_concentrations__i_to1 + (var_intracellular_ion_concentrations__i_NaK * (-2.0))) * var_intracellular_ion_concentrations__A_cap * 1.0) / (var_intracellular_ion_concentrations__V_myo * var_intracellular_ion_concentrations__F); // 'millimole per litre per second'
488  const NekDouble var_intracellular_ion_concentrations__d_Ca_ss_d_environment__time = var_intracellular_ion_concentrations__beta_SS * ((((var_intracellular_ion_concentrations__J_rel * var_intracellular_ion_concentrations__V_JSR) / var_intracellular_ion_concentrations__V_SS) - ((var_intracellular_ion_concentrations__J_xfer * var_intracellular_ion_concentrations__V_myo) / var_intracellular_ion_concentrations__V_SS)) - ((var_intracellular_ion_concentrations__i_Ca * var_intracellular_ion_concentrations__A_cap * 1.0) / (2.0 * var_intracellular_ion_concentrations__V_SS * var_intracellular_ion_concentrations__F))); // 'millimole per litre per second'
489  const NekDouble var_intracellular_ion_concentrations__d_Ca_JSR_d_environment__time = var_intracellular_ion_concentrations__beta_JSR * (var_intracellular_ion_concentrations__J_tr - var_intracellular_ion_concentrations__J_rel); // 'millimole per litre per second'
490  const NekDouble var_intracellular_ion_concentrations__d_Ca_NSR_d_environment__time = ((var_intracellular_ion_concentrations__J_up * var_intracellular_ion_concentrations__V_myo) / var_intracellular_ion_concentrations__V_NSR) - ((var_intracellular_ion_concentrations__J_tr * var_intracellular_ion_concentrations__V_JSR) / var_intracellular_ion_concentrations__V_NSR); // 'millimole per litre per second'
491  const NekDouble var_chaste_interface__fast_sodium_current_m_gate__d_m_d_environment__time_converter = var_fast_sodium_current__fast_sodium_current_m_gate__d_m_d_environment__time; // per_second
492  const NekDouble var_chaste_interface__fast_sodium_current_m_gate__d_m_d_environment__time = 0.001 * var_chaste_interface__fast_sodium_current_m_gate__d_m_d_environment__time_converter; // 'per millisecond'
493  const NekDouble var_chaste_interface__fast_sodium_current_h_gate__d_h_d_environment__time_converter = var_fast_sodium_current__fast_sodium_current_h_gate__d_h_d_environment__time; // per_second
494  const NekDouble var_chaste_interface__fast_sodium_current_h_gate__d_h_d_environment__time = 0.001 * var_chaste_interface__fast_sodium_current_h_gate__d_h_d_environment__time_converter; // 'per millisecond'
495  const NekDouble var_chaste_interface__fast_sodium_current_j_gate__d_j_d_environment__time_converter = var_fast_sodium_current__fast_sodium_current_j_gate__d_j_d_environment__time; // per_second
496  const NekDouble var_chaste_interface__fast_sodium_current_j_gate__d_j_d_environment__time = 0.001 * var_chaste_interface__fast_sodium_current_j_gate__d_j_d_environment__time_converter; // 'per millisecond'
497  const NekDouble var_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time_converter = var_rapid_activating_delayed_rectifiyer_K_current__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time; // per_second
498  const NekDouble var_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time = 0.001 * var_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time_converter; // 'per millisecond'
499  const NekDouble var_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time_converter = var_slow_activating_delayed_rectifiyer_K_current__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time; // per_second
500  const NekDouble var_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time = 0.001 * var_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time_converter; // 'per millisecond'
501  const NekDouble var_chaste_interface__transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time_converter = var_transient_outward_potassium_current__transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time; // per_second
502  const NekDouble var_chaste_interface__transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time = 0.001 * var_chaste_interface__transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time_converter; // 'per millisecond'
503  const NekDouble var_chaste_interface__transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time_converter = var_transient_outward_potassium_current__transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time; // per_second
504  const NekDouble var_chaste_interface__transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time = 0.001 * var_chaste_interface__transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time_converter; // 'per millisecond'
505  const NekDouble var_chaste_interface__L_type_Ca_current__d_O_d_environment__time_converter = var_L_type_Ca_current__d_O_d_environment__time; // per_second
506  const NekDouble var_chaste_interface__L_type_Ca_current__d_O_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_O_d_environment__time_converter; // 'per millisecond'
507  const NekDouble var_chaste_interface__L_type_Ca_current__d_O_Ca_d_environment__time_converter = var_L_type_Ca_current__d_O_Ca_d_environment__time; // per_second
508  const NekDouble var_chaste_interface__L_type_Ca_current__d_O_Ca_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_O_Ca_d_environment__time_converter; // 'per millisecond'
509  const NekDouble var_chaste_interface__L_type_Ca_current__d_C0_d_environment__time_converter = var_L_type_Ca_current__d_C0_d_environment__time; // per_second
510  const NekDouble var_chaste_interface__L_type_Ca_current__d_C0_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C0_d_environment__time_converter; // 'per millisecond'
511  const NekDouble var_chaste_interface__L_type_Ca_current__d_C1_d_environment__time_converter = var_L_type_Ca_current__d_C1_d_environment__time; // per_second
512  const NekDouble var_chaste_interface__L_type_Ca_current__d_C1_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C1_d_environment__time_converter; // 'per millisecond'
513  const NekDouble var_chaste_interface__L_type_Ca_current__d_C2_d_environment__time_converter = var_L_type_Ca_current__d_C2_d_environment__time; // per_second
514  const NekDouble var_chaste_interface__L_type_Ca_current__d_C2_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C2_d_environment__time_converter; // 'per millisecond'
515  const NekDouble var_chaste_interface__L_type_Ca_current__d_C3_d_environment__time_converter = var_L_type_Ca_current__d_C3_d_environment__time; // per_second
516  const NekDouble var_chaste_interface__L_type_Ca_current__d_C3_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C3_d_environment__time_converter; // 'per millisecond'
517  const NekDouble var_chaste_interface__L_type_Ca_current__d_C4_d_environment__time_converter = var_L_type_Ca_current__d_C4_d_environment__time; // per_second
518  const NekDouble var_chaste_interface__L_type_Ca_current__d_C4_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C4_d_environment__time_converter; // 'per millisecond'
519  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca0_d_environment__time_converter = var_L_type_Ca_current__d_C_Ca0_d_environment__time; // per_second
520  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca0_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C_Ca0_d_environment__time_converter; // 'per millisecond'
521  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca1_d_environment__time_converter = var_L_type_Ca_current__d_C_Ca1_d_environment__time; // per_second
522  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca1_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C_Ca1_d_environment__time_converter; // 'per millisecond'
523  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca2_d_environment__time_converter = var_L_type_Ca_current__d_C_Ca2_d_environment__time; // per_second
524  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca2_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C_Ca2_d_environment__time_converter; // 'per millisecond'
525  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca3_d_environment__time_converter = var_L_type_Ca_current__d_C_Ca3_d_environment__time; // per_second
526  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca3_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C_Ca3_d_environment__time_converter; // 'per millisecond'
527  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca4_d_environment__time_converter = var_L_type_Ca_current__d_C_Ca4_d_environment__time; // per_second
528  const NekDouble var_chaste_interface__L_type_Ca_current__d_C_Ca4_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current__d_C_Ca4_d_environment__time_converter; // 'per millisecond'
529  const NekDouble var_chaste_interface__L_type_Ca_current_y_gate__d_y_d_environment__time_converter = var_L_type_Ca_current__L_type_Ca_current_y_gate__d_y_d_environment__time; // per_second
530  const NekDouble var_chaste_interface__L_type_Ca_current_y_gate__d_y_d_environment__time = 0.001 * var_chaste_interface__L_type_Ca_current_y_gate__d_y_d_environment__time_converter; // 'per millisecond'
531  const NekDouble var_chaste_interface__RyR_channel__d_P_O1_d_environment__time_converter = var_RyR_channel__d_P_O1_d_environment__time; // per_second
532  const NekDouble var_chaste_interface__RyR_channel__d_P_O1_d_environment__time = 0.001 * var_chaste_interface__RyR_channel__d_P_O1_d_environment__time_converter; // 'per millisecond'
533  const NekDouble var_chaste_interface__RyR_channel__d_P_O2_d_environment__time_converter = var_RyR_channel__d_P_O2_d_environment__time; // per_second
534  const NekDouble var_chaste_interface__RyR_channel__d_P_O2_d_environment__time = 0.001 * var_chaste_interface__RyR_channel__d_P_O2_d_environment__time_converter; // 'per millisecond'
535  const NekDouble var_chaste_interface__RyR_channel__d_P_C1_d_environment__time_converter = var_RyR_channel__d_P_C1_d_environment__time; // per_second
536  const NekDouble var_chaste_interface__RyR_channel__d_P_C1_d_environment__time = 0.001 * var_chaste_interface__RyR_channel__d_P_C1_d_environment__time_converter; // 'per millisecond'
537  const NekDouble var_chaste_interface__RyR_channel__d_P_C2_d_environment__time_converter = var_RyR_channel__d_P_C2_d_environment__time; // per_second
538  const NekDouble var_chaste_interface__RyR_channel__d_P_C2_d_environment__time = 0.001 * var_chaste_interface__RyR_channel__d_P_C2_d_environment__time_converter; // 'per millisecond'
539  const NekDouble var_chaste_interface__intracellular_Ca_fluxes__d_HTRPNCa_d_environment__time_converter = var_intracellular_Ca_fluxes__d_HTRPNCa_d_environment__time; // ___units_89
540  const NekDouble var_chaste_interface__intracellular_Ca_fluxes__d_HTRPNCa_d_environment__time = 0.001 * var_chaste_interface__intracellular_Ca_fluxes__d_HTRPNCa_d_environment__time_converter; // 'millimolar per millisecond'
541  const NekDouble var_chaste_interface__intracellular_Ca_fluxes__d_LTRPNCa_d_environment__time_converter = var_intracellular_Ca_fluxes__d_LTRPNCa_d_environment__time; // millimole_per_litre_per_second
542  const NekDouble var_chaste_interface__intracellular_Ca_fluxes__d_LTRPNCa_d_environment__time = 0.001 * var_chaste_interface__intracellular_Ca_fluxes__d_LTRPNCa_d_environment__time_converter; // 'millimolar per millisecond'
543  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Nai_d_environment__time_converter = var_intracellular_ion_concentrations__d_Nai_d_environment__time; // millimole_per_litre_per_second
544  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Nai_d_environment__time = 0.001 * var_chaste_interface__intracellular_ion_concentrations__d_Nai_d_environment__time_converter; // 'millimolar per millisecond'
545  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Cai_d_environment__time_converter = var_intracellular_ion_concentrations__d_Cai_d_environment__time; // millimole_per_litre_per_second
546  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Cai_d_environment__time = 0.001 * var_chaste_interface__intracellular_ion_concentrations__d_Cai_d_environment__time_converter; // 'millimolar per millisecond'
547  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ki_d_environment__time_converter = var_intracellular_ion_concentrations__d_Ki_d_environment__time; // millimole_per_litre_per_second
548  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ki_d_environment__time = 0.001 * var_chaste_interface__intracellular_ion_concentrations__d_Ki_d_environment__time_converter; // 'millimolar per millisecond'
549  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ca_ss_d_environment__time_converter = var_intracellular_ion_concentrations__d_Ca_ss_d_environment__time; // millimole_per_litre_per_second
550  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ca_ss_d_environment__time = 0.001 * var_chaste_interface__intracellular_ion_concentrations__d_Ca_ss_d_environment__time_converter; // 'millimolar per millisecond'
551  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ca_JSR_d_environment__time_converter = var_intracellular_ion_concentrations__d_Ca_JSR_d_environment__time; // millimole_per_litre_per_second
552  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ca_JSR_d_environment__time = 0.001 * var_chaste_interface__intracellular_ion_concentrations__d_Ca_JSR_d_environment__time_converter; // 'millimolar per millisecond'
553  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ca_NSR_d_environment__time_converter = var_intracellular_ion_concentrations__d_Ca_NSR_d_environment__time; // millimole_per_litre_per_second
554  const NekDouble var_chaste_interface__intracellular_ion_concentrations__d_Ca_NSR_d_environment__time = 0.001 * var_chaste_interface__intracellular_ion_concentrations__d_Ca_NSR_d_environment__time_converter; // 'millimolar per millisecond'
555  const NekDouble d_dt_chaste_interface__fast_sodium_current_m_gate__m = var_chaste_interface__fast_sodium_current_m_gate__d_m_d_environment__time; // 'per millisecond'
556  const NekDouble d_dt_chaste_interface__fast_sodium_current_h_gate__h = var_chaste_interface__fast_sodium_current_h_gate__d_h_d_environment__time; // 'per millisecond'
557  const NekDouble d_dt_chaste_interface__fast_sodium_current_j_gate__j = var_chaste_interface__fast_sodium_current_j_gate__d_j_d_environment__time; // 'per millisecond'
558  const NekDouble d_dt_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr = var_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__d_X_kr_d_environment__time; // 'per millisecond'
559  const NekDouble d_dt_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks = var_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__d_X_ks_d_environment__time; // 'per millisecond'
560  const NekDouble d_dt_chaste_interface__transient_outward_potassium_current_X_to1_gate__X_to1 = var_chaste_interface__transient_outward_potassium_current_X_to1_gate__d_X_to1_d_environment__time; // 'per millisecond'
561  const NekDouble d_dt_chaste_interface__transient_outward_potassium_current_Y_to1_gate__Y_to1 = var_chaste_interface__transient_outward_potassium_current_Y_to1_gate__d_Y_to1_d_environment__time; // 'per millisecond'
562  const NekDouble d_dt_chaste_interface__L_type_Ca_current__O = var_chaste_interface__L_type_Ca_current__d_O_d_environment__time; // 'per millisecond'
563  const NekDouble d_dt_chaste_interface__L_type_Ca_current__O_Ca = var_chaste_interface__L_type_Ca_current__d_O_Ca_d_environment__time; // 'per millisecond'
564  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C0 = var_chaste_interface__L_type_Ca_current__d_C0_d_environment__time; // 'per millisecond'
565  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C1 = var_chaste_interface__L_type_Ca_current__d_C1_d_environment__time; // 'per millisecond'
566  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C2 = var_chaste_interface__L_type_Ca_current__d_C2_d_environment__time; // 'per millisecond'
567  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C3 = var_chaste_interface__L_type_Ca_current__d_C3_d_environment__time; // 'per millisecond'
568  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C4 = var_chaste_interface__L_type_Ca_current__d_C4_d_environment__time; // 'per millisecond'
569  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C_Ca0 = var_chaste_interface__L_type_Ca_current__d_C_Ca0_d_environment__time; // 'per millisecond'
570  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C_Ca1 = var_chaste_interface__L_type_Ca_current__d_C_Ca1_d_environment__time; // 'per millisecond'
571  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C_Ca2 = var_chaste_interface__L_type_Ca_current__d_C_Ca2_d_environment__time; // 'per millisecond'
572  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C_Ca3 = var_chaste_interface__L_type_Ca_current__d_C_Ca3_d_environment__time; // 'per millisecond'
573  const NekDouble d_dt_chaste_interface__L_type_Ca_current__C_Ca4 = var_chaste_interface__L_type_Ca_current__d_C_Ca4_d_environment__time; // 'per millisecond'
574  const NekDouble d_dt_chaste_interface__L_type_Ca_current_y_gate__y = var_chaste_interface__L_type_Ca_current_y_gate__d_y_d_environment__time; // 'per millisecond'
575  const NekDouble d_dt_chaste_interface__RyR_channel__P_O1 = var_chaste_interface__RyR_channel__d_P_O1_d_environment__time; // 'per millisecond'
576  const NekDouble d_dt_chaste_interface__RyR_channel__P_O2 = var_chaste_interface__RyR_channel__d_P_O2_d_environment__time; // 'per millisecond'
577  const NekDouble d_dt_chaste_interface__RyR_channel__P_C1 = var_chaste_interface__RyR_channel__d_P_C1_d_environment__time; // 'per millisecond'
578  const NekDouble d_dt_chaste_interface__RyR_channel__P_C2 = var_chaste_interface__RyR_channel__d_P_C2_d_environment__time; // 'per millisecond'
579  const NekDouble d_dt_chaste_interface__intracellular_Ca_fluxes__HTRPNCa = var_chaste_interface__intracellular_Ca_fluxes__d_HTRPNCa_d_environment__time; // 'millimole per litre per millisecond'
580  const NekDouble d_dt_chaste_interface__intracellular_Ca_fluxes__LTRPNCa = var_chaste_interface__intracellular_Ca_fluxes__d_LTRPNCa_d_environment__time; // 'millimole per litre per millisecond'
581  const NekDouble d_dt_chaste_interface__intracellular_ion_concentrations__Nai = var_chaste_interface__intracellular_ion_concentrations__d_Nai_d_environment__time; // 'millimole per litre per millisecond'
582  const NekDouble d_dt_chaste_interface__intracellular_ion_concentrations__Cai = var_chaste_interface__intracellular_ion_concentrations__d_Cai_d_environment__time; // 'millimole per litre per millisecond'
583  const NekDouble d_dt_chaste_interface__intracellular_ion_concentrations__Ki = var_chaste_interface__intracellular_ion_concentrations__d_Ki_d_environment__time; // 'millimole per litre per millisecond'
584  const NekDouble d_dt_chaste_interface__intracellular_ion_concentrations__Ca_ss = var_chaste_interface__intracellular_ion_concentrations__d_Ca_ss_d_environment__time; // 'millimole per litre per millisecond'
585  const NekDouble d_dt_chaste_interface__intracellular_ion_concentrations__Ca_JSR = var_chaste_interface__intracellular_ion_concentrations__d_Ca_JSR_d_environment__time; // 'millimole per litre per millisecond'
586  const NekDouble d_dt_chaste_interface__intracellular_ion_concentrations__Ca_NSR = var_chaste_interface__intracellular_ion_concentrations__d_Ca_NSR_d_environment__time; // 'millimole per litre per millisecond'
587 
588  const NekDouble var_membrane__C_sc = 0.001; // microF_per_cm2
589  const NekDouble var_membrane__i_Na = var_fast_sodium_current__i_Na; // microA_per_microF
590  const NekDouble var_membrane__i_Ca = var_L_type_Ca_current__i_Ca; // microA_per_microF
591  const NekDouble var_membrane__i_Ca_K = var_L_type_Ca_current__i_Ca_K; // microA_per_microF
592  const NekDouble var_membrane__i_Kr = var_rapid_activating_delayed_rectifiyer_K_current__i_Kr; // microA_per_microF
593  const NekDouble var_membrane__i_Ks = var_slow_activating_delayed_rectifiyer_K_current__i_Ks; // microA_per_microF
594  const NekDouble var_membrane__i_to1 = var_transient_outward_potassium_current__i_to1; // microA_per_microF
595  const NekDouble var_membrane__i_K1 = var_time_independent_potassium_current__i_K1; // microA_per_microF
596  const NekDouble var_membrane__i_Kp = var_plateau_potassium_current__i_Kp; // microA_per_microF
597  const NekDouble var_membrane__i_NaCa = var_Na_Ca_exchanger__i_NaCa; // microA_per_microF
598  const NekDouble var_membrane__i_NaK = var_sodium_potassium_pump__i_NaK; // microA_per_microF
599  const NekDouble var_membrane__i_p_Ca = var_sarcolemmal_calcium_pump__i_p_Ca; // microA_per_microF
600  const NekDouble var_membrane__i_Ca_b = var_calcium_background_current__i_Ca_b; // microA_per_microF
601  const NekDouble var_membrane__i_Na_b = var_sodium_background_current__i_Na_b; // microA_per_microF
602  const NekDouble var_chaste_interface__membrane__i_Stim = 0.0;
603  const NekDouble var_membrane__i_Stim_converter = var_chaste_interface__membrane__i_Stim; // uA_per_cm2
604  const NekDouble var_membrane__chaste_interface__chaste_membrane_capacitance = 1.0; // uF_per_cm2
605  const NekDouble var_membrane__i_Stim = var_membrane__i_Stim_converter / var_membrane__chaste_interface__chaste_membrane_capacitance; // microA_per_microF
606  const NekDouble var_membrane__d_V_d_environment__time = ((-1.0) * 1.0 * (var_membrane__i_Na + var_membrane__i_Ca + var_membrane__i_Ca_K + var_membrane__i_Kr + var_membrane__i_Ks + var_membrane__i_to1 + var_membrane__i_K1 + var_membrane__i_Kp + var_membrane__i_NaCa + var_membrane__i_NaK + var_membrane__i_p_Ca + var_membrane__i_Na_b + var_membrane__i_Ca_b + var_membrane__i_Stim)) / var_membrane__C_sc; // 'millivolt per second'
607  const NekDouble var_chaste_interface__membrane__d_V_d_environment__time_converter = var_membrane__d_V_d_environment__time; // ___units_1
608  const NekDouble var_chaste_interface__membrane__d_V_d_environment__time = 0.001 * var_chaste_interface__membrane__d_V_d_environment__time_converter; // 'millivolt per millisecond'
609  d_dt_chaste_interface__membrane__V = var_chaste_interface__membrane__d_V_d_environment__time; // 'millivolt per millisecond'
610  outarray[0][i] = d_dt_chaste_interface__membrane__V;
611  outarray[1][i] = d_dt_chaste_interface__fast_sodium_current_m_gate__m;
612  outarray[2][i] = d_dt_chaste_interface__fast_sodium_current_h_gate__h;
613  outarray[3][i] = d_dt_chaste_interface__fast_sodium_current_j_gate__j;
614  outarray[4][i] = d_dt_chaste_interface__rapid_activating_delayed_rectifiyer_K_current_X_kr_gate__X_kr;
615  outarray[5][i] = d_dt_chaste_interface__slow_activating_delayed_rectifiyer_K_current_X_ks_gate__X_ks;
616  outarray[6][i] = d_dt_chaste_interface__transient_outward_potassium_current_X_to1_gate__X_to1;
617  outarray[7][i] = d_dt_chaste_interface__transient_outward_potassium_current_Y_to1_gate__Y_to1;
618  outarray[8][i] = d_dt_chaste_interface__L_type_Ca_current__O;
619  outarray[9][i] = d_dt_chaste_interface__L_type_Ca_current__O_Ca;
620  outarray[10][i] = d_dt_chaste_interface__L_type_Ca_current__C0;
621  outarray[11][i] = d_dt_chaste_interface__L_type_Ca_current__C1;
622  outarray[12][i] = d_dt_chaste_interface__L_type_Ca_current__C2;
623  outarray[13][i] = d_dt_chaste_interface__L_type_Ca_current__C3;
624  outarray[14][i] = d_dt_chaste_interface__L_type_Ca_current__C4;
625  outarray[15][i] = d_dt_chaste_interface__L_type_Ca_current__C_Ca0;
626  outarray[16][i] = d_dt_chaste_interface__L_type_Ca_current__C_Ca1;
627  outarray[17][i] = d_dt_chaste_interface__L_type_Ca_current__C_Ca2;
628  outarray[18][i] = d_dt_chaste_interface__L_type_Ca_current__C_Ca3;
629  outarray[19][i] = d_dt_chaste_interface__L_type_Ca_current__C_Ca4;
630  outarray[20][i] = d_dt_chaste_interface__L_type_Ca_current_y_gate__y;
631  outarray[21][i] = d_dt_chaste_interface__RyR_channel__P_O1;
632  outarray[22][i] = d_dt_chaste_interface__RyR_channel__P_O2;
633  outarray[23][i] = d_dt_chaste_interface__RyR_channel__P_C1;
634  outarray[24][i] = d_dt_chaste_interface__RyR_channel__P_C2;
635  outarray[25][i] = d_dt_chaste_interface__intracellular_Ca_fluxes__HTRPNCa;
636  outarray[26][i] = d_dt_chaste_interface__intracellular_Ca_fluxes__LTRPNCa;
637  outarray[27][i] = d_dt_chaste_interface__intracellular_ion_concentrations__Nai;
638  outarray[28][i] = d_dt_chaste_interface__intracellular_ion_concentrations__Cai;
639  outarray[29][i] = d_dt_chaste_interface__intracellular_ion_concentrations__Ki;
640  outarray[30][i] = d_dt_chaste_interface__intracellular_ion_concentrations__Ca_ss;
641  outarray[31][i] = d_dt_chaste_interface__intracellular_ion_concentrations__Ca_JSR;
642  outarray[32][i] = d_dt_chaste_interface__intracellular_ion_concentrations__Ca_NSR;
643  }
644 
645  }
646 
647  /**
648  *
649  */
651  {
652  SolverUtils::AddSummaryItem(s, "Cell model", "Winslow99");
653  }
654 
655 
656  /**
657  *
658  */
660  {
661  Vmath::Fill(m_nq, -96.1638, m_cellSol[0], 1);
662  Vmath::Fill(m_nq, 0.0328302, m_cellSol[1], 1);
663  Vmath::Fill(m_nq, 0.988354, m_cellSol[2], 1);
664  Vmath::Fill(m_nq, 0.99254, m_cellSol[3], 1);
665  Vmath::Fill(m_nq, 0.51, m_cellSol[4], 1);
666  Vmath::Fill(m_nq, 0.264, m_cellSol[5], 1);
667  Vmath::Fill(m_nq, 2.63, m_cellSol[6], 1);
668  Vmath::Fill(m_nq, 0.99, m_cellSol[7], 1);
669  Vmath::Fill(m_nq, 9.84546e-21, m_cellSol[8], 1);
670  Vmath::Fill(m_nq, 0.0, m_cellSol[9], 1);
671  Vmath::Fill(m_nq, 0.997208, m_cellSol[10], 1);
672  Vmath::Fill(m_nq, 6.38897e-5, m_cellSol[11], 1);
673  Vmath::Fill(m_nq, 1.535e-9, m_cellSol[12], 1);
674  Vmath::Fill(m_nq, 1.63909e-14, m_cellSol[13], 1);
675  Vmath::Fill(m_nq, 6.56337e-20, m_cellSol[14], 1);
676  Vmath::Fill(m_nq, 0.00272826, m_cellSol[15], 1);
677  Vmath::Fill(m_nq, 6.99215e-7, m_cellSol[16], 1);
678  Vmath::Fill(m_nq, 6.71989e-11, m_cellSol[17], 1);
679  Vmath::Fill(m_nq, 2.87031e-15, m_cellSol[18], 1);
680  Vmath::Fill(m_nq, 4.59752e-20, m_cellSol[19], 1);
681  Vmath::Fill(m_nq, 0.798, m_cellSol[20], 1);
682  Vmath::Fill(m_nq, 0.0, m_cellSol[21], 1);
683  Vmath::Fill(m_nq, 0.0, m_cellSol[22], 1);
684  Vmath::Fill(m_nq, 0.47, m_cellSol[23], 1);
685  Vmath::Fill(m_nq, 0.53, m_cellSol[24], 1);
686  Vmath::Fill(m_nq, 0.98, m_cellSol[25], 1);
687  Vmath::Fill(m_nq, 0.078, m_cellSol[26], 1);
688  Vmath::Fill(m_nq, 10.0, m_cellSol[27], 1);
689  Vmath::Fill(m_nq, 0.00008, m_cellSol[28], 1);
690  Vmath::Fill(m_nq, 157.8, m_cellSol[29], 1);
691  Vmath::Fill(m_nq, 0.00011, m_cellSol[30], 1);
692  Vmath::Fill(m_nq, 0.257, m_cellSol[31], 1);
693  Vmath::Fill(m_nq, 0.257, m_cellSol[32], 1);
694  }
695 }
696 
Cell model base class.
Definition: CellModel.h:65
Array< OneD, Array< OneD, NekDouble > > m_cellSol
Cell model solution variables.
Definition: CellModel.h:125
int m_nq
Number of physical points.
Definition: CellModel.h:116
std::vector< int > m_concentrations
Indices of cell model variables which are concentrations.
Definition: CellModel.h:138
std::vector< int > m_gates
Indices of cell model variables which are gates.
Definition: CellModel.h:140
int m_nvar
Number of variables in cell model (inc. transmembrane voltage)
Definition: CellModel.h:118
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
virtual void v_Update(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)$.
Definition: Winslow99.cpp:96
virtual void v_GenerateSummary(SummaryList &s)
Prints a summary of the model parameters.
Definition: Winslow99.cpp:650
static CellModelSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const MultiRegions::ExpListSharedPtr &pField)
Creates an instance of this class.
Definition: Winslow99.h:47
static std::string className
Name of class.
Definition: Winslow99.h:55
virtual void v_SetInitialConditions()
Set initial conditions for cell model.
Definition: Winslow99.cpp:659
Winslow99(const LibUtilities::SessionReaderSharedPtr &pSession, const MultiRegions::ExpListSharedPtr &pField)
Constructor.
Definition: Winslow99.cpp:52
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
CellModelFactory & GetCellModelFactory()
Definition: CellModel.cpp:46
double NekDouble
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:278
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267