# Size of variable arrays: sizeAlgebraic = 80 sizeStates = 34 sizeConstants = 102 from math import * from numpy import * def createLegends(): legend_states = [""] * sizeStates legend_rates = [""] * sizeStates legend_algebraic = [""] * sizeAlgebraic legend_voi = "" legend_constants = [""] * sizeConstants legend_voi = "time in component environment (second)" legend_states[0] = "V in component membrane (millivolt)" legend_constants[0] = "R in component membrane (joule_per_mole_kelvin)" legend_constants[1] = "T in component membrane (kelvin)" legend_constants[2] = "F in component membrane (coulomb_per_millimole)" legend_constants[3] = "C_sc in component membrane (microF_per_cm2)" legend_algebraic[10] = "i_Stim in component membrane (microA_per_microF)" legend_algebraic[26] = "i_Na in component fast_sodium_current (microA_per_microF)" legend_algebraic[49] = "i_Ca in component L_type_Ca_current (microA_per_microF)" legend_algebraic[50] = "i_Ca_K in component L_type_Ca_current (microA_per_microF)" legend_algebraic[31] = "i_Kr in component rapid_activating_delayed_rectifiyer_K_current (microA_per_microF)" legend_algebraic[33] = "i_Ks in component slow_activating_delayed_rectifiyer_K_current (microA_per_microF)" legend_algebraic[34] = "i_to1 in component transient_outward_potassium_current (microA_per_microF)" legend_algebraic[36] = "i_K1 in component time_independent_potassium_current (microA_per_microF)" legend_algebraic[38] = "i_Kp in component plateau_potassium_current (microA_per_microF)" legend_algebraic[39] = "i_NaCa in component Na_Ca_exchanger (microA_per_microF)" legend_algebraic[41] = "i_NaK in component sodium_potassium_pump (microA_per_microF)" legend_algebraic[43] = "i_p_Ca in component sarcolemmal_calcium_pump (microA_per_microF)" legend_algebraic[45] = "i_Ca_b in component calcium_background_current (microA_per_microF)" legend_algebraic[46] = "i_Na_b in component sodium_background_current (microA_per_microF)" legend_algebraic[79] = "i_K_ATP in component ATP_sensitive_potassium_current (microA_per_microF)" legend_constants[4] = "stim_start in component membrane (second)" legend_constants[5] = "stim_end in component membrane (second)" legend_constants[6] = "stim_period in component membrane (second)" legend_constants[7] = "stim_duration in component membrane (second)" legend_constants[8] = "stim_amplitude in component membrane (microA_per_microF)" legend_algebraic[21] = "E_Na in component fast_sodium_current (millivolt)" legend_constants[9] = "g_Na in component fast_sodium_current (milliS_per_microF)" legend_constants[10] = "Na_o in component extracellular_ion_concentrations (millimolar)" legend_states[1] = "Na_i in component intracellular_ion_concentrations (millimolar)" legend_states[2] = "m in component fast_sodium_current_m_gate (dimensionless)" legend_states[3] = "h in component fast_sodium_current_h_gate (dimensionless)" legend_states[4] = "j in component fast_sodium_current_j_gate (dimensionless)" legend_algebraic[11] = "alpha_m in component fast_sodium_current_m_gate (per_second)" legend_algebraic[22] = "beta_m in component fast_sodium_current_m_gate (per_second)" legend_algebraic[0] = "E0_m in component fast_sodium_current_m_gate (millivolt)" legend_algebraic[1] = "alpha_h in component fast_sodium_current_h_gate (per_second)" legend_algebraic[12] = "beta_h in component fast_sodium_current_h_gate (per_second)" legend_algebraic[2] = "alpha_j in component fast_sodium_current_j_gate (per_second)" legend_algebraic[13] = "beta_j in component fast_sodium_current_j_gate (per_second)" legend_algebraic[29] = "E_K in component rapid_activating_delayed_rectifiyer_K_current (millivolt)" legend_constants[11] = "g_Kr in component rapid_activating_delayed_rectifiyer_K_current (milliS_per_microF)" legend_constants[99] = "f_K_o in component rapid_activating_delayed_rectifiyer_K_current (dimensionless)" legend_algebraic[30] = "R_V in component rapid_activating_delayed_rectifiyer_K_current (dimensionless)" legend_constants[12] = "K_o in component extracellular_ion_concentrations (millimolar)" legend_states[5] = "K_i in component intracellular_ion_concentrations (millimolar)" legend_states[6] = "X_kr in component rapid_activating_delayed_rectifiyer_K_current_X_kr_gate (dimensionless)" legend_algebraic[3] = "K12 in component rapid_activating_delayed_rectifiyer_K_current_X_kr_gate (dimensionless)" legend_algebraic[14] = "K21 in component rapid_activating_delayed_rectifiyer_K_current_X_kr_gate (dimensionless)" legend_algebraic[23] = "X_kr_inf in component rapid_activating_delayed_rectifiyer_K_current_X_kr_gate (dimensionless)" legend_algebraic[27] = "tau_X_kr in component rapid_activating_delayed_rectifiyer_K_current_X_kr_gate (second)" legend_constants[13] = "tau_factor in component rapid_activating_delayed_rectifiyer_K_current_X_kr_gate (dimensionless)" legend_constants[14] = "g_Ks in component slow_activating_delayed_rectifiyer_K_current (milliS_per_microF)" legend_algebraic[32] = "E_Ks in component slow_activating_delayed_rectifiyer_K_current (millivolt)" legend_states[7] = "X_ks in component slow_activating_delayed_rectifiyer_K_current_X_ks_gate (dimensionless)" legend_algebraic[15] = "tau_X_ks in component slow_activating_delayed_rectifiyer_K_current_X_ks_gate (second)" legend_algebraic[4] = "X_ks_infinity in component slow_activating_delayed_rectifiyer_K_current_X_ks_gate (dimensionless)" legend_constants[15] = "g_to1 in component transient_outward_potassium_current (milliS_per_microF)" legend_states[8] = "X_to1 in component transient_outward_potassium_current_X_to1_gate (dimensionless)" legend_states[9] = "Y_to1 in component transient_outward_potassium_current_Y_to1_gate (dimensionless)" legend_algebraic[5] = "alpha_X_to1 in component transient_outward_potassium_current_X_to1_gate (per_second)" legend_algebraic[16] = "beta_X_to1 in component transient_outward_potassium_current_X_to1_gate (per_second)" legend_algebraic[6] = "alpha_Y_to1 in component transient_outward_potassium_current_Y_to1_gate (per_second)" legend_algebraic[17] = "beta_Y_to1 in component transient_outward_potassium_current_Y_to1_gate (per_second)" legend_constants[16] = "g_K1 in component time_independent_potassium_current (milliS_per_microF)" legend_constants[17] = "K_mK1 in component time_independent_potassium_current (millimolar)" legend_algebraic[35] = "K1_infinity_V in component time_independent_potassium_current_K1_gate (dimensionless)" legend_constants[18] = "g_Kp in component plateau_potassium_current (milliS_per_microF)" legend_algebraic[37] = "Kp_V in component plateau_potassium_current_Kp_gate (dimensionless)" legend_constants[19] = "K_mCa in component Na_Ca_exchanger (millimolar)" legend_constants[20] = "K_mNa in component Na_Ca_exchanger (millimolar)" legend_constants[21] = "K_NaCa in component Na_Ca_exchanger (microA_per_microF)" legend_constants[22] = "K_sat in component Na_Ca_exchanger (dimensionless)" legend_constants[23] = "eta in component Na_Ca_exchanger (dimensionless)" legend_states[10] = "Ca_i in component intracellular_ion_concentrations (millimolar)" legend_constants[24] = "Ca_o in component extracellular_ion_concentrations (millimolar)" legend_constants[25] = "I_NaK in component sodium_potassium_pump (microA_per_microF)" legend_algebraic[40] = "f_NaK in component sodium_potassium_pump (dimensionless)" legend_constants[26] = "K_mNa_i in component sodium_potassium_pump (millimolar)" legend_constants[27] = "K_mK_o in component sodium_potassium_pump (millimolar)" legend_constants[100] = "sigma in component sodium_potassium_pump (dimensionless)" legend_algebraic[42] = "i_p_Ca_winslow in component sarcolemmal_calcium_pump (microA_per_microF)" legend_constants[28] = "K_mpCa in component sarcolemmal_calcium_pump (millimolar)" legend_constants[29] = "I_pCa in component sarcolemmal_calcium_pump (microA_per_microF)" legend_states[11] = "MgATP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_constants[30] = "MgATP_i0 in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_constants[31] = "g_Cab in component calcium_background_current (milliS_per_microF)" legend_algebraic[44] = "E_Ca in component calcium_background_current (millivolt)" legend_constants[32] = "g_Nab in component sodium_background_current (milliS_per_microF)" legend_constants[33] = "P_Ca in component L_type_Ca_current (cm_per_second)" legend_constants[34] = "p_prime_k in component L_type_Ca_current (cm_per_second)" legend_algebraic[48] = "f_Ca in component L_type_Ca_current (dimensionless)" legend_algebraic[47] = "P_rom in component L_type_Ca_current (dimensionless)" legend_constants[35] = "k_MgATP_ss in component L_type_Ca_current (millimolar)" legend_algebraic[7] = "alpha in component L_type_Ca_current (per_second)" legend_algebraic[18] = "beta in component L_type_Ca_current (per_second)" legend_algebraic[8] = "gamma in component L_type_Ca_current (per_second)" legend_algebraic[24] = "v_gamma in component L_type_Ca_current (per_second)" legend_algebraic[25] = "alpha_prime in component L_type_Ca_current (per_second)" legend_algebraic[28] = "beta_prime in component L_type_Ca_current (per_second)" legend_constants[36] = "a in component L_type_Ca_current (dimensionless)" legend_constants[37] = "b in component L_type_Ca_current (dimensionless)" legend_constants[38] = "g in component L_type_Ca_current (per_second)" legend_constants[39] = "f in component L_type_Ca_current (per_second)" legend_algebraic[19] = "v_omega in component L_type_Ca_current (per_second)" legend_constants[40] = "omega in component L_type_Ca_current (per_second)" legend_constants[101] = "x in component L_type_Ca_current (dimensionless)" legend_states[12] = "z in component L_type_Ca_current (dimensionless)" legend_states[13] = "v in component L_type_Ca_current (dimensionless)" legend_states[14] = "w in component L_type_Ca_current (dimensionless)" legend_states[15] = "Ca_ss in component intracellular_ion_concentrations (millimolar)" legend_states[16] = "MgATP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[17] = "y in component L_type_Ca_current_y_gate (dimensionless)" legend_algebraic[9] = "y_infinity in component L_type_Ca_current_y_gate (dimensionless)" legend_algebraic[20] = "tau_y in component L_type_Ca_current_y_gate (second)" legend_algebraic[75] = "g_K_ATP in component ATP_sensitive_potassium_current (milliS_per_microF)" legend_constants[41] = "G_K_ATP in component ATP_sensitive_potassium_current (milliS_per_microF)" legend_algebraic[73] = "f_K_ATP in component ATP_sensitive_potassium_current (dimensionless)" legend_algebraic[69] = "f_ATP in component ATP_sensitive_potassium_current (dimensionless)" legend_constants[42] = "k_ATP in component ATP_sensitive_potassium_current (millimolar)" legend_constants[43] = "k_MgADP in component ATP_sensitive_potassium_current (millimolar)" legend_constants[44] = "K_o_normal in component ATP_sensitive_potassium_current (millimolar)" legend_algebraic[51] = "f_MgADP_ in component ATP_sensitive_potassium_current (dimensionless)" legend_algebraic[53] = "f_MgADP in component ATP_sensitive_potassium_current (dimensionless)" legend_algebraic[52] = "c_MgADP in component ATP_sensitive_potassium_current (dimensionless)" legend_constants[45] = "go in component ATP_sensitive_potassium_current (dimensionless)" legend_constants[46] = "gd in component ATP_sensitive_potassium_current (dimensionless)" legend_algebraic[68] = "ATP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[18] = "MgADP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_algebraic[54] = "J_rel in component RyR_channel (millimolar_per_second)" legend_constants[47] = "v1 in component RyR_channel (per_second)" legend_constants[48] = "k_a_plus in component RyR_channel (millimolar4_per_second)" legend_constants[49] = "k_a_minus in component RyR_channel (per_second)" legend_constants[50] = "k_b_plus in component RyR_channel (millimolar3_per_second)" legend_constants[51] = "k_b_minus in component RyR_channel (per_second)" legend_constants[52] = "k_c_plus in component RyR_channel (per_second)" legend_constants[53] = "k_c_minus in component RyR_channel (per_second)" legend_states[19] = "P_O1 in component RyR_channel (dimensionless)" legend_states[20] = "P_O2 in component RyR_channel (dimensionless)" legend_states[21] = "P_C1 in component RyR_channel (dimensionless)" legend_states[22] = "P_C2 in component RyR_channel (dimensionless)" legend_constants[54] = "n in component RyR_channel (dimensionless)" legend_constants[55] = "m in component RyR_channel (dimensionless)" legend_states[23] = "Ca_JSR in component intracellular_ion_concentrations (millimolar)" legend_algebraic[58] = "J_up in component SERCA2a_pump (millimolar_per_second)" legend_algebraic[57] = "J_up_winslow in component SERCA2a_pump (millimolar_per_second)" legend_constants[56] = "K_fb in component SERCA2a_pump (millimolar)" legend_constants[57] = "K_rb in component SERCA2a_pump (millimolar)" legend_algebraic[55] = "fb in component SERCA2a_pump (dimensionless)" legend_algebraic[56] = "rb in component SERCA2a_pump (dimensionless)" legend_constants[58] = "Vmaxf in component SERCA2a_pump (millimolar_per_second)" legend_constants[59] = "Vmaxr in component SERCA2a_pump (millimolar_per_second)" legend_constants[60] = "K_SR in component SERCA2a_pump (dimensionless)" legend_constants[61] = "N_fb in component SERCA2a_pump (dimensionless)" legend_constants[62] = "N_rb in component SERCA2a_pump (dimensionless)" legend_states[24] = "Ca_NSR in component intracellular_ion_concentrations (millimolar)" legend_algebraic[60] = "J_tr in component intracellular_Ca_fluxes (millimolar_per_second)" legend_algebraic[59] = "J_xfer in component intracellular_Ca_fluxes (millimolar_per_second)" legend_algebraic[64] = "J_trpn in component intracellular_Ca_fluxes (millimolar_per_second)" legend_constants[63] = "tau_tr in component intracellular_Ca_fluxes (second)" legend_constants[64] = "tau_xfer in component intracellular_Ca_fluxes (second)" legend_states[25] = "HTRPNCa in component intracellular_Ca_fluxes (millimolar)" legend_states[26] = "LTRPNCa in component intracellular_Ca_fluxes (millimolar)" legend_algebraic[62] = "J_HTRPNCa in component intracellular_Ca_fluxes (millimolar_per_second)" legend_algebraic[63] = "J_LTRPNCa in component intracellular_Ca_fluxes (millimolar_per_second)" legend_constants[65] = "HTRPN_tot in component intracellular_Ca_fluxes (dimensionless)" legend_constants[66] = "LTRPN_tot in component intracellular_Ca_fluxes (dimensionless)" legend_constants[67] = "k_htrpn_plus in component intracellular_Ca_fluxes (per_millimolar_second)" legend_constants[68] = "k_htrpn_minus in component intracellular_Ca_fluxes (per_second)" legend_constants[69] = "k_ltrpn_plus in component intracellular_Ca_fluxes (per_millimolar_second)" legend_constants[70] = "k_ltrpn_minus in component intracellular_Ca_fluxes (per_second)" legend_constants[71] = "A_cap in component intracellular_ion_concentrations (cm2)" legend_constants[72] = "V_myo in component intracellular_ion_concentrations (microlitre)" legend_constants[73] = "V_JSR in component intracellular_ion_concentrations (microlitre)" legend_constants[74] = "V_NSR in component intracellular_ion_concentrations (microlitre)" legend_constants[75] = "V_ss in component intracellular_ion_concentrations (microlitre)" legend_constants[76] = "K_mCMDN in component intracellular_ion_concentrations (millimolar)" legend_constants[77] = "K_mEGTA in component intracellular_ion_concentrations (millimolar)" legend_constants[78] = "K_mCSQN in component intracellular_ion_concentrations (millimolar)" legend_constants[79] = "CMDN_tot in component intracellular_ion_concentrations (millimolar)" legend_constants[80] = "EGTA_tot in component intracellular_ion_concentrations (millimolar)" legend_constants[81] = "CSQN_tot in component intracellular_ion_concentrations (millimolar)" legend_algebraic[65] = "beta_i in component intracellular_ion_concentrations (dimensionless)" legend_algebraic[66] = "beta_SS in component intracellular_ion_concentrations (dimensionless)" legend_algebraic[61] = "beta_JSR in component intracellular_ion_concentrations (dimensionless)" legend_constants[82] = "k_plus_CaATP in component Ca_and_Mg_buffering_by_ATP (per_millimolar_second)" legend_constants[83] = "k_minus_CaATP in component Ca_and_Mg_buffering_by_ATP (per_second)" legend_constants[84] = "k_plus_CaADP in component Ca_and_Mg_buffering_by_ATP (per_millimolar_second)" legend_constants[85] = "k_minus_CaADP in component Ca_and_Mg_buffering_by_ATP (per_second)" legend_states[27] = "CaADP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[28] = "CaADP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[29] = "CaATP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_algebraic[74] = "ADP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_algebraic[70] = "ADP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_algebraic[67] = "ATP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[30] = "CaATP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[31] = "Mg_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[32] = "Mg_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_states[33] = "MgADP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_constants[86] = "ATP_tot in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_constants[87] = "k_plus_MgATP in component Ca_and_Mg_buffering_by_ATP (per_millimolar_second)" legend_constants[88] = "k_minus_MgATP in component Ca_and_Mg_buffering_by_ATP (per_second)" legend_algebraic[71] = "Jxfer_CaATP in component Ca_and_Mg_buffering_by_ATP (millimolar_per_second)" legend_algebraic[72] = "Jxfer_MgATP in component Ca_and_Mg_buffering_by_ATP (millimolar_per_second)" legend_algebraic[76] = "Jxfer_Mg in component Ca_and_Mg_buffering_by_ATP (millimolar_per_second)" legend_constants[89] = "tau_xfer_CaATP in component Ca_and_Mg_buffering_by_ATP (second)" legend_constants[90] = "tau_xfer_MgATP in component Ca_and_Mg_buffering_by_ATP (second)" legend_constants[91] = "tau_xfer_Mg in component Ca_and_Mg_buffering_by_ATP (second)" legend_constants[92] = "ADP_tot in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_constants[93] = "k_plus_MgADP in component Ca_and_Mg_buffering_by_ATP (per_millimolar_second)" legend_constants[94] = "k_minus_MgADP in component Ca_and_Mg_buffering_by_ATP (per_second)" legend_algebraic[77] = "Jxfer_CaADP in component Ca_and_Mg_buffering_by_ATP (millimolar_per_second)" legend_algebraic[78] = "Jxfer_MgADP in component Ca_and_Mg_buffering_by_ATP (millimolar_per_second)" legend_constants[95] = "tau_xfer_CaADP in component Ca_and_Mg_buffering_by_ATP (second)" legend_constants[96] = "tau_xfer_MgADP in component Ca_and_Mg_buffering_by_ATP (second)" legend_constants[97] = "V_myo in component model_parameters (microlitre)" legend_constants[98] = "V_ss in component model_parameters (microlitre)" legend_rates[0] = "d/dt V in component membrane (millivolt)" legend_rates[2] = "d/dt m in component fast_sodium_current_m_gate (dimensionless)" legend_rates[3] = "d/dt h in component fast_sodium_current_h_gate (dimensionless)" legend_rates[4] = "d/dt j in component fast_sodium_current_j_gate (dimensionless)" legend_rates[6] = "d/dt X_kr in component rapid_activating_delayed_rectifiyer_K_current_X_kr_gate (dimensionless)" legend_rates[7] = "d/dt X_ks in component slow_activating_delayed_rectifiyer_K_current_X_ks_gate (dimensionless)" legend_rates[8] = "d/dt X_to1 in component transient_outward_potassium_current_X_to1_gate (dimensionless)" legend_rates[9] = "d/dt Y_to1 in component transient_outward_potassium_current_Y_to1_gate (dimensionless)" legend_rates[13] = "d/dt v in component L_type_Ca_current (dimensionless)" legend_rates[14] = "d/dt w in component L_type_Ca_current (dimensionless)" legend_rates[12] = "d/dt z in component L_type_Ca_current (dimensionless)" legend_rates[17] = "d/dt y in component L_type_Ca_current_y_gate (dimensionless)" legend_rates[21] = "d/dt P_C1 in component RyR_channel (dimensionless)" legend_rates[19] = "d/dt P_O1 in component RyR_channel (dimensionless)" legend_rates[20] = "d/dt P_O2 in component RyR_channel (dimensionless)" legend_rates[22] = "d/dt P_C2 in component RyR_channel (dimensionless)" legend_rates[25] = "d/dt HTRPNCa in component intracellular_Ca_fluxes (millimolar)" legend_rates[26] = "d/dt LTRPNCa in component intracellular_Ca_fluxes (millimolar)" legend_rates[10] = "d/dt Ca_i in component intracellular_ion_concentrations (millimolar)" legend_rates[1] = "d/dt Na_i in component intracellular_ion_concentrations (millimolar)" legend_rates[5] = "d/dt K_i in component intracellular_ion_concentrations (millimolar)" legend_rates[15] = "d/dt Ca_ss in component intracellular_ion_concentrations (millimolar)" legend_rates[23] = "d/dt Ca_JSR in component intracellular_ion_concentrations (millimolar)" legend_rates[24] = "d/dt Ca_NSR in component intracellular_ion_concentrations (millimolar)" legend_rates[29] = "d/dt CaATP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[16] = "d/dt MgATP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[30] = "d/dt CaATP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[11] = "d/dt MgATP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[28] = "d/dt CaADP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[33] = "d/dt MgADP_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[27] = "d/dt CaADP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[18] = "d/dt MgADP_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[31] = "d/dt Mg_ss in component Ca_and_Mg_buffering_by_ATP (millimolar)" legend_rates[32] = "d/dt Mg_i in component Ca_and_Mg_buffering_by_ATP (millimolar)" return (legend_states, legend_algebraic, legend_voi, legend_constants) def initConsts(): constants = [0.0] * sizeConstants; states = [0.0] * sizeStates; states[0] = -96.1638 constants[0] = 8.314472 constants[1] = 310 constants[2] = 96.4853415 constants[3] = 0.001 constants[4] = 0.1 constants[5] = 100000000 constants[6] = 1 constants[7] = 0.0005 constants[8] = -100.0 constants[9] = 12.8 constants[10] = 138 states[1] = 10 states[2] = 0.0328302 states[3] = 0.988354 states[4] = 0.99254 constants[11] = 0.0034 constants[12] = 4 states[5] = 159.48 states[6] = 0.51 constants[13] = 1 constants[14] = 0.0027134 states[7] = 0.264 constants[15] = 0.23815 states[8] = 2.63 states[9] = 0.99 constants[16] = 2.8 constants[17] = 13 constants[18] = 0.002216 constants[19] = 1.38 constants[20] = 87.5 constants[21] = 0.3 constants[22] = 0.2 constants[23] = 0.35 states[10] = 8.464E-5 constants[24] = 2 constants[25] = 0.693 constants[26] = 10 constants[27] = 1.5 constants[28] = 0.00005 constants[29] = 0.05 states[11] = 6.4395 constants[30] = 2.888 constants[31] = 0.0003842 constants[32] = 0.0031 constants[33] = 3.594e-4 constants[34] = 6.658e-7 constants[35] = 1.4 constants[36] = 2 constants[37] = 2 constants[38] = 2000 constants[39] = 300 constants[40] = 10 states[12] = 0.1 states[13] = 0.1 states[14] = 0.1 states[15] = 1.315E-4 states[16] = 6.4395 states[17] = 0.798 constants[41] = 0.05 constants[42] = 0.6 constants[43] = 0.4 constants[44] = 4.5 constants[45] = 0.08 constants[46] = 0.89 states[18] = 0.298E-2 constants[47] = 1800 constants[48] = 1.215e13 constants[49] = 576 constants[50] = 4.05e9 constants[51] = 1930 constants[52] = 100 constants[53] = 0.8 states[19] = 0 states[20] = 0 states[21] = 0.47 states[22] = 0.53 constants[54] = 4 constants[55] = 3 states[23] = 0.2616 constants[56] = 0.000168 constants[57] = 3.29 constants[58] = 0.0813 constants[59] = 0.318 constants[60] = 1 constants[61] = 1.2 constants[62] = 1 states[24] = 0.2620 constants[63] = 0.0005747 constants[64] = 0.0267 states[25] = 0.98 states[26] = 0.078 constants[65] = 0.14 constants[66] = 0.07 constants[67] = 20000 constants[68] = 0.066 constants[69] = 40000 constants[70] = 40 constants[71] = 0.0001534 constants[72] = 0.00002584 constants[73] = 0.00000016 constants[74] = 0.0000021 constants[75] = 0.0000000012 constants[76] = 0.00238 constants[77] = 0.00015 constants[78] = 0.8 constants[79] = 0.05 constants[80] = 0 constants[81] = 15 constants[82] = 225000.0 constants[83] = 45000.0 constants[84] = 125000.0 constants[85] = 193500 states[27] = 0.11E-6 states[28] = 0.13E-6 states[29] = 0.25E-3 states[30] = 0.237E-3 states[31] = 1.0 states[32] = 1.0 states[33] = 0.298E-2 constants[86] = 7.0 constants[87] = 125000.0 constants[88] = 10875.0 constants[89] = 0.0534 constants[90] = 0.0534 constants[91] = 0.0267 constants[92] = 0.005 constants[93] = 125000.0 constants[94] = 84500.0 constants[95] = 0.0534 constants[96] = 0.0534 constants[97] = 0.00002584 constants[98] = 0.0000000012 constants[99] = power(constants[12]/4.00000, 1.0/2) constants[100] = (1.00000/7.00000)*(exp(constants[10]/67.3000)-1.00000) constants[101] = constants[39]/(constants[39]+constants[38]) return (states, constants) def computeRates(voi, states, constants): rates = [0.0] * sizeStates; algebraic = [0.0] * sizeAlgebraic rates[21] = -constants[48]*(power(states[15], constants[54]))*states[21]+constants[49]*states[19] rates[19] = (constants[48]*(power(states[15], constants[54]))*states[21]-(constants[49]*states[19]+constants[50]*(power(states[15], constants[55]))*states[19]+constants[52]*states[19]))+constants[51]*states[20]+constants[53]*states[22] rates[20] = constants[50]*(power(states[15], constants[55]))*states[19]-constants[51]*states[20] rates[22] = constants[52]*states[19]-constants[53]*states[22] algebraic[1] = custom_piecewise([less(states[0] , -40.0000), 135.000*exp((80.0000+states[0])/-6.80000) , True, 0.00000]) algebraic[12] = custom_piecewise([less(states[0] , -40.0000), 3560.00*exp(0.0790000*states[0])+310000.*exp(0.350000*states[0]) , True, 1000.00/(0.130000*(1.00000+exp((states[0]+10.6600)/-11.1000)))]) rates[3] = algebraic[1]*(1.00000-states[3])-algebraic[12]*states[3] algebraic[2] = custom_piecewise([less(states[0] , -40.0000), (1000.00*-(127140.*exp(0.244400*states[0])+3.47400e-05*exp(-0.0439100*states[0]))*(states[0]+37.7800))/(1.00000+exp(0.311000*(states[0]+79.2300))) , True, 0.00000]) algebraic[13] = custom_piecewise([less(states[0] , -40.0000), (121.200*exp(-0.0105200*states[0]))/(1.00000+exp(-0.137800*(states[0]+40.1400))) , True, (300.000*exp(-2.53500e-07*states[0]))/(1.00000+exp(-0.100000*(states[0]+32.0000)))]) rates[4] = algebraic[2]*(1.00000-states[4])-algebraic[13]*states[4] algebraic[15] = 0.00100000/((7.19000e-05*(states[0]-10.0000))/(1.00000-exp(-0.148000*(states[0]-10.0000)))+(0.000131000*(states[0]-10.0000))/(exp(0.0687000*(states[0]-10.0000))-1.00000)) algebraic[4] = 1.00000/(1.00000+exp(-(states[0]-24.7000)/13.6000)) rates[7] = (algebraic[4]-states[7])/algebraic[15] algebraic[5] = 45.1600*exp(0.0357700*states[0]) algebraic[16] = 98.9000*exp(-0.0623700*states[0]) rates[8] = algebraic[5]*(1.00000-states[8])-algebraic[16]*states[8] algebraic[6] = (5.41500*exp(-(states[0]+33.5000)/5.00000))/(1.00000+0.0513350*exp(-(states[0]+33.5000)/5.00000)) algebraic[17] = (5.41500*exp((states[0]+33.5000)/5.00000))/(1.00000+0.0513350*exp((states[0]+33.5000)/5.00000)) rates[9] = algebraic[6]*(1.00000-states[9])-algebraic[17]*states[9] algebraic[7] = 400.000*exp((states[0]+2.00000)/10.0000) algebraic[18] = 50.0000*exp(-(states[0]+2.00000)/13.0000) rates[13] = algebraic[7]*(1.00000-states[13])-algebraic[18]*states[13] algebraic[9] = 0.800000/(1.00000+exp((states[0]+12.5000)/5.00000))+0.200000 algebraic[20] = (20.0000+600.000/(1.00000+exp((states[0]+20.0000)/9.50000)))/1000.00 rates[17] = (algebraic[9]-states[17])/algebraic[20] algebraic[0] = states[0]+47.1300 algebraic[11] = custom_piecewise([less(fabs(algebraic[0]) , 1.00000e-05), 320.000/(0.100000-0.00500000*algebraic[0]) , True, (320.000*algebraic[0])/(1.00000-exp(-0.100000*algebraic[0]))]) algebraic[22] = 80.0000*exp(-states[0]/11.0000) rates[2] = custom_piecewise([greater_equal(states[0] , -90.0000), algebraic[11]*(1.00000-states[2])-algebraic[22]*states[2] , True, 0.00000]) algebraic[8] = (103.750*states[15])/1.00000 algebraic[24] = algebraic[8]*((power(constants[36], 0.00000))*(power(states[13], 0.00000))*(power(1.00000-states[13], 4.00000-0.00000))+(power(constants[36], 1.00000))*(power(states[13], 1.00000))*(power(1.00000-states[13], 4.00000-1.00000))+(power(constants[36], 2.00000))*(power(states[13], 2.00000))*(power(1.00000-states[13], 4.00000-2.00000))+(power(constants[36], 3.00000))*(power(states[13], 3.00000))*(power(1.00000-states[13], 4.00000-3.00000))+(power(constants[36], 4.00000))*(power(states[13], 4.00000))*(power(1.00000-states[13], 4.00000-4.00000))+-((power(constants[36], 4.00000))*(power(states[13], 4.00000))*constants[101])) algebraic[19] = constants[40]*((power(constants[37], -0.00000))*(power(states[14], 0.00000))*(power(1.00000-states[14], 4.00000-0.00000))+(power(constants[37], -1.00000))*(power(states[14], 1.00000))*(power(1.00000-states[14], 4.00000-1.00000))+(power(constants[37], -2.00000))*(power(states[14], 2.00000))*(power(1.00000-states[14], 4.00000-2.00000))+(power(constants[37], -3.00000))*(power(states[14], 3.00000))*(power(1.00000-states[14], 4.00000-3.00000))+(power(constants[37], -4.00000))*(power(states[14], 4.00000))*(power(1.00000-states[14], 4.00000-4.00000))) rates[12] = algebraic[19]*(1.00000-states[12])-algebraic[24]*states[12] algebraic[3] = exp(-5.49500+0.169100*states[0]) algebraic[14] = exp(-7.67700-0.0128000*states[0]) algebraic[23] = algebraic[3]/(algebraic[3]+algebraic[14]) algebraic[27] = 0.00100000/(algebraic[3]+algebraic[14])+constants[13]*0.0270000 rates[6] = (algebraic[23]-states[6])/algebraic[27] algebraic[25] = algebraic[7]*constants[36] algebraic[28] = algebraic[18]/constants[37] rates[14] = algebraic[25]*(1.00000-states[14])-algebraic[28]*states[14] algebraic[21] = ((constants[0]*constants[1])/constants[2])*log(constants[10]/states[1]) algebraic[26] = constants[9]*(power(states[2], 3.00000))*states[3]*states[4]*(states[0]-algebraic[21]) algebraic[39] = ((constants[21]*5000.00)/((power(constants[20], 3.00000)+power(constants[10], 3.00000))*(constants[19]+constants[24])*(1.00000+constants[22]*exp(((constants[23]-1.00000)*states[0]*constants[2])/(constants[0]*constants[1])))))*(exp((constants[23]*states[0]*constants[2])/(constants[0]*constants[1]))*(power(states[1], 3.00000))*constants[24]-exp(((constants[23]-1.00000)*states[0]*constants[2])/(constants[0]*constants[1]))*(power(constants[10], 3.00000))*states[10]) algebraic[40] = 1.00000/(1.00000+0.124500*exp((-0.100000*states[0]*constants[2])/(constants[0]*constants[1]))+0.0365000*constants[100]*exp((-states[0]*constants[2])/(constants[0]*constants[1]))) algebraic[41] = (((constants[25]*algebraic[40])/(1.00000+power(constants[26]/states[1], 1.50000)))*constants[12])/(constants[12]+constants[27]) algebraic[46] = constants[32]*(states[0]-algebraic[21]) rates[1] = (-0.00000*(algebraic[26]+algebraic[46]+algebraic[39]*3.00000+algebraic[41]*3.00000)*constants[71]*1.00000)/(constants[72]*constants[2]) rates[25] = constants[67]*states[10]*(1.00000-states[25])-constants[68]*states[25] algebraic[55] = power(states[10]/constants[56], constants[61]) algebraic[56] = power(states[24]/constants[57], constants[62]) algebraic[57] = (constants[60]*(constants[58]*algebraic[55]-constants[59]*algebraic[56]))/(1.00000+algebraic[55]+algebraic[56]) algebraic[58] = (states[11]/constants[30])*algebraic[57] algebraic[60] = (states[24]-states[23])/constants[63] rates[24] = (algebraic[58]*constants[72])/constants[74]-(algebraic[60]*constants[73])/constants[74] algebraic[54] = constants[47]*(states[19]+states[20])*(states[23]-states[15]) algebraic[61] = 1.00000/(1.00000+(constants[81]*constants[78])/(power(constants[78]+states[23], 2.00000))) rates[23] = algebraic[61]*(algebraic[60]-algebraic[54]) rates[26] = constants[69]*states[10]*(1.00000-states[26])-constants[70]*states[26] algebraic[48] = 1.00000/(1.00000+power(constants[35]/states[16], 2.60000)) algebraic[47] = (power(states[13], 4.00000))*states[17]*states[12]*(constants[39]/(constants[39]+constants[38])) algebraic[49] = algebraic[48]*(constants[33]/(1.00000*constants[3]))*((4.00000*states[0]*(power(constants[2], 2.00000)))/(constants[0]*constants[1]))*((0.00100000*exp((constants[2]*states[0])/(constants[0]*constants[1]))-0.340000*constants[24])/(exp((constants[2]*states[0])/(constants[0]*constants[1]))-1.00000))*algebraic[47] algebraic[59] = (states[15]-states[10])/constants[64] algebraic[66] = 1.00000/(1.00000+(constants[79]*constants[76])/(power(constants[76]+states[15], 2.00000))+(constants[80]*constants[77])/(power(constants[77]+states[15], 2.00000))) algebraic[70] = constants[92]-(states[28]+states[33]) algebraic[67] = constants[86]-(states[29]+states[16]) rates[15] = algebraic[66]*(((algebraic[54]*constants[73])/constants[75]+constants[83]*states[29]+constants[85]*states[28])-((algebraic[59]*constants[72])/constants[75]+algebraic[49]*((constants[71]*1.00000)/(2.00000*constants[75]*constants[2]))+constants[82]*states[15]*algebraic[67]+constants[84]*states[15]*algebraic[70])) algebraic[71] = (states[29]-states[30])/constants[89] rates[29] = constants[82]*states[15]*algebraic[67]-(algebraic[71]*(constants[97]/constants[98])+constants[83]*states[29]) algebraic[72] = (states[16]-states[11])/constants[90] rates[16] = constants[87]*states[31]*algebraic[67]-(algebraic[72]*(constants[97]/constants[98])+constants[88]*states[16]) algebraic[68] = constants[86]-(states[30]+states[11]) rates[30] = (algebraic[71]+constants[82]*states[10]*algebraic[68])-constants[83]*states[30] rates[11] = (algebraic[72]+constants[87]*states[32]*algebraic[68])-constants[88]*states[11] algebraic[42] = (constants[29]*states[10])/(constants[28]+states[10]) algebraic[43] = (states[11]/constants[30])*algebraic[42] algebraic[44] = ((constants[0]*constants[1])/(2.00000*constants[2]))*log(constants[24]/states[10]) algebraic[45] = constants[31]*(states[0]-algebraic[44]) algebraic[62] = rates[25] algebraic[63] = rates[26] algebraic[64] = constants[65]*algebraic[62]+constants[66]*algebraic[63] algebraic[65] = 1.00000/(1.00000+(constants[79]*constants[76])/(power(constants[76]+states[10], 2.00000))+(constants[80]*constants[77])/(power(constants[77]+states[10], 2.00000))) algebraic[74] = constants[92]-(states[27]+states[18]) rates[10] = algebraic[65]*((algebraic[59]-(algebraic[58]+algebraic[64]))+-(((algebraic[43]+algebraic[45])-2.00000*algebraic[39])*((constants[71]*constants[3])/(2.00000*constants[72]*constants[2])))+((constants[83]*states[30]+constants[85]*states[27])-(constants[82]*states[10]*algebraic[68]+constants[84]*states[10]*algebraic[74]))) algebraic[77] = (states[28]-states[27])/constants[95] rates[28] = constants[84]*states[15]*algebraic[70]-(algebraic[77]*(constants[97]/constants[98])+constants[85]*states[28]) algebraic[78] = (states[33]-states[18])/constants[96] rates[33] = constants[93]*states[31]*algebraic[70]-(algebraic[78]*(constants[97]/constants[98])+constants[94]*states[33]) rates[27] = (algebraic[77]+constants[84]*states[10]*algebraic[74])-constants[85]*states[27] rates[18] = (algebraic[78]+constants[93]*states[32]*algebraic[74])-constants[94]*states[18] algebraic[76] = (states[31]-states[32])/constants[91] rates[31] = (constants[88]*states[16]+constants[94]*states[33])-(constants[87]*states[31]*algebraic[67]+constants[93]*states[31]*algebraic[70]+algebraic[76]*(constants[97]/constants[98])) rates[32] = (algebraic[76]+constants[88]*states[11]+constants[94]*states[18])-(constants[87]*states[32]*algebraic[68]+constants[93]*states[32]*algebraic[74]) algebraic[10] = custom_piecewise([greater_equal(voi , constants[4]) & less_equal(voi , constants[5]) & less_equal((voi-constants[4])-floor((voi-constants[4])/constants[6])*constants[6] , constants[7]), constants[8] , True, 0.00000]) algebraic[50] = 1.00000*algebraic[48]*(constants[34]/constants[3])*algebraic[47] algebraic[29] = ((constants[0]*constants[1])/constants[2])*log(constants[12]/states[5]) algebraic[30] = 1.00000/(1.00000+1.49450*exp(0.0446000*states[0])) algebraic[31] = constants[11]*constants[99]*algebraic[30]*states[6]*(states[0]-algebraic[29]) algebraic[32] = ((constants[0]*constants[1])/constants[2])*log((constants[12]+0.0183300*constants[10])/(states[5]+0.0183300*states[1])) algebraic[33] = constants[14]*(power(states[7], 2.00000))*(states[0]-algebraic[32]) algebraic[34] = constants[15]*states[8]*states[9]*(states[0]-algebraic[29]) algebraic[35] = 1.00000/(2.00000+exp(((1.50000*constants[2])/(constants[0]*constants[1]))*(states[0]-algebraic[29]))) algebraic[36] = ((constants[16]*algebraic[35]*constants[12])/(constants[12]+constants[17]))*(states[0]-algebraic[29]) algebraic[37] = 1.00000/(1.00000+exp((7.48800-states[0])/5.98000)) algebraic[38] = constants[18]*algebraic[37]*(states[0]-algebraic[29]) algebraic[69] = power(1.00000-algebraic[68]/(algebraic[68]+constants[42]), 4.00000) algebraic[51] = power(states[18]/(states[18]+constants[43]), 2.00000) algebraic[52] = 1.00000-algebraic[51] algebraic[53] = 1.00000-power(algebraic[52], 4.00000) algebraic[73] = constants[45]*algebraic[69]*(1.00000-algebraic[53])+constants[46]*algebraic[69]*algebraic[53] algebraic[75] = constants[41]*algebraic[73]*(power(constants[12]/constants[44], 0.240000)) algebraic[79] = algebraic[75]*(states[0]-algebraic[29]) rates[0] = (-1.00000*1.00000*(algebraic[26]+algebraic[49]+algebraic[50]+algebraic[31]+algebraic[33]+algebraic[34]+algebraic[36]+algebraic[38]+algebraic[39]+algebraic[41]+algebraic[43]+algebraic[46]+algebraic[45]+algebraic[79]+algebraic[10]))/constants[3] rates[5] = (-0.00000*(algebraic[50]+algebraic[31]+algebraic[33]+algebraic[36]+algebraic[38]+algebraic[34]+algebraic[79]+algebraic[41]*-2.00000)*constants[71]*1.00000)/(constants[72]*constants[2]) return(rates) def computeAlgebraic(constants, states, voi): algebraic = array([[0.0] * len(voi)] * sizeAlgebraic) states = array(states) voi = array(voi) algebraic[1] = custom_piecewise([less(states[0] , -40.0000), 135.000*exp((80.0000+states[0])/-6.80000) , True, 0.00000]) algebraic[12] = custom_piecewise([less(states[0] , -40.0000), 3560.00*exp(0.0790000*states[0])+310000.*exp(0.350000*states[0]) , True, 1000.00/(0.130000*(1.00000+exp((states[0]+10.6600)/-11.1000)))]) algebraic[2] = custom_piecewise([less(states[0] , -40.0000), (1000.00*-(127140.*exp(0.244400*states[0])+3.47400e-05*exp(-0.0439100*states[0]))*(states[0]+37.7800))/(1.00000+exp(0.311000*(states[0]+79.2300))) , True, 0.00000]) algebraic[13] = custom_piecewise([less(states[0] , -40.0000), (121.200*exp(-0.0105200*states[0]))/(1.00000+exp(-0.137800*(states[0]+40.1400))) , True, (300.000*exp(-2.53500e-07*states[0]))/(1.00000+exp(-0.100000*(states[0]+32.0000)))]) algebraic[15] = 0.00100000/((7.19000e-05*(states[0]-10.0000))/(1.00000-exp(-0.148000*(states[0]-10.0000)))+(0.000131000*(states[0]-10.0000))/(exp(0.0687000*(states[0]-10.0000))-1.00000)) algebraic[4] = 1.00000/(1.00000+exp(-(states[0]-24.7000)/13.6000)) algebraic[5] = 45.1600*exp(0.0357700*states[0]) algebraic[16] = 98.9000*exp(-0.0623700*states[0]) algebraic[6] = (5.41500*exp(-(states[0]+33.5000)/5.00000))/(1.00000+0.0513350*exp(-(states[0]+33.5000)/5.00000)) algebraic[17] = (5.41500*exp((states[0]+33.5000)/5.00000))/(1.00000+0.0513350*exp((states[0]+33.5000)/5.00000)) algebraic[7] = 400.000*exp((states[0]+2.00000)/10.0000) algebraic[18] = 50.0000*exp(-(states[0]+2.00000)/13.0000) algebraic[9] = 0.800000/(1.00000+exp((states[0]+12.5000)/5.00000))+0.200000 algebraic[20] = (20.0000+600.000/(1.00000+exp((states[0]+20.0000)/9.50000)))/1000.00 algebraic[0] = states[0]+47.1300 algebraic[11] = custom_piecewise([less(fabs(algebraic[0]) , 1.00000e-05), 320.000/(0.100000-0.00500000*algebraic[0]) , True, (320.000*algebraic[0])/(1.00000-exp(-0.100000*algebraic[0]))]) algebraic[22] = 80.0000*exp(-states[0]/11.0000) algebraic[8] = (103.750*states[15])/1.00000 algebraic[24] = algebraic[8]*((power(constants[36], 0.00000))*(power(states[13], 0.00000))*(power(1.00000-states[13], 4.00000-0.00000))+(power(constants[36], 1.00000))*(power(states[13], 1.00000))*(power(1.00000-states[13], 4.00000-1.00000))+(power(constants[36], 2.00000))*(power(states[13], 2.00000))*(power(1.00000-states[13], 4.00000-2.00000))+(power(constants[36], 3.00000))*(power(states[13], 3.00000))*(power(1.00000-states[13], 4.00000-3.00000))+(power(constants[36], 4.00000))*(power(states[13], 4.00000))*(power(1.00000-states[13], 4.00000-4.00000))+-((power(constants[36], 4.00000))*(power(states[13], 4.00000))*constants[101])) algebraic[19] = constants[40]*((power(constants[37], -0.00000))*(power(states[14], 0.00000))*(power(1.00000-states[14], 4.00000-0.00000))+(power(constants[37], -1.00000))*(power(states[14], 1.00000))*(power(1.00000-states[14], 4.00000-1.00000))+(power(constants[37], -2.00000))*(power(states[14], 2.00000))*(power(1.00000-states[14], 4.00000-2.00000))+(power(constants[37], -3.00000))*(power(states[14], 3.00000))*(power(1.00000-states[14], 4.00000-3.00000))+(power(constants[37], -4.00000))*(power(states[14], 4.00000))*(power(1.00000-states[14], 4.00000-4.00000))) algebraic[3] = exp(-5.49500+0.169100*states[0]) algebraic[14] = exp(-7.67700-0.0128000*states[0]) algebraic[23] = algebraic[3]/(algebraic[3]+algebraic[14]) algebraic[27] = 0.00100000/(algebraic[3]+algebraic[14])+constants[13]*0.0270000 algebraic[25] = algebraic[7]*constants[36] algebraic[28] = algebraic[18]/constants[37] algebraic[21] = ((constants[0]*constants[1])/constants[2])*log(constants[10]/states[1]) algebraic[26] = constants[9]*(power(states[2], 3.00000))*states[3]*states[4]*(states[0]-algebraic[21]) algebraic[39] = ((constants[21]*5000.00)/((power(constants[20], 3.00000)+power(constants[10], 3.00000))*(constants[19]+constants[24])*(1.00000+constants[22]*exp(((constants[23]-1.00000)*states[0]*constants[2])/(constants[0]*constants[1])))))*(exp((constants[23]*states[0]*constants[2])/(constants[0]*constants[1]))*(power(states[1], 3.00000))*constants[24]-exp(((constants[23]-1.00000)*states[0]*constants[2])/(constants[0]*constants[1]))*(power(constants[10], 3.00000))*states[10]) algebraic[40] = 1.00000/(1.00000+0.124500*exp((-0.100000*states[0]*constants[2])/(constants[0]*constants[1]))+0.0365000*constants[100]*exp((-states[0]*constants[2])/(constants[0]*constants[1]))) algebraic[41] = (((constants[25]*algebraic[40])/(1.00000+power(constants[26]/states[1], 1.50000)))*constants[12])/(constants[12]+constants[27]) algebraic[46] = constants[32]*(states[0]-algebraic[21]) algebraic[55] = power(states[10]/constants[56], constants[61]) algebraic[56] = power(states[24]/constants[57], constants[62]) algebraic[57] = (constants[60]*(constants[58]*algebraic[55]-constants[59]*algebraic[56]))/(1.00000+algebraic[55]+algebraic[56]) algebraic[58] = (states[11]/constants[30])*algebraic[57] algebraic[60] = (states[24]-states[23])/constants[63] algebraic[54] = constants[47]*(states[19]+states[20])*(states[23]-states[15]) algebraic[61] = 1.00000/(1.00000+(constants[81]*constants[78])/(power(constants[78]+states[23], 2.00000))) algebraic[48] = 1.00000/(1.00000+power(constants[35]/states[16], 2.60000)) algebraic[47] = (power(states[13], 4.00000))*states[17]*states[12]*(constants[39]/(constants[39]+constants[38])) algebraic[49] = algebraic[48]*(constants[33]/(1.00000*constants[3]))*((4.00000*states[0]*(power(constants[2], 2.00000)))/(constants[0]*constants[1]))*((0.00100000*exp((constants[2]*states[0])/(constants[0]*constants[1]))-0.340000*constants[24])/(exp((constants[2]*states[0])/(constants[0]*constants[1]))-1.00000))*algebraic[47] algebraic[59] = (states[15]-states[10])/constants[64] algebraic[66] = 1.00000/(1.00000+(constants[79]*constants[76])/(power(constants[76]+states[15], 2.00000))+(constants[80]*constants[77])/(power(constants[77]+states[15], 2.00000))) algebraic[70] = constants[92]-(states[28]+states[33]) algebraic[67] = constants[86]-(states[29]+states[16]) algebraic[71] = (states[29]-states[30])/constants[89] algebraic[72] = (states[16]-states[11])/constants[90] algebraic[68] = constants[86]-(states[30]+states[11]) algebraic[42] = (constants[29]*states[10])/(constants[28]+states[10]) algebraic[43] = (states[11]/constants[30])*algebraic[42] algebraic[44] = ((constants[0]*constants[1])/(2.00000*constants[2]))*log(constants[24]/states[10]) algebraic[45] = constants[31]*(states[0]-algebraic[44]) algebraic[62] = rates[25] algebraic[63] = rates[26] algebraic[64] = constants[65]*algebraic[62]+constants[66]*algebraic[63] algebraic[65] = 1.00000/(1.00000+(constants[79]*constants[76])/(power(constants[76]+states[10], 2.00000))+(constants[80]*constants[77])/(power(constants[77]+states[10], 2.00000))) algebraic[74] = constants[92]-(states[27]+states[18]) algebraic[77] = (states[28]-states[27])/constants[95] algebraic[78] = (states[33]-states[18])/constants[96] algebraic[76] = (states[31]-states[32])/constants[91] algebraic[10] = custom_piecewise([greater_equal(voi , constants[4]) & less_equal(voi , constants[5]) & less_equal((voi-constants[4])-floor((voi-constants[4])/constants[6])*constants[6] , constants[7]), constants[8] , True, 0.00000]) algebraic[50] = 1.00000*algebraic[48]*(constants[34]/constants[3])*algebraic[47] algebraic[29] = ((constants[0]*constants[1])/constants[2])*log(constants[12]/states[5]) algebraic[30] = 1.00000/(1.00000+1.49450*exp(0.0446000*states[0])) algebraic[31] = constants[11]*constants[99]*algebraic[30]*states[6]*(states[0]-algebraic[29]) algebraic[32] = ((constants[0]*constants[1])/constants[2])*log((constants[12]+0.0183300*constants[10])/(states[5]+0.0183300*states[1])) algebraic[33] = constants[14]*(power(states[7], 2.00000))*(states[0]-algebraic[32]) algebraic[34] = constants[15]*states[8]*states[9]*(states[0]-algebraic[29]) algebraic[35] = 1.00000/(2.00000+exp(((1.50000*constants[2])/(constants[0]*constants[1]))*(states[0]-algebraic[29]))) algebraic[36] = ((constants[16]*algebraic[35]*constants[12])/(constants[12]+constants[17]))*(states[0]-algebraic[29]) algebraic[37] = 1.00000/(1.00000+exp((7.48800-states[0])/5.98000)) algebraic[38] = constants[18]*algebraic[37]*(states[0]-algebraic[29]) algebraic[69] = power(1.00000-algebraic[68]/(algebraic[68]+constants[42]), 4.00000) algebraic[51] = power(states[18]/(states[18]+constants[43]), 2.00000) algebraic[52] = 1.00000-algebraic[51] algebraic[53] = 1.00000-power(algebraic[52], 4.00000) algebraic[73] = constants[45]*algebraic[69]*(1.00000-algebraic[53])+constants[46]*algebraic[69]*algebraic[53] algebraic[75] = constants[41]*algebraic[73]*(power(constants[12]/constants[44], 0.240000)) algebraic[79] = algebraic[75]*(states[0]-algebraic[29]) return algebraic def custom_piecewise(cases): """Compute result of a piecewise function""" return select(cases[0::2],cases[1::2]) def solve_model(): """Solve model with ODE solver""" from scipy.integrate import ode # Initialise constants and state variables (init_states, constants) = initConsts() # Set timespan to solve over voi = linspace(0, 10, 500) # Construct ODE object to solve r = ode(computeRates) r.set_integrator('vode', method='bdf', atol=1e-06, rtol=1e-06, max_step=1) r.set_initial_value(init_states, voi[0]) r.set_f_params(constants) # Solve model states = array([[0.0] * len(voi)] * sizeStates) states[:,0] = init_states for (i,t) in enumerate(voi[1:]): if r.successful(): r.integrate(t) states[:,i+1] = r.y else: break # Compute algebraic variables algebraic = computeAlgebraic(constants, states, voi) return (voi, states, algebraic) def plot_model(voi, states, algebraic): """Plot variables against variable of integration""" import pylab (legend_states, legend_algebraic, legend_voi, legend_constants) = createLegends() pylab.figure(1) pylab.plot(voi,vstack((states,algebraic)).T) pylab.xlabel(legend_voi) pylab.legend(legend_states + legend_algebraic, loc='best') pylab.show() if __name__ == "__main__": (voi, states, algebraic) = solve_model() plot_model(voi, states, algebraic)