00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #include "HodgkinHuxleySquidAxon1952OriginalOdeSystem.hpp"
00029 #include "AbstractOdeSystem.hpp"
00030 #include "OdeSystemInformation.hpp"
00031 #include <cmath>
00032
00033
00037 HodgkinHuxleySquidAxon1952OriginalOdeSystem::HodgkinHuxleySquidAxon1952OriginalOdeSystem(
00038 AbstractIvpOdeSolver *pOdeSolver,
00039 AbstractStimulusFunction *pIntracellularStimulus)
00040 : AbstractCardiacCell(pOdeSolver, 4, 0, pIntracellularStimulus)
00041 {
00042 mpSystemInfo = OdeSystemInformation<HodgkinHuxleySquidAxon1952OriginalOdeSystem>::Instance();
00043
00044 AbstractCardiacCell::Init();
00045 }
00046
00050 HodgkinHuxleySquidAxon1952OriginalOdeSystem::~HodgkinHuxleySquidAxon1952OriginalOdeSystem(void)
00051 {
00052 }
00053
00060 void HodgkinHuxleySquidAxon1952OriginalOdeSystem::EvaluateYDerivatives(double time, const std::vector<double> &rY, std::vector<double>& rDY)
00061 {
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071 double membrane_V = rY[0];
00072 double potassium_channel_n_gate_n = rY[1];
00073 double sodium_channel_h_gate_h = rY[2];
00074 double sodium_channel_m_gate_m = rY[3];
00075
00076
00077
00078
00079
00080 double leakage_current_E_L = membrane_E_R+10.613;
00081 double leakage_current_i_L = leakage_current_g_L*(membrane_V-leakage_current_E_L);
00082
00083 double membrane_i_Stim = GetStimulus(time);
00084
00085 double sodium_channel_E_Na = membrane_E_R+115.0;
00086 double sodium_channel_i_Na = sodium_channel_g_Na*sodium_channel_m_gate_m*sodium_channel_m_gate_m*sodium_channel_m_gate_m*sodium_channel_h_gate_h*(membrane_V-sodium_channel_E_Na);
00087 double potassium_channel_E_K = membrane_E_R-12.0;
00088 double potassium_channel_i_K = potassium_channel_g_K*potassium_channel_n_gate_n*potassium_channel_n_gate_n*potassium_channel_n_gate_n*potassium_channel_n_gate_n*(membrane_V-potassium_channel_E_K);
00089
00090 double membrane_V_prime = -(-membrane_i_Stim+sodium_channel_i_Na+potassium_channel_i_K+leakage_current_i_L)/membrane_Cm;
00091
00092 if (mSetVoltageDerivativeToZero)
00093 {
00094 membrane_V_prime = 0;
00095 }
00096
00097 double potassium_channel_n_gate_alpha_n;
00098 if (-65.0001<membrane_V && membrane_V<-64.9999)
00099 {
00100 potassium_channel_n_gate_alpha_n = 0.1;
00101 }
00102 else
00103 {
00104 potassium_channel_n_gate_alpha_n = -0.01*(membrane_V+65.0)/(exp(-(membrane_V+65.0)/10.0)-1.0);
00105 }
00106
00107 double potassium_channel_n_gate_beta_n = 0.125*exp((membrane_V+75.0)/80.0);
00108 double potassium_channel_n_gate_n_prime = potassium_channel_n_gate_alpha_n*(1.0-potassium_channel_n_gate_n)-potassium_channel_n_gate_beta_n*potassium_channel_n_gate_n;
00109 double sodium_channel_h_gate_alpha_h = 0.07*exp(-(membrane_V+75.0)/20.0);
00110 double sodium_channel_h_gate_beta_h = 1.0/(exp(-(membrane_V+45.0)/10.0)+1.0);
00111 double sodium_channel_h_gate_h_prime = sodium_channel_h_gate_alpha_h*(1.0-sodium_channel_h_gate_h)-sodium_channel_h_gate_beta_h*sodium_channel_h_gate_h;
00112
00113 double sodium_channel_m_gate_alpha_m;
00114 if (-50.0001<membrane_V && membrane_V<-49.9999)
00115 {
00116 sodium_channel_m_gate_alpha_m = 1;
00117 }
00118 else
00119 {
00120 sodium_channel_m_gate_alpha_m = -0.1*(membrane_V+50.0)/(exp(-(membrane_V+50.0)/10.0)-1.0);
00121 }
00122 double sodium_channel_m_gate_beta_m = 4.0*exp(-(membrane_V+75.0)/18.0);
00123 double sodium_channel_m_gate_m_prime = sodium_channel_m_gate_alpha_m*(1.0-sodium_channel_m_gate_m)-sodium_channel_m_gate_beta_m*sodium_channel_m_gate_m;
00124
00125 rDY[0] = membrane_V_prime;
00126 rDY[1] = potassium_channel_n_gate_n_prime;
00127 rDY[2] = sodium_channel_h_gate_h_prime;
00128 rDY[3] = sodium_channel_m_gate_m_prime;
00129 }
00130
00131
00132 double HodgkinHuxleySquidAxon1952OriginalOdeSystem::GetIIonic()
00133 {
00134 double membrane_V = mStateVariables[mVoltageIndex];
00135 double potassium_channel_n_gate_n = mStateVariables[1];
00136 double sodium_channel_h_gate_h = mStateVariables[2];
00137 double sodium_channel_m_gate_m = mStateVariables[3];
00138
00139
00140
00141
00142
00143 double leakage_current_E_L = membrane_E_R+10.613;
00144 double leakage_current_i_L = leakage_current_g_L*(membrane_V-leakage_current_E_L);
00145
00146 double sodium_channel_E_Na = membrane_E_R+115.0;
00147 double sodium_channel_i_Na = sodium_channel_g_Na*sodium_channel_m_gate_m*sodium_channel_m_gate_m*sodium_channel_m_gate_m*sodium_channel_h_gate_h*(membrane_V-sodium_channel_E_Na);
00148 double potassium_channel_E_K = membrane_E_R-12.0;
00149 double potassium_channel_i_K = potassium_channel_g_K*potassium_channel_n_gate_n*potassium_channel_n_gate_n*potassium_channel_n_gate_n*potassium_channel_n_gate_n*(membrane_V-potassium_channel_E_K);
00150 double i_ionic = sodium_channel_i_Na+potassium_channel_i_K+leakage_current_i_L;
00151 return i_ionic;
00152 }
00153
00154
00155 template<>
00156 void OdeSystemInformation<HodgkinHuxleySquidAxon1952OriginalOdeSystem>::Initialise(void)
00157 {
00158
00159
00160
00161 this->mVariableNames.push_back("V");
00162 this->mVariableUnits.push_back("mV");
00163 this->mInitialConditions.push_back(-75.0);
00164
00165 this->mVariableNames.push_back("n");
00166 this->mVariableUnits.push_back("");
00167 this->mInitialConditions.push_back(0.325);
00168
00169 this->mVariableNames.push_back("h");
00170 this->mVariableUnits.push_back("");
00171 this->mInitialConditions.push_back(0.6);
00172
00173 this->mVariableNames.push_back("m");
00174 this->mVariableUnits.push_back("");
00175 this->mInitialConditions.push_back(0.05);
00176
00177 this->mInitialised = true;
00178 }