39 #include <lifev/electrophysiology/solver/IonicModels/IonicMinimalModel.hpp> 40 #include <lifev/eta/fem/ETFESpace.hpp> 41 #include <lifev/eta/expression/Integrate.hpp> 42 #include <boost/typeof/typeof.hpp> 50 IonicMinimalModel::IonicMinimalModel() :
81 M_restingConditions.at (0) = 0.0;
82 M_restingConditions.at (1) = 1.0;
83 M_restingConditions.at (2) = 1.0;
84 M_restingConditions.at (3) = 0.021553043080281;
88 IonicMinimalModel::IonicMinimalModel ( Teuchos::ParameterList& parameterList ) :
91 M_uo = parameterList.get (
"uo", 1000.0 );
92 M_uu = parameterList.get (
"uu", 1.61 );
93 M_tetav = parameterList.get (
"tetav", 0.30 );
94 M_tetaw = parameterList.get (
"tetaw", 0.130 );
95 M_tetavm = parameterList.get (
"tetavm", 0.10 );
96 M_tetao = parameterList.get (
"tetao", 0.005 );
97 M_tauv1 = parameterList.get (
"tauv1", 80.0 );
98 M_tauv2 = parameterList.get (
"tauv2", 1.4506 );
99 M_tauvp = parameterList.get (
"tauvp", 1.4506 );
100 M_tauw1 = parameterList.get (
"tauw1", 70.0 );
101 M_tauw2 = parameterList.get (
"tauw2", 8.0 );
102 M_kw = parameterList.get (
"kw", 200.0 );
103 M_uw = parameterList.get (
"uw", 0.016 );
104 M_tauwp = parameterList.get (
"tauwp", 280.0 );
105 M_taufi = parameterList.get (
"taufi", 0.078 );
106 M_tauo1 = parameterList.get (
"tauo1", 410.0 );
107 M_tauo2 = parameterList.get (
"tauo2", 7.0 );
108 M_tauso1 = parameterList.get (
"tauso1", 91.0 );
109 M_tauso2 = parameterList.get (
"tauso2", 0.8 );
110 M_kso = parameterList.get (
"kso", 2.1 );
111 M_uso = parameterList.get (
"uso", 0.6 );
112 M_taus1 = parameterList.get (
"taus1", 2.7342 );
113 M_taus2 = parameterList.get (
"taus2", 4.0 );
114 M_ks = parameterList.get (
"ks", 2.0994 );
115 M_us = parameterList.get (
"us", 0.9087 );
116 M_tausi = parameterList.get (
"tausi", 3.3849 );
117 M_tauwinf = parameterList.get (
"tauwinf", 0.01 );
118 M_winfstar = parameterList.get (
"winfstar", 0.5 );
120 M_restingConditions.at (0) = 0.0;
121 M_restingConditions.at (1) = 1.0;
122 M_restingConditions.at (2) = 1.0;
123 M_restingConditions.at (3) = 0.021553043080281;
126 IonicMinimalModel::IonicMinimalModel (
const IonicMinimalModel& model )
131 M_tetav = model.M_tetav;
132 M_tetaw = model.M_tetaw;
133 M_tetavm = model.M_tetavm;
134 M_tetao = model.M_tetao;
135 M_tauv1 = model.M_tauv1;
136 M_tauv2 = model.M_tauv2;
137 M_tauvp = model.M_tauvp;
138 M_tauw1 = model.M_tauw1;
139 M_tauw2 = model.M_tauw2;
142 M_tauwp = model.M_tauwp;
143 M_taufi = model.M_taufi;
144 M_tauo1 = model.M_tauo1;
145 M_tauo2 = model.M_tauo2;
146 M_tauso1 = model.M_tauso1;
147 M_tauso2 = model.M_tauso2;
150 M_taus1 = model.M_taus1;
151 M_taus2 = model.M_taus2;
154 M_tausi = model.M_tausi;
155 M_tauwinf = model.M_tauwinf;
156 M_winfstar = model.M_winfstar;
158 M_numberOfEquations = model.M_numberOfEquations;
159 M_restingConditions = model.M_restingConditions;
166 IonicMinimalModel& IonicMinimalModel::operator= (
const IonicMinimalModel& model )
170 M_tetav = model.M_tetav;
171 M_tetaw = model.M_tetaw;
172 M_tetavm = model.M_tetavm;
173 M_tetao = model.M_tetao;
174 M_tauv1 = model.M_tauv1;
175 M_tauv2 = model.M_tauv2;
176 M_tauvp = model.M_tauvp;
177 M_tauw1 = model.M_tauw1;
178 M_tauw2 = model.M_tauw2;
181 M_tauwp = model.M_tauwp;
182 M_taufi = model.M_taufi;
183 M_tauo1 = model.M_tauo1;
184 M_tauo2 = model.M_tauo2;
185 M_tauso1 = model.M_tauso1;
186 M_tauso2 = model.M_tauso2;
189 M_taus1 = model.M_taus1;
190 M_taus2 = model.M_taus2;
193 M_tausi = model.M_tausi;
194 M_tauwinf = model.M_tauwinf;
195 M_winfstar = model.M_winfstar;
197 M_numberOfEquations = model.M_numberOfEquations;
198 M_restingConditions = model.M_restingConditions;
208 void IonicMinimalModel::computeGatingRhs (
const std::vector<Real>& v,
209 std::vector<Real>& rhs )
217 Real tauvm = ( 1.0 - Heaviside ( U - M_tetavm ) ) * M_tauv1 + Heaviside ( U - M_tetavm ) * M_tauv2;
218 Real tauwm = M_tauw1 + ( M_tauw2 - M_tauw1 ) * ( 1.0 + std::tanh ( M_kw * ( U - M_uw ) ) ) / 2.0;
219 Real taus = ( 1.0 - Heaviside ( U - M_tetaw ) ) * M_taus1 + Heaviside ( U - M_tetaw ) * M_taus2;
221 Real vinf = Heaviside ( M_tetavm - U );
222 Real winf = ( 1.0 - Heaviside ( U - M_tetao ) ) * ( 1.0 - U / M_tauwinf ) + Heaviside ( U - M_tetao ) * M_winfstar;
224 rhs[0] = ( 1.0 - Heaviside ( U - M_tetav ) ) * ( vinf - V ) / tauvm - Heaviside ( U - M_tetav ) * V / M_tauvp;
225 rhs[1] = ( 1.0 - Heaviside ( U - M_tetaw ) ) * ( winf - W ) / tauwm - Heaviside ( U - M_tetaw ) * W / M_tauwp;
226 rhs[2] = ( ( 1.0 + std::tanh ( M_ks * ( U - M_us ) ) ) / 2.0 - S ) / taus;
230 void IonicMinimalModel::computeRhs (
const std::vector<Real>& v,
231 std::vector<Real>& rhs )
239 Real tauvm = ( 1.0 - Heaviside ( U - M_tetavm ) ) * M_tauv1 + Heaviside ( U - M_tetavm ) * M_tauv2;
240 Real tauwm = M_tauw1 + ( M_tauw2 - M_tauw1 ) * ( 1.0 + std::tanh ( M_kw * ( U - M_uw ) ) ) / 2.0;
241 Real tauso = M_tauso1 + ( M_tauso2 - M_tauso1 ) * ( 1.0 + std::tanh ( M_kso * ( U - M_uso ) ) ) / 2.0;
242 Real taus = ( 1.0 - Heaviside ( U - M_tetaw ) ) * M_taus1 + Heaviside ( U - M_tetaw ) * M_taus2;
243 Real tauo = ( 1.0 - Heaviside ( U - M_tetao ) ) * M_tauo1 + Heaviside ( U - M_tetao ) * M_tauo2;
245 Real vinf = Heaviside ( M_tetavm - U );
246 Real winf = ( 1.0 - Heaviside ( U - M_tetao ) ) * ( 1.0 - U / M_tauwinf ) + Heaviside ( U - M_tetao ) * M_winfstar;
248 Real Jfi = - V * Heaviside ( U - M_tetav ) * ( U - M_tetav ) * ( M_uu - U ) / M_taufi;
249 Real Jso = ( U - M_uo ) * ( 1.0 - Heaviside ( U - M_tetaw ) ) / tauo + Heaviside ( U - M_tetaw ) / tauso;
250 Real Jsi = - Heaviside ( U - M_tetaw ) * W * S / M_tausi;
252 rhs[0] = - ( Jfi + Jso + Jsi );
253 rhs[1] = ( 1.0 - Heaviside ( U - M_tetav ) ) * ( vinf - V ) / tauvm - Heaviside ( U - M_tetav ) * V / M_tauvp;
254 rhs[2] = ( 1.0 - Heaviside ( U - M_tetaw ) ) * ( winf - W ) / tauwm - Heaviside ( U - M_tetaw ) * W / M_tauwp;
255 rhs[3] = ( ( 1.0 + std::tanh ( M_ks * ( U - M_us ) ) ) / 2.0 - S ) / taus;
260 void IonicMinimalModel::computeGatingVariablesWithRushLarsen ( std::vector<Real>& v,
const Real dt )
262 std::cout <<
"\n\nRush Larsen method, for minimal model not implemented!!!\n";
263 std::cout <<
"\n\nI will use Forward Euler!!!\n";
264 std::vector<Real> rhs (3, 0.0);
265 computeGatingRhs ( v, rhs );
272 Real IonicMinimalModel::computeLocalPotentialRhs (
const std::vector<Real>& v )
274 Real dPotential (0.0);
281 Real tauso = M_tauso1 + ( M_tauso2 - M_tauso1 ) * ( 1.0 + std::tanh ( M_kso * ( U - M_uso ) ) ) / 2.0;
282 Real tauo = ( 1.0 - Heaviside ( U - M_tetao ) ) * M_tauo1 + Heaviside ( U - M_tetao ) * M_tauo2;
284 Real Jfi = - V * Heaviside ( U - M_tetav ) * ( U - M_tetav ) * ( M_uu - U ) / M_taufi;
285 Real Jso = ( U - M_uo ) * ( 1.0 - Heaviside ( U - M_tetaw ) ) / tauo + Heaviside ( U - M_tetaw ) / tauso;
286 Real Jsi = - Heaviside ( U - M_tetaw ) * W * S / M_tausi;
288 dPotential = - ( Jfi + Jso + Jsi );
295 void IonicMinimalModel::showMe()
298 std::cout <<
"\n\tHi, I'm the minimal model\n\t See you soon\n\n";
300 std::cout <<
"\nuo " << M_uo ;
301 std::cout <<
"\nuu " << M_uu ;
302 std::cout <<
"\ntetav " << M_tetav ;
303 std::cout <<
"\ntetaw " << M_tetaw ;
304 std::cout <<
"\ntetavm " << M_tetavm ;
305 std::cout <<
"\ntetao " << M_tetao ;
306 std::cout <<
"\ntauv1 " << M_tauv1 ;
307 std::cout <<
"\ntauv2 " << M_tauv2 ;
308 std::cout <<
"\ntauvp " << M_tauvp ;
309 std::cout <<
"\ntauw1 " << M_tauw1 ;
310 std::cout <<
"\ntauw2 " << M_tauw2 ;
311 std::cout <<
"\nkw " << M_kw ;
312 std::cout <<
"\nuw " << M_uw ;
313 std::cout <<
"\ntauwp " << M_tauwp ;
314 std::cout <<
"\ntaufi " << M_taufi ;
315 std::cout <<
"\ntauo1 " << M_tauo1 ;
316 std::cout <<
"\ntauo2 " << M_tauo2 ;
317 std::cout <<
"\ntauso1 " << M_tauso1 ;
318 std::cout <<
"\ntauso2 " << M_tauso2 ;
319 std::cout <<
"\nkso " << M_kso ;
320 std::cout <<
"\nuso " << M_uso ;
321 std::cout <<
"\ntaus1 " << M_taus1 ;
322 std::cout <<
"\ntaus2 " << M_taus2 ;
323 std::cout <<
"\nks " << M_ks ;
324 std::cout <<
"\nus " << M_us ;
325 std::cout <<
"\ntausi " << M_tausi ;
326 std::cout <<
"\ntauwinf " << M_tauwinf ;
327 std::cout <<
"\nwinfstar " << M_winfstar <<
"\n" ;
330 void IonicMinimalModel::computePotentialRhsSVI (
const std::vector<vectorPtr_Type>& v,
331 std::vector<vectorPtr_Type>& rhs,
332 FESpace<mesh_Type, MapEpetra>& uFESpace)
334 typedef ETFESpace<mesh_Type, MapEpetra, 3, 1> ETFESpace_Type;
335 typedef std::shared_ptr<ETFESpace_Type> ETFESpacePtr_Type;
339 if (uFESpace.mapPtr() -> commPtr() -> MyPID() == 0)
341 std::cout <<
"\nMinimal Model: Assembling SVI using ETA!\n";
344 ETFESpacePtr_Type spaceScalar (
345 new ETFESpace_Type (uFESpace.mesh(), &feTetraP1, uFESpace.mapPtr() -> commPtr() ) );
348 using namespace ExpressionAssembly;
350 std::shared_ptr<MMTanhFunctor> tanh (
new MMTanhFunctor);
351 std::shared_ptr<MMHFunctor> H (
new MMHFunctor);
352 std::shared_ptr<MMSV> sv (
new MMSV);
354 BOOST_AUTO_TPL (U, value (spaceScalar, * (v[0] ) ) );
355 BOOST_AUTO_TPL (V, value (spaceScalar, * (v[1] ) ) );
356 BOOST_AUTO_TPL (W, value (spaceScalar, * (v[2] ) ) );
357 BOOST_AUTO_TPL (S, value (spaceScalar, * (v[3] ) ) );
358 BOOST_AUTO_TPL (Iapp, value (spaceScalar, *M_appliedCurrentPtr ) );
361 BOOST_AUTO_TPL (tauso, M_tauso1 + ( M_tauso2 - M_tauso1 ) * ( 1.0 + eval (tanh, M_kso * ( U - M_uso ) ) ) / 2.0);
362 BOOST_AUTO_TPL (tauo, ( 1.0 - eval (H, U - M_tetao ) ) * M_tauo1 + eval (H, U - M_tetao ) * M_tauo2);
363 BOOST_AUTO_TPL (Jfi, value (-1.0) * V * eval (H, U - M_tetav ) * ( U - M_tetav ) * ( M_uu - U ) / M_taufi);
364 BOOST_AUTO_TPL (Jso, ( U - M_uo ) * ( 1.0 - eval (H, U - M_tetaw ) ) / tauo + eval (H, U - M_tetaw ) / tauso);
365 BOOST_AUTO_TPL (Jsi, value (-1.0) * eval (H, U - M_tetaw ) * W * S / M_tausi);
367 integrate ( elements ( uFESpace.mesh() ),
370 ( Iapp - ( Jfi + Jso + Jsi ) ) * phi_i ) >> rhs.at (0);
373 rhs.at (0) -> globalAssemble();