module mechanism implicit none INTEGER, PARAMETER :: WP = SELECTED_REAL_KIND(14,300) REAL(KIND=WP) :: working_precision_variables integer, parameter :: nspec = 22 integer, parameter :: nreac = 173 integer, parameter :: isc_T = 1 integer, parameter :: neq = nspec + 1 ! QSS variables ! Number of QSS species integer, parameter :: nqss = 12 ! QSS species integer, dimension(nspec + nqss) :: iqss real(WP), dimension(nspec + nqss) :: W_sp,Cp_sp,h_sp,dh_sp character(len=15), dimension(nspec + nqss) :: gname character(len=5), dimension(nreac) :: rname ! Post processing variables ! Number of groups integer, parameter :: ng = 35 ! Max number of species in groups integer, parameter :: maxppn = 1 ! Number of species in each group integer, dimension(ng) :: ppn ! Species in each group integer, dimension(ng,maxppn) :: pp ! Name of species in each group character(len=30), dimension(ng) :: ppname ! Actual expression of each reaction character(len=65), dimension(nreac) :: reacexp ! Link between backward and forward rates integer, dimension(nreac) :: fofb ! Index of species integer, parameter :: sN2 = 1 integer, parameter :: sO = 2 integer, parameter :: sO2 = 3 integer, parameter :: sH = 4 integer, parameter :: sOH = 5 integer, parameter :: sH2 = 6 integer, parameter :: sH2O = 7 integer, parameter :: sH2O2 = 8 integer, parameter :: sHO2 = 9 integer, parameter :: sCO = 10 integer, parameter :: sCH2O = 11 integer, parameter :: sCH3 = 12 integer, parameter :: sC2H6 = 13 integer, parameter :: sCH4 = 14 integer, parameter :: sC2H4 = 15 integer, parameter :: sCO2 = 16 integer, parameter :: sCH3O2 = 17 integer, parameter :: sCH3O2H = 18 integer, parameter :: sC2H2 = 19 integer, parameter :: sC3H6 = 20 integer, parameter :: sC3H5O = 21 integer, parameter :: sC3H8 = 22 integer, parameter :: sCH2GSGXCH2 = 23 integer, parameter :: sCH3O = 24 integer, parameter :: sC2H5 = 25 integer, parameter :: sHCO = 26 integer, parameter :: sHCCO = 27 integer, parameter :: sC2H3 = 28 integer, parameter :: sCH2CHO = 29 integer, parameter :: sC3H5XAXC3H5 = 30 integer, parameter :: sIXC3H7 = 31 integer, parameter :: sNXC3H7 = 32 integer, parameter :: sIXC3H7O2 = 33 integer, parameter :: sNXC3H7O2 = 34 ! Index of reactions integer, parameter :: rH1 = 1 integer, parameter :: rH2 = 2 integer, parameter :: rH3f = 3 integer, parameter :: rH3b = 4 integer, parameter :: rH4 = 5 integer, parameter :: rH5f = 6 integer, parameter :: rH6f = 7 integer, parameter :: rH6b = 8 integer, parameter :: rH7f = 9 integer, parameter :: rH7b = 10 integer, parameter :: rH8f = 11 integer, parameter :: rH9f = 12 integer, parameter :: rH9b = 13 integer, parameter :: rH10f = 14 integer, parameter :: rH10b = 15 integer, parameter :: rH11f = 16 integer, parameter :: rH11b = 17 integer, parameter :: rH12 = 18 integer, parameter :: rH13 = 19 integer, parameter :: rH14 = 20 integer, parameter :: rH15 = 21 integer, parameter :: rH16f = 22 integer, parameter :: rH16b = 23 integer, parameter :: rH17 = 24 integer, parameter :: rH18 = 25 integer, parameter :: rH19f = 26 integer, parameter :: rH19b = 27 integer, parameter :: rH20 = 28 integer, parameter :: rH21 = 29 integer, parameter :: rH22 = 30 integer, parameter :: rH23 = 31 integer, parameter :: rH24f = 32 integer, parameter :: rH24b = 33 integer, parameter :: rH25f = 34 integer, parameter :: rH26f = 35 integer, parameter :: rH26b = 36 integer, parameter :: rH27 = 37 integer, parameter :: rH28f = 38 integer, parameter :: rH28b = 39 integer, parameter :: rH29 = 40 integer, parameter :: rH30f = 41 integer, parameter :: rH31 = 42 integer, parameter :: rH32 = 43 integer, parameter :: rH33 = 44 integer, parameter :: rH34 = 45 integer, parameter :: rH35 = 46 integer, parameter :: rH36 = 47 integer, parameter :: rH37f = 48 integer, parameter :: rH37b = 49 integer, parameter :: rH38f = 50 integer, parameter :: rH38b = 51 integer, parameter :: rH39f = 52 integer, parameter :: rH39b = 53 integer, parameter :: rH40f = 54 integer, parameter :: rH40b = 55 integer, parameter :: rH41 = 56 integer, parameter :: rH42 = 57 integer, parameter :: rH43f = 58 integer, parameter :: rH43b = 59 integer, parameter :: rH44f = 60 integer, parameter :: rH44b = 61 integer, parameter :: rH45 = 62 integer, parameter :: rH46 = 63 integer, parameter :: rH47 = 64 integer, parameter :: rH48f = 65 integer, parameter :: rH48b = 66 integer, parameter :: rH49 = 67 integer, parameter :: rH50 = 68 integer, parameter :: rH51 = 69 integer, parameter :: rH52 = 70 integer, parameter :: rH53 = 71 integer, parameter :: rH54 = 72 integer, parameter :: rH55 = 73 integer, parameter :: rH56 = 74 integer, parameter :: rH57 = 75 integer, parameter :: rH58 = 76 integer, parameter :: rH59 = 77 integer, parameter :: rH60 = 78 integer, parameter :: rH61 = 79 integer, parameter :: rH62 = 80 integer, parameter :: rH63f = 81 integer, parameter :: rH63b = 82 integer, parameter :: rH64 = 83 integer, parameter :: rH65 = 84 integer, parameter :: rH66f = 85 integer, parameter :: rH66b = 86 integer, parameter :: rH67 = 87 integer, parameter :: rH68 = 88 integer, parameter :: rH69 = 89 integer, parameter :: rH70 = 90 integer, parameter :: rH71 = 91 integer, parameter :: rH72 = 92 integer, parameter :: rH73 = 93 integer, parameter :: rH74 = 94 integer, parameter :: rH75f = 95 integer, parameter :: rH75b = 96 integer, parameter :: rH76 = 97 integer, parameter :: rH77f = 98 integer, parameter :: rH77b = 99 integer, parameter :: rH78f = 100 integer, parameter :: rH78b = 101 integer, parameter :: rH79 = 102 integer, parameter :: rH80f = 103 integer, parameter :: rH81 = 104 integer, parameter :: rH82 = 105 integer, parameter :: rH83 = 106 integer, parameter :: rH84 = 107 integer, parameter :: rH85f = 108 integer, parameter :: rH85b = 109 integer, parameter :: rH86f = 110 integer, parameter :: rH86b = 111 integer, parameter :: rH87 = 112 integer, parameter :: rH88 = 113 integer, parameter :: rH89 = 114 integer, parameter :: rH90 = 115 integer, parameter :: rH91 = 116 integer, parameter :: rH92 = 117 integer, parameter :: rH93f = 118 integer, parameter :: rH93b = 119 integer, parameter :: rH94 = 120 integer, parameter :: rH95 = 121 integer, parameter :: rH96 = 122 integer, parameter :: rH97 = 123 integer, parameter :: rH98 = 124 integer, parameter :: rH99 = 125 integer, parameter :: rH100 = 126 integer, parameter :: rH101 = 127 integer, parameter :: rH102 = 128 integer, parameter :: rH103f = 129 integer, parameter :: rH103b = 130 integer, parameter :: rH104 = 131 integer, parameter :: rH105 = 132 integer, parameter :: rH106f = 133 integer, parameter :: rH106b = 134 integer, parameter :: rH107f = 135 integer, parameter :: rH107b = 136 integer, parameter :: rH108 = 137 integer, parameter :: rH109f = 138 integer, parameter :: rH109b = 139 integer, parameter :: rH110 = 140 integer, parameter :: rH111 = 141 integer, parameter :: rH112f = 142 integer, parameter :: rH112b = 143 integer, parameter :: rH113 = 144 integer, parameter :: rH114f = 145 integer, parameter :: rH114b = 146 integer, parameter :: rH115f = 147 integer, parameter :: rH115b = 148 integer, parameter :: rH116 = 149 integer, parameter :: rH117 = 150 integer, parameter :: rH118 = 151 integer, parameter :: rH119 = 152 integer, parameter :: rH120 = 153 integer, parameter :: rH121 = 154 integer, parameter :: rH122 = 155 integer, parameter :: rH123 = 156 integer, parameter :: rH124f = 157 integer, parameter :: rH124b = 158 integer, parameter :: rH125 = 159 integer, parameter :: rH126f = 160 integer, parameter :: rH127f = 161 integer, parameter :: rH127b = 162 integer, parameter :: rH128 = 163 integer, parameter :: rH129f = 164 integer, parameter :: rH129b = 165 integer, parameter :: rH130f = 166 integer, parameter :: rH130b = 167 integer, parameter :: rH5b = 168 integer, parameter :: rH8b = 169 integer, parameter :: rH25b = 170 integer, parameter :: rH30b = 171 integer, parameter :: rH80b = 172 integer, parameter :: rH126b = 173 ! Index of third bodies integer, parameter :: mM11 = 1 integer, parameter :: mM5 = 2 integer, parameter :: mM6 = 3 integer, parameter :: mM4 = 4 integer, parameter :: mM9 = 5 integer, parameter :: mM12 = 6 integer, parameter :: mM2 = 7 integer, parameter :: mM1 = 8 integer, parameter :: mM8 = 9 integer, parameter :: mM7 = 10 integer, parameter :: mM3 = 11 integer, parameter :: mM0 = 12 contains ! Name of mechanism ! subroutine which_mechanism ! use parallel !implicit none ! if (irank.eq.iroot) print*, 'Using mechanism Spec_34_Reac_173_QSS_12.mech' ! return ! end subroutine which_mechanism ! Subroutine to define groups for post-processing subroutine pp_data implicit none ! Number of species in each group ppn(1) = 1 ppn(2) = 1 ppn(3) = 1 ppn(4) = 1 ppn(5) = 1 ppn(6) = 1 ppn(7) = 1 ppn(8) = 1 ppn(9) = 1 ppn(10) = 1 ppn(11) = 1 ppn(12) = 1 ppn(13) = 1 ppn(14) = 1 ppn(15) = 1 ppn(16) = 1 ppn(17) = 1 ppn(18) = 1 ppn(19) = 1 ppn(20) = 1 ppn(21) = 1 ppn(22) = 1 ppn(23) = 1 ppn(24) = 1 ppn(25) = 1 ppn(26) = 1 ppn(27) = 1 ppn(28) = 1 ppn(29) = 1 ppn(30) = 1 ppn(31) = 1 ppn(32) = 1 ppn(33) = 1 ppn(34) = 1 ppn(35) = 1 ! Indices of species in each group pp(1,1) = 1 pp(2,1) = 2 pp(3,1) = 3 pp(4,1) = 4 pp(5,1) = 5 pp(6,1) = 6 pp(7,1) = 7 pp(8,1) = 8 pp(9,1) = 9 pp(10,1) = 10 pp(11,1) = 11 pp(12,1) = 12 pp(13,1) = 13 pp(14,1) = 14 pp(15,1) = 15 pp(16,1) = 16 pp(17,1) = 17 pp(18,1) = 18 pp(19,1) = 19 pp(20,1) = 20 pp(21,1) = 21 pp(22,1) = 22 pp(23,1) = 23 pp(24,1) = 24 pp(25,1) = 25 pp(26,1) = 26 pp(27,1) = 27 pp(28,1) = 28 pp(29,1) = 29 pp(30,1) = 30 pp(31,1) = 31 pp(32,1) = 32 pp(33,1) = 33 pp(34,1) = 34 pp(35,1) = 1 ! Name of group of species ppname(1) = trim(gname(sN2)) ppname(2) = trim(gname(sO)) ppname(3) = trim(gname(sO2)) ppname(4) = trim(gname(sH)) ppname(5) = trim(gname(sOH)) ppname(6) = trim(gname(sH2)) ppname(7) = trim(gname(sH2O)) ppname(8) = trim(gname(sH2O2)) ppname(9) = trim(gname(sHO2)) ppname(10) = trim(gname(sCH2GSGXCH2)) ppname(11) = trim(gname(sCO)) ppname(12) = trim(gname(sCH2O)) ppname(13) = trim(gname(sCH3)) ppname(14) = trim(gname(sC2H6)) ppname(15) = trim(gname(sCH3O)) ppname(16) = trim(gname(sCH4)) ppname(17) = trim(gname(sC2H5)) ppname(18) = trim(gname(sC2H4)) ppname(19) = trim(gname(sCO2)) ppname(20) = trim(gname(sHCO)) ppname(21) = trim(gname(sCH3O2)) ppname(22) = trim(gname(sCH3O2H)) ppname(23) = trim(gname(sC2H2)) ppname(24) = trim(gname(sHCCO)) ppname(25) = trim(gname(sC2H3)) ppname(26) = trim(gname(sCH2CHO)) ppname(27) = trim(gname(sC3H5XAXC3H5)) ppname(28) = trim(gname(sC3H6)) ppname(29) = trim(gname(sC3H5O)) ppname(30) = trim(gname(sIXC3H7)) ppname(31) = trim(gname(sNXC3H7)) ppname(32) = trim(gname(sC3H8)) ppname(33) = trim(gname(sIXC3H7O2)) ppname(34) = trim(gname(sNXC3H7O2)) ppname(35) = 'N2X' return end subroutine pp_data ! Molar mass subroutine molar_mass implicit none W_sp(sN2) = 0.02802_WP W_sp(sO) = 0.016_WP W_sp(sO2) = 0.032_WP W_sp(sH) = 0.001008_WP W_sp(sOH) = 0.017008_WP W_sp(sH2) = 0.002016_WP W_sp(sH2O) = 0.018016_WP W_sp(sH2O2) = 0.034016_WP W_sp(sHO2) = 0.033008_WP W_sp(sCH2GSGXCH2) = 0.014026_WP W_sp(sCO) = 0.02801_WP W_sp(sCH2O) = 0.030026_WP W_sp(sCH3) = 0.015034_WP W_sp(sC2H6) = 0.030068_WP W_sp(sCH3O) = 0.031034_WP W_sp(sCH4) = 0.016042_WP W_sp(sC2H5) = 0.02906_WP W_sp(sC2H4) = 0.028052_WP W_sp(sCO2) = 0.04401_WP W_sp(sHCO) = 0.029018_WP W_sp(sCH3O2) = 0.047034_WP W_sp(sCH3O2H) = 0.048042_WP W_sp(sC2H2) = 0.026036_WP W_sp(sHCCO) = 0.041028_WP W_sp(sC2H3) = 0.027044_WP W_sp(sCH2CHO) = 0.043044_WP W_sp(sC3H5XAXC3H5) = 0.04107_WP W_sp(sC3H6) = 0.042078_WP W_sp(sC3H5O) = 0.05707_WP W_sp(sIXC3H7) = 0.043086_WP W_sp(sNXC3H7) = 0.043086_WP W_sp(sC3H8) = 0.044094_WP W_sp(sIXC3H7O2) = 0.075086_WP W_sp(sNXC3H7O2) = 0.075086_WP return end subroutine molar_mass ! Species names subroutine species_name implicit none gname(sN2) = 'N2' gname(sO) = 'O' gname(sO2) = 'O2' gname(sH) = 'H' gname(sOH) = 'OH' gname(sH2) = 'H2' gname(sH2O) = 'H2O' gname(sH2O2) = 'H2O2' gname(sHO2) = 'HO2' gname(sCH2GSGXCH2) = 'CH2GSG-CH2' gname(sCO) = 'CO' gname(sCH2O) = 'CH2O' gname(sCH3) = 'CH3' gname(sC2H6) = 'C2H6' gname(sCH3O) = 'CH3O' gname(sCH4) = 'CH4' gname(sC2H5) = 'C2H5' gname(sC2H4) = 'C2H4' gname(sCO2) = 'CO2' gname(sHCO) = 'HCO' gname(sCH3O2) = 'CH3O2' gname(sCH3O2H) = 'CH3O2H' gname(sC2H2) = 'C2H2' gname(sHCCO) = 'HCCO' gname(sC2H3) = 'C2H3' gname(sCH2CHO) = 'CH2CHO' gname(sC3H5XAXC3H5) = 'C3H5-A-C3H5' gname(sC3H6) = 'C3H6' gname(sC3H5O) = 'C3H5O' gname(sIXC3H7) = 'I-C3H7' gname(sNXC3H7) = 'N-C3H7' gname(sC3H8) = 'C3H8' gname(sIXC3H7O2) = 'I-C3H7O2' gname(sNXC3H7O2) = 'N-C3H7O2' return end subroutine species_name ! Reaction names subroutine reaction_name implicit none rname(rH1) = 'H1' rname(rH2) = 'H2' rname(rH3f) = 'H3f' rname(rH3b) = 'H3b' rname(rH4) = 'H4' rname(rH5f) = 'H5f' rname(rH6f) = 'H6f' rname(rH6b) = 'H6b' rname(rH7f) = 'H7f' rname(rH7b) = 'H7b' rname(rH8f) = 'H8f' rname(rH9f) = 'H9f' rname(rH9b) = 'H9b' rname(rH10f) = 'H10f' rname(rH10b) = 'H10b' rname(rH11f) = 'H11f' rname(rH11b) = 'H11b' rname(rH12) = 'H12' rname(rH13) = 'H13' rname(rH14) = 'H14' rname(rH15) = 'H15' rname(rH16f) = 'H16f' rname(rH16b) = 'H16b' rname(rH17) = 'H17' rname(rH18) = 'H18' rname(rH19f) = 'H19f' rname(rH19b) = 'H19b' rname(rH20) = 'H20' rname(rH21) = 'H21' rname(rH22) = 'H22' rname(rH23) = 'H23' rname(rH24f) = 'H24f' rname(rH24b) = 'H24b' rname(rH25f) = 'H25f' rname(rH26f) = 'H26f' rname(rH26b) = 'H26b' rname(rH27) = 'H27' rname(rH28f) = 'H28f' rname(rH28b) = 'H28b' rname(rH29) = 'H29' rname(rH30f) = 'H30f' rname(rH31) = 'H31' rname(rH32) = 'H32' rname(rH33) = 'H33' rname(rH34) = 'H34' rname(rH35) = 'H35' rname(rH36) = 'H36' rname(rH37f) = 'H37f' rname(rH37b) = 'H37b' rname(rH38f) = 'H38f' rname(rH38b) = 'H38b' rname(rH39f) = 'H39f' rname(rH39b) = 'H39b' rname(rH40f) = 'H40f' rname(rH40b) = 'H40b' rname(rH41) = 'H41' rname(rH42) = 'H42' rname(rH43f) = 'H43f' rname(rH43b) = 'H43b' rname(rH44f) = 'H44f' rname(rH44b) = 'H44b' rname(rH45) = 'H45' rname(rH46) = 'H46' rname(rH47) = 'H47' rname(rH48f) = 'H48f' rname(rH48b) = 'H48b' rname(rH49) = 'H49' rname(rH50) = 'H50' rname(rH51) = 'H51' rname(rH52) = 'H52' rname(rH53) = 'H53' rname(rH54) = 'H54' rname(rH55) = 'H55' rname(rH56) = 'H56' rname(rH57) = 'H57' rname(rH58) = 'H58' rname(rH59) = 'H59' rname(rH60) = 'H60' rname(rH61) = 'H61' rname(rH62) = 'H62' rname(rH63f) = 'H63f' rname(rH63b) = 'H63b' rname(rH64) = 'H64' rname(rH65) = 'H65' rname(rH66f) = 'H66f' rname(rH66b) = 'H66b' rname(rH67) = 'H67' rname(rH68) = 'H68' rname(rH69) = 'H69' rname(rH70) = 'H70' rname(rH71) = 'H71' rname(rH72) = 'H72' rname(rH73) = 'H73' rname(rH74) = 'H74' rname(rH75f) = 'H75f' rname(rH75b) = 'H75b' rname(rH76) = 'H76' rname(rH77f) = 'H77f' rname(rH77b) = 'H77b' rname(rH78f) = 'H78f' rname(rH78b) = 'H78b' rname(rH79) = 'H79' rname(rH80f) = 'H80f' rname(rH81) = 'H81' rname(rH82) = 'H82' rname(rH83) = 'H83' rname(rH84) = 'H84' rname(rH85f) = 'H85f' rname(rH85b) = 'H85b' rname(rH86f) = 'H86f' rname(rH86b) = 'H86b' rname(rH87) = 'H87' rname(rH88) = 'H88' rname(rH89) = 'H89' rname(rH90) = 'H90' rname(rH91) = 'H91' rname(rH92) = 'H92' rname(rH93f) = 'H93f' rname(rH93b) = 'H93b' rname(rH94) = 'H94' rname(rH95) = 'H95' rname(rH96) = 'H96' rname(rH97) = 'H97' rname(rH98) = 'H98' rname(rH99) = 'H99' rname(rH100) = 'H100' rname(rH101) = 'H101' rname(rH102) = 'H102' rname(rH103f) = 'H103f' rname(rH103b) = 'H103b' rname(rH104) = 'H104' rname(rH105) = 'H105' rname(rH106f) = 'H106f' rname(rH106b) = 'H106b' rname(rH107f) = 'H107f' rname(rH107b) = 'H107b' rname(rH108) = 'H108' rname(rH109f) = 'H109f' rname(rH109b) = 'H109b' rname(rH110) = 'H110' rname(rH111) = 'H111' rname(rH112f) = 'H112f' rname(rH112b) = 'H112b' rname(rH113) = 'H113' rname(rH114f) = 'H114f' rname(rH114b) = 'H114b' rname(rH115f) = 'H115f' rname(rH115b) = 'H115b' rname(rH116) = 'H116' rname(rH117) = 'H117' rname(rH118) = 'H118' rname(rH119) = 'H119' rname(rH120) = 'H120' rname(rH121) = 'H121' rname(rH122) = 'H122' rname(rH123) = 'H123' rname(rH124f) = 'H124f' rname(rH124b) = 'H124b' rname(rH125) = 'H125' rname(rH126f) = 'H126f' rname(rH127f) = 'H127f' rname(rH127b) = 'H127b' rname(rH128) = 'H128' rname(rH129f) = 'H129f' rname(rH129b) = 'H129b' rname(rH130f) = 'H130f' rname(rH130b) = 'H130b' rname(rH5b) = 'H5b' rname(rH8b) = 'H8b' rname(rH25b) = 'H25b' rname(rH30b) = 'H30b' rname(rH80b) = 'H80b' rname(rH126b) = 'H126b' return end subroutine reaction_name ! List of QSS species subroutine QSS_list implicit none iqss(sN2) = 0 iqss(sO) = 0 iqss(sO2) = 0 iqss(sH) = 0 iqss(sOH) = 0 iqss(sH2) = 0 iqss(sH2O) = 0 iqss(sH2O2) = 0 iqss(sHO2) = 0 iqss(sCH2GSGXCH2) = 1 iqss(sCO) = 0 iqss(sCH2O) = 0 iqss(sCH3) = 0 iqss(sC2H6) = 0 iqss(sCH3O) = 1 iqss(sCH4) = 0 iqss(sC2H5) = 1 iqss(sC2H4) = 0 iqss(sCO2) = 0 iqss(sHCO) = 1 iqss(sCH3O2) = 0 iqss(sCH3O2H) = 0 iqss(sC2H2) = 0 iqss(sHCCO) = 1 iqss(sC2H3) = 1 iqss(sCH2CHO) = 1 iqss(sC3H5XAXC3H5) = 1 iqss(sC3H6) = 0 iqss(sC3H5O) = 0 iqss(sIXC3H7) = 1 iqss(sNXC3H7) = 1 iqss(sC3H8) = 0 iqss(sIXC3H7O2) = 1 iqss(sNXC3H7O2) = 1 return end subroutine QSS_list ! Subroutine for pressure dependent coefficients real(WP) function getlindratecoeff (Tloc,k0,kinf,fc,concin,Ploc) implicit none real(WP) :: Tloc,k0,kinf,fc,Ploc real(WP), parameter :: R = 8.31434_WP real(WP) :: ntmp,ccoeff,dcoeff,lgknull real(WP) :: f real(WP) :: conc, concin if (concin.gt.0.0_WP) then conc = concin else conc = Ploc / ( R * Tloc ) end if ntmp = 0.75_WP - 1.27_WP * dlog10( fc ) ccoeff = - 0.4_WP - 0.67_WP * dlog10( fc ) dcoeff = 0.14_WP k0 = k0 * conc / max(kinf, 1.0e-60_WP) lgknull = dlog10(k0) f = (lgknull+ccoeff)/(ntmp-dcoeff*(lgknull+ccoeff)) f = fc**(1.0_WP / ( f * f + 1.0_WP )) getlindratecoeff = kinf * f * k0 / ( 1.0_WP + k0 ) end function getlindratecoeff end module mechanism subroutine TEST(P, T, Y, WDOT) use mechanism implicit none real(WP), dimension(nspec) :: Y,c,WDOT real(WP), dimension(nspec + nqss) :: H,dH,Cp real(WP), dimension(nqss) :: cqss real(WP), dimension(nreac) :: w,k real(WP), dimension(13) :: m real(WP) :: T,P CALL compute_thermodata(H,Cp,dH,T) CALL YtoC(c,P,T,Y) CALL get_thirdbodies(M,c) CALL get_rate_coefficients(k,M,T,P) CALL get_QSS(cqss,c,k,M) CALL get_reaction_rates(w,k,M,c,cqss) CALL get_production_rates(WDOT,w) return end subroutine TEST subroutine YtoC(c,P,T,Y) use mechanism implicit none real(WP), dimension(nspec) :: c,Y real(WP) :: bla, P, T integer :: K call molar_mass c(sN2) = Y(sN2) / W_sp(sN2) c(sO) = Y(sO) / W_sp(sO) c(sO2) = Y(sO2) / W_sp(sO2) c(sH) = Y(sH) / W_sp(sH) c(sOH) = Y(sOH) / W_sp(sOH) c(sH2) = Y(sH2) / W_sp(sH2) c(sH2O) = Y(sH2O) / W_sp(sH2O) c(sH2O2) = Y(sH2O2) / W_sp(sH2O2) c(sHO2) = Y(sHO2) / W_sp(sHO2) !c(sCH2GSGXCH2) = Y(sCH2GSGXCH2) / W_sp(sCH2GSGXCH2) c(sCO) = Y(sCO) / W_sp(sCO) c(sCH2O) = Y(sCH2O) / W_sp(sCH2O) c(sCH3) = Y(sCH3) / W_sp(sCH3) c(sC2H6) = Y(sC2H6) / W_sp(sC2H6) !c(sCH3O) = Y(sCH3O) / W_sp(sCH3O) c(sCH4) = Y(sCH4) / W_sp(sCH4) !c(sC2H5) = Y(sC2H5) / W_sp(sC2H5) c(sC2H4) = Y(sC2H4) / W_sp(sC2H4) c(sCO2) = Y(sCO2) / W_sp(sCO2) !c(sHCO) = Y(sHCO) / W_sp(sHCO) c(sCH3O2) = Y(sCH3O2) / W_sp(sCH3O2) c(sCH3O2H) = Y(sCH3O2H) / W_sp(sCH3O2H) c(sC2H2) = Y(sC2H2) / W_sp(sC2H2) !c(sHCCO) = Y(sHCCO) / W_sp(sHCCO) !c(sC2H3) = Y(sC2H3) / W_sp(sC2H3) !c(sCH2CHO) = Y(sCH2CHO) / W_sp(sCH2CHO) !c(sC3H5XAXC3H5) = Y(sC3H5XAXC3H5) / W_sp(sC3H5XAXC3H5) c(sC3H6) = Y(sC3H6) / W_sp(sC3H6) c(sC3H5O) = Y(sC3H5O) / W_sp(sC3H5O) !c(sIXC3H7) = Y(sIXC3H7) / W_sp(sIXC3H7) !c(sNXC3H7) = Y(sNXC3H7) / W_sp(sNXC3H7) c(sC3H8) = Y(sC3H8) / W_sp(sC3H8) !c(sIXC3H7O2) = Y(sIXC3H7O2) / W_sp(sIXC3H7O2) !c(sNXC3H7O2) = Y(sNXC3H7O2) / W_sp(sNXC3H7O2) bla = 0.0_WP DO K = 1, nspec bla = bla + c(K) ENDDO bla = P/(bla*T*8.31451_WP) DO K = 1, nspec c(K) = max(c(K),1.0e-60_WP) * bla ENDDO return end subroutine YtoC ! Subroutine for Cp and H computations subroutine compute_thermodata(H,Cp,dH,T) use mechanism implicit none real(WP), dimension(nspec + nqss) :: H,Cp,dH real(WP) :: T if (T.gt.1000_WP) then H(sN2) = 296.716630977873 * ( T * ((2.92664000e+00_WP) & + T *((0.0007439885_WP) & + T *((-1.89492033333333e-07_WP) & + T *((2.52426e-11_WP) & + T *((-1.3506702e-15_WP)))))) + (-9.22797700e+02_WP)) Cp(sN2) = 296.716630977873 * ((2.92664000e+00_WP) & + T *((1.48797700e-03_WP) & + T *((-5.68476100e-07_WP) & + T *((1.00970400e-10_WP) & + T *(-6.75335100e-15_WP))))) dH(sN2) = 296.716630977873 * ( (2.92664000e+00_WP) & + T *((1.48797700e-03_WP) & + T *((-5.68476100e-07_WP) & + T *((1.00970400e-10_WP) & + T *((-6.75335100e-15_WP))))) ) H(sO) = 519.625 * ( T * ((2.54206000e+00_WP) & + T *((-1.377531e-05_WP) & + T *((-1.03426766666667e-09_WP) & + T *((1.13776675e-12_WP) & + T *((-8.736104e-17_WP)))))) + (2.92308000e+04_WP)) Cp(sO) = 519.625 * ((2.54206000e+00_WP) & + T *((-2.75506200e-05_WP) & + T *((-3.10280300e-09_WP) & + T *((4.55106700e-12_WP) & + T *(-4.36805200e-16_WP))))) dH(sO) = 519.625 * ( (2.54206000e+00_WP) & + T *((-2.75506200e-05_WP) & + T *((-3.10280300e-09_WP) & + T *((4.55106700e-12_WP) & + T *((-4.36805200e-16_WP))))) ) H(sO2) = 259.8125 * ( T * ((3.69757800e+00_WP) & + T *((0.00030675985_WP) & + T *((-4.19614e-08_WP) & + T *((4.4382025e-12_WP) & + T *((-2.27287e-16_WP)))))) + (-1.23393000e+03_WP)) Cp(sO2) = 259.8125 * ((3.69757800e+00_WP) & + T *((6.13519700e-04_WP) & + T *((-1.25884200e-07_WP) & + T *((1.77528100e-11_WP) & + T *(-1.13643500e-15_WP))))) dH(sO2) = 259.8125 * ( (3.69757800e+00_WP) & + T *((6.13519700e-04_WP) & + T *((-1.25884200e-07_WP) & + T *((1.77528100e-11_WP) & + T *((-1.13643500e-15_WP))))) ) H(sH) = 8248.01587301587 * ( T * ((2.50104422e+00_WP) & + T *((0_WP) & + T *((0_WP) & + T *((0_WP) & + T *((0_WP)))))) + (2.54747466e+04_WP)) Cp(sH) = 8248.01587301587 * ((2.50104422e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *(0.00000000e+00_WP))))) dH(sH) = 8248.01587301587 * ( (2.50104422e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP))))) ) H(sOH) = 488.828786453434 * ( T * ((3.69039275e+00_WP) & + T *((-0.000290401159_WP) & + T *((2.54434418666667e-07_WP) & + T *((-4.506551475e-11_WP) & + T *((0_WP)))))) + (3.67872158e+03_WP)) Cp(sOH) = 488.828786453434 * ((3.69039275e+00_WP) & + T *((-5.80802318e-04_WP) & + T *((7.63303256e-07_WP) & + T *((-1.80262059e-10_WP) & + T *(0.00000000e+00_WP))))) dH(sOH) = 488.828786453434 * ( (3.69039275e+00_WP) & + T *((-5.80802318e-04_WP) & + T *((7.63303256e-07_WP) & + T *((-1.80262059e-10_WP) & + T *((0.00000000e+00_WP))))) ) H(sH2) = 4124.00793650794 * ( T * ((2.99142300e+00_WP) & + T *((0.0003500322_WP) & + T *((-1.877943e-08_WP) & + T *((-2.3078945e-12_WP) & + T *((3.165504e-16_WP)))))) + (-8.35034000e+02_WP)) Cp(sH2) = 4124.00793650794 * ((2.99142300e+00_WP) & + T *((7.00064400e-04_WP) & + T *((-5.63382900e-08_WP) & + T *((-9.23157800e-12_WP) & + T *(1.58275200e-15_WP))))) dH(sH2) = 4124.00793650794 * ( (2.99142300e+00_WP) & + T *((7.00064400e-04_WP) & + T *((-5.63382900e-08_WP) & + T *((-9.23157800e-12_WP) & + T *((1.58275200e-15_WP))))) ) H(sH2O) = 461.478685612789 * ( T * ((2.67214600e+00_WP) & + T *((0.0015281465_WP) & + T *((-2.91008666666667e-07_WP) & + T *((3.00249e-11_WP) & + T *((-1.2783236e-15_WP)))))) + (-2.98992100e+04_WP)) Cp(sH2O) = 461.478685612789 * ((2.67214600e+00_WP) & + T *((3.05629300e-03_WP) & + T *((-8.73026000e-07_WP) & + T *((1.20099600e-10_WP) & + T *(-6.39161800e-15_WP))))) dH(sH2O) = 461.478685612789 * ( (2.67214600e+00_WP) & + T *((3.05629300e-03_WP) & + T *((-8.73026000e-07_WP) & + T *((1.20099600e-10_WP) & + T *((-6.39161800e-15_WP))))) ) H(sH2O2) = 244.414393226717 * ( T * ((4.57316700e+00_WP) & + T *((0.002168068_WP) & + T *((-4.91563e-07_WP) & + T *((5.87226e-11_WP) & + T *((-2.863308e-15_WP)))))) + (-1.80069600e+04_WP)) Cp(sH2O2) = 244.414393226717 * ((4.57316700e+00_WP) & + T *((4.33613600e-03_WP) & + T *((-1.47468900e-06_WP) & + T *((2.34890400e-10_WP) & + T *(-1.43165400e-14_WP))))) dH(sH2O2) = 244.414393226717 * ( (4.57316700e+00_WP) & + T *((4.33613600e-03_WP) & + T *((-1.47468900e-06_WP) & + T *((2.34890400e-10_WP) & + T *((-1.43165400e-14_WP))))) ) H(sHO2) = 251.878332525448 * ( T * ((3.00201579e+00_WP) & + T *((0.00219276212_WP) & + T *((-6.58647396666667e-07_WP) & + T *((8.302252425e-11_WP) & + T *((0_WP)))))) + (8.44944001e+02_WP)) Cp(sHO2) = 251.878332525448 * ((3.00201579e+00_WP) & + T *((4.38552424e-03_WP) & + T *((-1.97594219e-06_WP) & + T *((3.32090097e-10_WP) & + T *(0.00000000e+00_WP))))) dH(sHO2) = 251.878332525448 * ( (3.00201579e+00_WP) & + T *((4.38552424e-03_WP) & + T *((-1.97594219e-06_WP) & + T *((3.32090097e-10_WP) & + T *((0.00000000e+00_WP))))) ) H(sCH2GSGXCH2) = 592.756309710538 * ( T * ((3.55288900e+00_WP) & + T *((0.001033394_WP) & + T *((-6.38038666666667e-08_WP) & + T *((-2.7616825e-11_WP) & + T *((4.0427e-15_WP)))))) + (4.98497500e+04_WP)) Cp(sCH2GSGXCH2) = 592.756309710538 * ((3.55288900e+00_WP) & + T *((2.06678800e-03_WP) & + T *((-1.91411600e-07_WP) & + T *((-1.10467300e-10_WP) & + T *(2.02135000e-14_WP))))) dH(sCH2GSGXCH2) = 592.756309710538 * ( (3.55288900e+00_WP) & + T *((2.06678800e-03_WP) & + T *((-1.91411600e-07_WP) & + T *((-1.10467300e-10_WP) & + T *((2.02135000e-14_WP))))) ) H(sCO) = 296.822563370225 * ( T * ((3.02507800e+00_WP) & + T *((0.0007213445_WP) & + T *((-1.87694266666667e-07_WP) & + T *((2.5464525e-11_WP) & + T *((-1.3821904e-15_WP)))))) + (-1.42683500e+04_WP)) Cp(sCO) = 296.822563370225 * ((3.02507800e+00_WP) & + T *((1.44268900e-03_WP) & + T *((-5.63082800e-07_WP) & + T *((1.01858100e-10_WP) & + T *(-6.91095200e-15_WP))))) dH(sCO) = 296.822563370225 * ( (3.02507800e+00_WP) & + T *((1.44268900e-03_WP) & + T *((-5.63082800e-07_WP) & + T *((1.01858100e-10_WP) & + T *((-6.91095200e-15_WP))))) ) H(sCH2O) = 276.89335908879 * ( T * ((2.99560600e+00_WP) & + T *((0.0033406605_WP) & + T *((-8.76318333333333e-07_WP) & + T *((1.18428825e-10_WP) & + T *((-6.425034e-15_WP)))))) + (-1.53203700e+04_WP)) Cp(sCH2O) = 276.89335908879 * ((2.99560600e+00_WP) & + T *((6.68132100e-03_WP) & + T *((-2.62895500e-06_WP) & + T *((4.73715300e-10_WP) & + T *(-3.21251700e-14_WP))))) dH(sCH2O) = 276.89335908879 * ( (2.99560600e+00_WP) & + T *((6.68132100e-03_WP) & + T *((-2.62895500e-06_WP) & + T *((4.73715300e-10_WP) & + T *((-3.21251700e-14_WP))))) ) H(sCH3) = 553.013170147665 * ( T * ((2.84405200e+00_WP) & + T *((0.003068987_WP) & + T *((-7.43448333333333e-07_WP) & + T *((9.4629025e-11_WP) & + T *((-4.904318e-15_WP)))))) + (1.64378100e+04_WP)) Cp(sCH3) = 553.013170147665 * ((2.84405200e+00_WP) & + T *((6.13797400e-03_WP) & + T *((-2.23034500e-06_WP) & + T *((3.78516100e-10_WP) & + T *(-2.45215900e-14_WP))))) dH(sCH3) = 553.013170147665 * ( (2.84405200e+00_WP) & + T *((6.13797400e-03_WP) & + T *((-2.23034500e-06_WP) & + T *((3.78516100e-10_WP) & + T *((-2.45215900e-14_WP))))) ) H(sC2H6) = 276.506585073833 * ( T * ((2.78247515e-01_WP) & + T *((0.0115844838_WP) & + T *((-3.41064786666667e-06_WP) & + T *((4.1924835e-10_WP) & + T *((0_WP)))))) + (-1.13124605e+04_WP)) Cp(sC2H6) = 276.506585073833 * ((2.78247515e-01_WP) & + T *((2.31689676e-02_WP) & + T *((-1.02319436e-05_WP) & + T *((1.67699340e-09_WP) & + T *(0.00000000e+00_WP))))) dH(sC2H6) = 276.506585073833 * ( (2.78247515e-01_WP) & + T *((2.31689676e-02_WP) & + T *((-1.02319436e-05_WP) & + T *((1.67699340e-09_WP) & + T *((0.00000000e+00_WP))))) ) H(sCH3O) = 267.899722884578 * ( T * ((3.77080000e+00_WP) & + T *((0.0039357485_WP) & + T *((-8.85461333333333e-07_WP) & + T *((9.8610775e-11_WP) & + T *((-4.225232e-15_WP)))))) + (1.27832500e+02_WP)) Cp(sCH3O) = 267.899722884578 * ((3.77080000e+00_WP) & + T *((7.87149700e-03_WP) & + T *((-2.65638400e-06_WP) & + T *((3.94443100e-10_WP) & + T *(-2.11261600e-14_WP))))) dH(sCH3O) = 267.899722884578 * ( (3.77080000e+00_WP) & + T *((7.87149700e-03_WP) & + T *((-2.65638400e-06_WP) & + T *((3.94443100e-10_WP) & + T *((-2.11261600e-14_WP))))) ) H(sCH4) = 518.264555541703 * ( T * ((1.68347900e+00_WP) & + T *((0.00511862_WP) & + T *((-1.29170966666667e-06_WP) & + T *((1.69639625e-10_WP) & + T *((-9.006846e-15_WP)))))) + (-1.00807900e+04_WP)) Cp(sCH4) = 518.264555541703 * ((1.68347900e+00_WP) & + T *((1.02372400e-02_WP) & + T *((-3.87512900e-06_WP) & + T *((6.78558500e-10_WP) & + T *(-4.50342300e-14_WP))))) dH(sCH4) = 518.264555541703 * ( (1.68347900e+00_WP) & + T *((1.02372400e-02_WP) & + T *((-3.87512900e-06_WP) & + T *((6.78558500e-10_WP) & + T *((-4.50342300e-14_WP))))) ) H(sC2H5) = 286.097728836889 * ( T * ((7.19048000e+00_WP) & + T *((0.0032420385_WP) & + T *((-2.14268833333333e-07_WP) & + T *((-5.8696975e-11_WP) & + T *((7.761754e-15_WP)))))) + (1.06745500e+04_WP)) Cp(sC2H5) = 286.097728836889 * ((7.19048000e+00_WP) & + T *((6.48407700e-03_WP) & + T *((-6.42806500e-07_WP) & + T *((-2.34787900e-10_WP) & + T *(3.88087700e-14_WP))))) dH(sC2H5) = 286.097728836889 * ( (7.19048000e+00_WP) & + T *((6.48407700e-03_WP) & + T *((-6.42806500e-07_WP) & + T *((-2.34787900e-10_WP) & + T *((3.88087700e-14_WP))))) ) H(sC2H4) = 296.378154855269 * ( T * ((3.52841900e+00_WP) & + T *((0.00574259_WP) & + T *((-1.472795e-06_WP) & + T *((1.96115025e-10_WP) & + T *((-1.0533696e-14_WP)))))) + (4.42828900e+03_WP)) Cp(sC2H4) = 296.378154855269 * ((3.52841900e+00_WP) & + T *((1.14851800e-02_WP) & + T *((-4.41838500e-06_WP) & + T *((7.84460100e-10_WP) & + T *(-5.26684800e-14_WP))))) dH(sC2H4) = 296.378154855269 * ( (3.52841900e+00_WP) & + T *((1.14851800e-02_WP) & + T *((-4.41838500e-06_WP) & + T *((7.84460100e-10_WP) & + T *((-5.26684800e-14_WP))))) ) H(sCO2) = 188.911610997501 * ( T * ((4.45362300e+00_WP) & + T *((0.0015700845_WP) & + T *((-4.26137e-07_WP) & + T *((5.9849925e-11_WP) & + T *((-3.338066e-15_WP)))))) + (-4.89669600e+04_WP)) Cp(sCO2) = 188.911610997501 * ((4.45362300e+00_WP) & + T *((3.14016900e-03_WP) & + T *((-1.27841100e-06_WP) & + T *((2.39399700e-10_WP) & + T *(-1.66903300e-14_WP))))) dH(sCO2) = 188.911610997501 * ( (4.45362300e+00_WP) & + T *((3.14016900e-03_WP) & + T *((-1.27841100e-06_WP) & + T *((2.39399700e-10_WP) & + T *((-1.66903300e-14_WP))))) ) H(sHCO) = 286.5118202495 * ( T * ((3.55727100e+00_WP) & + T *((0.0016727865_WP) & + T *((-4.45002e-07_WP) & + T *((6.1764325e-11_WP) & + T *((-3.427702e-15_WP)))))) + (3.91632400e+03_WP)) Cp(sHCO) = 286.5118202495 * ((3.55727100e+00_WP) & + T *((3.34557300e-03_WP) & + T *((-1.33500600e-06_WP) & + T *((2.47057300e-10_WP) & + T *(-1.71385100e-14_WP))))) dH(sHCO) = 286.5118202495 * ( (3.55727100e+00_WP) & + T *((3.34557300e-03_WP) & + T *((-1.33500600e-06_WP) & + T *((2.47057300e-10_WP) & + T *((-1.71385100e-14_WP))))) ) H(sCH3O2) = 176.765743929923 * ( T * ((3.72224006e+00_WP) & + T *((0.00577152465_WP) & + T *((-1.51686254333333e-06_WP) & + T *((1.67577132e-10_WP) & + T *((0_WP)))))) + (-5.34692750e+02_WP)) Cp(sCH3O2) = 176.765743929923 * ((3.72224006e+00_WP) & + T *((1.15430493e-02_WP) & + T *((-4.55058763e-06_WP) & + T *((6.70308528e-10_WP) & + T *(0.00000000e+00_WP))))) dH(sCH3O2) = 176.765743929923 * ( (3.72224006e+00_WP) & + T *((1.15430493e-02_WP) & + T *((-4.55058763e-06_WP) & + T *((6.70308528e-10_WP) & + T *((0.00000000e+00_WP))))) ) H(sCH3O2H) = 173.056908538362 * ( T * ((4.46640989e+00_WP) & + T *((0.0075510663_WP) & + T *((-2.26701845e-06_WP) & + T *((2.83131345e-10_WP) & + T *((0_WP)))))) + (-1.80236075e+04_WP)) Cp(sCH3O2H) = 173.056908538362 * ((4.46640989e+00_WP) & + T *((1.51021326e-02_WP) & + T *((-6.80105535e-06_WP) & + T *((1.13252538e-09_WP) & + T *(0.00000000e+00_WP))))) dH(sCH3O2H) = 173.056908538362 * ( (4.46640989e+00_WP) & + T *((1.51021326e-02_WP) & + T *((-6.80105535e-06_WP) & + T *((1.13252538e-09_WP) & + T *((0.00000000e+00_WP))))) ) H(sC2H2) = 319.327085573821 * ( T * ((4.43677000e+00_WP) & + T *((0.0026880195_WP) & + T *((-6.37605666666667e-07_WP) & + T *((8.2159475e-11_WP) & + T *((-4.31342e-15_WP)))))) + (2.56676600e+04_WP)) Cp(sC2H2) = 319.327085573821 * ((4.43677000e+00_WP) & + T *((5.37603900e-03_WP) & + T *((-1.91281700e-06_WP) & + T *((3.28637900e-10_WP) & + T *(-2.15671000e-14_WP))))) dH(sC2H2) = 319.327085573821 * ( (4.43677000e+00_WP) & + T *((5.37603900e-03_WP) & + T *((-1.91281700e-06_WP) & + T *((3.28637900e-10_WP) & + T *((-2.15671000e-14_WP))))) ) H(sHCCO) = 202.64209807936 * ( T * ((6.75807300e+00_WP) & + T *((0.0010002_WP) & + T *((-6.75869e-08_WP) & + T *((-2.60283e-11_WP) & + T *((3.93033e-15_WP)))))) + (1.90151300e+04_WP)) Cp(sHCCO) = 202.64209807936 * ((6.75807300e+00_WP) & + T *((2.00040000e-03_WP) & + T *((-2.02760700e-07_WP) & + T *((-1.04113200e-10_WP) & + T *(1.96516500e-14_WP))))) dH(sHCCO) = 202.64209807936 * ( (6.75807300e+00_WP) & + T *((2.00040000e-03_WP) & + T *((-2.02760700e-07_WP) & + T *((-1.04113200e-10_WP) & + T *((1.96516500e-14_WP))))) ) H(sC2H3) = 307.424937139476 * ( T * ((5.93346800e+00_WP) & + T *((0.002008873_WP) & + T *((-1.32224666666667e-07_WP) & + T *((-3.6031675e-11_WP) & + T *((4.757288e-15_WP)))))) + (3.18543500e+04_WP)) Cp(sC2H3) = 307.424937139476 * ((5.93346800e+00_WP) & + T *((4.01774600e-03_WP) & + T *((-3.96674000e-07_WP) & + T *((-1.44126700e-10_WP) & + T *(2.37864400e-14_WP))))) dH(sC2H3) = 307.424937139476 * ( (5.93346800e+00_WP) & + T *((4.01774600e-03_WP) & + T *((-3.96674000e-07_WP) & + T *((-1.44126700e-10_WP) & + T *((2.37864400e-14_WP))))) ) H(sCH2CHO) = 193.15119412694 * ( T * ((5.97567000e+00_WP) & + T *((0.0040652955_WP) & + T *((-9.14541333333333e-07_WP) & + T *((1.017576e-10_WP) & + T *((-4.352034e-15_WP)))))) + (4.90321800e+02_WP)) Cp(sCH2CHO) = 193.15119412694 * ((5.97567000e+00_WP) & + T *((8.13059100e-03_WP) & + T *((-2.74362400e-06_WP) & + T *((4.07030400e-10_WP) & + T *(-2.17601700e-14_WP))))) dH(sCH2CHO) = 193.15119412694 * ( (5.97567000e+00_WP) & + T *((8.13059100e-03_WP) & + T *((-2.74362400e-06_WP) & + T *((4.07030400e-10_WP) & + T *((-2.17601700e-14_WP))))) ) H(sC3H5XAXC3H5) = 202.434867299732 * ( T * ((3.66851913e+00_WP) & + T *((0.0098101658_WP) & + T *((-2.84310523e-06_WP) & + T *((3.45283085e-10_WP) & + T *((0_WP)))))) + (1.83798091e+04_WP)) Cp(sC3H5XAXC3H5) = 202.434867299732 * ((3.66851913e+00_WP) & + T *((1.96203316e-02_WP) & + T *((-8.52931569e-06_WP) & + T *((1.38113234e-09_WP) & + T *(0.00000000e+00_WP))))) dH(sC3H5XAXC3H5) = 202.434867299732 * ( (3.66851913e+00_WP) & + T *((1.96203316e-02_WP) & + T *((-8.52931569e-06_WP) & + T *((1.38113234e-09_WP) & + T *((0.00000000e+00_WP))))) ) H(sC3H6) = 197.585436570179 * ( T * ((1.68134957e+00_WP) & + T *((0.0124845359_WP) & + T *((-3.70920136666667e-06_WP) & + T *((4.61033815e-10_WP) & + T *((0_WP)))))) + (7.40726095e+02_WP)) Cp(sC3H6) = 197.585436570179 * ((1.68134957e+00_WP) & + T *((2.49690718e-02_WP) & + T *((-1.11276041e-05_WP) & + T *((1.84413526e-09_WP) & + T *(0.00000000e+00_WP))))) dH(sC3H6) = 197.585436570179 * ( (1.68134957e+00_WP) & + T *((2.49690718e-02_WP) & + T *((-1.11276041e-05_WP) & + T *((1.84413526e-09_WP) & + T *((0.00000000e+00_WP))))) ) H(sC3H5O) = 145.680742947258 * ( T * ((3.39074577e+00_WP) & + T *((0.012065081_WP) & + T *((-3.78836313333333e-06_WP) & + T *((4.94752345e-10_WP) & + T *((0_WP)))))) + (9.00757452e+03_WP)) Cp(sC3H5O) = 145.680742947258 * ((3.39074577e+00_WP) & + T *((2.41301620e-02_WP) & + T *((-1.13650894e-05_WP) & + T *((1.97900938e-09_WP) & + T *(0.00000000e+00_WP))))) dH(sC3H5O) = 145.680742947258 * ( (3.39074577e+00_WP) & + T *((2.41301620e-02_WP) & + T *((-1.13650894e-05_WP) & + T *((1.97900938e-09_WP) & + T *((0.00000000e+00_WP))))) ) H(sIXC3H7) = 192.962911386529 * ( T * ((8.06336900e+00_WP) & + T *((0.00787244_WP) & + T *((-1.727464e-06_WP) & + T *((1.86931125e-10_WP) & + T *((-7.708844e-15_WP)))))) + (5.31387100e+03_WP)) Cp(sIXC3H7) = 192.962911386529 * ((8.06336900e+00_WP) & + T *((1.57448800e-02_WP) & + T *((-5.18239200e-06_WP) & + T *((7.47724500e-10_WP) & + T *(-3.85442200e-14_WP))))) dH(sIXC3H7) = 192.962911386529 * ( (8.06336900e+00_WP) & + T *((1.57448800e-02_WP) & + T *((-5.18239200e-06_WP) & + T *((7.47724500e-10_WP) & + T *((-3.85442200e-14_WP))))) ) H(sNXC3H7) = 192.962911386529 * ( T * ((7.97829100e+00_WP) & + T *((0.007880565_WP) & + T *((-1.72441433333333e-06_WP) & + T *((1.860973e-10_WP) & + T *((-7.649956e-15_WP)))))) + (7.57940200e+03_WP)) Cp(sNXC3H7) = 192.962911386529 * ((7.97829100e+00_WP) & + T *((1.57611300e-02_WP) & + T *((-5.17324300e-06_WP) & + T *((7.44389200e-10_WP) & + T *(-3.82497800e-14_WP))))) dH(sNXC3H7) = 192.962911386529 * ( (7.97829100e+00_WP) & + T *((1.57611300e-02_WP) & + T *((-5.17324300e-06_WP) & + T *((7.44389200e-10_WP) & + T *((-3.82497800e-14_WP))))) ) H(sC3H8) = 188.551730394158 * ( T * ((7.52521700e+00_WP) & + T *((0.00944517_WP) & + T *((-2.09464133333333e-06_WP) & + T *((2.29484325e-10_WP) & + T *((-9.62482e-15_WP)))))) + (-1.64645500e+04_WP)) Cp(sC3H8) = 188.551730394158 * ((7.52521700e+00_WP) & + T *((1.88903400e-02_WP) & + T *((-6.28392400e-06_WP) & + T *((9.17937300e-10_WP) & + T *(-4.81241000e-14_WP))))) dH(sC3H8) = 188.551730394158 * ( (7.52521700e+00_WP) & + T *((1.88903400e-02_WP) & + T *((-6.28392400e-06_WP) & + T *((9.17937300e-10_WP) & + T *((-4.81241000e-14_WP))))) ) H(sIXC3H7O2) = 110.72636709906 * ( T * ((6.12023736e+00_WP) & + T *((0.0142615975_WP) & + T *((-4.10011726666667e-06_WP) & + T *((4.895046625e-10_WP) & + T *((0_WP)))))) + (-1.13282549e+04_WP)) Cp(sIXC3H7O2) = 110.72636709906 * ((6.12023736e+00_WP) & + T *((2.85231950e-02_WP) & + T *((-1.23003518e-05_WP) & + T *((1.95801865e-09_WP) & + T *(0.00000000e+00_WP))))) dH(sIXC3H7O2) = 110.72636709906 * ( (6.12023736e+00_WP) & + T *((2.85231950e-02_WP) & + T *((-1.23003518e-05_WP) & + T *((1.95801865e-09_WP) & + T *((0.00000000e+00_WP))))) ) H(sNXC3H7O2) = 110.72636709906 * ( T * ((4.71102090e+00_WP) & + T *((0.0153071024_WP) & + T *((-4.48109563333333e-06_WP) & + T *((5.437036375e-10_WP) & + T *((0_WP)))))) + (-8.53286429e+03_WP)) Cp(sNXC3H7O2) = 110.72636709906 * ((4.71102090e+00_WP) & + T *((3.06142048e-02_WP) & + T *((-1.34432869e-05_WP) & + T *((2.17481455e-09_WP) & + T *(0.00000000e+00_WP))))) dH(sNXC3H7O2) = 110.72636709906 * ( (4.71102090e+00_WP) & + T *((3.06142048e-02_WP) & + T *((-1.34432869e-05_WP) & + T *((2.17481455e-09_WP) & + T *((0.00000000e+00_WP))))) ) else H(sN2) = 296.716630977873 * ( T * ((3.29867700e+00_WP) & + T *((0.00070412_WP) & + T *((-1.321074e-06_WP) & + T *((1.41037875e-09_WP) & + T *((-4.88971e-13_WP)))))) + (-1.02090000e+03_WP)) Cp(sN2) = 296.716630977873 * ((3.29867700e+00_WP) & + T *((1.40824000e-03_WP) & + T *((-3.96322200e-06_WP) & + T *((5.64151500e-09_WP) & + T *(-2.44485500e-12_WP))))) dH(sN2) = 296.716630977873 * ( (3.29867700e+00_WP) & + T *((1.40824000e-03_WP) & + T *((-3.96322200e-06_WP) & + T *((5.64151500e-09_WP) & + T *((-2.44485500e-12_WP))))) ) H(sO) = 519.625 * ( T * ((2.94642900e+00_WP) & + T *((-0.000819083_WP) & + T *((8.07010666666667e-07_WP) & + T *((-4.0071075e-10_WP) & + T *((7.781392e-14_WP)))))) + (2.91476400e+04_WP)) Cp(sO) = 519.625 * ((2.94642900e+00_WP) & + T *((-1.63816600e-03_WP) & + T *((2.42103200e-06_WP) & + T *((-1.60284300e-09_WP) & + T *(3.89069600e-13_WP))))) dH(sO) = 519.625 * ( (2.94642900e+00_WP) & + T *((-1.63816600e-03_WP) & + T *((2.42103200e-06_WP) & + T *((-1.60284300e-09_WP) & + T *((3.89069600e-13_WP))))) ) H(sO2) = 259.8125 * ( T * ((3.21293600e+00_WP) & + T *((0.000563743_WP) & + T *((-1.91871666666667e-07_WP) & + T *((3.2846925e-10_WP) & + T *((-1.7537108e-13_WP)))))) + (-1.00524900e+03_WP)) Cp(sO2) = 259.8125 * ((3.21293600e+00_WP) & + T *((1.12748600e-03_WP) & + T *((-5.75615000e-07_WP) & + T *((1.31387700e-09_WP) & + T *(-8.76855400e-13_WP))))) dH(sO2) = 259.8125 * ( (3.21293600e+00_WP) & + T *((1.12748600e-03_WP) & + T *((-5.75615000e-07_WP) & + T *((1.31387700e-09_WP) & + T *((-8.76855400e-13_WP))))) ) H(sH) = 8248.01587301587 * ( T * ((2.50104422e+00_WP) & + T *((0_WP) & + T *((0_WP) & + T *((0_WP) & + T *((0_WP)))))) + (2.54747466e+04_WP)) Cp(sH) = 8248.01587301587 * ((2.50104422e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *(0.00000000e+00_WP))))) dH(sH) = 8248.01587301587 * ( (2.50104422e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP) & + T *((0.00000000e+00_WP))))) ) H(sOH) = 488.828786453434 * ( T * ((3.43586219e+00_WP) & + T *((0.000101117902_WP) & + T *((-3.7848804e-08_WP) & + T *((6.061128725e-11_WP) & + T *((-1.487302062e-14_WP)))))) + (3.74321252e+03_WP)) Cp(sOH) = 488.828786453434 * ((3.43586219e+00_WP) & + T *((2.02235804e-04_WP) & + T *((-1.13546412e-07_WP) & + T *((2.42445149e-10_WP) & + T *(-7.43651031e-14_WP))))) dH(sOH) = 488.828786453434 * ( (3.43586219e+00_WP) & + T *((2.02235804e-04_WP) & + T *((-1.13546412e-07_WP) & + T *((2.42445149e-10_WP) & + T *((-7.43651031e-14_WP))))) ) H(sH2) = 4124.00793650794 * ( T * ((3.29812400e+00_WP) & + T *((0.0004124721_WP) & + T *((-2.71433833333333e-07_WP) & + T *((-2.3688585e-11_WP) & + T *((8.269744e-14_WP)))))) + (-1.01252100e+03_WP)) Cp(sH2) = 4124.00793650794 * ((3.29812400e+00_WP) & + T *((8.24944200e-04_WP) & + T *((-8.14301500e-07_WP) & + T *((-9.47543400e-11_WP) & + T *(4.13487200e-13_WP))))) dH(sH2) = 4124.00793650794 * ( (3.29812400e+00_WP) & + T *((8.24944200e-04_WP) & + T *((-8.14301500e-07_WP) & + T *((-9.47543400e-11_WP) & + T *((4.13487200e-13_WP))))) ) H(sH2O) = 461.478685612789 * ( T * ((3.38684200e+00_WP) & + T *((0.001737491_WP) & + T *((-2.118232e-06_WP) & + T *((1.74214525e-09_WP) & + T *((-5.013176e-13_WP)))))) + (-3.02081100e+04_WP)) Cp(sH2O) = 461.478685612789 * ((3.38684200e+00_WP) & + T *((3.47498200e-03_WP) & + T *((-6.35469600e-06_WP) & + T *((6.96858100e-09_WP) & + T *(-2.50658800e-12_WP))))) dH(sH2O) = 461.478685612789 * ( (3.38684200e+00_WP) & + T *((3.47498200e-03_WP) & + T *((-6.35469600e-06_WP) & + T *((6.96858100e-09_WP) & + T *((-2.50658800e-12_WP))))) ) H(sH2O2) = 244.414393226717 * ( T * ((3.38875400e+00_WP) & + T *((0.003284613_WP) & + T *((-4.95004333333333e-08_WP) & + T *((-1.1564515e-09_WP) & + T *((4.94303e-13_WP)))))) + (-1.76631500e+04_WP)) Cp(sH2O2) = 244.414393226717 * ((3.38875400e+00_WP) & + T *((6.56922600e-03_WP) & + T *((-1.48501300e-07_WP) & + T *((-4.62580600e-09_WP) & + T *(2.47151500e-12_WP))))) dH(sH2O2) = 244.414393226717 * ( (3.38875400e+00_WP) & + T *((6.56922600e-03_WP) & + T *((-1.48501300e-07_WP) & + T *((-4.62580600e-09_WP) & + T *((2.47151500e-12_WP))))) ) H(sHO2) = 251.878332525448 * ( T * ((3.18310656e+00_WP) & + T *((0.00183383975_WP) & + T *((-3.10795040666667e-07_WP) & + T *((-8.146322975e-11_WP) & + T *((3.02279824e-14_WP)))))) + (8.09181013e+02_WP)) Cp(sHO2) = 251.878332525448 * ((3.18310656e+00_WP) & + T *((3.66767950e-03_WP) & + T *((-9.32385122e-07_WP) & + T *((-3.25852919e-10_WP) & + T *(1.51139912e-13_WP))))) dH(sHO2) = 251.878332525448 * ( (3.18310656e+00_WP) & + T *((3.66767950e-03_WP) & + T *((-9.32385122e-07_WP) & + T *((-3.25852919e-10_WP) & + T *((1.51139912e-13_WP))))) ) H(sCH2GSGXCH2) = 592.756309710538 * ( T * ((3.97126500e+00_WP) & + T *((-8.495445e-05_WP) & + T *((3.41789666666667e-07_WP) & + T *((6.2313775e-10_WP) & + T *((-3.962532e-13_WP)))))) + (4.98936800e+04_WP)) Cp(sCH2GSGXCH2) = 592.756309710538 * ((3.97126500e+00_WP) & + T *((-1.69908900e-04_WP) & + T *((1.02536900e-06_WP) & + T *((2.49255100e-09_WP) & + T *(-1.98126600e-12_WP))))) dH(sCH2GSGXCH2) = 592.756309710538 * ( (3.97126500e+00_WP) & + T *((-1.69908900e-04_WP) & + T *((1.02536900e-06_WP) & + T *((2.49255100e-09_WP) & + T *((-1.98126600e-12_WP))))) ) H(sCO) = 296.822563370225 * ( T * ((3.26245200e+00_WP) & + T *((0.0007559705_WP) & + T *((-1.29391833333333e-06_WP) & + T *((1.395486e-09_WP) & + T *((-4.949902e-13_WP)))))) + (-1.43105400e+04_WP)) Cp(sCO) = 296.822563370225 * ((3.26245200e+00_WP) & + T *((1.51194100e-03_WP) & + T *((-3.88175500e-06_WP) & + T *((5.58194400e-09_WP) & + T *(-2.47495100e-12_WP))))) dH(sCO) = 296.822563370225 * ( (3.26245200e+00_WP) & + T *((1.51194100e-03_WP) & + T *((-3.88175500e-06_WP) & + T *((5.58194400e-09_WP) & + T *((-2.47495100e-12_WP))))) ) H(sCH2O) = 276.89335908879 * ( T * ((1.65273100e+00_WP) & + T *((0.00631572_WP) & + T *((-6.29389333333333e-06_WP) & + T *((5.1250775e-09_WP) & + T *((-1.6826474e-12_WP)))))) + (-1.48654000e+04_WP)) Cp(sCH2O) = 276.89335908879 * ((1.65273100e+00_WP) & + T *((1.26314400e-02_WP) & + T *((-1.88816800e-05_WP) & + T *((2.05003100e-08_WP) & + T *(-8.41323700e-12_WP))))) dH(sCH2O) = 276.89335908879 * ( (1.65273100e+00_WP) & + T *((1.26314400e-02_WP) & + T *((-1.88816800e-05_WP) & + T *((2.05003100e-08_WP) & + T *((-8.41323700e-12_WP))))) ) H(sCH3) = 553.013170147665 * ( T * ((2.43044300e+00_WP) & + T *((0.00556205_WP) & + T *((-5.60073333333333e-06_WP) & + T *((4.0545725e-09_WP) & + T *((-1.1729906e-12_WP)))))) + (1.64237800e+04_WP)) Cp(sCH3) = 553.013170147665 * ((2.43044300e+00_WP) & + T *((1.11241000e-02_WP) & + T *((-1.68022000e-05_WP) & + T *((1.62182900e-08_WP) & + T *(-5.86495300e-12_WP))))) dH(sCH3) = 553.013170147665 * ( (2.43044300e+00_WP) & + T *((1.11241000e-02_WP) & + T *((-1.68022000e-05_WP) & + T *((1.62182900e-08_WP) & + T *((-5.86495300e-12_WP))))) ) H(sC2H6) = 276.506585073833 * ( T * ((-2.52854344e-02_WP) & + T *((0.0120382377_WP) & + T *((-3.7297824e-06_WP) & + T *((5.208522525e-10_WP) & + T *((-1.059737232e-14_WP)))))) + (-1.12345534e+04_WP)) Cp(sC2H6) = 276.506585073833 * ((-2.52854344e-02_WP) & + T *((2.40764754e-02_WP) & + T *((-1.11893472e-05_WP) & + T *((2.08340901e-09_WP) & + T *(-5.29868616e-14_WP))))) dH(sC2H6) = 276.506585073833 * ( (-2.52854344e-02_WP) & + T *((2.40764754e-02_WP) & + T *((-1.11893472e-05_WP) & + T *((2.08340901e-09_WP) & + T *((-5.29868616e-14_WP))))) ) H(sCH3O) = 267.899722884578 * ( T * ((2.10620400e+00_WP) & + T *((0.0036082975_WP) & + T *((1.77949066666667e-06_WP) & + T *((-1.844409e-09_WP) & + T *((4.151222e-13_WP)))))) + (9.78601100e+02_WP)) Cp(sCH3O) = 267.899722884578 * ((2.10620400e+00_WP) & + T *((7.21659500e-03_WP) & + T *((5.33847200e-06_WP) & + T *((-7.37763600e-09_WP) & + T *(2.07561100e-12_WP))))) dH(sCH3O) = 267.899722884578 * ( (2.10620400e+00_WP) & + T *((7.21659500e-03_WP) & + T *((5.33847200e-06_WP) & + T *((-7.37763600e-09_WP) & + T *((2.07561100e-12_WP))))) ) H(sCH4) = 518.264555541703 * ( T * ((7.78741500e-01_WP) & + T *((0.00873834_WP) & + T *((-9.27803e-06_WP) & + T *((7.62427e-09_WP) & + T *((-2.447862e-12_WP)))))) + (-9.82522900e+03_WP)) Cp(sCH4) = 518.264555541703 * ((7.78741500e-01_WP) & + T *((1.74766800e-02_WP) & + T *((-2.78340900e-05_WP) & + T *((3.04970800e-08_WP) & + T *(-1.22393100e-11_WP))))) dH(sCH4) = 518.264555541703 * ( (7.78741500e-01_WP) & + T *((1.74766800e-02_WP) & + T *((-2.78340900e-05_WP) & + T *((3.04970800e-08_WP) & + T *((-1.22393100e-11_WP))))) ) H(sC2H5) = 286.097728836889 * ( T * ((2.69070200e+00_WP) & + T *((0.0043595665_WP) & + T *((1.47327966666667e-06_WP) & + T *((2.33467575e-10_WP) & + T *((-7.855546e-13_WP)))))) + (1.28704000e+04_WP)) Cp(sC2H5) = 286.097728836889 * ((2.69070200e+00_WP) & + T *((8.71913300e-03_WP) & + T *((4.41983900e-06_WP) & + T *((9.33870300e-10_WP) & + T *(-3.92777300e-12_WP))))) dH(sC2H5) = 286.097728836889 * ( (2.69070200e+00_WP) & + T *((8.71913300e-03_WP) & + T *((4.41983900e-06_WP) & + T *((9.33870300e-10_WP) & + T *((-3.92777300e-12_WP))))) ) H(sC2H4) = 296.378154855269 * ( T * ((-8.61488000e-01_WP) & + T *((0.013980815_WP) & + T *((-1.129559e-05_WP) & + T *((6.96288e-09_WP) & + T *((-1.9475758e-12_WP)))))) + (5.57304600e+03_WP)) Cp(sC2H4) = 296.378154855269 * ((-8.61488000e-01_WP) & + T *((2.79616300e-02_WP) & + T *((-3.38867700e-05_WP) & + T *((2.78515200e-08_WP) & + T *(-9.73787900e-12_WP))))) dH(sC2H4) = 296.378154855269 * ( (-8.61488000e-01_WP) & + T *((2.79616300e-02_WP) & + T *((-3.38867700e-05_WP) & + T *((2.78515200e-08_WP) & + T *((-9.73787900e-12_WP))))) ) H(sCO2) = 188.911610997501 * ( T * ((2.27572500e+00_WP) & + T *((0.004961036_WP) & + T *((-3.46970333333333e-06_WP) & + T *((1.71667175e-09_WP) & + T *((-4.23456e-13_WP)))))) + (-4.83731400e+04_WP)) Cp(sCO2) = 188.911610997501 * ((2.27572500e+00_WP) & + T *((9.92207200e-03_WP) & + T *((-1.04091100e-05_WP) & + T *((6.86668700e-09_WP) & + T *(-2.11728000e-12_WP))))) dH(sCO2) = 188.911610997501 * ( (2.27572500e+00_WP) & + T *((9.92207200e-03_WP) & + T *((-1.04091100e-05_WP) & + T *((6.86668700e-09_WP) & + T *((-2.11728000e-12_WP))))) ) H(sHCO) = 286.5118202495 * ( T * ((2.89833000e+00_WP) & + T *((0.0030995735_WP) & + T *((-3.20769466666667e-06_WP) & + T *((2.7245625e-09_WP) & + T *((-9.14977e-13_WP)))))) + (4.15992200e+03_WP)) Cp(sHCO) = 286.5118202495 * ((2.89833000e+00_WP) & + T *((6.19914700e-03_WP) & + T *((-9.62308400e-06_WP) & + T *((1.08982500e-08_WP) & + T *(-4.57488500e-12_WP))))) dH(sHCO) = 286.5118202495 * ( (2.89833000e+00_WP) & + T *((6.19914700e-03_WP) & + T *((-9.62308400e-06_WP) & + T *((1.08982500e-08_WP) & + T *((-4.57488500e-12_WP))))) ) H(sCH3O2) = 176.765743929923 * ( T * ((4.26146906e+00_WP) & + T *((0.00504367995_WP) & + T *((-1.07168728e-06_WP) & + T *((5.235231675e-11_WP) & + T *((8.36678206e-15_WP)))))) + (-6.84394259e+02_WP)) Cp(sCH3O2) = 176.765743929923 * ((4.26146906e+00_WP) & + T *((1.00873599e-02_WP) & + T *((-3.21506184e-06_WP) & + T *((2.09409267e-10_WP) & + T *(4.18339103e-14_WP))))) dH(sCH3O2) = 176.765743929923 * ( (4.26146906e+00_WP) & + T *((1.00873599e-02_WP) & + T *((-3.21506184e-06_WP) & + T *((2.09409267e-10_WP) & + T *((4.18339103e-14_WP))))) ) H(sCH3O2H) = 173.056908538362 * ( T * ((3.23442817e+00_WP) & + T *((0.00950648835_WP) & + T *((-3.7795429e-06_WP) & + T *((8.507666325e-10_WP) & + T *((-8.23660444e-14_WP)))))) + (-1.77197926e+04_WP)) Cp(sCH3O2H) = 173.056908538362 * ((3.23442817e+00_WP) & + T *((1.90129767e-02_WP) & + T *((-1.13386287e-05_WP) & + T *((3.40306653e-09_WP) & + T *(-4.11830222e-13_WP))))) dH(sCH3O2H) = 173.056908538362 * ( (3.23442817e+00_WP) & + T *((1.90129767e-02_WP) & + T *((-1.13386287e-05_WP) & + T *((3.40306653e-09_WP) & + T *((-4.11830222e-13_WP))))) ) H(sC2H2) = 319.327085573821 * ( T * ((2.01356200e+00_WP) & + T *((0.007595225_WP) & + T *((-5.38773e-06_WP) & + T *((2.269748e-09_WP) & + T *((-3.825492e-13_WP)))))) + (2.61244400e+04_WP)) Cp(sC2H2) = 319.327085573821 * ((2.01356200e+00_WP) & + T *((1.51904500e-02_WP) & + T *((-1.61631900e-05_WP) & + T *((9.07899200e-09_WP) & + T *(-1.91274600e-12_WP))))) dH(sC2H2) = 319.327085573821 * ( (2.01356200e+00_WP) & + T *((1.51904500e-02_WP) & + T *((-1.61631900e-05_WP) & + T *((9.07899200e-09_WP) & + T *((-1.91274600e-12_WP))))) ) H(sHCCO) = 202.64209807936 * ( T * ((5.04796500e+00_WP) & + T *((0.002226739_WP) & + T *((7.56094333333333e-08_WP) & + T *((-3.7052375e-10_WP) & + T *((4.501484e-14_WP)))))) + (1.96589200e+04_WP)) Cp(sHCCO) = 202.64209807936 * ((5.04796500e+00_WP) & + T *((4.45347800e-03_WP) & + T *((2.26828300e-07_WP) & + T *((-1.48209500e-09_WP) & + T *(2.25074200e-13_WP))))) dH(sHCCO) = 202.64209807936 * ( (5.04796500e+00_WP) & + T *((4.45347800e-03_WP) & + T *((2.26828300e-07_WP) & + T *((-1.48209500e-09_WP) & + T *((2.25074200e-13_WP))))) ) H(sC2H3) = 307.424937139476 * ( T * ((2.45927600e+00_WP) & + T *((0.003685738_WP) & + T *((7.03291e-07_WP) & + T *((-3.304105e-10_WP) & + T *((-2.369568e-13_WP)))))) + (3.33522500e+04_WP)) Cp(sC2H3) = 307.424937139476 * ((2.45927600e+00_WP) & + T *((7.37147600e-03_WP) & + T *((2.10987300e-06_WP) & + T *((-1.32164200e-09_WP) & + T *(-1.18478400e-12_WP))))) dH(sC2H3) = 307.424937139476 * ( (2.45927600e+00_WP) & + T *((7.37147600e-03_WP) & + T *((2.10987300e-06_WP) & + T *((-1.32164200e-09_WP) & + T *((-1.18478400e-12_WP))))) ) H(sCH2CHO) = 193.15119412694 * ( T * ((3.40906200e+00_WP) & + T *((0.005369285_WP) & + T *((6.30497333333333e-07_WP) & + T *((-1.78964575e-09_WP) & + T *((5.73477e-13_WP)))))) + (1.52147700e+03_WP)) Cp(sCH2CHO) = 193.15119412694 * ((3.40906200e+00_WP) & + T *((1.07385700e-02_WP) & + T *((1.89149200e-06_WP) & + T *((-7.15858300e-09_WP) & + T *(2.86738500e-12_WP))))) dH(sCH2CHO) = 193.15119412694 * ( (3.40906200e+00_WP) & + T *((1.07385700e-02_WP) & + T *((1.89149200e-06_WP) & + T *((-7.15858300e-09_WP) & + T *((2.86738500e-12_WP))))) ) H(sC3H5XAXC3H5) = 202.434867299732 * ( T * ((-5.29131958e-01_WP) & + T *((0.016727955_WP) & + T *((-8.4467009e-06_WP) & + T *((2.57164385e-09_WP) & + T *((-3.4651668e-13_WP)))))) + (1.93834226e+04_WP)) Cp(sC3H5XAXC3H5) = 202.434867299732 * ((-5.29131958e-01_WP) & + T *((3.34559100e-02_WP) & + T *((-2.53401027e-05_WP) & + T *((1.02865754e-08_WP) & + T *(-1.73258340e-12_WP))))) dH(sC3H5XAXC3H5) = 202.434867299732 * ( (-5.29131958e-01_WP) & + T *((3.34559100e-02_WP) & + T *((-2.53401027e-05_WP) & + T *((1.02865754e-08_WP) & + T *((-1.73258340e-12_WP))))) ) H(sC3H6) = 197.585436570179 * ( T * ((3.94615444e-01_WP) & + T *((0.0144553831_WP) & + T *((-5.1628936e-06_WP) & + T *((9.720355225e-10_WP) & + T *((-6.75780704e-14_WP)))))) + (1.06688164e+03_WP)) Cp(sC3H6) = 197.585436570179 * ((3.94615444e-01_WP) & + T *((2.89107662e-02_WP) & + T *((-1.54886808e-05_WP) & + T *((3.88814209e-09_WP) & + T *(-3.37890352e-13_WP))))) dH(sC3H6) = 197.585436570179 * ( (3.94615444e-01_WP) & + T *((2.89107662e-02_WP) & + T *((-1.54886808e-05_WP) & + T *((3.88814209e-09_WP) & + T *((-3.37890352e-13_WP))))) ) H(sC3H5O) = 145.680742947258 * ( T * ((1.19822582e+00_WP) & + T *((0.01527899185_WP) & + T *((-6.0210092e-06_WP) & + T *((1.2153750825e-09_WP) & + T *((-8.39709124e-14_WP)))))) + (9.58217784e+03_WP)) Cp(sC3H5O) = 145.680742947258 * ((1.19822582e+00_WP) & + T *((3.05579837e-02_WP) & + T *((-1.80630276e-05_WP) & + T *((4.86150033e-09_WP) & + T *(-4.19854562e-13_WP))))) dH(sC3H5O) = 145.680742947258 * ( (1.19822582e+00_WP) & + T *((3.05579837e-02_WP) & + T *((-1.80630276e-05_WP) & + T *((4.86150033e-09_WP) & + T *((-4.19854562e-13_WP))))) ) H(sIXC3H7) = 192.962911386529 * ( T * ((1.71330000e+00_WP) & + T *((0.01271308_WP) & + T *((5.26936e-07_WP) & + T *((-4.553215e-09_WP) & + T *((1.765542e-12_WP)))))) + (7.53580900e+03_WP)) Cp(sIXC3H7) = 192.962911386529 * ((1.71330000e+00_WP) & + T *((2.54261600e-02_WP) & + T *((1.58080800e-06_WP) & + T *((-1.82128600e-08_WP) & + T *(8.82771000e-12_WP))))) dH(sIXC3H7) = 192.962911386529 * ( (1.71330000e+00_WP) & + T *((2.54261600e-02_WP) & + T *((1.58080800e-06_WP) & + T *((-1.82128600e-08_WP) & + T *((8.82771000e-12_WP))))) ) H(sNXC3H7) = 192.962911386529 * ( T * ((1.92253700e+00_WP) & + T *((0.012394635_WP) & + T *((6.03416333333333e-07_WP) & + T *((-4.458165e-09_WP) & + T *((1.7165992e-12_WP)))))) + (9.71328100e+03_WP)) Cp(sNXC3H7) = 192.962911386529 * ((1.92253700e+00_WP) & + T *((2.47892700e-02_WP) & + T *((1.81024900e-06_WP) & + T *((-1.78326600e-08_WP) & + T *(8.58299600e-12_WP))))) dH(sNXC3H7) = 192.962911386529 * ( (1.92253700e+00_WP) & + T *((2.47892700e-02_WP) & + T *((1.81024900e-06_WP) & + T *((-1.78326600e-08_WP) & + T *((8.58299600e-12_WP))))) ) H(sC3H8) = 188.551730394158 * ( T * ((8.96920800e-01_WP) & + T *((0.01334493_WP) & + T *((1.810475e-06_WP) & + T *((-5.3150025e-09_WP) & + T *((1.848666e-12_WP)))))) + (-1.39549200e+04_WP)) Cp(sC3H8) = 188.551730394158 * ((8.96920800e-01_WP) & + T *((2.66898600e-02_WP) & + T *((5.43142500e-06_WP) & + T *((-2.12600100e-08_WP) & + T *(9.24333000e-12_WP))))) dH(sC3H8) = 188.551730394158 * ( (8.96920800e-01_WP) & + T *((2.66898600e-02_WP) & + T *((5.43142500e-06_WP) & + T *((-2.12600100e-08_WP) & + T *((9.24333000e-12_WP))))) ) H(sIXC3H7O2) = 110.72636709906 * ( T * ((1.49941639e+00_WP) & + T *((0.02215406025_WP) & + T *((-1.07471485333333e-05_WP) & + T *((3.2421784e-09_WP) & + T *((-4.46741138e-13_WP)))))) + (-1.02587980e+04_WP)) Cp(sIXC3H7O2) = 110.72636709906 * ((1.49941639e+00_WP) & + T *((4.43081205e-02_WP) & + T *((-3.22414456e-05_WP) & + T *((1.29687136e-08_WP) & + T *(-2.23370569e-12_WP))))) dH(sIXC3H7O2) = 110.72636709906 * ( (1.49941639e+00_WP) & + T *((4.43081205e-02_WP) & + T *((-3.22414456e-05_WP) & + T *((1.29687136e-08_WP) & + T *((-2.23370569e-12_WP))))) ) H(sNXC3H7O2) = 110.72636709906 * ( T * ((2.10731492e+00_WP) & + T *((0.0198082493_WP) & + T *((-8.31638663333333e-06_WP) & + T *((2.14862575e-09_WP) & + T *((-2.6248066e-13_WP)))))) + (-7.93745567e+03_WP)) Cp(sNXC3H7O2) = 110.72636709906 * ((2.10731492e+00_WP) & + T *((3.96164986e-02_WP) & + T *((-2.49491599e-05_WP) & + T *((8.59450300e-09_WP) & + T *(-1.31240330e-12_WP))))) dH(sNXC3H7O2) = 110.72636709906 * ( (2.10731492e+00_WP) & + T *((3.96164986e-02_WP) & + T *((-2.49491599e-05_WP) & + T *((8.59450300e-09_WP) & + T *((-1.31240330e-12_WP))))) ) end if return end subroutine compute_thermodata ! --- Thirdbodies --- ! subroutine get_thirdbodies(M,c) use mechanism implicit none real(WP), dimension(nspec) :: c real(WP), dimension(13) :: M M(mM11) = (1_WP)*c(sH2) & + (4_WP)*c(sH2O) & + (2_WP)*c(sCO2) & + (1_WP)*c(sCO) & + sum(c) M(mM5) = sum(c) M(mM6) = (1_WP)*c(sCO) & + (2_WP)*c(sCO2) & + (4_WP)*c(sH2O) & + (1_WP)*c(sH2) & + sum(c) M(mM4) = (1.5_WP)*c(sH2) & + (11_WP)*c(sH2O) & + (0.9_WP)*c(sCO) & + (2.8_WP)*c(sCO2) & + sum(c) M(mM9) = (5_WP)*c(sH2O) & + (2.8_WP)*c(sCO2) & + (0.9_WP)*c(sCO) & + (1.5_WP)*c(sH2) & + sum(c) M(mM12) = sum(c) M(mM2) = (1.5_WP)*c(sH2) & + (2.8_WP)*c(sCO2) & + (0.9_WP)*c(sCO) & + (11_WP)*c(sH2O) & + sum(c) M(mM1) = (1.5_WP)*c(sH2) & + (11_WP)*c(sH2O) & + (0.9_WP)*c(sCO) & + (2.8_WP)*c(sCO2) & + sum(c) M(mM8) = (1.5_WP)*c(sH2) & + (2.8_WP)*c(sCO2) & + (0.9_WP)*c(sCO) & + (11_WP)*c(sH2O) & + sum(c) M(mM7) = (1_WP)*c(sH2) & + (4_WP)*c(sH2O) & + (2_WP)*c(sCO2) & + (1_WP)*c(sCO) & + sum(c) M(mM3) = (1.5_WP)*c(sH2) & + (2.8_WP)*c(sCO2) & + (0.9_WP)*c(sCO) & + (11_WP)*c(sH2O) & + sum(c) M(mM0) = (11_WP)*c(sH2O) & + (2.8_WP)*c(sCO2) & + (0.9_WP)*c(sCO) & + (1.5_WP)*c(sH2) & + sum(c) return end subroutine get_thirdbodies ! --- Rate coefficients --- ! subroutine get_rate_coefficients(k,M,Tloc,Ploc) use mechanism implicit none real(WP), dimension(nreac) :: k real(WP), dimension(13) :: M real(WP) :: Tloc,Ploc real(WP) :: kH5f_0, kH5f_inf, FCH5f real(WP) :: kH8f_0, kH8f_inf, FCH8f real(WP) :: kH25f_0, kH25f_inf, FCH25f real(WP) :: kH29_0, kH29_inf, FCH29 real(WP) :: kH42_0, kH42_inf, FCH42 real(WP) :: kH59_0, kH59_inf, FCH59 real(WP) :: kH73_0, kH73_inf, FCH73 real(WP) :: kH74_0, kH74_inf, FCH74 real(WP) :: kH80f_0, kH80f_inf, FCH80f real(WP) :: kH5b_0, kH5b_inf, FCH5b real(WP) :: kH8b_0, kH8b_inf, FCH8b real(WP) :: kH25b_0, kH25b_inf, FCH25b real(WP) :: kH80b_0, kH80b_inf, FCH80b ! Rate coefficients k(rH1) = (6.17000000e+03_WP)*Tloc**(-0.500_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH2) = (4.72000000e+06_WP)*Tloc**(-1.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH3f) = (5.08000000e-02_WP)*Tloc**(2.670_WP)* & exp(-(2.633e+04_WP)/(8.314_WP*Tloc)) k(rH3b) = (2.23100000e-02_WP)*Tloc**(2.670_WP)* & exp(-(1.756e+04_WP)/(8.314_WP*Tloc)) k(rH4) = (2.25000000e+10_WP)*Tloc**(-2.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) kH5f_0 = (3.04100000e+18_WP)*Tloc**(-4.630_WP)* & exp(-(8.573e+03_WP)/(8.314_WP*Tloc)) kH5f_inf = (1.23600000e+08_WP)*Tloc**(-0.370_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) FCH5f = (5.300e-01_WP)* & exp(-Tloc/(1.000e+02_WP)) + (4.700e-01_WP)* & exp(-Tloc/(2.000e+03_WP)) + (1.000e+00_WP)* & exp(-(1.000e+15_WP)/Tloc) k(rH5f) = & getlindratecoeff & (Tloc, kH5f_0, kH5f_inf, FCH5f, M(mM2), Ploc ) k(rH6f) = (2.16000000e+02_WP)*Tloc**(1.510_WP)* & exp(-(1.435e+04_WP)/(8.314_WP*Tloc)) k(rH6b) = (9.35200000e+02_WP)*Tloc**(1.510_WP)* & exp(-(7.774e+04_WP)/(8.314_WP*Tloc)) k(rH7f) = (2.97000000e+00_WP)*Tloc**(2.020_WP)* & exp(-(5.607e+04_WP)/(8.314_WP*Tloc)) k(rH7b) = (3.01300000e-01_WP)*Tloc**(2.020_WP)* & exp(-(-1.611e+04_WP)/(8.314_WP*Tloc)) kH8f_0 = (3.50000000e+04_WP)*Tloc**(-0.410_WP)* & exp(-(-4.669e+03_WP)/(8.314_WP*Tloc)) kH8f_inf = (1.47500000e+06_WP)*Tloc**(0.600_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) FCH8f = (5.000e-01_WP)* & exp(-Tloc/(1.000e-30_WP)) + (5.000e-01_WP)* & exp(-Tloc/(1.000e+30_WP)) + (1.000e+00_WP)* & exp(-(1.000e+100_WP)/Tloc) k(rH8f) = & getlindratecoeff & (Tloc, kH8f_0, kH8f_inf, FCH8f, M(mM4), Ploc ) k(rH9f) = (1.97000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(6.920e+04_WP)/(8.314_WP*Tloc)) k(rH9b) = (1.55500000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.778e+03_WP)/(8.314_WP*Tloc)) k(rH10f) = (2.89000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(-2.092e+03_WP)/(8.314_WP*Tloc)) k(rH10b) = (6.88800000e+09_WP)*Tloc**(-0.330_WP)* & exp(-(3.018e+05_WP)/(8.314_WP*Tloc)) k(rH11f) = (3.25000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH11b) = (7.85700000e+08_WP)*Tloc**(-0.330_WP)* & exp(-(2.318e+05_WP)/(8.314_WP*Tloc)) k(rH12) = (1.66000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(3.431e+03_WP)/(8.314_WP*Tloc)) k(rH13) = (1.30000000e+05_WP)*Tloc**(0.000_WP)* & exp(-(-6.816e+03_WP)/(8.314_WP*Tloc)) k(rH14) = (7.08000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.255e+03_WP)/(8.314_WP*Tloc)) k(rH15) = (4.20000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(5.012e+04_WP)/(8.314_WP*Tloc)) k(rH16f) = (4.82000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(3.326e+04_WP)/(8.314_WP*Tloc)) k(rH16b) = (1.87500000e+06_WP)*Tloc**(0.330_WP)* & exp(-(1.015e+05_WP)/(8.314_WP*Tloc)) k(rH17) = (9.55000000e+00_WP)*Tloc**(2.000_WP)* & exp(-(1.661e+04_WP)/(8.314_WP*Tloc)) k(rH18) = (1.00000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH19f) = (5.80000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(4.000e+04_WP)/(8.314_WP*Tloc)) k(rH19b) = (9.77100000e+07_WP)*Tloc**(0.330_WP)* & exp(-(1.716e+05_WP)/(8.314_WP*Tloc)) k(rH20) = (2.41000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.661e+04_WP)/(8.314_WP*Tloc)) k(rH21) = (3.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH22) = (7.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH23) = (3.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH24f) = (7.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH24b) = (2.48200000e+11_WP)*Tloc**(-0.890_WP)* & exp(-(6.749e+04_WP)/(8.314_WP*Tloc)) kH25f_0 = (1.13500000e+24_WP)*Tloc**(-5.246_WP)* & exp(-(7.134e+03_WP)/(8.314_WP*Tloc)) kH25f_inf = (9.21400000e+10_WP)*Tloc**(-1.170_WP)* & exp(-(2.660e+03_WP)/(8.314_WP*Tloc)) FCH25f = (5.950e-01_WP)* & exp(-Tloc/(1.120e+03_WP)) + (4.050e-01_WP)* & exp(-Tloc/(6.960e+01_WP)) + (1.000e+00_WP)* & exp(-(1.000e+15_WP)/Tloc) k(rH25f) = & getlindratecoeff & (Tloc, kH25f_0, kH25f_inf, FCH25f, M(mM6), Ploc ) k(rH26f) = (1.99500000e+12_WP)*Tloc**(-1.570_WP)* & exp(-(1.222e+05_WP)/(8.314_WP*Tloc)) k(rH26b) = (3.58500000e+12_WP)*Tloc**(-1.590_WP)* & exp(-(-6.824e+03_WP)/(8.314_WP*Tloc)) k(rH27) = (1.10000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH28f) = (2.65000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(9.146e+03_WP)/(8.314_WP*Tloc)) k(rH28b) = (3.23600000e+04_WP)*Tloc**(0.890_WP)* & exp(-(5.067e+03_WP)/(8.314_WP*Tloc)) kH29_0 = (3.31000000e+18_WP)*Tloc**(-4.000_WP)* & exp(-(8.820e+03_WP)/(8.314_WP*Tloc)) kH29_inf = (2.13800000e+09_WP)*Tloc**(-0.400_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) FCH29 = (1.000e+00_WP)* & exp(-Tloc/(1.000e-15_WP)) + (0.000e+00_WP)* & exp(-Tloc/(1.000e-15_WP)) + (1.000e+00_WP)* & exp(-(4.000e+01_WP)/Tloc) k(rH29) = & getlindratecoeff & (Tloc, kH29_0, kH29_inf, FCH29, M(mM7), Ploc ) k(rH30f) = (6.84000000e+06_WP)*Tloc**(0.100_WP)* & exp(-(4.435e+04_WP)/(8.314_WP*Tloc)) k(rH31) = (7.47000000e+05_WP)*Tloc**(0.000_WP)* & exp(-(5.962e+04_WP)/(8.314_WP*Tloc)) k(rH32) = (8.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH33) = (3.60000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH34) = (2.25000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.799e+04_WP)/(8.314_WP*Tloc)) k(rH35) = (2.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH36) = (3.36500000e+05_WP)*Tloc**(-0.330_WP)* & exp(-(1.047e+04_WP)/(8.314_WP*Tloc)) k(rH37f) = (1.93000000e-01_WP)*Tloc**(2.400_WP)* & exp(-(8.812e+03_WP)/(8.314_WP*Tloc)) k(rH37b) = (3.19900000e-02_WP)*Tloc**(2.400_WP)* & exp(-(7.021e+04_WP)/(8.314_WP*Tloc)) k(rH38f) = (1.72700000e-02_WP)*Tloc**(3.000_WP)* & exp(-(3.441e+04_WP)/(8.314_WP*Tloc)) k(rH38b) = (6.61000000e-04_WP)*Tloc**(3.000_WP)* & exp(-(3.240e+04_WP)/(8.314_WP*Tloc)) k(rH39f) = (4.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH39b) = (5.42900000e+09_WP)*Tloc**(-0.890_WP)* & exp(-(6.548e+04_WP)/(8.314_WP*Tloc)) k(rH40f) = (3.15000000e+06_WP)*Tloc**(0.500_WP)* & exp(-(4.305e+04_WP)/(8.314_WP*Tloc)) k(rH40b) = (5.29600000e+04_WP)*Tloc**(0.500_WP)* & exp(-(3.228e+04_WP)/(8.314_WP*Tloc)) k(rH41) = (3.01000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(9.623e+04_WP)/(8.314_WP*Tloc)) kH42_0 = (1.35000000e+12_WP)*Tloc**(-2.788_WP)* & exp(-(1.754e+04_WP)/(8.314_WP*Tloc)) kH42_inf = (1.80000000e+04_WP)*Tloc**(0.000_WP)* & exp(-(9.975e+03_WP)/(8.314_WP*Tloc)) FCH42 = + (1.000e+00_WP)* & exp(-(0.000e+00_WP)/Tloc) k(rH42) = & getlindratecoeff & (Tloc, kH42_0, kH42_inf, FCH42, M(mM8), Ploc ) k(rH43f) = (1.40000000e-01_WP)*Tloc**(1.950_WP)* & exp(-(-5.636e+03_WP)/(8.314_WP*Tloc)) k(rH43b) = (1.56800000e+01_WP)*Tloc**(1.950_WP)* & exp(-(8.782e+04_WP)/(8.314_WP*Tloc)) k(rH44f) = (1.62000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.996e+05_WP)/(8.314_WP*Tloc)) k(rH44b) = (1.43300000e+08_WP)*Tloc**(0.000_WP)* & exp(-(2.256e+05_WP)/(8.314_WP*Tloc)) k(rH45) = (3.02000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH46) = (3.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH47) = (7.34000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH48f) = (1.86000000e+11_WP)*Tloc**(-1.000_WP)* & exp(-(7.113e+04_WP)/(8.314_WP*Tloc)) k(rH48b) = (6.46700000e+01_WP)*Tloc**(0.000_WP)* & exp(-(-1.849e+03_WP)/(8.314_WP*Tloc)) k(rH49) = (1.02000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH50) = (1.21000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH51) = (7.58000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(1.715e+03_WP)/(8.314_WP*Tloc)) k(rH52) = (5.82000000e-09_WP)*Tloc**(4.530_WP)* & exp(-(2.743e+04_WP)/(8.314_WP*Tloc)) k(rH53) = (3.43000000e+03_WP)*Tloc**(1.180_WP)* & exp(-(-1.870e+03_WP)/(8.314_WP*Tloc)) k(rH54) = (9.33400000e+02_WP)*Tloc**(1.500_WP)* & exp(-(1.245e+04_WP)/(8.314_WP*Tloc)) k(rH55) = (3.63600000e-12_WP)*Tloc**(5.420_WP)* & exp(-(4.176e+03_WP)/(8.314_WP*Tloc)) k(rH56) = (4.16000000e+05_WP)*Tloc**(0.570_WP)* & exp(-(1.156e+04_WP)/(8.314_WP*Tloc)) k(rH57) = (2.05000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.630e+05_WP)/(8.314_WP*Tloc)) k(rH58) = (5.50000000e+04_WP)*Tloc**(0.000_WP)* & exp(-(1.014e+04_WP)/(8.314_WP*Tloc)) kH59_0 = (2.34400000e+19_WP)*Tloc**(-2.700_WP)* & exp(-(1.280e+05_WP)/(8.314_WP*Tloc)) kH59_inf = (5.45000000e+13_WP)*Tloc**(0.000_WP)* & exp(-(5.648e+04_WP)/(8.314_WP*Tloc)) FCH59 = + (1.000e+00_WP)* & exp(-(0.000e+00_WP)/Tloc) k(rH59) = & getlindratecoeff & (Tloc, kH59_0, kH59_inf, FCH59, M(mM5), Ploc ) k(rH60) = (3.00000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH61) = (1.75000000e+04_WP)*Tloc**(0.000_WP)* & exp(-(-1.370e+04_WP)/(8.314_WP*Tloc)) k(rH62) = (1.40000000e+10_WP)*Tloc**(-1.610_WP)* & exp(-(7.782e+03_WP)/(8.314_WP*Tloc)) k(rH63f) = (4.34300000e+21_WP)*Tloc**(-3.420_WP)* & exp(-(1.275e+05_WP)/(8.314_WP*Tloc)) k(rH63b) = (5.44000000e+13_WP)*Tloc**(-3.300_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH64) = (1.99000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(4.878e+04_WP)/(8.314_WP*Tloc)) k(rH65) = (7.00000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(-4.184e+03_WP)/(8.314_WP*Tloc)) k(rH66f) = (6.31000000e+14_WP)*Tloc**(0.000_WP)* & exp(-(1.770e+05_WP)/(8.314_WP*Tloc)) k(rH66b) = (1.16600000e+05_WP)*Tloc**(0.600_WP)* & exp(-(-7.410e+03_WP)/(8.314_WP*Tloc)) k(rH67) = (1.43000000e+01_WP)*Tloc**(2.000_WP)* & exp(-(7.950e+03_WP)/(8.314_WP*Tloc)) k(rH68) = (2.00000000e+02_WP)*Tloc**(1.500_WP)* & exp(-(1.259e+05_WP)/(8.314_WP*Tloc)) k(rH69) = (2.12000000e-12_WP)*Tloc**(6.000_WP)* & exp(-(3.968e+04_WP)/(8.314_WP*Tloc)) k(rH70) = (3.50000000e+08_WP)*Tloc**(-0.610_WP)* & exp(-(2.201e+04_WP)/(8.314_WP*Tloc)) k(rH71) = (1.70000000e+23_WP)*Tloc**(-5.310_WP)* & exp(-(2.720e+04_WP)/(8.314_WP*Tloc)) k(rH72) = (2.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.046e+04_WP)/(8.314_WP*Tloc)) kH73_0 = (1.16400000e+33_WP)*Tloc**(-6.821_WP)* & exp(-(1.862e+05_WP)/(8.314_WP*Tloc)) kH73_inf = (1.60600000e+10_WP)*Tloc**(1.028_WP)* & exp(-(1.695e+05_WP)/(8.314_WP*Tloc)) FCH73 = (0.000e+00_WP)* & exp(-Tloc/(1.000e-15_WP)) + (1.000e+00_WP)* & exp(-Tloc/(6.750e+02_WP)) + (1.000e+00_WP)* & exp(-(1.000e+15_WP)/Tloc) k(rH73) = & getlindratecoeff & (Tloc, kH73_0, kH73_inf, FCH73, M(mM11), Ploc ) kH74_0 = (1.50000000e+09_WP)*Tloc**(0.000_WP)* & exp(-(2.320e+05_WP)/(8.314_WP*Tloc)) kH74_inf = (1.80000000e+13_WP)*Tloc**(0.000_WP)* & exp(-(3.180e+05_WP)/(8.314_WP*Tloc)) FCH74 = + (1.000e+00_WP)* & exp(-(0.000e+00_WP)/Tloc) k(rH74) = & getlindratecoeff & (Tloc, kH74_0, kH74_inf, FCH74, M(mM5), Ploc ) k(rH75f) = (8.42000000e-09_WP)*Tloc**(4.620_WP)* & exp(-(1.081e+04_WP)/(8.314_WP*Tloc)) k(rH75b) = (5.72300000e-07_WP)*Tloc**(3.790_WP)* & exp(-(1.353e+04_WP)/(8.314_WP*Tloc)) k(rH76) = (2.23000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(7.192e+04_WP)/(8.314_WP*Tloc)) k(rH77f) = (2.05000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(2.490e+04_WP)/(8.314_WP*Tloc)) k(rH77b) = (6.03300000e+09_WP)*Tloc**(-0.830_WP)* & exp(-(9.104e+04_WP)/(8.314_WP*Tloc)) k(rH78f) = (6.62000000e-06_WP)*Tloc**(3.700_WP)* & exp(-(3.975e+04_WP)/(8.314_WP*Tloc)) k(rH78b) = (1.44000000e-06_WP)*Tloc**(4.020_WP)* & exp(-(2.290e+04_WP)/(8.314_WP*Tloc)) k(rH79) = (1.02000000e+01_WP)*Tloc**(1.880_WP)* & exp(-(7.490e+02_WP)/(8.314_WP*Tloc)) kH80f_0 = (1.11200000e+22_WP)*Tloc**(-5.000_WP)* & exp(-(1.861e+04_WP)/(8.314_WP*Tloc)) kH80f_inf = (1.08100000e+06_WP)*Tloc**(0.450_WP)* & exp(-(7.623e+03_WP)/(8.314_WP*Tloc)) FCH80f = (0.000e+00_WP)* & exp(-Tloc/(1.000e-15_WP)) + (1.000e+00_WP)* & exp(-Tloc/(9.500e+01_WP)) + (1.000e+00_WP)* & exp(-(2.000e+02_WP)/Tloc) k(rH80f) = & getlindratecoeff & (Tloc, kH80f_0, kH80f_inf, FCH80f, M(mM12), Ploc ) k(rH81) = (3.39000000e+00_WP)*Tloc**(1.880_WP)* & exp(-(7.490e+02_WP)/(8.314_WP*Tloc)) k(rH82) = (2.67900000e+02_WP)*Tloc**(0.890_WP)* & exp(-(-8.042e+03_WP)/(8.314_WP*Tloc)) k(rH83) = (1.95000000e+07_WP)*Tloc**(-0.500_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH84) = (5.83100000e+05_WP)*Tloc**(0.599_WP)* & exp(-(-1.219e+04_WP)/(8.314_WP*Tloc)) k(rH85f) = (1.22000000e+24_WP)*Tloc**(-5.760_WP)* & exp(-(4.226e+04_WP)/(8.314_WP*Tloc)) k(rH85b) = (1.25900000e+24_WP)*Tloc**(-5.630_WP)* & exp(-(9.330e+04_WP)/(8.314_WP*Tloc)) k(rH86f) = (5.54000000e-04_WP)*Tloc**(3.500_WP)* & exp(-(2.162e+04_WP)/(8.314_WP*Tloc)) k(rH86b) = (1.35500000e-07_WP)*Tloc**(4.060_WP)* & exp(-(3.706e+04_WP)/(8.314_WP*Tloc)) k(rH87) = (5.80000000e+01_WP)*Tloc**(1.730_WP)* & exp(-(4.853e+03_WP)/(8.314_WP*Tloc)) k(rH88) = (1.20000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH89) = (1.51000000e-13_WP)*Tloc**(6.000_WP)* & exp(-(2.530e+04_WP)/(8.314_WP*Tloc)) k(rH90) = (1.30000000e+01_WP)*Tloc**(2.130_WP)* & exp(-(2.172e+04_WP)/(8.314_WP*Tloc)) k(rH91) = (1.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH92) = (8.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH93f) = (1.10000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH93b) = (2.04600000e+06_WP)*Tloc**(0.890_WP)* & exp(-(1.164e+05_WP)/(8.314_WP*Tloc)) k(rH94) = (2.40000000e+05_WP)*Tloc**(0.000_WP)* & exp(-(-3.573e+03_WP)/(8.314_WP*Tloc)) k(rH95) = (2.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.757e+04_WP)/(8.314_WP*Tloc)) k(rH96) = (4.88700000e+50_WP)*Tloc**(-12.250_WP)* & exp(-(1.175e+05_WP)/(8.314_WP*Tloc)) k(rH97) = (7.00000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(-4.184e+03_WP)/(8.314_WP*Tloc)) k(rH98) = (2.39700000e+48_WP)*Tloc**(-9.900_WP)* & exp(-(3.434e+05_WP)/(8.314_WP*Tloc)) k(rH99) = (7.00000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(-4.184e+03_WP)/(8.314_WP*Tloc)) k(rH100) = (9.72000000e+23_WP)*Tloc**(-5.710_WP)* & exp(-(8.975e+04_WP)/(8.314_WP*Tloc)) k(rH101) = (7.14000000e+09_WP)*Tloc**(-1.210_WP)* & exp(-(8.807e+04_WP)/(8.314_WP*Tloc)) k(rH102) = (3.33200000e+04_WP)*Tloc**(0.340_WP)* & exp(-(-2.326e+03_WP)/(8.314_WP*Tloc)) k(rH103f) = (6.30000000e+02_WP)*Tloc**(1.900_WP)* & exp(-(7.611e+04_WP)/(8.314_WP*Tloc)) k(rH103b) = (1.09700000e+02_WP)*Tloc**(1.890_WP)* & exp(-(6.628e+04_WP)/(8.314_WP*Tloc)) k(rH104) = (3.12000000e+00_WP)*Tloc**(2.000_WP)* & exp(-(-1.247e+03_WP)/(8.314_WP*Tloc)) k(rH105) = (5.24000000e+05_WP)*Tloc**(0.700_WP)* & exp(-(2.462e+04_WP)/(8.314_WP*Tloc)) k(rH106f) = (1.73000000e-01_WP)*Tloc**(2.500_WP)* & exp(-(1.043e+04_WP)/(8.314_WP*Tloc)) k(rH106b) = (7.93300000e-02_WP)*Tloc**(2.510_WP)* & exp(-(8.167e+04_WP)/(8.314_WP*Tloc)) k(rH107f) = (4.83000000e+27_WP)*Tloc**(-5.810_WP)* & exp(-(7.740e+04_WP)/(8.314_WP*Tloc)) k(rH107b) = (2.31300000e+27_WP)*Tloc**(-5.900_WP)* & exp(-(1.323e+05_WP)/(8.314_WP*Tloc)) k(rH108) = (1.58000000e+01_WP)*Tloc**(1.760_WP)* & exp(-(-5.088e+03_WP)/(8.314_WP*Tloc)) k(rH109f) = (2.73000000e+62_WP)*Tloc**(-13.280_WP)* & exp(-(5.155e+05_WP)/(8.314_WP*Tloc)) k(rH109b) = (4.71200000e+53_WP)*Tloc**(-13.190_WP)* & exp(-(1.236e+05_WP)/(8.314_WP*Tloc)) k(rH110) = (2.21000000e-06_WP)*Tloc**(3.500_WP)* & exp(-(2.374e+04_WP)/(8.314_WP*Tloc)) k(rH111) = (4.50000000e+05_WP)*Tloc**(0.000_WP)* & exp(-(2.100e+04_WP)/(8.314_WP*Tloc)) k(rH112f) = (8.56900000e+18_WP)*Tloc**(-1.570_WP)* & exp(-(1.688e+05_WP)/(8.314_WP*Tloc)) k(rH112b) = (1.30000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(6.527e+03_WP)/(8.314_WP*Tloc)) k(rH113) = (3.00000000e+05_WP)*Tloc**(0.000_WP)* & exp(-(1.255e+04_WP)/(8.314_WP*Tloc)) k(rH114f) = (2.66700000e+15_WP)*Tloc**(-0.640_WP)* & exp(-(1.541e+05_WP)/(8.314_WP*Tloc)) k(rH114b) = (1.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.046e+04_WP)/(8.314_WP*Tloc)) k(rH115f) = (2.28400000e+14_WP)*Tloc**(-0.550_WP)* & exp(-(1.188e+05_WP)/(8.314_WP*Tloc)) k(rH115b) = (4.10000000e+05_WP)*Tloc**(0.000_WP)* & exp(-(3.014e+04_WP)/(8.314_WP*Tloc)) k(rH116) = (2.08000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH117) = (2.81000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(2.176e+04_WP)/(8.314_WP*Tloc)) k(rH118) = (4.67000000e+01_WP)*Tloc**(1.610_WP)* & exp(-(-1.460e+02_WP)/(8.314_WP*Tloc)) k(rH119) = (1.29000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(4.853e+04_WP)/(8.314_WP*Tloc)) k(rH120) = (1.68000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(8.548e+04_WP)/(8.314_WP*Tloc)) k(rH121) = (5.60000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(7.406e+04_WP)/(8.314_WP*Tloc)) k(rH122) = (1.05400000e+04_WP)*Tloc**(0.970_WP)* & exp(-(6.636e+03_WP)/(8.314_WP*Tloc)) k(rH123) = (3.98000000e+05_WP)*Tloc**(0.000_WP)* & exp(-(3.975e+04_WP)/(8.314_WP*Tloc)) k(rH124f) = (4.00000000e+07_WP)*Tloc**(0.000_WP)* & exp(-(1.987e+05_WP)/(8.314_WP*Tloc)) k(rH124b) = (2.08000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH125) = (1.13000000e+08_WP)*Tloc**(0.000_WP)* & exp(-(3.284e+04_WP)/(8.314_WP*Tloc)) k(rH126f) = (3.97200000e+00_WP)*Tloc**(2.750_WP)* & exp(-(2.827e+04_WP)/(8.314_WP*Tloc)) k(rH127f) = (1.30000000e+00_WP)*Tloc**(2.400_WP)* & exp(-(1.871e+04_WP)/(8.314_WP*Tloc)) k(rH127b) = (4.70900000e-01_WP)*Tloc**(2.150_WP)* & exp(-(5.096e+04_WP)/(8.314_WP*Tloc)) k(rH128) = (2.02800000e+12_WP)*Tloc**(0.090_WP)* & exp(-(9.858e+04_WP)/(8.314_WP*Tloc)) k(rH129f) = (2.80300000e+17_WP)*Tloc**(-0.620_WP)* & exp(-(1.508e+05_WP)/(8.314_WP*Tloc)) k(rH129b) = (7.54000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) k(rH130f) = (3.36400000e+19_WP)*Tloc**(-1.320_WP)* & exp(-(1.496e+05_WP)/(8.314_WP*Tloc)) k(rH130b) = (4.52000000e+06_WP)*Tloc**(0.000_WP)* & exp(-(0.000e+00_WP)/(8.314_WP*Tloc)) kH5b_0 = (4.07335643e+30_WP)*Tloc**(-5.83_WP)* & exp(-(2.2815e+05_WP)/(8.314_WP*Tloc)) kH5b_inf = (1.65559637e+20_WP)*Tloc**(-1.57_WP)* & exp(-(2.1958e+05_WP)/(8.314_WP*Tloc)) FCH5b = (5.300e-01_WP)* & exp(-Tloc/(1.000e+02_WP)) + (4.700e-01_WP)* & exp(-Tloc/(2.000e+03_WP)) + (1.000e+00_WP)* & exp(-(1.000e+15_WP)/Tloc) k(rH5b) = & getlindratecoeff & (Tloc, kH5b_0, kH5b_inf, FCH5b, M(mM2), Ploc ) kH8b_0 = (5.93678616e+10_WP)*Tloc**(-0.45_WP)* & exp(-(1.9683e+05_WP)/(8.314_WP*Tloc)) kH8b_inf = (2.50193131e+12_WP)*Tloc**(0.56_WP)* & exp(-(2.0150e+05_WP)/(8.314_WP*Tloc)) FCH8b = (5.000e-01_WP)* & exp(-Tloc/(1.000e-30_WP)) + (5.000e-01_WP)* & exp(-Tloc/(1.000e+30_WP)) + (1.000e+00_WP)* & exp(-(1.000e+100_WP)/Tloc) k(rH8b) = & getlindratecoeff & (Tloc, kH8b_0, kH8b_inf, FCH8b, M(mM4), Ploc ) kH25b_0 = (7.03201592e+38_WP)*Tloc**(-6.81_WP)* & exp(-(3.9194e+05_WP)/(8.314_WP*Tloc)) kH25b_inf = (5.70863389e+25_WP)*Tloc**(-2.73_WP)* & exp(-(3.8747e+05_WP)/(8.314_WP*Tloc)) FCH25b = (5.950e-01_WP)* & exp(-Tloc/(1.120e+03_WP)) + (4.050e-01_WP)* & exp(-Tloc/(6.960e+01_WP)) + (1.000e+00_WP)* & exp(-(1.000e+15_WP)/Tloc) k(rH25b) = & getlindratecoeff & (Tloc, kH25b_0, kH25b_inf, FCH25b, M(mM6), Ploc ) k(rH30b) = (1.35039738e+11_WP)*Tloc**(-0.81_WP)* & exp(-(6.7200e+03_WP)/(8.314_WP*Tloc)) kH80b_0 = (3.30091021e+27_WP)*Tloc**(-4.90_WP)* & exp(-(1.7163e+05_WP)/(8.314_WP*Tloc)) kH80b_inf = (3.20888843e+11_WP)*Tloc**(0.55_WP)* & exp(-(1.6065e+05_WP)/(8.314_WP*Tloc)) FCH80b = (0.000e+00_WP)* & exp(-Tloc/(1.000e-15_WP)) + (1.000e+00_WP)* & exp(-Tloc/(9.500e+01_WP)) + (1.000e+00_WP)* & exp(-(2.000e+02_WP)/Tloc) k(rH80b) = & getlindratecoeff & (Tloc, kH80b_0, kH80b_inf, FCH80b, M(mM12), Ploc ) k(rH126b) = (9.85582206e-03_WP)*Tloc**(3.24_WP)* & exp(-(4.3191e+04_WP)/(8.314_WP*Tloc)) return end subroutine get_rate_coefficients ! --- Reaction rates --- ! subroutine get_reaction_rates(w,k,m,c,cqss) use mechanism implicit none real(WP), dimension(nspec) :: c real(WP), dimension(nqss) :: cqss real(WP), dimension(nreac) :: w,k real(WP), dimension(13) :: m w(rH1) = k(rH1) * c(sO)**2_WP * m(mM0) w(rH2) = k(rH2) * c(sO) * c(sH) * m(mM1) w(rH3f) = k(rH3f) * c(sO) * c(sH2) w(rH3b) = k(rH3b) * c(sH) * c(sOH) w(rH4) = k(rH4) * c(sH) * c(sOH) * m(mM3) w(rH5f) = k(rH5f) * c(sOH)**2_WP w(rH6f) = k(rH6f) * c(sOH) * c(sH2) w(rH6b) = k(rH6b) * c(sH) * c(sH2O) w(rH7f) = k(rH7f) * c(sO) * c(sH2O) w(rH7b) = k(rH7b) * c(sOH)**2_WP w(rH8f) = k(rH8f) * c(sH) * c(sO2) w(rH9f) = k(rH9f) * c(sH) * c(sO2) w(rH9b) = k(rH9b) * c(sO) * c(sOH) w(rH10f) = k(rH10f) * c(sHO2) * c(sOH) w(rH10b) = k(rH10b) * c(sH2O) * c(sO2) w(rH11f) = k(rH11f) * c(sHO2) * c(sO) w(rH11b) = k(rH11b) * c(sOH) * c(sO2) w(rH12) = k(rH12) * c(sHO2) * c(sH) w(rH13) = k(rH13) * c(sHO2)**2_WP w(rH14) = k(rH14) * c(sHO2) * c(sH) w(rH15) = k(rH15) * c(sHO2)**2_WP w(rH16f) = k(rH16f) * c(sH2O2) * c(sH) w(rH16b) = k(rH16b) * c(sH2) * c(sHO2) w(rH17) = k(rH17) * c(sH2O2) * c(sO) w(rH18) = k(rH18) * c(sH2O2) * c(sOH) w(rH19f) = k(rH19f) * c(sH2O2) * c(sOH) w(rH19b) = k(rH19b) * c(sH2O) * c(sHO2) w(rH20) = k(rH20) * c(sH2O2) * c(sH) w(rH21) = k(rH21) * cqss(sCH2GSGXCH2-nspec) * c(sO) w(rH22) = k(rH22) * cqss(sCH2GSGXCH2-nspec) * c(sO2) w(rH23) = k(rH23) * cqss(sCH2GSGXCH2-nspec) * c(sOH) w(rH24f) = k(rH24f) * cqss(sCH2GSGXCH2-nspec) * c(sH2) w(rH24b) = k(rH24b) * c(sCH3) * c(sH) w(rH25f) = k(rH25f) * c(sCH3)**2_WP w(rH26f) = k(rH26f) * c(sCH3) * c(sO2) w(rH26b) = k(rH26b) * cqss(sCH3O-nspec) * c(sO) w(rH27) = k(rH27) * c(sCH3) * c(sHO2) w(rH28f) = k(rH28f) * c(sCH3) * c(sOH) w(rH28b) = k(rH28b) * cqss(sCH2GSGXCH2-nspec) * c(sH2O) w(rH29) = k(rH29) * c(sCH3) * c(sH) w(rH30f) = k(rH30f) * c(sCH3)**2_WP w(rH31) = k(rH31) * c(sCH3) * c(sO2) w(rH32) = k(rH32) * c(sCH3) * c(sO) w(rH33) = k(rH33) * c(sCH3) * c(sHO2) w(rH34) = k(rH34) * c(sCH3) * c(sOH) w(rH35) = k(rH35) * cqss(sCH2GSGXCH2-nspec) * c(sCH3) w(rH36) = k(rH36) * c(sCH3) * c(sH2O2) w(rH37f) = k(rH37f) * c(sCH4) * c(sOH) w(rH37b) = k(rH37b) * c(sCH3) * c(sH2O) w(rH38f) = k(rH38f) * c(sCH4) * c(sH) w(rH38b) = k(rH38b) * c(sCH3) * c(sH2) w(rH39f) = k(rH39f) * cqss(sCH2GSGXCH2-nspec) * c(sCH4) w(rH39b) = k(rH39b) * c(sCH3)**2_WP w(rH40f) = k(rH40f) * c(sCH4) * c(sO) w(rH40b) = k(rH40b) * c(sCH3) * c(sOH) w(rH41) = k(rH41) * c(sCO) * c(sHO2) w(rH42) = k(rH42) * c(sCO) * c(sO) w(rH43f) = k(rH43f) * c(sCO) * c(sOH) w(rH43b) = k(rH43b) * c(sCO2) * c(sH) w(rH44f) = k(rH44f) * c(sCO) * c(sO2) w(rH44b) = k(rH44b) * c(sCO2) * c(sO) w(rH45) = k(rH45) * cqss(sHCO-nspec) * c(sO) w(rH46) = k(rH46) * cqss(sHCO-nspec) * c(sO) w(rH47) = k(rH47) * cqss(sHCO-nspec) * c(sH) w(rH48f) = k(rH48f) * cqss(sHCO-nspec) * m(mM9) w(rH48b) = k(rH48b) * c(sH) * c(sCO) * m(mM9) w(rH49) = k(rH49) * cqss(sHCO-nspec) * c(sOH) w(rH50) = k(rH50) * cqss(sHCO-nspec) * c(sCH3) w(rH51) = k(rH51) * cqss(sHCO-nspec) * c(sO2) w(rH52) = k(rH52) * c(sCH2O) * c(sHO2) w(rH53) = k(rH53) * c(sCH2O) * c(sOH) w(rH54) = k(rH54) * c(sCH2O) * c(sH) w(rH55) = k(rH55) * c(sCH2O) * c(sCH3) w(rH56) = k(rH56) * c(sCH2O) * c(sO) w(rH57) = k(rH57) * c(sCH2O) * c(sO2) w(rH58) = k(rH58) * cqss(sCH3O-nspec) * c(sO2) w(rH59) = k(rH59) * cqss(sCH3O-nspec) w(rH60) = k(rH60) * cqss(sCH2GSGXCH2-nspec) * c(sCO2) w(rH61) = k(rH61) * c(sCH3O2) * c(sHO2) w(rH62) = k(rH62) * c(sCH3O2)**2_WP w(rH63f) = k(rH63f) * c(sCH3O2) * m(mM5) w(rH63b) = k(rH63b) * c(sCH3) * c(sO2) * m(mM5) w(rH64) = k(rH64) * c(sCH3O2) * c(sCH2O) w(rH65) = k(rH65) * c(sCH3O2) * c(sCH3) w(rH66f) = k(rH66f) * c(sCH3O2H) w(rH66b) = k(rH66b) * cqss(sCH3O-nspec) * c(sOH) w(rH67) = k(rH67) * c(sC2H2) * c(sO) w(rH68) = k(rH68) * c(sC2H2) * c(sO2) w(rH69) = k(rH69) * cqss(sC2H3-nspec) * c(sO2) w(rH70) = k(rH70) * cqss(sC2H3-nspec) * c(sO2) w(rH71) = k(rH71) * cqss(sC2H3-nspec) * c(sO2) w(rH72) = k(rH72) * cqss(sC2H3-nspec) * c(sH) w(rH73) = k(rH73) * cqss(sC2H3-nspec) w(rH74) = k(rH74) * c(sC2H4) w(rH75f) = k(rH75f) * c(sC2H4) * c(sH) w(rH75b) = k(rH75b) * cqss(sC2H3-nspec) * c(sH2) w(rH76) = k(rH76) * c(sC2H4) * c(sCH3O2) w(rH77f) = k(rH77f) * c(sC2H4) * c(sOH) w(rH77b) = k(rH77b) * cqss(sC2H3-nspec) * c(sH2O) w(rH78f) = k(rH78f) * c(sC2H4) * c(sCH3) w(rH78b) = k(rH78b) * cqss(sC2H3-nspec) * c(sCH4) w(rH79) = k(rH79) * c(sC2H4) * c(sO) w(rH80f) = k(rH80f) * c(sH) * c(sC2H4) w(rH81) = k(rH81) * c(sC2H4) * c(sO) w(rH82) = k(rH82) * cqss(sC2H5-nspec) * c(sHO2) w(rH83) = k(rH83) * c(sCH3) * cqss(sC2H5-nspec) w(rH84) = k(rH84) * c(sH) * cqss(sC2H5-nspec) w(rH85f) = k(rH85f) * cqss(sC2H5-nspec) * c(sO2) w(rH85b) = k(rH85b) * c(sC2H4) * c(sHO2) w(rH86f) = k(rH86f) * c(sC2H6) * c(sH) w(rH86b) = k(rH86b) * cqss(sC2H5-nspec) * c(sH2) w(rH87) = k(rH87) * c(sC2H6) * c(sOH) w(rH88) = k(rH88) * cqss(sCH2GSGXCH2-nspec) * c(sC2H6) w(rH89) = k(rH89) * c(sC2H6) * c(sCH3) w(rH90) = k(rH90) * c(sC2H6) * c(sO) w(rH91) = k(rH91) * cqss(sHCCO-nspec) * c(sOH) w(rH92) = k(rH92) * cqss(sHCCO-nspec) * c(sO) w(rH93f) = k(rH93f) * cqss(sHCCO-nspec) * c(sH) w(rH93b) = k(rH93b) * cqss(sCH2GSGXCH2-nspec) * c(sCO) w(rH94) = k(rH94) * cqss(sHCCO-nspec) * c(sO2) w(rH95) = k(rH95) * cqss(sCH2CHO-nspec) * c(sO2) w(rH96) = k(rH96) * cqss(sC3H5XAXC3H5-nspec) * c(sH) w(rH97) = k(rH97) * cqss(sC3H5XAXC3H5-nspec) * c(sHO2) w(rH98) = k(rH98) * cqss(sC3H5XAXC3H5-nspec) w(rH99) = k(rH99) * cqss(sC3H5XAXC3H5-nspec) * c(sCH3O2) w(rH100) = k(rH100) * cqss(sC3H5XAXC3H5-nspec) * c(sO2) w(rH101) = k(rH101) * cqss(sC3H5XAXC3H5-nspec) * c(sO2) w(rH102) = k(rH102) * cqss(sC3H5XAXC3H5-nspec) * c(sHO2) w(rH103f) = k(rH103f) * cqss(sC3H5XAXC3H5-nspec) * c(sCH2O) w(rH103b) = k(rH103b) * c(sC3H6) * cqss(sHCO-nspec) w(rH104) = k(rH104) * c(sC3H6) * c(sOH) w(rH105) = k(rH105) * c(sC3H6) * c(sO) w(rH106f) = k(rH106f) * c(sC3H6) * c(sH) w(rH106b) = k(rH106b) * cqss(sC3H5XAXC3H5-nspec) * c(sH2) w(rH107f) = k(rH107f) * c(sC3H6) * c(sH) w(rH107b) = k(rH107b) * c(sC2H4) * c(sCH3) w(rH108) = k(rH108) * c(sC3H6) * c(sO) w(rH109f) = k(rH109f) * c(sC3H6) w(rH109b) = k(rH109b) * cqss(sC2H3-nspec) * c(sCH3) w(rH110) = k(rH110) * c(sC3H6) * c(sCH3) w(rH111) = k(rH111) * cqss(sIXC3H7-nspec) * c(sO2) w(rH112f) = k(rH112f) * cqss(sIXC3H7-nspec) w(rH112b) = k(rH112b) * c(sH) * c(sC3H6) w(rH113) = k(rH113) * cqss(sNXC3H7-nspec) * c(sO2) w(rH114f) = k(rH114f) * cqss(sNXC3H7-nspec) w(rH114b) = k(rH114b) * c(sH) * c(sC3H6) w(rH115f) = k(rH115f) * cqss(sNXC3H7-nspec) w(rH115b) = k(rH115b) * c(sCH3) * c(sC2H4) w(rH116) = k(rH116) * cqss(sNXC3H7-nspec) * c(sHO2) w(rH117) = k(rH117) * c(sC3H8) * c(sO) w(rH118) = k(rH118) * c(sC3H8) * c(sOH) w(rH119) = k(rH119) * c(sCH3) * c(sC3H8) w(rH120) = k(rH120) * c(sC3H8) * c(sHO2) w(rH121) = k(rH121) * c(sC3H8) * c(sHO2) w(rH122) = k(rH122) * c(sC3H8) * c(sOH) w(rH123) = k(rH123) * c(sCH3) * c(sC3H8) w(rH124f) = k(rH124f) * c(sC3H8) * c(sO2) w(rH124b) = k(rH124b) * cqss(sIXC3H7-nspec) * c(sHO2) w(rH125) = k(rH125) * c(sC3H8) * c(sO) w(rH126f) = k(rH126f) * c(sH) * c(sC3H8) w(rH127f) = k(rH127f) * c(sH) * c(sC3H8) w(rH127b) = k(rH127b) * c(sH2) * cqss(sIXC3H7-nspec) w(rH128) = k(rH128) * c(sC3H5O) w(rH129f) = k(rH129f) * cqss(sIXC3H7O2-nspec) w(rH129b) = k(rH129b) * cqss(sIXC3H7-nspec) * c(sO2) w(rH130f) = k(rH130f) * cqss(sNXC3H7O2-nspec) w(rH130b) = k(rH130b) * cqss(sNXC3H7-nspec) * c(sO2) w(rH5b) = k(rH5b) * c(sH2O2) w(rH8b) = k(rH8b) * c(sHO2) w(rH25b) = k(rH25b) * c(sC2H6) w(rH30b) = k(rH30b) * c(sH) * cqss(sC2H5-nspec) w(rH80b) = k(rH80b) * cqss(sC2H5-nspec) w(rH126b) = k(rH126b) * c(sH2) * cqss(sNXC3H7-nspec) return end subroutine get_reaction_rates ! --- Production rates --- ! subroutine get_production_rates(cdot,w) use mechanism implicit none real(WP), dimension(nspec) :: cdot real(WP), dimension(nreac) :: w cdot(sN2) = 0.0_WP cdot(sO) = 0.0_WP - 2_WP * & w(rH1) - & w(rH2) - & w(rH3f) + & w(rH3b) - & w(rH7f) + & w(rH7b) + & w(rH9f) - & w(rH9b) - & w(rH11f) + & w(rH11b) - & w(rH17) - & w(rH21) + & w(rH26f) - & w(rH26b) - & w(rH32) - & w(rH40f) + & w(rH40b) - & w(rH42) + & w(rH44f) - & w(rH44b) - & w(rH45) - & w(rH46) - & w(rH56) - & w(rH67) + & w(rH70) - & w(rH79) - & w(rH81) - & w(rH90) - & w(rH92) - & w(rH105) - & w(rH108) - & w(rH117) - & w(rH125) cdot(sO2) = 0.0_WP + & w(rH1) - & w(rH8f) - & w(rH9f) + & w(rH9b) + & w(rH10f) - & w(rH10b) + & w(rH11f) - & w(rH11b) + & w(rH12) + & w(rH13) + & w(rH15) - & w(rH22) - & w(rH26f) + & w(rH26b) - & w(rH31) + & w(rH33) - & w(rH44f) + & w(rH44b) - & w(rH51) - & w(rH57) - & w(rH58) + & w(rH61) + & w(rH62) + & w(rH63f) - & w(rH63b) - & w(rH68) - & w(rH69) - & w(rH70) - & w(rH71) + & w(rH82) - & w(rH85f) + & w(rH85b) - & w(rH94) - & w(rH95) - & w(rH100) - & w(rH101) + & w(rH102) - & w(rH111) - & w(rH113) + & w(rH116) - & w(rH124f) + & w(rH124b) + & w(rH129f) - & w(rH129b) + & w(rH130f) - & w(rH130b) + & w(rH8b) cdot(sH) = 0.0_WP - & w(rH2) + & w(rH3f) - & w(rH3b) - & w(rH4) + & w(rH6f) - & w(rH6b) - & w(rH8f) - & w(rH9f) + & w(rH9b) - & w(rH12) - & w(rH14) - & w(rH16f) + & w(rH16b) - & w(rH20) + 2_WP * & w(rH21) + & w(rH22) + & w(rH23) + & w(rH24f) - & w(rH24b) - & w(rH29) + & w(rH30f) + & w(rH32) + & w(rH35) - & w(rH38f) + & w(rH38b) + & w(rH43f) - & w(rH43b) + & w(rH46) - & w(rH47) + & w(rH48f) - & w(rH48b) - & w(rH54) + & w(rH59) + & w(rH67) - & w(rH72) + & w(rH73) - & w(rH75f) + & w(rH75b) - & w(rH80f) + & w(rH81) - & w(rH84) - & w(rH86f) + & w(rH86b) + & w(rH92) - & w(rH93f) + & w(rH93b) - & w(rH96) - & w(rH106f) + & w(rH106b) - & w(rH107f) cdot(sH) = cdot(sH) + & w(rH107b) + & w(rH112f) - & w(rH112b) + & w(rH114f) - & w(rH114b) - & w(rH126f) - & w(rH127f) + & w(rH127b) + & w(rH8b) - & w(rH30b) + & w(rH80b) + & w(rH126b) cdot(sOH) = 0.0_WP + & w(rH2) + & w(rH3f) - & w(rH3b) - & w(rH4) - 2_WP * & w(rH5f) - & w(rH6f) + & w(rH6b) + 2_WP * & w(rH7f) - 2_WP * & w(rH7b) + & w(rH9f) - & w(rH9b) - & w(rH10f) + & w(rH10b) + & w(rH11f) - & w(rH11b) + 2_WP * & w(rH14) + & w(rH17) - & w(rH18) - & w(rH19f) + & w(rH19b) + & w(rH20) + & w(rH22) - & w(rH23) + & w(rH27) - & w(rH28f) + & w(rH28b) + & w(rH31) - & w(rH34) - & w(rH37f) + & w(rH37b) + & w(rH40f) - & w(rH40b) + & w(rH41) - & w(rH43f) + & w(rH43b) + & w(rH45) - & w(rH49) - & w(rH53) + & w(rH56) + & w(rH66f) - & w(rH66b) + & w(rH68) - & w(rH77f) + & w(rH77b) - & w(rH87) + & w(rH90) - & w(rH91) + & w(rH95) + & w(rH97) + & w(rH100) cdot(sOH) = cdot(sOH) - & w(rH104) + & w(rH105) + & w(rH117) - & w(rH118) - & w(rH122) + & w(rH125) + 2_WP * & w(rH5b) cdot(sH2) = 0.0_WP - & w(rH3f) + & w(rH3b) - & w(rH6f) + & w(rH6b) + & w(rH12) + & w(rH16f) - & w(rH16b) - & w(rH24f) + & w(rH24b) + & w(rH34) + & w(rH38f) - & w(rH38b) + & w(rH47) + & w(rH54) + & w(rH72) + & w(rH74) + & w(rH75f) - & w(rH75b) + & w(rH86f) - & w(rH86b) + & w(rH106f) - & w(rH106b) + & w(rH126f) + & w(rH127f) - & w(rH127b) - & w(rH126b) cdot(sH2O) = 0.0_WP + & w(rH4) + & w(rH6f) - & w(rH6b) - & w(rH7f) + & w(rH7b) + & w(rH10f) - & w(rH10b) + & w(rH18) + & w(rH19f) - & w(rH19b) + & w(rH20) + & w(rH28f) - & w(rH28b) + & w(rH37f) - & w(rH37b) + & w(rH49) + & w(rH53) + & w(rH77f) - & w(rH77b) + & w(rH87) + & w(rH104) + & w(rH118) + & w(rH122) cdot(sH2O2) = 0.0_WP + & w(rH5f) + & w(rH13) + & w(rH15) - & w(rH16f) + & w(rH16b) - & w(rH17) - & w(rH18) - & w(rH19f) + & w(rH19b) - & w(rH20) - & w(rH36) + & w(rH52) + & w(rH120) + & w(rH121) - & w(rH5b) cdot(sHO2) = 0.0_WP + & w(rH8f) - & w(rH10f) + & w(rH10b) - & w(rH11f) + & w(rH11b) - & w(rH12) - 2_WP * & w(rH13) - & w(rH14) - 2_WP * & w(rH15) + & w(rH16f) - & w(rH16b) + & w(rH17) + & w(rH18) + & w(rH19f) - & w(rH19b) - & w(rH27) - & w(rH33) + & w(rH36) - & w(rH41) + & w(rH51) - & w(rH52) + & w(rH57) + & w(rH58) - & w(rH61) + & w(rH69) - & w(rH82) + & w(rH85f) - & w(rH85b) - & w(rH97) - & w(rH102) + & w(rH111) + & w(rH113) - & w(rH116) - & w(rH120) - & w(rH121) + & w(rH124f) - & w(rH124b) - & w(rH8b) !cdot(sCH2GSGXCH2) = 0.0_WP cdot(sCO) = 0.0_WP + & w(rH21) + & w(rH22) - & w(rH41) - & w(rH42) - & w(rH43f) + & w(rH43b) - & w(rH44f) + & w(rH44b) + & w(rH45) + & w(rH47) + & w(rH48f) - & w(rH48b) + & w(rH49) + & w(rH50) + & w(rH51) + & w(rH60) + 2_WP * & w(rH92) + & w(rH93f) - & w(rH93b) + & w(rH95) cdot(sCH2O) = 0.0_WP + & w(rH23) + & w(rH31) + & w(rH32) + & w(rH34) - & w(rH52) - & w(rH53) - & w(rH54) - & w(rH55) - & w(rH56) - & w(rH57) + & w(rH58) + & w(rH59) + & w(rH60) - & w(rH64) + & w(rH71) + & w(rH95) + & w(rH100) + & w(rH101) - & w(rH103f) + & w(rH103b) + & w(rH128) cdot(sCH3) = 0.0_WP + & w(rH24f) - & w(rH24b) - 2_WP * & w(rH25f) - & w(rH26f) + & w(rH26b) - & w(rH27) - & w(rH28f) + & w(rH28b) - & w(rH29) - 2_WP * & w(rH30f) - & w(rH31) - & w(rH32) - & w(rH33) - & w(rH34) - & w(rH35) - & w(rH36) + & w(rH37f) - & w(rH37b) + & w(rH38f) - & w(rH38b) + 2_WP * & w(rH39f) - 2_WP * & w(rH39b) + & w(rH40f) - & w(rH40b) - & w(rH50) - & w(rH55) + & w(rH63f) - & w(rH63b) - & w(rH65) - & w(rH78f) + & w(rH78b) + & w(rH79) - & w(rH83) + & w(rH88) - & w(rH89) + & w(rH98) + & w(rH107f) - & w(rH107b) + & w(rH109f) - & w(rH109b) - & w(rH110) + & w(rH115f) - & w(rH115b) - & w(rH119) - & w(rH123) + 2_WP * & w(rH25b) + 2_WP * & w(rH30b) cdot(sC2H6) = 0.0_WP + & w(rH25f) + & w(rH82) + & w(rH84) - & w(rH86f) + & w(rH86b) - & w(rH87) - & w(rH88) - & w(rH89) - & w(rH90) - & w(rH25b) !cdot(sCH3O) = 0.0_WP cdot(sCH4) = 0.0_WP + & w(rH29) + & w(rH33) + & w(rH36) - & w(rH37f) + & w(rH37b) - & w(rH38f) + & w(rH38b) - & w(rH39f) + & w(rH39b) - & w(rH40f) + & w(rH40b) + & w(rH50) + & w(rH55) + & w(rH78f) - & w(rH78b) + & w(rH83) + & w(rH89) + & w(rH110) + & w(rH119) + & w(rH123) !cdot(sC2H5) = 0.0_WP cdot(sC2H4) = 0.0_WP + & w(rH35) - & w(rH74) - & w(rH75f) + & w(rH75b) - & w(rH76) - & w(rH77f) + & w(rH77b) - & w(rH78f) + & w(rH78b) - & w(rH79) - & w(rH80f) - & w(rH81) + & w(rH83) + & w(rH85f) - & w(rH85b) + & w(rH107f) - & w(rH107b) + & w(rH115f) - & w(rH115b) + & w(rH80b) cdot(sCO2) = 0.0_WP + & w(rH41) + & w(rH42) + & w(rH43f) - & w(rH43b) + & w(rH44f) - & w(rH44b) + & w(rH46) - & w(rH60) + & w(rH94) !cdot(sHCO) = 0.0_WP cdot(sCH3O2) = 0.0_WP - & w(rH61) - 2_WP * & w(rH62) - & w(rH63f) + & w(rH63b) - & w(rH64) - & w(rH65) - & w(rH76) - & w(rH99) cdot(sCH3O2H) = 0.0_WP + & w(rH61) + & w(rH64) - & w(rH66f) + & w(rH66b) + & w(rH76) cdot(sC2H2) = 0.0_WP - & w(rH67) - & w(rH68) + & w(rH69) + & w(rH72) + & w(rH73) + & w(rH74) + & w(rH98) + & w(rH100) !cdot(sHCCO) = 0.0_WP !cdot(sC2H3) = 0.0_WP !cdot(sCH2CHO) = 0.0_WP !cdot(sC3H5XAXC3H5) = 0.0_WP cdot(sC3H6) = 0.0_WP + & w(rH96) + & w(rH102) + & w(rH103f) - & w(rH103b) - & w(rH104) - & w(rH105) - & w(rH106f) + & w(rH106b) - & w(rH107f) + & w(rH107b) - & w(rH108) - & w(rH109f) + & w(rH109b) - & w(rH110) + & w(rH111) + & w(rH112f) - & w(rH112b) + & w(rH113) + & w(rH114f) - & w(rH114b) cdot(sC3H5O) = 0.0_WP + & w(rH97) + & w(rH99) - & w(rH128) !cdot(sIXC3H7) = 0.0_WP !cdot(sNXC3H7) = 0.0_WP cdot(sC3H8) = 0.0_WP + & w(rH116) - & w(rH117) - & w(rH118) - & w(rH119) - & w(rH120) - & w(rH121) - & w(rH122) - & w(rH123) - & w(rH124f) + & w(rH124b) - & w(rH125) - & w(rH126f) - & w(rH127f) + & w(rH127b) + & w(rH126b) !cdot(sIXC3H7O2) = 0.0_WP !cdot(sNXC3H7O2) = 0.0_WP return end subroutine get_production_rates ! --- Actual reactions --- ! subroutine reaction_expressions use mechanism implicit none reacexp(1) = '2 O + M0 -> O2 + M0' reacexp(2) = 'O + H + M1 -> OH + M1' reacexp(3) = 'O + H2 -> H + OH' reacexp(4) = 'H + OH -> O + H2' reacexp(5) = 'H + OH + M3 -> H2O + M3' reacexp(6) = '2 OH + M2 -> H2O2 + M2' reacexp(7) = 'OH + H2 -> H + H2O' reacexp(8) = 'H + H2O -> OH + H2' reacexp(9) = 'O + H2O -> 2 OH' reacexp(10) = '2 OH -> O + H2O' reacexp(11) = 'H + O2 + M4 -> HO2 + M4' reacexp(12) = 'H + O2 -> O + OH' reacexp(13) = 'O + OH -> H + O2' reacexp(14) = 'HO2 + OH -> H2O + O2' reacexp(15) = 'H2O + O2 -> HO2 + OH' reacexp(16) = 'HO2 + O -> OH + O2' reacexp(17) = 'OH + O2 -> HO2 + O' reacexp(18) = 'HO2 + H -> H2 + O2' reacexp(19) = '2 HO2 -> H2O2 + O2' reacexp(20) = 'HO2 + H -> 2 OH' reacexp(21) = '2 HO2 -> H2O2 + O2' reacexp(22) = 'H2O2 + H -> H2 + HO2' reacexp(23) = 'H2 + HO2 -> H2O2 + H' reacexp(24) = 'H2O2 + O -> OH + HO2' reacexp(25) = 'H2O2 + OH -> H2O + HO2' reacexp(26) = 'H2O2 + OH -> H2O + HO2' reacexp(27) = 'H2O + HO2 -> H2O2 + OH' reacexp(28) = 'H2O2 + H -> H2O + OH' reacexp(29) = 'CH2GSG-CH2 + O -> CO + 2 H' reacexp(30) = 'CH2GSG-CH2 + O2 -> CO + OH + H' reacexp(31) = 'CH2GSG-CH2 + OH -> CH2O + H' reacexp(32) = 'CH2GSG-CH2 + H2 -> CH3 + H' reacexp(33) = 'CH3 + H -> CH2GSG-CH2 + H2' reacexp(34) = '2 CH3 + M6 -> C2H6 + M6' reacexp(35) = 'CH3 + O2 -> CH3O + O' reacexp(36) = 'CH3O + O -> CH3 + O2' reacexp(37) = 'CH3 + HO2 -> CH3O + OH' reacexp(38) = 'CH3 + OH -> CH2GSG-CH2 + H2O' reacexp(39) = 'CH2GSG-CH2 + H2O -> CH3 + OH' reacexp(40) = 'CH3 + H + M7 -> CH4 + M7' reacexp(41) = '2 CH3 -> H + C2H5' reacexp(42) = 'CH3 + O2 -> CH2O + OH' reacexp(43) = 'CH3 + O -> CH2O + H' reacexp(44) = 'CH3 + HO2 -> CH4 + O2' reacexp(45) = 'CH3 + OH -> CH2O + H2' reacexp(46) = 'CH2GSG-CH2 + CH3 -> C2H4 + H' reacexp(47) = 'CH3 + H2O2 -> CH4 + HO2' reacexp(48) = 'CH4 + OH -> CH3 + H2O' reacexp(49) = 'CH3 + H2O -> CH4 + OH' reacexp(50) = 'CH4 + H -> CH3 + H2' reacexp(51) = 'CH3 + H2 -> CH4 + H' reacexp(52) = 'CH2GSG-CH2 + CH4 -> 2 CH3' reacexp(53) = '2 CH3 -> CH2GSG-CH2 + CH4' reacexp(54) = 'CH4 + O -> CH3 + OH' reacexp(55) = 'CH3 + OH -> CH4 + O' reacexp(56) = 'CO + HO2 -> CO2 + OH' reacexp(57) = 'CO + O + M8 -> CO2 + M8' reacexp(58) = 'CO + OH -> CO2 + H' reacexp(59) = 'CO2 + H -> CO + OH' reacexp(60) = 'CO + O2 -> CO2 + O' reacexp(61) = 'CO2 + O -> CO + O2' reacexp(62) = 'HCO + O -> CO + OH' reacexp(63) = 'HCO + O -> CO2 + H' reacexp(64) = 'HCO + H -> CO + H2' reacexp(65) = 'HCO + M9 -> H + CO + M9' reacexp(66) = 'H + CO + M9 -> HCO + M9' reacexp(67) = 'HCO + OH -> CO + H2O' reacexp(68) = 'HCO + CH3 -> CH4 + CO' reacexp(69) = 'HCO + O2 -> CO + HO2' reacexp(70) = 'CH2O + HO2 -> HCO + H2O2' reacexp(71) = 'CH2O + OH -> HCO + H2O' reacexp(72) = 'CH2O + H -> HCO + H2' reacexp(73) = 'CH2O + CH3 -> HCO + CH4' reacexp(74) = 'CH2O + O -> HCO + OH' reacexp(75) = 'CH2O + O2 -> HCO + HO2' reacexp(76) = 'CH3O + O2 -> CH2O + HO2' reacexp(77) = 'CH3O + M5 -> CH2O + H + M5' reacexp(78) = 'CH2GSG-CH2 + CO2 -> CH2O + CO' reacexp(79) = 'CH3O2 + HO2 -> CH3O2H + O2' reacexp(80) = '2 CH3O2 -> O2 + 2 CH3O' reacexp(81) = 'CH3O2 + M5 -> CH3 + O2 + M5' reacexp(82) = 'CH3 + O2 + M5 -> CH3O2 + M5' reacexp(83) = 'CH3O2 + CH2O -> CH3O2H + HCO' reacexp(84) = 'CH3O2 + CH3 -> 2 CH3O' reacexp(85) = 'CH3O2H -> CH3O + OH' reacexp(86) = 'CH3O + OH -> CH3O2H' reacexp(87) = 'C2H2 + O -> HCCO + H' reacexp(88) = 'C2H2 + O2 -> HCCO + OH' reacexp(89) = 'C2H3 + O2 -> C2H2 + HO2' reacexp(90) = 'C2H3 + O2 -> CH2CHO + O' reacexp(91) = 'C2H3 + O2 -> CH2O + HCO' reacexp(92) = 'C2H3 + H -> C2H2 + H2' reacexp(93) = 'C2H3 + M11 -> H + C2H2 + M11' reacexp(94) = 'C2H4 + M5 -> C2H2 + H2 + M5' reacexp(95) = 'C2H4 + H -> C2H3 + H2' reacexp(96) = 'C2H3 + H2 -> C2H4 + H' reacexp(97) = 'C2H4 + CH3O2 -> C2H3 + CH3O2H' reacexp(98) = 'C2H4 + OH -> C2H3 + H2O' reacexp(99) = 'C2H3 + H2O -> C2H4 + OH' reacexp(100) = 'C2H4 + CH3 -> C2H3 + CH4' reacexp(101) = 'C2H3 + CH4 -> C2H4 + CH3' reacexp(102) = 'C2H4 + O -> CH3 + HCO' reacexp(103) = 'H + C2H4 + M12 -> C2H5 + M12' reacexp(104) = 'C2H4 + O -> CH2CHO + H' reacexp(105) = 'C2H5 + HO2 -> C2H6 + O2' reacexp(106) = 'CH3 + C2H5 -> CH4 + C2H4' reacexp(107) = 'H + C2H5 -> C2H6' reacexp(108) = 'C2H5 + O2 -> C2H4 + HO2' reacexp(109) = 'C2H4 + HO2 -> C2H5 + O2' reacexp(110) = 'C2H6 + H -> C2H5 + H2' reacexp(111) = 'C2H5 + H2 -> C2H6 + H' reacexp(112) = 'C2H6 + OH -> C2H5 + H2O' reacexp(113) = 'CH2GSG-CH2 + C2H6 -> CH3 + C2H5' reacexp(114) = 'C2H6 + CH3 -> C2H5 + CH4' reacexp(115) = 'C2H6 + O -> C2H5 + OH' reacexp(116) = 'HCCO + OH -> 2 HCO' reacexp(117) = 'HCCO + O -> H + 2 CO' reacexp(118) = 'HCCO + H -> CH2GSG-CH2 + CO' reacexp(119) = 'CH2GSG-CH2 + CO -> HCCO + H' reacexp(120) = 'HCCO + O2 -> CO2 + HCO' reacexp(121) = 'CH2CHO + O2 -> CH2O + CO + OH' reacexp(122) = 'C3H5-A-C3H5 + H -> C3H6' reacexp(123) = 'C3H5-A-C3H5 + HO2 -> C3H5O + OH' reacexp(124) = 'C3H5-A-C3H5 -> C2H2 + CH3' reacexp(125) = 'C3H5-A-C3H5 + CH3O2 -> C3H5O + CH3O' reacexp(126) = 'C3H5-A-C3H5 + O2 -> C2H2 + CH2O + OH' reacexp(127) = 'C3H5-A-C3H5 + O2 -> CH2CHO + CH2O' reacexp(128) = 'C3H5-A-C3H5 + HO2 -> C3H6 + O2' reacexp(129) = 'C3H5-A-C3H5 + CH2O -> C3H6 + HCO' reacexp(130) = 'C3H6 + HCO -> C3H5-A-C3H5 + CH2O' reacexp(131) = 'C3H6 + OH -> C3H5-A-C3H5 + H2O' reacexp(132) = 'C3H6 + O -> C3H5-A-C3H5 + OH' reacexp(133) = 'C3H6 + H -> C3H5-A-C3H5 + H2' reacexp(134) = 'C3H5-A-C3H5 + H2 -> C3H6 + H' reacexp(135) = 'C3H6 + H -> C2H4 + CH3' reacexp(136) = 'C2H4 + CH3 -> C3H6 + H' reacexp(137) = 'C3H6 + O -> C2H5 + HCO' reacexp(138) = 'C3H6 -> C2H3 + CH3' reacexp(139) = 'C2H3 + CH3 -> C3H6' reacexp(140) = 'C3H6 + CH3 -> C3H5-A-C3H5 + CH4' reacexp(141) = 'I-C3H7 + O2 -> C3H6 + HO2' reacexp(142) = 'I-C3H7 -> H + C3H6' reacexp(143) = 'H + C3H6 -> I-C3H7' reacexp(144) = 'N-C3H7 + O2 -> C3H6 + HO2' reacexp(145) = 'N-C3H7 -> H + C3H6' reacexp(146) = 'H + C3H6 -> N-C3H7' reacexp(147) = 'N-C3H7 -> CH3 + C2H4' reacexp(148) = 'CH3 + C2H4 -> N-C3H7' reacexp(149) = 'N-C3H7 + HO2 -> C3H8 + O2' reacexp(150) = 'C3H8 + O -> I-C3H7 + OH' reacexp(151) = 'C3H8 + OH -> I-C3H7 + H2O' reacexp(152) = 'CH3 + C3H8 -> CH4 + N-C3H7' reacexp(153) = 'C3H8 + HO2 -> N-C3H7 + H2O2' reacexp(154) = 'C3H8 + HO2 -> I-C3H7 + H2O2' reacexp(155) = 'C3H8 + OH -> N-C3H7 + H2O' reacexp(156) = 'CH3 + C3H8 -> CH4 + I-C3H7' reacexp(157) = 'C3H8 + O2 -> I-C3H7 + HO2' reacexp(158) = 'I-C3H7 + HO2 -> C3H8 + O2' reacexp(159) = 'C3H8 + O -> N-C3H7 + OH' reacexp(160) = 'H + C3H8 -> H2 + N-C3H7' reacexp(161) = 'H + C3H8 -> H2 + I-C3H7' reacexp(162) = 'H2 + I-C3H7 -> H + C3H8' reacexp(163) = 'C3H5O -> C2H3 + CH2O' reacexp(164) = 'I-C3H7O2 -> I-C3H7 + O2' reacexp(165) = 'I-C3H7 + O2 -> I-C3H7O2' reacexp(166) = 'N-C3H7O2 -> N-C3H7 + O2' reacexp(167) = 'N-C3H7 + O2 -> N-C3H7O2' reacexp(168) = 'Reverse of 2 OH + M2 -> H2O2 + M2' reacexp(169) = 'Reverse of H + O2 + M4 -> HO2 + M4' reacexp(170) = 'Reverse of 2 CH3 + M6 -> C2H6 + M6' reacexp(171) = 'Reverse of 2 CH3 -> H + C2H5' reacexp(172) = 'Reverse of H + C2H4 + M12 -> C2H5 + M12' reacexp(173) = 'Reverse of H + C3H8 -> H2 + N-C3H7' return end subroutine reaction_expressions ! --- Forward/Backward link --- ! subroutine reverse_reactions use mechanism implicit none fofb = 0 ! Attach corresponding forward reaction to each backward reaction fofb(4) = 3 fofb(8) = 7 fofb(10) = 9 fofb(13) = 12 fofb(15) = 14 fofb(17) = 16 fofb(23) = 22 fofb(27) = 26 fofb(33) = 32 fofb(36) = 35 fofb(39) = 38 fofb(49) = 48 fofb(51) = 50 fofb(53) = 52 fofb(55) = 54 fofb(59) = 58 fofb(61) = 60 fofb(66) = 65 fofb(82) = 81 fofb(86) = 85 fofb(96) = 95 fofb(99) = 98 fofb(101) = 100 fofb(109) = 108 fofb(111) = 110 fofb(119) = 118 fofb(130) = 129 fofb(134) = 133 fofb(136) = 135 fofb(139) = 138 fofb(143) = 142 fofb(146) = 145 fofb(148) = 147 fofb(158) = 157 fofb(162) = 161 fofb(165) = 164 fofb(167) = 166 fofb(168) = 6 fofb(169) = 11 fofb(170) = 34 fofb(171) = 41 fofb(172) = 103 fofb(173) = 160 return end subroutine reverse_reactions ! --- Evaluation of QSS concentrations --- ! subroutine get_QSS(cqss,c,k,M) use mechanism implicit none real(WP), dimension(nqss) :: cqss real(WP), dimension(nspec) :: c real(WP), dimension(nreac) :: k real(WP), dimension(13) :: M real(WP) :: NXC3H7O2_denom1 & , NXC3H7O2_ct1 & , NXC3H7O2_denom2 & , NXC3H7O2_ct2 real(WP) :: HCCO_denom1 & , HCCO_ct1 & , HCCO_denom2 & , HCCO_ct2 & , HCCO_CH2GSGXCH2 & , HCCO_CH2GSGXCH2_coeff real(WP) :: CH2GSGXCH2_denom1 & , CH2GSGXCH2_ct1 & , CH2GSGXCH2_denom2 & , CH2GSGXCH2_ct2 real(WP) :: IXC3H7O2_denom1 & , IXC3H7O2_ct1 & , IXC3H7O2_denom2 & , IXC3H7O2_ct2 & , IXC3H7O2_IXC3H7 & , IXC3H7O2_IXC3H7_coeff real(WP) :: IXC3H7_denom1 & , IXC3H7_ct1 & , IXC3H7_denom2 & , IXC3H7_ct2 real(WP) :: C3H5XAXC3H5_denom1 & , C3H5XAXC3H5_ct1 & , C3H5XAXC3H5_denom2 & , C3H5XAXC3H5_ct2 real(WP) :: NXC3H7_denom1 & , NXC3H7_ct1 & , NXC3H7_denom2 & , NXC3H7_ct2 & , NXC3H7_NXC3H7O2 & , NXC3H7_NXC3H7O2_coeff real(WP) :: HCO_denom1 & , HCO_ct1 & , HCO_denom2 & , HCO_ct2 & , HCO_C3H5XAXC3H5 & , HCO_C3H5XAXC3H5_coeff ! c(sIXC3H7O2) c(sIXC3H7) (coupled) -------------------- ! Primary denominators----------------------- IXC3H7O2_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH129f) & ) IXC3H7_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH111) * c(sO2) & + k(rH112f) & + k(rH124b) * c(sHO2) & + k(rH127b) * c(sH2) & + k(rH129b) * c(sO2) & ) ! Primary constant parts ----------------------- IXC3H7O2_ct1 = ( 0.0_WP & ) IXC3H7_ct1 = ( 0.0_WP & + k(rH112b) * c(sH) * c(sC3H6) & + k(rH117) * c(sC3H8) * c(sO) & + k(rH118) * c(sC3H8) * c(sOH) & + k(rH121) * c(sC3H8) * c(sHO2) & + k(rH123) * c(sCH3) * c(sC3H8) & + k(rH124f) * c(sC3H8) * c(sO2) & + k(rH127f) * c(sH) * c(sC3H8) & ) ! IXC3H7O2 --------------------------------------- IXC3H7O2_denom2 = tiny(1.0_WP) + ( IXC3H7O2_denom1 & ) IXC3H7O2_ct2 = ( IXC3H7O2_ct1 & ) / IXC3H7O2_denom2 IXC3H7O2_IXC3H7 = ( 0.0_WP & + k(rH129b) * c(sO2) & ) / IXC3H7O2_denom2 IXC3H7O2_IXC3H7_coeff = ( 0.0_WP & + k(rH129f) & ) ! IXC3H7 --------------------------------------- IXC3H7_denom2 = tiny(1.0_WP) + ( IXC3H7_denom1 & - IXC3H7O2_IXC3H7_coeff * IXC3H7O2_IXC3H7 & ) IXC3H7_ct2 = ( IXC3H7_ct1 & + IXC3H7O2_IXC3H7_coeff * IXC3H7O2_ct2 & ) / IXC3H7_denom2 ! Reconstruction ------------------------------------ cqss(sIXC3H7-nspec) = ( IXC3H7_ct2 & ) cqss(sIXC3H7O2-nspec) = ( IXC3H7O2_ct2 & + IXC3H7O2_IXC3H7 * cqss(sIXC3H7-nspec) & ) ! c(sNXC3H7) c(sNXC3H7O2) (coupled) -------------------- ! Primary denominators----------------------- NXC3H7_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH113) * c(sO2) & + k(rH114f) & + k(rH115f) & + k(rH116) * c(sHO2) & + k(rH126b) * c(sH2) & + k(rH130b) * c(sO2) & ) NXC3H7O2_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH130f) & ) ! Primary constant parts ----------------------- NXC3H7_ct1 = ( 0.0_WP & + k(rH114b) * c(sH) * c(sC3H6) & + k(rH115b) * c(sCH3) * c(sC2H4) & + k(rH119) * c(sCH3) * c(sC3H8) & + k(rH120) * c(sC3H8) * c(sHO2) & + k(rH122) * c(sC3H8) * c(sOH) & + k(rH125) * c(sC3H8) * c(sO) & + k(rH126f) * c(sH) * c(sC3H8) & ) NXC3H7O2_ct1 = ( 0.0_WP & ) ! NXC3H7 --------------------------------------- NXC3H7_denom2 = tiny(1.0_WP) + ( NXC3H7_denom1 & ) NXC3H7_ct2 = ( NXC3H7_ct1 & ) / NXC3H7_denom2 NXC3H7_NXC3H7O2 = ( 0.0_WP & + k(rH130f) & ) / NXC3H7_denom2 NXC3H7_NXC3H7O2_coeff = ( 0.0_WP & + k(rH130b) * c(sO2) & ) ! NXC3H7O2 --------------------------------------- NXC3H7O2_denom2 = tiny(1.0_WP) + ( NXC3H7O2_denom1 & - NXC3H7_NXC3H7O2_coeff * NXC3H7_NXC3H7O2 & ) NXC3H7O2_ct2 = ( NXC3H7O2_ct1 & + NXC3H7_NXC3H7O2_coeff * NXC3H7_ct2 & ) / NXC3H7O2_denom2 ! Reconstruction ------------------------------------ cqss(sNXC3H7O2-nspec) = ( NXC3H7O2_ct2 & ) cqss(sNXC3H7-nspec) = ( NXC3H7_ct2 & + NXC3H7_NXC3H7O2 * cqss(sNXC3H7O2-nspec) & ) ! cqss(sC2H3) (uncoupled) -------------------- cqss(sC2H3-nspec) = ( 0.0_WP & + k(rH75f) * c(sC2H4) * c(sH) & + k(rH76) * c(sC2H4) * c(sCH3O2) & + k(rH77f) * c(sC2H4) * c(sOH) & + k(rH78f) * c(sC2H4) * c(sCH3) & + k(rH109f) * c(sC3H6) & + k(rH128) * c(sC3H5O) & ) / ( tiny(1.0_WP) + ( & + k(rH69) * c(sO2) & + k(rH70) * c(sO2) & + k(rH71) * c(sO2) & + k(rH72) * c(sH) & + k(rH73) & + k(rH75b) * c(sH2) & + k(rH77b) * c(sH2O) & + k(rH78b) * c(sCH4) & + k(rH109b) * c(sCH3) & ) ) ! c(sHCCO) c(sCH2GSGXCH2) (coupled) -------------------- ! Primary denominators----------------------- HCCO_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH91) * c(sOH) & + k(rH92) * c(sO) & + k(rH93f) * c(sH) & + k(rH94) * c(sO2) & ) CH2GSGXCH2_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH21) * c(sO) & + k(rH22) * c(sO2) & + k(rH23) * c(sOH) & + k(rH24f) * c(sH2) & + k(rH28b) * c(sH2O) & + k(rH35) * c(sCH3) & + k(rH39f) * c(sCH4) & + k(rH60) * c(sCO2) & + k(rH88) * c(sC2H6) & + k(rH93b) * c(sCO) & ) ! Primary constant parts ----------------------- HCCO_ct1 = ( 0.0_WP & + k(rH67) * c(sC2H2) * c(sO) & + k(rH68) * c(sC2H2) * c(sO2) & ) CH2GSGXCH2_ct1 = ( 0.0_WP & + k(rH24b) * c(sCH3) * c(sH) & + k(rH28f) * c(sCH3) * c(sOH) & + k(rH39b) * c(sCH3) * c(sCH3) & ) ! HCCO --------------------------------------- HCCO_denom2 = tiny(1.0_WP) + ( HCCO_denom1 & ) HCCO_ct2 = ( HCCO_ct1 & ) / HCCO_denom2 HCCO_CH2GSGXCH2 = ( 0.0_WP & + k(rH93b) * c(sCO) & ) / HCCO_denom2 HCCO_CH2GSGXCH2_coeff = ( 0.0_WP & + k(rH93f) * c(sH) & ) ! CH2GSGXCH2 --------------------------------------- CH2GSGXCH2_denom2 = tiny(1.0_WP) + ( CH2GSGXCH2_denom1 & - HCCO_CH2GSGXCH2_coeff * HCCO_CH2GSGXCH2 & ) CH2GSGXCH2_ct2 = ( CH2GSGXCH2_ct1 & + HCCO_CH2GSGXCH2_coeff * HCCO_ct2 & ) / CH2GSGXCH2_denom2 ! Reconstruction ------------------------------------ cqss(sCH2GSGXCH2-nspec) = ( CH2GSGXCH2_ct2 & ) cqss(sHCCO-nspec) = ( HCCO_ct2 & + HCCO_CH2GSGXCH2 * cqss(sCH2GSGXCH2-nspec) & ) ! c(sHCO) c(sC3H5XAXC3H5) (coupled) -------------------- ! Primary denominators----------------------- HCO_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH45) * c(sO) & + k(rH46) * c(sO) & + k(rH47) * c(sH) & + k(rH48f) * M(mM9) & + k(rH49) * c(sOH) & + k(rH50) * c(sCH3) & + k(rH51) * c(sO2) & + k(rH103b) * c(sC3H6) & ) C3H5XAXC3H5_denom1 = tiny(1.0_WP) + ( 0.0_WP & + k(rH96) * c(sH) & + k(rH97) * c(sHO2) & + k(rH98) & + k(rH99) * c(sCH3O2) & + k(rH100) * c(sO2) & + k(rH101) * c(sO2) & + k(rH102) * c(sHO2) & + k(rH103f) * c(sCH2O) & + k(rH106b) * c(sH2) & ) ! Primary constant parts ----------------------- HCO_ct1 = ( 0.0_WP & + k(rH48b) * c(sH) * c(sCO) * M(mM9) & + k(rH52) * c(sCH2O) * c(sHO2) & + k(rH53) * c(sCH2O) * c(sOH) & + k(rH54) * c(sCH2O) * c(sH) & + k(rH55) * c(sCH2O) * c(sCH3) & + k(rH56) * c(sCH2O) * c(sO) & + k(rH57) * c(sCH2O) * c(sO2) & + k(rH64) * c(sCH3O2) * c(sCH2O) & + k(rH71) * cqss(sC2H3-nspec) *c(sO2) & + k(rH79) * c(sC2H4) * c(sO) & + k(rH91) * 2.0 * cqss(sHCCO-nspec) *c(sOH) & + k(rH94) * cqss(sHCCO-nspec) *c(sO2) & + k(rH108) * c(sC3H6) * c(sO) & ) C3H5XAXC3H5_ct1 = ( 0.0_WP & + k(rH104) * c(sC3H6) * c(sOH) & + k(rH105) * c(sC3H6) * c(sO) & + k(rH106f) * c(sC3H6) * c(sH) & + k(rH110) * c(sC3H6) * c(sCH3) & ) ! HCO --------------------------------------- HCO_denom2 = tiny(1.0_WP) + ( HCO_denom1 & ) HCO_ct2 = ( HCO_ct1 & ) / HCO_denom2 HCO_C3H5XAXC3H5 = ( 0.0_WP & + k(rH103f) * c(sCH2O) & ) / HCO_denom2 HCO_C3H5XAXC3H5_coeff = ( 0.0_WP & + k(rH103b) * c(sC3H6) & ) ! C3H5XAXC3H5 --------------------------------------- C3H5XAXC3H5_denom2 = tiny(1.0_WP) + ( C3H5XAXC3H5_denom1 & - HCO_C3H5XAXC3H5_coeff * HCO_C3H5XAXC3H5 & ) C3H5XAXC3H5_ct2 = ( C3H5XAXC3H5_ct1 & + HCO_C3H5XAXC3H5_coeff * HCO_ct2 & ) / C3H5XAXC3H5_denom2 ! Reconstruction ------------------------------------ cqss(sC3H5XAXC3H5-nspec) = ( C3H5XAXC3H5_ct2 & ) cqss(sHCO-nspec) = ( HCO_ct2 & + HCO_C3H5XAXC3H5 * cqss(sC3H5XAXC3H5-nspec) & ) ! cqss(sCH3O) (uncoupled) -------------------- cqss(sCH3O-nspec) = ( 0.0_WP & + k(rH26f) * c(sCH3) * c(sO2) & + k(rH27) * c(sCH3) * c(sHO2) & + k(rH62) * 2.0 * c(sCH3O2) * c(sCH3O2) & + k(rH65) * 2.0 * c(sCH3O2) * c(sCH3) & + k(rH66f) * c(sCH3O2H) & + k(rH99) * cqss(sC3H5XAXC3H5-nspec) *c(sCH3O2) & ) / ( tiny(1.0_WP) + ( & + k(rH26b) * c(sO) & + k(rH58) * c(sO2) & + k(rH59) & + k(rH66b) * c(sOH) & ) ) ! cqss(sCH2CHO) (uncoupled) -------------------- cqss(sCH2CHO-nspec) = ( 0.0_WP & + k(rH70) * cqss(sC2H3-nspec) *c(sO2) & + k(rH81) * c(sC2H4) * c(sO) & + k(rH101) * cqss(sC3H5XAXC3H5-nspec) *c(sO2) & ) / ( tiny(1.0_WP) + ( & + k(rH95) * c(sO2) & ) ) ! cqss(sC2H5) (uncoupled) -------------------- cqss(sC2H5-nspec) = ( 0.0_WP & + k(rH30f) * c(sCH3) * c(sCH3) & + k(rH80f) * c(sH) * c(sC2H4) & + k(rH85b) * c(sC2H4) * c(sHO2) & + k(rH86f) * c(sC2H6) * c(sH) & + k(rH87) * c(sC2H6) * c(sOH) & + k(rH88) * cqss(sCH2GSGXCH2-nspec) *c(sC2H6) & + k(rH89) * c(sC2H6) * c(sCH3) & + k(rH90) * c(sC2H6) * c(sO) & + k(rH108) * c(sC3H6) * c(sO) & ) / ( tiny(1.0_WP) + ( & + k(rH30b) * c(sH) & + k(rH80b) & + k(rH82) * c(sHO2) & + k(rH83) * c(sCH3) & + k(rH84) * c(sH) & + k(rH85f) * c(sO2) & + k(rH86b) * c(sH2) & ) ) return end subroutine get_QSS