luo_rudy_1991Opt.cpp

Go to the documentation of this file.
00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013 #include "luo_rudy_1991Opt.hpp"
00014 #include <cmath>
00015 #include <cassert>
00016 #include "Exception.hpp"
00017 #include "OdeSystemInformation.hpp"
00018 
00019 class Cellluo_rudy_1991FromCellMLOpt_LookupTables
00020 {
00021 public:
00022     static Cellluo_rudy_1991FromCellMLOpt_LookupTables* Instance()
00023     {
00024         if (mpInstance == NULL)
00025         {
00026             mpInstance = new Cellluo_rudy_1991FromCellMLOpt_LookupTables;
00027         }
00028         return mpInstance;
00029     }
00030     
00031     // Row lookup methods
00032     // using linear interpolation
00033     double* _lookup_0_row(unsigned i, double factor)
00034     {
00035         for (unsigned j=0; j<16; j++)
00036         {
00037             double y1 = _lookup_table_0[i][j];
00038             double y2 = _lookup_table_0[i+1][j];
00039             _lookup_table_0_row[j] = y1 + (y2-y1)*factor;
00040         }
00041         return _lookup_table_0_row;
00042     }
00043     
00044     
00045 protected:
00046     Cellluo_rudy_1991FromCellMLOpt_LookupTables(const Cellluo_rudy_1991FromCellMLOpt_LookupTables&);
00047     Cellluo_rudy_1991FromCellMLOpt_LookupTables& operator= (const Cellluo_rudy_1991FromCellMLOpt_LookupTables&);
00048     Cellluo_rudy_1991FromCellMLOpt_LookupTables()
00049     {
00050         assert(mpInstance == NULL);
00051         for (int i=0 ; i<16001; i++)
00052         {
00053             double var_membrane__V = -100.0001 + i*0.01;
00054             _lookup_table_0[i][0] = (0.32 * (var_membrane__V + 47.13)) / (1.0 - exp( -0.1 * (var_membrane__V + 47.13)));
00055         }
00056         
00057         for (int i=0 ; i<16001; i++)
00058         {
00059             double var_membrane__V = -100.0001 + i*0.01;
00060             _lookup_table_0[i][1] = 0.08 * exp((-var_membrane__V) * 0.0909090909091);
00061         }
00062         
00063         for (int i=0 ; i<16001; i++)
00064         {
00065             double var_membrane__V = -100.0001 + i*0.01;
00066             _lookup_table_0[i][2] = (var_membrane__V <  -40.0) ? (0.135 * exp((80.0 + var_membrane__V) *  -0.147058823529)) : 0.0;
00067         }
00068         
00069         for (int i=0 ; i<16001; i++)
00070         {
00071             double var_membrane__V = -100.0001 + i*0.01;
00072             _lookup_table_0[i][3] = (var_membrane__V <  -40.0) ? ((3.56 * exp(0.079 * var_membrane__V)) + (310000.0 * exp(0.35 * var_membrane__V))) : (1.0 / (0.13 * (1.0 + exp((var_membrane__V + 10.66) *  -0.0900900900901))));
00073         }
00074         
00075         for (int i=0 ; i<16001; i++)
00076         {
00077             double var_membrane__V = -100.0001 + i*0.01;
00078             _lookup_table_0[i][4] = (var_membrane__V <  -40.0) ? (((( -127140.0 * exp(0.2444 * var_membrane__V)) - (3.474e-05 * exp( -0.04391 * var_membrane__V))) * (var_membrane__V + 37.78)) / (1.0 + exp(0.311 * (var_membrane__V + 79.23)))) : 0.0;
00079         }
00080         
00081         for (int i=0 ; i<16001; i++)
00082         {
00083             double var_membrane__V = -100.0001 + i*0.01;
00084             _lookup_table_0[i][5] = (var_membrane__V <  -40.0) ? ((0.1212 * exp( -0.01052 * var_membrane__V)) / (1.0 + exp( -0.1378 * (var_membrane__V + 40.14)))) : ((0.3 * exp( -2.535e-07 * var_membrane__V)) / (1.0 + exp( -0.1 * (var_membrane__V + 32.0))));
00085         }
00086         
00087         for (int i=0 ; i<16001; i++)
00088         {
00089             double var_membrane__V = -100.0001 + i*0.01;
00090             _lookup_table_0[i][6] = (0.095 * exp( -0.01 * (var_membrane__V - 5.0))) / (1.0 + exp( -0.072 * (var_membrane__V - 5.0)));
00091         }
00092         
00093         for (int i=0 ; i<16001; i++)
00094         {
00095             double var_membrane__V = -100.0001 + i*0.01;
00096             _lookup_table_0[i][7] = (0.07 * exp( -0.017 * (var_membrane__V + 44.0))) / (1.0 + exp(0.05 * (var_membrane__V + 44.0)));
00097         }
00098         
00099         for (int i=0 ; i<16001; i++)
00100         {
00101             double var_membrane__V = -100.0001 + i*0.01;
00102             _lookup_table_0[i][8] = (0.012 * exp( -0.008 * (var_membrane__V + 28.0))) / (1.0 + exp(0.15 * (var_membrane__V + 28.0)));
00103         }
00104         
00105         for (int i=0 ; i<16001; i++)
00106         {
00107             double var_membrane__V = -100.0001 + i*0.01;
00108             _lookup_table_0[i][9] = (0.0065 * exp( -0.02 * (var_membrane__V + 30.0))) / (1.0 + exp( -0.2 * (var_membrane__V + 30.0)));
00109         }
00110         
00111         for (int i=0 ; i<16001; i++)
00112         {
00113             double var_membrane__V = -100.0001 + i*0.01;
00114             _lookup_table_0[i][10] = (var_membrane__V >  -100.0) ? ((2.837 * (exp(0.04 * (var_membrane__V + 77.0)) - 1.0)) / ((var_membrane__V + 77.0) * exp(0.04 * (var_membrane__V + 35.0)))) : 1.0;
00115         }
00116         
00117         for (int i=0 ; i<16001; i++)
00118         {
00119             double var_membrane__V = -100.0001 + i*0.01;
00120             _lookup_table_0[i][11] = (0.0005 * exp(0.083 * (var_membrane__V + 50.0))) / (1.0 + exp(0.057 * (var_membrane__V + 50.0)));
00121         }
00122         
00123         for (int i=0 ; i<16001; i++)
00124         {
00125             double var_membrane__V = -100.0001 + i*0.01;
00126             _lookup_table_0[i][12] = (0.0013 * exp( -0.06 * (var_membrane__V + 20.0))) / (1.0 + exp( -0.04 * (var_membrane__V + 20.0)));
00127         }
00128         
00129         for (int i=0 ; i<16001; i++)
00130         {
00131             double var_membrane__V = -100.0001 + i*0.01;
00132             _lookup_table_0[i][13] = ((0.49124 * exp(0.08032 * ((var_membrane__V + 5.476) -  -87.8929017138))) + (1.0 * exp(0.06175 * (var_membrane__V - 506.417098286)))) / (1.0 + exp( -0.5143 * ((var_membrane__V -  -87.8929017138) + 4.753)));
00133         }
00134         
00135         for (int i=0 ; i<16001; i++)
00136         {
00137             double var_membrane__V = -100.0001 + i*0.01;
00138             _lookup_table_0[i][14] = 1.02 / (1.0 + exp(0.2385 * ((var_membrane__V -  -87.8929017138) - 59.215)));
00139         }
00140         
00141         for (int i=0 ; i<16001; i++)
00142         {
00143             double var_membrane__V = -100.0001 + i*0.01;
00144             _lookup_table_0[i][15] = 0.0183 * (1.0 / (1.0 + exp((7.488 - var_membrane__V) * 0.167224080268))) * (var_membrane__V -  -87.8929017138);
00145         }
00146         
00147     }
00148     
00149 private:
00151     static Cellluo_rudy_1991FromCellMLOpt_LookupTables *mpInstance;
00152 
00153     // Row lookup methods memory
00154     double _lookup_table_0_row[16];
00155     
00156     // Lookup tables
00157     double _lookup_table_0[16001][16];
00158     
00159 };
00160 
00161 Cellluo_rudy_1991FromCellMLOpt_LookupTables* Cellluo_rudy_1991FromCellMLOpt_LookupTables::mpInstance = NULL;
00162 
00163     double Cellluo_rudy_1991FromCellMLOpt::Get_membrane__I_stim()
00164     {
00165         return var_membrane__I_stim;
00166     }
00167     
00168     double Cellluo_rudy_1991FromCellMLOpt::Get_membrane__i_Na()
00169     {
00170         return var_membrane__i_Na;
00171     }
00172     
00173     double Cellluo_rudy_1991FromCellMLOpt::Get_membrane__i_si()
00174     {
00175         return var_membrane__i_si;
00176     }
00177     
00178     double Cellluo_rudy_1991FromCellMLOpt::Get_membrane__i_K()
00179     {
00180         return var_membrane__i_K;
00181     }
00182     
00183     double Cellluo_rudy_1991FromCellMLOpt::Get_membrane__i_K1()
00184     {
00185         return var_membrane__i_K1;
00186     }
00187     
00188     double Cellluo_rudy_1991FromCellMLOpt::Get_membrane__i_Kp()
00189     {
00190         return var_membrane__i_Kp;
00191     }
00192     
00193     double Cellluo_rudy_1991FromCellMLOpt::Get_membrane__i_b()
00194     {
00195         return var_membrane__i_b;
00196     }
00197     
00198     Cellluo_rudy_1991FromCellMLOpt::Cellluo_rudy_1991FromCellMLOpt(boost::shared_ptr<AbstractIvpOdeSolver> pSolver, boost::shared_ptr<AbstractStimulusFunction> pIntracellularStimulus)
00199         : AbstractCardiacCell(
00200                 pSolver,
00201                 8,
00202                 0,
00203                 pIntracellularStimulus)
00204     {
00205         // Time units: millisecond
00206         // 
00207         mpSystemInfo = OdeSystemInformation<Cellluo_rudy_1991FromCellMLOpt>::Instance();
00208         Init();
00209 
00210     }
00211     
00212     Cellluo_rudy_1991FromCellMLOpt::~Cellluo_rudy_1991FromCellMLOpt()
00213     {
00214     }
00215     
00216     void Cellluo_rudy_1991FromCellMLOpt::VerifyStateVariables()
00217     {}
00218 
00219     double Cellluo_rudy_1991FromCellMLOpt::GetIIonic()
00220     {
00221         std::vector<double>& rY = rGetStateVariables();
00222         double var_membrane__V = rY[0];
00223         // Units: millivolt; Initial value: -83.853
00224         double var_fast_sodium_current_m_gate__m = rY[1];
00225         // Units: dimensionless; Initial value: 0.00187018
00226         double var_fast_sodium_current_h_gate__h = rY[2];
00227         // Units: dimensionless; Initial value: 0.9804713
00228         double var_fast_sodium_current_j_gate__j = rY[3];
00229         // Units: dimensionless; Initial value: 0.98767124
00230         double var_slow_inward_current_d_gate__d = rY[4];
00231         // Units: dimensionless; Initial value: 0.00316354
00232         double var_slow_inward_current_f_gate__f = rY[5];
00233         // Units: dimensionless; Initial value: 0.99427859
00234         double var_time_dependent_potassium_current_X_gate__X = rY[6];
00235         // Units: dimensionless; Initial value: 0.16647703
00236         double var_intracellular_calcium_concentration__Cai = rY[7];
00237         // Units: millimolar; Initial value: 0.0002
00238         
00239         // Lookup table indexing
00240         if (var_membrane__V>59.9999 || var_membrane__V<-100.0001)
00241         {
00242 #define COVERAGE_IGNORE
00243             EXCEPTION(DumpState("V outside lookup table range", rY));
00244 #undef COVERAGE_IGNORE
00245         }
00246         
00247         double _offset_0 = var_membrane__V - -100.0001;
00248         double _offset_0_over_table_step = _offset_0 * 100.0;
00249         unsigned _table_index_0 = (unsigned)(_offset_0_over_table_step);
00250         double _factor_0 = _offset_0_over_table_step - _table_index_0;
00251         double* _lt_0_row = Cellluo_rudy_1991FromCellMLOpt_LookupTables::Instance()->_lookup_0_row(_table_index_0, _factor_0);
00252         
00253         var_membrane__i_Na = 23.0 * pow(var_fast_sodium_current_m_gate__m, 3.0) * var_fast_sodium_current_h_gate__h * var_fast_sodium_current_j_gate__j * (var_membrane__V - 54.7944639351);
00254         const double var_slow_inward_current__i_si = 0.09 * var_slow_inward_current_d_gate__d * var_slow_inward_current_f_gate__f * (var_membrane__V - (7.7 - (13.0287 * log(var_intracellular_calcium_concentration__Cai * 1.0))));
00255         var_membrane__i_si = var_slow_inward_current__i_si;
00256         var_membrane__i_K = 0.282 * var_time_dependent_potassium_current_X_gate__X * _lt_0_row[10] * (var_membrane__V -  -77.5675843853);
00257         const double var_time_independent_potassium_current_K1_gate__alpha_K1 = _lt_0_row[14];
00258         var_membrane__i_K1 = 0.6047 * (var_time_independent_potassium_current_K1_gate__alpha_K1 / (var_time_independent_potassium_current_K1_gate__alpha_K1 + _lt_0_row[13])) * (var_membrane__V -  -87.8929017138);
00259         var_membrane__i_Kp = _lt_0_row[15];
00260         var_membrane__i_b = 0.03921 * (var_membrane__V -  -59.87);
00261         
00262         return (var_membrane__i_Na+var_membrane__i_si+var_membrane__i_K+var_membrane__i_K1+var_membrane__i_Kp+var_membrane__i_b);
00263     }
00264     
00265     void Cellluo_rudy_1991FromCellMLOpt::EvaluateYDerivatives(double var_environment__time, const std::vector<double>& rY, std::vector<double>& rDY)
00266     {
00267         // Inputs:
00268         // Time units: millisecond
00269         var_environment__time *= 1.0;
00270         double var_membrane__V = rY[0];
00271         // Units: millivolt; Initial value: -83.853
00272         double var_fast_sodium_current_m_gate__m = rY[1];
00273         // Units: dimensionless; Initial value: 0.00187018
00274         double var_fast_sodium_current_h_gate__h = rY[2];
00275         // Units: dimensionless; Initial value: 0.9804713
00276         double var_fast_sodium_current_j_gate__j = rY[3];
00277         // Units: dimensionless; Initial value: 0.98767124
00278         double var_slow_inward_current_d_gate__d = rY[4];
00279         // Units: dimensionless; Initial value: 0.00316354
00280         double var_slow_inward_current_f_gate__f = rY[5];
00281         // Units: dimensionless; Initial value: 0.99427859
00282         double var_time_dependent_potassium_current_X_gate__X = rY[6];
00283         // Units: dimensionless; Initial value: 0.16647703
00284         double var_intracellular_calcium_concentration__Cai = rY[7];
00285         // Units: millimolar; Initial value: 0.0002
00286         
00287         
00288         // Lookup table indexing
00289         if (var_membrane__V>59.9999 || var_membrane__V<-100.0001)
00290         {
00291 #define COVERAGE_IGNORE
00292             EXCEPTION(DumpState("V outside lookup table range", rY));
00293 #undef COVERAGE_IGNORE
00294         }
00295         
00296         double _offset_0 = var_membrane__V - -100.0001;
00297         double _offset_0_over_table_step = _offset_0 * 100.0;
00298         unsigned _table_index_0 = (unsigned)(_offset_0_over_table_step);
00299         double _factor_0 = _offset_0_over_table_step - _table_index_0;
00300         double* _lt_0_row = Cellluo_rudy_1991FromCellMLOpt_LookupTables::Instance()->_lookup_0_row(_table_index_0, _factor_0);
00301         
00302         // Mathematics
00303         var_membrane__I_stim = GetStimulus((1.0/1)*var_environment__time);
00304         var_membrane__i_Na = 23.0 * pow(var_fast_sodium_current_m_gate__m, 3.0) * var_fast_sodium_current_h_gate__h * var_fast_sodium_current_j_gate__j * (var_membrane__V - 54.7944639351);
00305         const double var_slow_inward_current__i_si = 0.09 * var_slow_inward_current_d_gate__d * var_slow_inward_current_f_gate__f * (var_membrane__V - (7.7 - (13.0287 * log(var_intracellular_calcium_concentration__Cai * 1.0))));
00306         var_membrane__i_si = var_slow_inward_current__i_si;
00307         var_membrane__i_K = 0.282 * var_time_dependent_potassium_current_X_gate__X * _lt_0_row[10] * (var_membrane__V -  -77.5675843853);
00308         const double var_time_independent_potassium_current_K1_gate__alpha_K1 = _lt_0_row[14];
00309         var_membrane__i_K1 = 0.6047 * (var_time_independent_potassium_current_K1_gate__alpha_K1 / (var_time_independent_potassium_current_K1_gate__alpha_K1 + _lt_0_row[13])) * (var_membrane__V -  -87.8929017138);
00310         var_membrane__i_Kp = _lt_0_row[15];
00311         var_membrane__i_b = 0.03921 * (var_membrane__V -  -59.87);
00312         
00313         double d_dt_membrane__V;
00314         if (mSetVoltageDerivativeToZero)
00315         {
00316             d_dt_membrane__V = 0.0;
00317         }
00318         else
00319         {
00320             d_dt_membrane__V =  -1.0 * (var_membrane__I_stim + var_membrane__i_Na + var_membrane__i_si + var_membrane__i_K + var_membrane__i_K1 + var_membrane__i_Kp + var_membrane__i_b);
00321         }
00322         
00323         const double d_dt_fast_sodium_current_m_gate__m = (_lt_0_row[0] * (1.0 - var_fast_sodium_current_m_gate__m)) - (_lt_0_row[1] * var_fast_sodium_current_m_gate__m);
00324         const double d_dt_fast_sodium_current_h_gate__h = (_lt_0_row[2] * (1.0 - var_fast_sodium_current_h_gate__h)) - (_lt_0_row[3] * var_fast_sodium_current_h_gate__h);
00325         const double d_dt_fast_sodium_current_j_gate__j = (_lt_0_row[4] * (1.0 - var_fast_sodium_current_j_gate__j)) - (_lt_0_row[5] * var_fast_sodium_current_j_gate__j);
00326         const double d_dt_slow_inward_current_d_gate__d = (_lt_0_row[6] * (1.0 - var_slow_inward_current_d_gate__d)) - (_lt_0_row[7] * var_slow_inward_current_d_gate__d);
00327         const double d_dt_slow_inward_current_f_gate__f = (_lt_0_row[8] * (1.0 - var_slow_inward_current_f_gate__f)) - (_lt_0_row[9] * var_slow_inward_current_f_gate__f);
00328         const double d_dt_time_dependent_potassium_current_X_gate__X = (_lt_0_row[11] * (1.0 - var_time_dependent_potassium_current_X_gate__X)) - (_lt_0_row[12] * var_time_dependent_potassium_current_X_gate__X);
00329         const double d_dt_intracellular_calcium_concentration__Cai = ( -0.0001 * var_slow_inward_current__i_si) + (0.07 * (0.0001 - var_intracellular_calcium_concentration__Cai));
00330         
00331         rDY[0] = 1.0*d_dt_membrane__V;
00332         rDY[1] = 1.0*d_dt_fast_sodium_current_m_gate__m;
00333         rDY[2] = 1.0*d_dt_fast_sodium_current_h_gate__h;
00334         rDY[3] = 1.0*d_dt_fast_sodium_current_j_gate__j;
00335         rDY[4] = 1.0*d_dt_slow_inward_current_d_gate__d;
00336         rDY[5] = 1.0*d_dt_slow_inward_current_f_gate__f;
00337         rDY[6] = 1.0*d_dt_time_dependent_potassium_current_X_gate__X;
00338         rDY[7] = 1.0*d_dt_intracellular_calcium_concentration__Cai;
00339     }
00340     
00341 template<>
00342 void OdeSystemInformation<Cellluo_rudy_1991FromCellMLOpt>::Initialise(void)
00343 {
00344     // Time units: millisecond
00345     // 
00346     this->mVariableNames.push_back("membrane__V");
00347     this->mVariableUnits.push_back("millivolt");
00348     this->mInitialConditions.push_back(-83.853);
00349 
00350     this->mVariableNames.push_back("fast_sodium_current_m_gate__m");
00351     this->mVariableUnits.push_back("dimensionless");
00352     this->mInitialConditions.push_back(0.00187018);
00353 
00354     this->mVariableNames.push_back("fast_sodium_current_h_gate__h");
00355     this->mVariableUnits.push_back("dimensionless");
00356     this->mInitialConditions.push_back(0.9804713);
00357 
00358     this->mVariableNames.push_back("fast_sodium_current_j_gate__j");
00359     this->mVariableUnits.push_back("dimensionless");
00360     this->mInitialConditions.push_back(0.98767124);
00361 
00362     this->mVariableNames.push_back("slow_inward_current_d_gate__d");
00363     this->mVariableUnits.push_back("dimensionless");
00364     this->mInitialConditions.push_back(0.00316354);
00365 
00366     this->mVariableNames.push_back("slow_inward_current_f_gate__f");
00367     this->mVariableUnits.push_back("dimensionless");
00368     this->mInitialConditions.push_back(0.99427859);
00369 
00370     this->mVariableNames.push_back("time_dependent_potassium_current_X_gate__X");
00371     this->mVariableUnits.push_back("dimensionless");
00372     this->mInitialConditions.push_back(0.16647703);
00373 
00374     this->mVariableNames.push_back("intracellular_calcium_concentration__Cai");
00375     this->mVariableUnits.push_back("millimolar");
00376     this->mInitialConditions.push_back(0.0002);
00377 
00378     this->mInitialised = true;
00379 }
00380 
00381 
00382 // Serialization for Boost >= 1.36
00383 #include "SerializationExportWrapperForCpp.hpp"
00384 CHASTE_CLASS_EXPORT(Cellluo_rudy_1991FromCellMLOpt)

Generated by  doxygen 1.6.2