[57] | 1 | c***************************************************************** |
---|
| 2 | c |
---|
| 3 | c Photochemical routines for Scheme B |
---|
| 4 | c |
---|
| 5 | c Author: Franck Lefevre |
---|
| 6 | c ------ |
---|
| 7 | c |
---|
| 8 | c update 12/06/2003: compatibility with dust and iceparty |
---|
| 9 | c (S. Lebonnois) |
---|
| 10 | c |
---|
| 11 | c***************************************************************** |
---|
| 12 | c |
---|
| 13 | subroutine photochemist_B(lswitch,zycol, sza, ptimestep, press, |
---|
| 14 | c $ temp, dens, dist_sol) |
---|
| 15 | $ temp, dens, dist_sol, jo3) |
---|
| 16 | c |
---|
| 17 | c |
---|
| 18 | implicit none |
---|
| 19 | c |
---|
| 20 | #include "dimensions.h" |
---|
| 21 | #include "dimphys.h" |
---|
| 22 | #include "chimiedata.h" |
---|
| 23 | #include "callkeys.h" |
---|
| 24 | c |
---|
| 25 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 26 | c input/output: |
---|
| 27 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 28 | c |
---|
| 29 | real zycol(nlayermx,nqmx) ! chemical species volume mixing ratio |
---|
| 30 | c if iceparty, zycol(l,nqmx-1) is ice surface area (in microns^2/cm^3) |
---|
| 31 | c |
---|
| 32 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 33 | c inputs: |
---|
| 34 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 35 | c |
---|
| 36 | real sza ! solar zenith angle (deg) |
---|
| 37 | real ptimestep ! physics timestep (s) |
---|
| 38 | real press(nlayermx) ! pressure (hpa) |
---|
| 39 | real temp(nlayermx) ! temperature (k) |
---|
| 40 | real dens(nlayermx) ! density (cm-3) |
---|
| 41 | real dist_sol ! sun distance (au) |
---|
| 42 | c |
---|
| 43 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 44 | c output: |
---|
| 45 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 46 | c |
---|
| 47 | real jo3(nlayermx) ! photodissociation rate O3 -> O1D |
---|
| 48 | c |
---|
| 49 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 50 | c local: |
---|
| 51 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 52 | c |
---|
| 53 | integer phychemrat ! ratio physics timestep/chemistry timestep |
---|
| 54 | integer istep |
---|
| 55 | integer i_o3,lev |
---|
| 56 | integer nesp |
---|
| 57 | integer lswitch |
---|
| 58 | c |
---|
| 59 | parameter (nesp = 15) ! number of species in the chemistry |
---|
| 60 | c |
---|
| 61 | real ctimestep ! chemistry timestep (s) |
---|
| 62 | real rm(nlayermx,nesp) ! species volume mixing ratio |
---|
| 63 | real j(nlayermx,nd) ! interpolated photolysis rates (s-1) |
---|
| 64 | real icesurf(nlayermx) ! ice surface area (cm^2/cm^3) |
---|
| 65 | real rmo3(nlayermx) ! ozone mixing ratio |
---|
| 66 | c |
---|
| 67 | c reaction rates |
---|
| 68 | c |
---|
| 69 | real a001(nlayermx), a002(nlayermx), a003(nlayermx) |
---|
| 70 | real b001(nlayermx), b002(nlayermx), b003(nlayermx), |
---|
| 71 | $ b004(nlayermx), b005(nlayermx), b006(nlayermx) |
---|
| 72 | real c001(nlayermx), c002(nlayermx), c003(nlayermx), |
---|
| 73 | $ c004(nlayermx), c005(nlayermx), c006(nlayermx), |
---|
| 74 | $ c007(nlayermx), c008(nlayermx), c009(nlayermx), |
---|
| 75 | $ c010(nlayermx), c011(nlayermx), c012(nlayermx), |
---|
| 76 | $ c013(nlayermx), c014(nlayermx), c015(nlayermx), |
---|
| 77 | $ c016(nlayermx), c017(nlayermx), c018(nlayermx) |
---|
| 78 | real d001(nlayermx), d002(nlayermx), d003(nlayermx) |
---|
| 79 | real e001(nlayermx), e002(nlayermx) |
---|
| 80 | real h001(nlayermx), h002(nlayermx) |
---|
| 81 | c |
---|
| 82 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 83 | c ctimestep : chemistry timestep (s) c |
---|
| 84 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 85 | c |
---|
| 86 | ctimestep = 10.*60. |
---|
| 87 | c |
---|
| 88 | phychemrat = int(ptimestep/ctimestep) |
---|
| 89 | c print*,"PHYCHEMRAT = ",phychemrat |
---|
| 90 | c |
---|
| 91 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 92 | c initialisation of chemical species and families c |
---|
| 93 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 94 | c |
---|
| 95 | call gcmtochim(zycol, lswitch, nesp, rm) |
---|
| 96 | c |
---|
| 97 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 98 | c ice surface area for heterogenous chemistry |
---|
| 99 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 100 | c |
---|
| 101 | if (iceparty) then |
---|
| 102 | do lev=1,lswitch-1 |
---|
| 103 | icesurf(lev)= zycol(lev,nqmx-1)*1.e-8 ! microns^2 -> cm^2 |
---|
| 104 | enddo |
---|
| 105 | else |
---|
| 106 | do lev=1,lswitch-1 |
---|
| 107 | icesurf(lev)=0.0e0 |
---|
| 108 | enddo |
---|
| 109 | endif |
---|
| 110 | c |
---|
| 111 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 112 | c compute reaction rates c |
---|
| 113 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 114 | c |
---|
| 115 | call chemrates(lswitch, dens, press, temp, icesurf, |
---|
| 116 | $ a001, a002, a003, |
---|
| 117 | $ b001, b002, b003, b004, b005, b006, |
---|
| 118 | $ c001, c002, c003, c004, c005, c006, |
---|
| 119 | $ c007, c008, c009, c010, c011, c012, |
---|
| 120 | $ c013, c014, c015, c016, c017, c018, |
---|
| 121 | $ d001, d002, d003, |
---|
| 122 | $ e001, e002, |
---|
| 123 | $ h001, h002) |
---|
| 124 | c |
---|
| 125 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 126 | c interpolation of photolysis rates in the lookup table c |
---|
| 127 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 128 | c |
---|
| 129 | i_o3 = 5 |
---|
| 130 | c |
---|
| 131 | do lev=1,lswitch-1 |
---|
| 132 | rmo3(lev)=rm(lev,i_o3) |
---|
| 133 | enddo |
---|
| 134 | |
---|
| 135 | call phot(lswitch, press, temp, sza, dist_sol, rmo3, j) |
---|
| 136 | c |
---|
| 137 | do lev=1,lswitch-1 |
---|
| 138 | jo3(lev) = j(lev,5) |
---|
| 139 | enddo |
---|
| 140 | |
---|
| 141 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 142 | c loop over time within the physical timestep c |
---|
| 143 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 144 | c |
---|
| 145 | do istep = 1,phychemrat |
---|
| 146 | call chimie_B(lswitch,nesp, rm, j, dens, ctimestep, |
---|
| 147 | $ press, temp, sza, dist_sol, |
---|
| 148 | $ a001, a002, a003, |
---|
| 149 | $ b001, b002, b003, b004, b005, b006, |
---|
| 150 | $ c001, c002, c003, c004, c005, c006, |
---|
| 151 | $ c007, c008, c009, c010, c011, c012, |
---|
| 152 | $ c013, c014, c015, c016, c017, c018, |
---|
| 153 | $ d001, d002, d003, |
---|
| 154 | $ e001, e002, |
---|
| 155 | $ h001, h002) |
---|
| 156 | c |
---|
| 157 | c print*,'appel chimie ok iteration ', istep |
---|
| 158 | c |
---|
| 159 | end do |
---|
| 160 | c |
---|
| 161 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 162 | c save chemical species for the gcm c |
---|
| 163 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 164 | c |
---|
| 165 | call chimtogcm(zycol, lswitch, nesp, rm) |
---|
| 166 | c |
---|
| 167 | return |
---|
| 168 | end |
---|
| 169 | c |
---|
| 170 | c***************************************************************** |
---|
| 171 | c |
---|
| 172 | subroutine chimie_B(lswitch, nesp, rm, j, dens, dt, |
---|
| 173 | $ press, t, sza, dist_sol, |
---|
| 174 | $ a001, a002, a003, |
---|
| 175 | $ b001, b002, b003, b004, b005, b006, |
---|
| 176 | $ c001, c002, c003, c004, c005, c006, |
---|
| 177 | $ c007, c008, c009, c010, c011, c012, |
---|
| 178 | $ c013, c014, c015, c016, c017, c018, |
---|
| 179 | $ d001, d002, d003, |
---|
| 180 | $ e001, e002, |
---|
| 181 | $ h001, h002) |
---|
| 182 | c |
---|
| 183 | c***************************************************************** |
---|
| 184 | c |
---|
| 185 | implicit none |
---|
| 186 | c |
---|
| 187 | #include "dimensions.h" |
---|
| 188 | #include "dimphys.h" |
---|
| 189 | #include "chimiedata.h" |
---|
| 190 | #include "callkeys.h" |
---|
| 191 | c |
---|
| 192 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 193 | c input/output: |
---|
| 194 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 195 | c |
---|
| 196 | integer lswitch ! interface level between chemistries |
---|
| 197 | integer nesp ! number of species |
---|
| 198 | c |
---|
| 199 | real rm(nlayermx,nesp) ! volume mixing ratios |
---|
| 200 | c |
---|
| 201 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 202 | c inputs: |
---|
| 203 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 204 | c |
---|
| 205 | real dens(nlayermx) ! density (cm-3) |
---|
| 206 | real dt ! chemistry timestep (s) |
---|
| 207 | real j(nlayermx,nd) ! interpolated photolysis rates (s-1) |
---|
| 208 | real press(nlayermx) ! pressure (hpa) |
---|
| 209 | real t(nlayermx) ! temperature (k) |
---|
| 210 | real sza ! solar zenith angle (deg) |
---|
| 211 | real dist_sol ! sun distance (au) |
---|
| 212 | c |
---|
| 213 | c reaction rates |
---|
| 214 | c |
---|
| 215 | real a001(nlayermx), a002(nlayermx), a003(nlayermx) |
---|
| 216 | real b001(nlayermx), b002(nlayermx), b003(nlayermx), |
---|
| 217 | $ b004(nlayermx), b005(nlayermx), b006(nlayermx) |
---|
| 218 | real c001(nlayermx), c002(nlayermx), c003(nlayermx), |
---|
| 219 | $ c004(nlayermx), c005(nlayermx), c006(nlayermx), |
---|
| 220 | $ c007(nlayermx), c008(nlayermx), c009(nlayermx), |
---|
| 221 | $ c010(nlayermx), c011(nlayermx), c012(nlayermx), |
---|
| 222 | $ c013(nlayermx), c014(nlayermx), c015(nlayermx), |
---|
| 223 | $ c016(nlayermx), c017(nlayermx), c018(nlayermx) |
---|
| 224 | real d001(nlayermx), d002(nlayermx), d003(nlayermx) |
---|
| 225 | real e001(nlayermx), e002(nlayermx) |
---|
| 226 | real h001(nlayermx), h002(nlayermx) |
---|
| 227 | c |
---|
| 228 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 229 | c local: |
---|
| 230 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 231 | c |
---|
| 232 | integer HETERO |
---|
| 233 | parameter (HETERO=0) |
---|
| 234 | |
---|
| 235 | integer iesp, iter, l, niter |
---|
| 236 | integer i_co2, i_co, i_o2, i_h2, i_h2o, i_h2o2, i_hox, i_ox, |
---|
| 237 | $ i_o1d, i_o, i_o3, i_h, i_oh, i_ho2, i_n2 |
---|
| 238 | integer j_o2_o, j_o2_o1d, j_co2_o, j_co2_o1d, |
---|
| 239 | $ j_o3_o1d, j_o3_o, j_h2o, j_hdo, j_h2o2, |
---|
| 240 | $ j_ho2, j_no2 |
---|
| 241 | c |
---|
| 242 | parameter (niter = 5) ! iterations in the chemical scheme |
---|
| 243 | c |
---|
| 244 | real*8 cc0(nlayermx,nesp) ! initial number density (cm-3) |
---|
| 245 | real*8 cc(nlayermx,nesp) ! final number density (cm-3) |
---|
| 246 | real*8 co2(nlayermx) ! co2 number density (cm-3) |
---|
| 247 | real*8 nox(nlayermx) ! nox number density (cm-3) |
---|
| 248 | real*8 no(nlayermx) ! no number density (cm-3) |
---|
| 249 | real*8 no2(nlayermx) ! no2 number density (cm-3) |
---|
| 250 | real*8 production(nlayermx,nesp) ! production rate |
---|
| 251 | real*8 loss(nlayermx,nesp) ! loss rate |
---|
| 252 | c |
---|
| 253 | real c007h(nlayermx) ! c007 rate modified |
---|
| 254 | c ! by heterogenous reactions h001 and h002 |
---|
| 255 | c |
---|
| 256 | real ro_o3, rh_ho2, roh_ho2, rno2_no |
---|
| 257 | c |
---|
| 258 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 259 | c tracer numbering in the chemistry |
---|
| 260 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 261 | c |
---|
| 262 | i_co2 = 1 |
---|
| 263 | i_co = 2 |
---|
| 264 | i_o = 3 |
---|
| 265 | i_o1d = 4 |
---|
| 266 | i_o2 = 5 |
---|
| 267 | i_o3 = 6 |
---|
| 268 | i_h = 7 |
---|
| 269 | i_h2 = 8 |
---|
| 270 | i_oh = 9 |
---|
| 271 | i_ho2 = 10 |
---|
| 272 | i_h2o2 = 11 |
---|
| 273 | i_h2o = 12 |
---|
| 274 | i_n2 = 13 |
---|
| 275 | i_hox = 14 |
---|
| 276 | i_ox = 15 |
---|
| 277 | c |
---|
| 278 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 279 | c numbering of photolysis rates |
---|
| 280 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 281 | c |
---|
| 282 | j_o2_o = 1 ! o2 + hv -> o + o |
---|
| 283 | j_o2_o1d = 2 ! o2 + hv -> o + o(1d) |
---|
| 284 | j_co2_o = 3 ! co2 + hv -> co + o |
---|
| 285 | j_co2_o1d = 4 ! co2 + hv -> co + o(1d) |
---|
| 286 | j_o3_o1d = 5 ! o3 + hv -> o2 + o(1d) |
---|
| 287 | j_o3_o = 6 ! o3 + hv -> o2 + o |
---|
| 288 | j_h2o = 7 ! h2o + hv -> h + oh |
---|
| 289 | j_hdo = 8 ! hdo + hv -> d + oh |
---|
| 290 | j_h2o2 = 9 ! h2o2 + hv -> oh + oh |
---|
| 291 | j_ho2 = 10 ! ho2 + hv -> oh + o |
---|
| 292 | j_no2 = 11 ! no2 + hv -> no + o |
---|
| 293 | c |
---|
| 294 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 295 | c volume mixing ratio -> density conversion |
---|
| 296 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 297 | c |
---|
| 298 | do iesp = 1,nesp |
---|
| 299 | do l = 1,lswitch-1 |
---|
| 300 | cc0(l,iesp) = rm(l,iesp)*dens(l) |
---|
| 301 | cc(l,iesp) = cc0(l,iesp) |
---|
| 302 | end do |
---|
| 303 | end do |
---|
| 304 | c |
---|
| 305 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 306 | c co2 and nox number densities (cm-3) |
---|
| 307 | c |
---|
| 308 | c co2 mixing ratio: 0.953 |
---|
| 309 | c nox mixing ratio: 6.e-10 (yung and demore, 1999) |
---|
| 310 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 311 | c |
---|
| 312 | do l = 1,lswitch-1 |
---|
| 313 | co2(l) = cc(l,i_co2) |
---|
| 314 | nox(l) = 6.e-10*dens(l) |
---|
| 315 | end do |
---|
| 316 | c |
---|
| 317 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 318 | c loop over iterations in the chemical scheme |
---|
| 319 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 320 | c |
---|
| 321 | do iter = 1,niter |
---|
| 322 | c |
---|
| 323 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 324 | c modification of c007 |
---|
| 325 | c by heterogenous reactions h001 and h002 |
---|
| 326 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 327 | c |
---|
| 328 | if (HETERO.eq.1) then |
---|
| 329 | do l = 1,lswitch-1 |
---|
| 330 | c007h(l) = c007(l) |
---|
| 331 | $ + h001(l)/max(cc(l,i_oh), 1.e-30*dens(l)) |
---|
| 332 | $ + h002(l)/max(cc(l,i_ho2),1.e-30*dens(l)) |
---|
| 333 | c write(*,*) "C007 = ",c007(l)," / C007H = ",c007h(l) |
---|
| 334 | end do |
---|
| 335 | else |
---|
| 336 | do l = 1,lswitch-1 |
---|
| 337 | c007h(l) = c007(l) |
---|
| 338 | end do |
---|
| 339 | endif |
---|
| 340 | c |
---|
| 341 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 342 | c nox species |
---|
| 343 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 344 | c no2/no |
---|
| 345 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 346 | c |
---|
| 347 | do l = 1,lswitch-1 |
---|
| 348 | c |
---|
| 349 | rno2_no = (d002(l)*cc(l,i_o3) + d003(l)*cc(l,i_ho2)) |
---|
| 350 | $ /(j(l,j_no2) + |
---|
| 351 | $ d001(l)*max(cc(l,i_o),1.e-30*dens(l))) |
---|
| 352 | c |
---|
| 353 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 354 | c no |
---|
| 355 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 356 | c |
---|
| 357 | no(l) = nox(l)/(1. + rno2_no) |
---|
| 358 | c |
---|
| 359 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 360 | c no2 |
---|
| 361 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 362 | c |
---|
| 363 | no2(l) = no(l)*rno2_no |
---|
| 364 | c |
---|
| 365 | end do |
---|
| 366 | c |
---|
| 367 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 368 | c hox species |
---|
| 369 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 370 | c photochemical equilibrium for oh and ho2 |
---|
| 371 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 372 | c h/ho2 |
---|
| 373 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 374 | c |
---|
| 375 | do l = 1,lswitch-1 |
---|
| 376 | c |
---|
| 377 | rh_ho2 = (c001(l)*cc(l,i_o) |
---|
| 378 | $ + c004(l)*cc(l,i_h) |
---|
| 379 | $ + c005(l)*cc(l,i_h) |
---|
| 380 | $ + c006(l)*cc(l,i_h) |
---|
| 381 | $ + c007h(l)*cc(l,i_oh) |
---|
| 382 | $ + 2.*c008(l)*cc(l,i_ho2) |
---|
| 383 | $ + c015(l)*cc(l,i_o3) |
---|
| 384 | $ + 2.*c016(l)*cc(l,i_ho2) |
---|
| 385 | $ + j(l,j_ho2)) |
---|
| 386 | $ /(c011(l)*cc(l,i_o2)) |
---|
| 387 | c |
---|
| 388 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 389 | c oh/ho2 |
---|
| 390 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 391 | c |
---|
| 392 | roh_ho2 = (c001(l)*cc(l,i_o) |
---|
| 393 | $ + c003(l)*cc(l,i_o3)*rh_ho2 |
---|
| 394 | $ + 2.*c004(l)*cc(l,i_h) |
---|
| 395 | $ + 2.*c008(l)*cc(l,i_ho2) |
---|
| 396 | $ + c015(l)*cc(l,i_o3) |
---|
| 397 | $ + d003(l)*no(l) |
---|
| 398 | $ + j(l,j_ho2) |
---|
| 399 | $ + j(l,j_h2o)*cc(l,i_h2o) |
---|
| 400 | $ /max(cc(l,i_ho2),dens(l)*1.e-30)) |
---|
| 401 | $ /(c002(l)*cc(l,i_o) |
---|
| 402 | $ + c007h(l)*cc(l,i_ho2) |
---|
| 403 | $ + e001(l)*cc(l,i_co)) |
---|
| 404 | c |
---|
| 405 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 406 | c h |
---|
| 407 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 408 | c |
---|
| 409 | cc(l,i_h) = cc(l,i_hox) |
---|
| 410 | $ /(1. + (1. + roh_ho2)/rh_ho2) |
---|
| 411 | c |
---|
| 412 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 413 | c ho2 |
---|
| 414 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 415 | c |
---|
| 416 | cc(l,i_ho2) = cc(l,i_h)/rh_ho2 |
---|
| 417 | c |
---|
| 418 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 419 | c oh |
---|
| 420 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 421 | c |
---|
| 422 | cc(l,i_oh) = cc(l,i_ho2)*roh_ho2 |
---|
| 423 | c |
---|
| 424 | end do |
---|
| 425 | c |
---|
| 426 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 427 | c modification of c007 |
---|
| 428 | c by heterogenous reactions h001 and h002 |
---|
| 429 | c after update of cc(l,i_oh) and cc(l,i_ho2) |
---|
| 430 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 431 | c |
---|
| 432 | if (HETERO.eq.1) then |
---|
| 433 | do l = 1,lswitch-1 |
---|
| 434 | c007h(l) = c007(l) |
---|
| 435 | $ + h001(l)/max(cc(l,i_oh), 1.e-30*dens(l)) |
---|
| 436 | $ + h002(l)/max(cc(l,i_ho2),1.e-30*dens(l)) |
---|
| 437 | c write(*,*) "C007 = ",c007(l)," / C007H = ",c007h(l) |
---|
| 438 | end do |
---|
| 439 | endif |
---|
| 440 | c |
---|
| 441 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 442 | c ox species |
---|
| 443 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 444 | c day: |
---|
| 445 | c - o1d at photochemical equilibrium |
---|
| 446 | c - o3 at photochemical equilibrium |
---|
| 447 | c - continuity equation for ox |
---|
| 448 | c night: |
---|
| 449 | c - o1d = 0 |
---|
| 450 | c - continuity equation for o3 |
---|
| 451 | c - continuity equation for o |
---|
| 452 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 453 | c |
---|
| 454 | if (sza .le. 95.) then |
---|
| 455 | c |
---|
| 456 | do l = 1,lswitch-1 |
---|
| 457 | c |
---|
| 458 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 459 | c o(1d) |
---|
| 460 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 461 | c |
---|
| 462 | cc(l,i_o1d) = (j(l,j_co2_o1d)*co2(l) |
---|
| 463 | $ + j(l,j_o2_o1d)*cc(l,i_o2) |
---|
| 464 | $ + j(l,j_o3_o1d)*cc(l,i_o3)) |
---|
| 465 | $ /( b001(l)*co2(l) |
---|
| 466 | $ + b002(l)*cc(l,i_h2o) |
---|
| 467 | $ + b003(l)*cc(l,i_h2) |
---|
| 468 | $ + b004(l)*cc(l,i_o2) |
---|
| 469 | $ + b005(l)*cc(l,i_o3) |
---|
| 470 | $ + b006(l)*cc(l,i_o3)) |
---|
| 471 | c |
---|
| 472 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 473 | c o/o3 |
---|
| 474 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 475 | c |
---|
| 476 | ro_o3 = (j(l,j_o3_o1d) + j(l,j_o3_o) |
---|
| 477 | $ + a003(l)*cc(l,i_o) |
---|
| 478 | $ + c003(l)*cc(l,i_h) |
---|
| 479 | $ + c014(l)*cc(l,i_oh) |
---|
| 480 | $ + c015(l)*cc(l,i_ho2)) |
---|
| 481 | $ /(a001(l)*cc(l,i_o2)) |
---|
| 482 | c |
---|
| 483 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 484 | c o3 |
---|
| 485 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 486 | c |
---|
| 487 | cc(l,i_o3) = cc(l,i_ox)/(1. + ro_o3) |
---|
| 488 | c |
---|
| 489 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 490 | c o |
---|
| 491 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 492 | c |
---|
| 493 | cc(l,i_o) = cc(l,i_o3)*ro_o3 |
---|
| 494 | c |
---|
| 495 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 496 | c ox = o + o3 |
---|
| 497 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 498 | c |
---|
| 499 | production(l,i_ox) = |
---|
| 500 | $ + j(l,j_co2_o)*co2(l) |
---|
| 501 | $ + j(l,j_co2_o1d)*co2(l) |
---|
| 502 | $ + j(l,j_ho2)*cc(l,i_ho2) |
---|
| 503 | $ + 2.*j(l,j_o2_o)*cc(l,i_o2) |
---|
| 504 | $ + 2.*j(l,j_o2_o1d)*cc(l,i_o2) |
---|
| 505 | $ + c006(l)*cc(l,i_h)*cc(l,i_ho2) |
---|
| 506 | $ + c013(l)*cc(l,i_oh)*cc(l,i_oh) |
---|
| 507 | $ + d003(l)*cc(l,i_ho2)*no(l) |
---|
| 508 | c |
---|
| 509 | loss(l,i_ox) = 2.*a002(l)*cc(l,i_o)*cc(l,i_o) |
---|
| 510 | $ + 2.*a003(l)*cc(l,i_o)*cc(l,i_o3) |
---|
| 511 | $ + c001(l)*cc(l,i_ho2)*cc(l,i_o) |
---|
| 512 | $ + c002(l)*cc(l,i_oh)*cc(l,i_o) |
---|
| 513 | $ + c003(l)*cc(l,i_h)*cc(l,i_o3) |
---|
| 514 | $ + c012(l)*cc(l,i_o)*cc(l,i_h2o2) |
---|
| 515 | $ + c014(l)*cc(l,i_o3)*cc(l,i_oh) |
---|
| 516 | $ + c015(l)*cc(l,i_o3)*cc(l,i_ho2) |
---|
| 517 | $ + d001(l)*cc(l,i_o)*no2(l) |
---|
| 518 | $ + e002(l)*cc(l,i_o)*cc(l,i_co) |
---|
| 519 | c |
---|
| 520 | loss(l,i_ox) = loss(l,i_ox)/cc(l,i_ox) |
---|
| 521 | c |
---|
| 522 | end do |
---|
| 523 | c |
---|
| 524 | else |
---|
| 525 | c |
---|
| 526 | do l = 1,lswitch-1 |
---|
| 527 | c |
---|
| 528 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 529 | c o(1d) |
---|
| 530 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 531 | c |
---|
| 532 | cc(l,i_o1d) = 0. |
---|
| 533 | c |
---|
| 534 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 535 | c o3 |
---|
| 536 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 537 | c |
---|
| 538 | production(l,i_o3) = a001(l)*cc(l,i_o2)*cc(l,i_o) |
---|
| 539 | c |
---|
| 540 | loss(l,i_o3) = a003(l)*cc(l,i_o) |
---|
| 541 | $ + c003(l)*cc(l,i_h) |
---|
| 542 | $ + c014(l)*cc(l,i_oh) |
---|
| 543 | $ + c015(l)*cc(l,i_ho2) |
---|
| 544 | c |
---|
| 545 | |
---|
| 546 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 547 | c o |
---|
| 548 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 549 | c |
---|
| 550 | production(l,i_o) = c006(l)*cc(l,i_h)*cc(l,i_ho2) |
---|
| 551 | $ + c013(l)*cc(l,i_oh)*cc(l,i_oh) |
---|
| 552 | c |
---|
| 553 | loss(l,i_o) = a001(l)*cc(l,i_o2) |
---|
| 554 | $ + 2.*a002(l)*cc(l,i_o) |
---|
| 555 | $ + a003(l)*cc(l,i_o3) |
---|
| 556 | $ + c001(l)*cc(l,i_ho2) |
---|
| 557 | $ + c002(l)*cc(l,i_oh) |
---|
| 558 | $ + c012(l)*cc(l,i_h2o2) |
---|
| 559 | $ + e002(l)*cc(l,i_co) |
---|
| 560 | c |
---|
| 561 | end do |
---|
| 562 | end if |
---|
| 563 | c |
---|
| 564 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 565 | c other species |
---|
| 566 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 567 | c |
---|
| 568 | do l = 1,lswitch-1 |
---|
| 569 | c |
---|
| 570 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 571 | c co |
---|
| 572 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 573 | c |
---|
| 574 | production(l,i_co) = j(l,j_co2_o)*co2(l) |
---|
| 575 | $ + j(l,j_co2_o1d)*co2(l) |
---|
| 576 | c |
---|
| 577 | loss(l,i_co) = e001(l)*cc(l,i_oh) |
---|
| 578 | $ + e002(l)*cc(l,i_o) |
---|
| 579 | c |
---|
| 580 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 581 | c o2 |
---|
| 582 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 583 | c |
---|
| 584 | production(l,i_o2) = |
---|
| 585 | $ j(l,j_o3_o)*cc(l,i_o3) |
---|
| 586 | $ + j(l,j_o3_o1d)*cc(l,i_o3) |
---|
| 587 | $ + a002(l)*cc(l,i_o)*cc(l,i_o) |
---|
| 588 | $ + 2.*a003(l)*cc(l,i_o)*cc(l,i_o3) |
---|
| 589 | $ + 2.*b005(l)*cc(l,i_o1d)*cc(l,i_o3) |
---|
| 590 | $ + b006(l)*cc(l,i_o1d)*cc(l,i_o3) |
---|
| 591 | $ + c001(l)*cc(l,i_o)*cc(l,i_ho2) |
---|
| 592 | $ + c002(l)*cc(l,i_o)*cc(l,i_oh) |
---|
| 593 | $ + c003(l)*cc(l,i_h)*cc(l,i_o3) |
---|
| 594 | $ + c005(l)*cc(l,i_h)*cc(l,i_ho2) |
---|
| 595 | $ + c007h(l)*cc(l,i_oh)*cc(l,i_ho2) |
---|
| 596 | $ + c008(l)*cc(l,i_ho2)*cc(l,i_ho2) |
---|
| 597 | $ + c014(l)*cc(l,i_o3)*cc(l,i_oh) |
---|
| 598 | $ + 2.*c015(l)*cc(l,i_o3)*cc(l,i_ho2) |
---|
| 599 | $ + c016(l)*cc(l,i_ho2)*cc(l,i_ho2) |
---|
| 600 | $ + d001(l)*cc(l,i_o)*no2(l) |
---|
| 601 | c |
---|
| 602 | loss(l,i_o2) = j(l,j_o2_o) |
---|
| 603 | $ + j(l,j_o2_o1d) |
---|
| 604 | $ + a001(l)*cc(l,i_o) |
---|
| 605 | $ + c011(l)*cc(l,i_h) |
---|
| 606 | c |
---|
| 607 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 608 | c h2 |
---|
| 609 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 610 | c |
---|
| 611 | production(l,i_h2) = c005(l)*cc(l,i_h)*cc(l,i_ho2) |
---|
| 612 | $ + c018(l)*cc(l,i_h)*cc(l,i_h) |
---|
| 613 | c |
---|
| 614 | loss(l,i_h2) = b003(l)*cc(l,i_o1d) |
---|
| 615 | $ + c010(l)*cc(l,i_oh) |
---|
| 616 | c |
---|
| 617 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 618 | c h2o |
---|
| 619 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 620 | c |
---|
| 621 | production(l,i_h2o) = |
---|
| 622 | $ c006(l)*cc(l,i_h)*cc(l,i_ho2) |
---|
| 623 | $ + c007h(l)*cc(l,i_oh)*cc(l,i_ho2) |
---|
| 624 | $ + c009(l)*cc(l,i_oh)*cc(l,i_h2o2) |
---|
| 625 | $ + c010(l)*cc(l,i_oh)*cc(l,i_h2) |
---|
| 626 | $ + c013(l)*cc(l,i_oh)*cc(l,i_oh) |
---|
| 627 | c |
---|
| 628 | loss(l,i_h2o) = j(l,j_h2o) |
---|
| 629 | $ + b002(l)*cc(l,i_o1d) |
---|
| 630 | c |
---|
| 631 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 632 | c h2o2 |
---|
| 633 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 634 | c |
---|
| 635 | production(l,i_h2o2) = |
---|
| 636 | $ c008(l)*cc(l,i_ho2)*cc(l,i_ho2) |
---|
| 637 | $ + c016(l)*cc(l,i_ho2)*cc(l,i_ho2) |
---|
| 638 | $ + c017(l)*cc(l,i_oh)*cc(l,i_oh) |
---|
| 639 | c |
---|
| 640 | loss(l,i_h2o2) = j(l,j_h2o2) |
---|
| 641 | $ + c009(l)*cc(l,i_oh) |
---|
| 642 | $ + c012(l)*cc(l,i_o) |
---|
| 643 | c |
---|
| 644 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 645 | c hox = h + oh + ho2 |
---|
| 646 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 647 | c |
---|
| 648 | production(l,i_hox) = |
---|
| 649 | $ 2.*j(l,j_h2o)*cc(l,i_h2o) |
---|
| 650 | $ + 2.*j(l,j_h2o2)*cc(l,i_h2o2) |
---|
| 651 | $ + 2.*b002(l)*cc(l,i_o1d)*cc(l,i_h2o) |
---|
| 652 | $ + 2.*b003(l)*cc(l,i_o1d)*cc(l,i_h2) |
---|
| 653 | $ + 2.*c012(l)*cc(l,i_o)*cc(l,i_h2o2) |
---|
| 654 | c |
---|
| 655 | loss(l,i_hox) = 2.*c005(l)*cc(l,i_h)*cc(l,i_ho2) |
---|
| 656 | $ + 2.*c006(l)*cc(l,i_h)*cc(l,i_ho2) |
---|
| 657 | $ + 2.*c007h(l)*cc(l,i_oh)*cc(l,i_ho2) |
---|
| 658 | $ + 2.*c008(l)*cc(l,i_ho2)*cc(l,i_ho2) |
---|
| 659 | $ + 2.*c013(l)*cc(l,i_oh)*cc(l,i_oh) |
---|
| 660 | $ + 2.*c016(l)*cc(l,i_ho2)*cc(l,i_ho2) |
---|
| 661 | $ + 2.*c017(l)*cc(l,i_oh)*cc(l,i_oh) |
---|
| 662 | $ + 2.*c018(l)*cc(l,i_h)*cc(l,i_h) |
---|
| 663 | c |
---|
| 664 | loss(l,i_hox) = loss(l,i_hox)/cc(l,i_hox) |
---|
| 665 | c |
---|
| 666 | end do |
---|
| 667 | c |
---|
| 668 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 669 | c update number densities |
---|
| 670 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 671 | c |
---|
| 672 | c long-lived species |
---|
| 673 | c |
---|
| 674 | do l = 1,lswitch-1 |
---|
| 675 | cc(l,i_co) = (cc0(l,i_co) + production(l,i_co)*dt) |
---|
| 676 | $ /(1. + loss(l,i_co)*dt) |
---|
| 677 | cc(l,i_o2) = (cc0(l,i_o2) + production(l,i_o2)*dt) |
---|
| 678 | $ /(1. + loss(l,i_o2)*dt) |
---|
| 679 | cc(l,i_h2) = (cc0(l,i_h2) + production(l,i_h2)*dt) |
---|
| 680 | $ /(1. + loss(l,i_h2)*dt) |
---|
| 681 | cc(l,i_h2o2)= (cc0(l,i_h2o2) + production(l,i_h2o2)*dt) |
---|
| 682 | $ /(1. + loss(l,i_h2o2)*dt) |
---|
| 683 | cc(l,i_h2o) = (cc0(l,i_h2o) + production(l,i_h2o)*dt) |
---|
| 684 | $ /(1. + loss(l,i_h2o)*dt) |
---|
| 685 | cc(l,i_hox) = (cc0(l,i_hox) + production(l,i_hox)*dt) |
---|
| 686 | $ /(1. + loss(l,i_hox)*dt) |
---|
| 687 | end do |
---|
| 688 | c |
---|
| 689 | c ox species |
---|
| 690 | c |
---|
| 691 | if (sza .le. 95.) then |
---|
| 692 | do l = 1,lswitch-1 |
---|
| 693 | cc(l,i_ox) = (cc0(l,i_ox) + production(l,i_ox)*dt) |
---|
| 694 | $ /(1. + loss(l,i_ox)*dt) |
---|
| 695 | end do |
---|
| 696 | else |
---|
| 697 | do l = 1,lswitch-1 |
---|
| 698 | cc(l,i_o) = (cc0(l,i_o) + production(l,i_o)*dt) |
---|
| 699 | $ /(1. + loss(l,i_o)*dt) |
---|
| 700 | cc(l,i_o3) = (cc0(l,i_o3) + production(l,i_o3)*dt) |
---|
| 701 | $ /(1. + loss(l,i_o3)*dt) |
---|
| 702 | cc(l,i_ox) = cc(l,i_o) + cc(l,i_o3) |
---|
| 703 | end do |
---|
| 704 | end if |
---|
| 705 | c |
---|
| 706 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 707 | c end of loop over iterations |
---|
| 708 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 709 | c |
---|
| 710 | end do |
---|
| 711 | c |
---|
| 712 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 713 | c density -> volume mixing ratio conversion |
---|
| 714 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 715 | c |
---|
| 716 | do iesp = 1,nesp |
---|
| 717 | do l = 1,lswitch-1 |
---|
| 718 | rm(l,iesp) = max(cc(l,iesp)/dens(l), 1.e-30) |
---|
| 719 | end do |
---|
| 720 | end do |
---|
| 721 | c |
---|
| 722 | return |
---|
| 723 | end |
---|
| 724 | c |
---|
| 725 | c***************************************************************** |
---|
| 726 | c |
---|
| 727 | subroutine phot(lswitch, press, temp, sza, dist_sol, rmo3, j) |
---|
| 728 | c |
---|
| 729 | c***************************************************************** |
---|
| 730 | c |
---|
| 731 | implicit none |
---|
| 732 | c |
---|
| 733 | #include "dimensions.h" |
---|
| 734 | #include "dimphys.h" |
---|
| 735 | #include "chimiedata.h" |
---|
| 736 | c |
---|
| 737 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 738 | c inputs: |
---|
| 739 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 740 | c |
---|
| 741 | integer lswitch ! interface level between chemistries |
---|
| 742 | real press(nlayermx) ! pressure (hPa) |
---|
| 743 | real temp(nlayermx) ! temperature (K) |
---|
| 744 | real sza ! solar zenith angle (deg) |
---|
| 745 | real dist_sol ! sun distance (AU) |
---|
| 746 | real rmo3(nlayermx) ! ozone mixing ratio |
---|
| 747 | c |
---|
| 748 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 749 | c output: |
---|
| 750 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 751 | c |
---|
| 752 | real j(nlayermx,nd) ! interpolated photolysis rates (s-1) |
---|
| 753 | c |
---|
| 754 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 755 | c local: |
---|
| 756 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 757 | c |
---|
| 758 | integer icol, ij, indsza, indcol, indozo, indtemp, |
---|
| 759 | $ iozo, isza, it, iztop, l |
---|
| 760 | c |
---|
| 761 | real col(nlayermx) ! overhead air column (molecule cm-2) |
---|
| 762 | real colo3(nlayermx) ! overhead ozone column (molecule cm-2) |
---|
| 763 | real poids(2,2,2,2) ! 4D interpolation weights |
---|
| 764 | real tref ! temperature at 1.9 hPa in the gcm (K) |
---|
| 765 | real table_temp(ntemp) ! temperatures at 1.9 hPa in jmars (K) |
---|
| 766 | real cinf, csup, cicol, |
---|
| 767 | $ ciozo, cisza, citemp |
---|
| 768 | real colo3min, dp, scaleh |
---|
| 769 | real ratio_o3(nlayermx) |
---|
| 770 | c |
---|
| 771 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 772 | c day/night criterion |
---|
| 773 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 774 | c |
---|
| 775 | if (sza .le. 95.) then |
---|
| 776 | c |
---|
| 777 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 778 | c temperatures at 1.9 hPa in jmars |
---|
| 779 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 780 | c |
---|
| 781 | table_temp(1) = 226.2 |
---|
| 782 | table_temp(2) = 206.2 |
---|
| 783 | table_temp(3) = 186.2 |
---|
| 784 | table_temp(4) = 169.8 |
---|
| 785 | c |
---|
| 786 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 787 | c interpolation in solar zenith angle |
---|
| 788 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 789 | c |
---|
| 790 | indsza = nsza - 1 |
---|
| 791 | do isza = 1,nsza |
---|
| 792 | if (szatab(isza) .ge. sza) then |
---|
| 793 | indsza = min(indsza,isza - 1) |
---|
| 794 | indsza = max(indsza, 1) |
---|
| 795 | end if |
---|
| 796 | end do |
---|
| 797 | cisza = (sza - szatab(indsza)) |
---|
| 798 | $ /(szatab(indsza + 1) - szatab(indsza)) |
---|
| 799 | c |
---|
| 800 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 801 | c air and ozone columns |
---|
| 802 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 803 | c |
---|
| 804 | c air column at model top: assign standard value |
---|
| 805 | c |
---|
| 806 | scaleh = 10. |
---|
| 807 | iztop = nint(scaleh*log(press(1)/press(lswitch-1))) |
---|
| 808 | iztop = min(iztop,200) |
---|
| 809 | col(lswitch-1) = colairtab(iztop) |
---|
| 810 | c |
---|
| 811 | c ozone column at model top |
---|
| 812 | c |
---|
| 813 | colo3(lswitch-1) = 0. |
---|
| 814 | c |
---|
| 815 | c air and ozone columns for other levels |
---|
| 816 | c |
---|
| 817 | do l = lswitch-2,1,-1 |
---|
| 818 | dp = (press(l) - press(l+1))*100. |
---|
| 819 | col(l) = col(l+1) + factor*dp |
---|
| 820 | col(l) = min(col(l), colairtab(0)) |
---|
| 821 | colo3(l) = colo3(l+1) |
---|
| 822 | $ + factor*(rmo3(l+1) + rmo3(l))*0.5*dp |
---|
| 823 | end do |
---|
| 824 | c |
---|
| 825 | c ratio ozone column/minimal theoretical column (ro3 = 7.171e-10) |
---|
| 826 | c |
---|
| 827 | do l = 1,lswitch-1 |
---|
| 828 | colo3min = col(l)*7.171e-10 |
---|
| 829 | ratio_o3(l) = colo3(l)/colo3min |
---|
| 830 | ratio_o3(l) = min(ratio_o3(l), table_ozo(nozo)*10.) |
---|
| 831 | ratio_o3(l) = max(ratio_o3(l), table_ozo(1)*10.) |
---|
| 832 | end do |
---|
| 833 | c |
---|
| 834 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 835 | c temperature dependence |
---|
| 836 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 837 | c |
---|
| 838 | c 1) search for temperature at 1.9 hPa (tref): vertical interpolation |
---|
| 839 | c |
---|
| 840 | c print*,'****** recherche de T a 1.9 hPa ******' |
---|
| 841 | tref = temp(1) |
---|
| 842 | do l = (lswitch-1)-1,1,-1 |
---|
| 843 | if (press(l) .gt. 1.9) then |
---|
| 844 | cinf = (press(l) - 1.9) |
---|
| 845 | $ /(press(l) - press(l+1)) |
---|
| 846 | csup = 1. - cinf |
---|
| 847 | tref = cinf*temp(l+1) + csup*temp(l) |
---|
| 848 | c |
---|
| 849 | c print*, press(l+1), temp(l+1) |
---|
| 850 | c print*, press(l ), temp(l ) |
---|
| 851 | c print*, 'tref = ', tref |
---|
| 852 | go to 10 |
---|
| 853 | end if |
---|
| 854 | end do |
---|
| 855 | 10 continue |
---|
| 856 | c |
---|
| 857 | c 2) interpolation in temperature |
---|
| 858 | c |
---|
| 859 | tref = min(tref,table_temp(1)) |
---|
| 860 | tref = max(tref,table_temp(ntemp)) |
---|
| 861 | c |
---|
| 862 | do it = 2, ntemp |
---|
| 863 | if (table_temp(it) .le. tref) then |
---|
| 864 | citemp = (log(tref) - log(table_temp(it))) |
---|
| 865 | $ /(log(table_temp(it-1)) - log(table_temp(it))) |
---|
| 866 | indtemp = it - 1 |
---|
| 867 | c print*,'interpolation entre ', table_temp(it), ' et ', |
---|
| 868 | c $ table_temp(it-1) |
---|
| 869 | c print*,'indtemp = ', indtemp, ' citemp = ', citemp |
---|
| 870 | goto 20 |
---|
| 871 | end if |
---|
| 872 | end do |
---|
| 873 | 20 continue |
---|
| 874 | c print*,'****** fin interpolation en T ********' |
---|
| 875 | c |
---|
| 876 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 877 | c loop over vertical levels |
---|
| 878 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 879 | c |
---|
| 880 | do l = 1,lswitch-1 |
---|
| 881 | c |
---|
| 882 | c interpolation in air column |
---|
| 883 | c |
---|
| 884 | do icol = 0,200 |
---|
| 885 | if (colairtab(icol) .lt. col(l)) then |
---|
| 886 | cicol = (log(col(l)) - log(colairtab(icol))) |
---|
| 887 | $ /(log(colairtab(icol-1)) - log(colairtab(icol))) |
---|
| 888 | indcol = icol - 1 |
---|
| 889 | goto 30 |
---|
| 890 | end if |
---|
| 891 | end do |
---|
| 892 | 30 continue |
---|
| 893 | c |
---|
| 894 | cc interpolation in ozone column |
---|
| 895 | c |
---|
| 896 | indozo = nozo - 1 |
---|
| 897 | do iozo = 1,nozo |
---|
| 898 | if (table_ozo(iozo)*10. .ge. ratio_o3(l)) then |
---|
| 899 | indozo = min(indozo, iozo - 1) |
---|
| 900 | indozo = max(indozo, 1) |
---|
| 901 | end if |
---|
| 902 | end do |
---|
| 903 | ciozo = (ratio_o3(l) - table_ozo(indozo)*10.) |
---|
| 904 | $ /(table_ozo(indozo + 1)*10. - table_ozo(indozo)*10.) |
---|
| 905 | c |
---|
| 906 | c write(6,'(i2,4e13.4,5f10.4)') |
---|
| 907 | c $ l, press(l), col(l), |
---|
| 908 | c $ colairtab(indcol), colairtab(indcol+1), |
---|
| 909 | c $ colo3(l)/2.687e15, |
---|
| 910 | c $ ratio_o3(l), table_ozo(indozo), table_ozo(indozo+1), |
---|
| 911 | c $ ciozo |
---|
| 912 | c |
---|
| 913 | cc 4-dimensional interpolation weights |
---|
| 914 | c |
---|
| 915 | poids(1,1,1,1) = citemp*(1.-cisza)*cicol*(1.-ciozo) |
---|
| 916 | poids(1,1,1,2) = citemp*(1.-cisza)*cicol*ciozo |
---|
| 917 | poids(1,1,2,1) = citemp*(1.-cisza)*(1.-cicol)*(1.-ciozo) |
---|
| 918 | poids(1,1,2,2) = citemp*(1.-cisza)*(1.-cicol)*ciozo |
---|
| 919 | poids(1,2,1,1) = citemp*cisza*cicol*(1.-ciozo) |
---|
| 920 | poids(1,2,1,2) = citemp*cisza*cicol*ciozo |
---|
| 921 | poids(1,2,2,1) = citemp*cisza*(1.-cicol)*(1.-ciozo) |
---|
| 922 | poids(1,2,2,2) = citemp*cisza*(1.-cicol)*ciozo |
---|
| 923 | poids(2,1,1,1) = (1.-citemp)*(1.-cisza)*cicol*(1.-ciozo) |
---|
| 924 | poids(2,1,1,2) = (1.-citemp)*(1.-cisza)*cicol*ciozo |
---|
| 925 | poids(2,1,2,1) = (1.-citemp)*(1.-cisza)*(1.-cicol)*(1.-ciozo) |
---|
| 926 | poids(2,1,2,2) = (1.-citemp)*(1.-cisza)*(1.-cicol)*ciozo |
---|
| 927 | poids(2,2,1,1) = (1.-citemp)*cisza*cicol*(1.-ciozo) |
---|
| 928 | poids(2,2,1,2) = (1.-citemp)*cisza*cicol*ciozo |
---|
| 929 | poids(2,2,2,1) = (1.-citemp)*cisza*(1.-cicol)*(1.-ciozo) |
---|
| 930 | poids(2,2,2,2) = (1.-citemp)*cisza*(1.-cicol)*ciozo |
---|
| 931 | c |
---|
| 932 | cc 4-dimensional interpolation in the lookup table |
---|
| 933 | c |
---|
| 934 | do ij = 1,nd |
---|
| 935 | j(l,ij) = |
---|
| 936 | $ poids(1,1,1,1)*jphot(indtemp,indsza,indcol,indozo,ij) |
---|
| 937 | $ + poids(1,1,1,2)*jphot(indtemp,indsza,indcol,indozo+1,ij) |
---|
| 938 | $ + poids(1,1,2,1)*jphot(indtemp,indsza,indcol+1,indozo,ij) |
---|
| 939 | $ + poids(1,1,2,2)*jphot(indtemp,indsza,indcol+1,indozo+1,ij) |
---|
| 940 | $ + poids(1,2,1,1)*jphot(indtemp,indsza+1,indcol,indozo,ij) |
---|
| 941 | $ + poids(1,2,1,2)*jphot(indtemp,indsza+1,indcol,indozo+1,ij) |
---|
| 942 | $ + poids(1,2,2,1)*jphot(indtemp,indsza+1,indcol+1,indozo,ij) |
---|
| 943 | $ + poids(1,2,2,2)*jphot(indtemp,indsza+1,indcol+1,indozo+1,ij) |
---|
| 944 | $ + poids(2,1,1,1)*jphot(indtemp+1,indsza,indcol,indozo,ij) |
---|
| 945 | $ + poids(2,1,1,2)*jphot(indtemp+1,indsza,indcol,indozo+1,ij) |
---|
| 946 | $ + poids(2,1,2,1)*jphot(indtemp+1,indsza,indcol+1,indozo,ij) |
---|
| 947 | $ + poids(2,1,2,2)*jphot(indtemp+1,indsza,indcol+1,indozo+1,ij) |
---|
| 948 | $ + poids(2,2,1,1)*jphot(indtemp+1,indsza+1,indcol,indozo,ij) |
---|
| 949 | $ + poids(2,2,1,2)*jphot(indtemp+1,indsza+1,indcol,indozo+1,ij) |
---|
| 950 | $ + poids(2,2,2,1)*jphot(indtemp+1,indsza+1,indcol+1,indozo,ij) |
---|
| 951 | $ + poids(2,2,2,2)*jphot(indtemp+1,indsza+1,indcol+1,indozo+1,ij) |
---|
| 952 | end do |
---|
| 953 | c |
---|
| 954 | cc correction for sun distance |
---|
| 955 | c |
---|
| 956 | do ij = 1,nd |
---|
| 957 | j(l,ij) = j(l,ij)*(1.52/dist_sol)**2. |
---|
| 958 | end do |
---|
| 959 | c |
---|
| 960 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 961 | c end of loop over vertical levels |
---|
| 962 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 963 | c |
---|
| 964 | end do |
---|
| 965 | c |
---|
| 966 | else |
---|
| 967 | c |
---|
| 968 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 969 | c night |
---|
| 970 | ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 971 | c |
---|
| 972 | do ij = 1,nd |
---|
| 973 | do l = 1,lswitch-1 |
---|
| 974 | j(l,ij) = 0. |
---|
| 975 | end do |
---|
| 976 | end do |
---|
| 977 | c |
---|
| 978 | end if |
---|
| 979 | c |
---|
| 980 | return |
---|
| 981 | end |
---|
| 982 | c |
---|
| 983 | c***************************************************************** |
---|
| 984 | c |
---|
| 985 | subroutine gcmtochim(zycol, lswitch, nesp, rm) |
---|
| 986 | c |
---|
| 987 | c***************************************************************** |
---|
| 988 | c |
---|
| 989 | implicit none |
---|
| 990 | c |
---|
| 991 | #include "dimensions.h" |
---|
| 992 | #include "dimphys.h" |
---|
| 993 | #include "callkeys.h" |
---|
| 994 | c |
---|
| 995 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 996 | c inputs: |
---|
| 997 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 998 | c |
---|
| 999 | real zycol(nlayermx,nqmx)! species volume mixing ratio in the gcm |
---|
| 1000 | c |
---|
| 1001 | integer nesp ! number of species in the chemistry |
---|
| 1002 | integer lswitch ! interface level between chemistries |
---|
| 1003 | c |
---|
| 1004 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1005 | c outputs: |
---|
| 1006 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1007 | c |
---|
| 1008 | real rm(nlayermx,nesp) ! species volume mixing ratio |
---|
| 1009 | c |
---|
| 1010 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1011 | c local: |
---|
| 1012 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1013 | c |
---|
| 1014 | integer l |
---|
| 1015 | integer i_co2, i_co, i_o2, i_h2, i_h2o, i_h2o2, i_hox, i_ox, |
---|
| 1016 | $ i_o1d, i_o, i_o3, i_h, i_oh, i_ho2, i_n2 |
---|
| 1017 | c |
---|
| 1018 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1019 | c tracers numbering in the gcm |
---|
| 1020 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1021 | c |
---|
| 1022 | c co2 = nqchem_min |
---|
| 1023 | c co = nqchem_min + 1 |
---|
| 1024 | c o = nqchem_min + 2 |
---|
| 1025 | c o(1d) = nqchem_min + 3 |
---|
| 1026 | c o2 = nqchem_min + 4 |
---|
| 1027 | c o3 = nqchem_min + 5 |
---|
| 1028 | c h = nqchem_min + 6 |
---|
| 1029 | c h2 = nqchem_min + 7 |
---|
| 1030 | c oh = nqchem_min + 8 |
---|
| 1031 | c ho2 = nqchem_min + 9 |
---|
| 1032 | c h2o2 = nqchem_min + 10 |
---|
| 1033 | c n2 = nqchem_min + 11 |
---|
| 1034 | c ar = nqchem_min + 12 |
---|
| 1035 | c h2o = nqmx |
---|
| 1036 | c |
---|
| 1037 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1038 | c tracers numbering in the chemistry |
---|
| 1039 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1040 | c |
---|
| 1041 | i_co2 = 1 |
---|
| 1042 | i_co = 2 |
---|
| 1043 | i_o = 3 |
---|
| 1044 | i_o1d = 4 |
---|
| 1045 | i_o2 = 5 |
---|
| 1046 | i_o3 = 6 |
---|
| 1047 | i_h = 7 |
---|
| 1048 | i_h2 = 8 |
---|
| 1049 | i_oh = 9 |
---|
| 1050 | i_ho2 = 10 |
---|
| 1051 | i_h2o2 = 11 |
---|
| 1052 | i_h2o = 12 |
---|
| 1053 | i_n2 = 13 |
---|
| 1054 | i_hox = 14 |
---|
| 1055 | i_ox = 15 |
---|
| 1056 | c |
---|
| 1057 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1058 | c initialise chemical species |
---|
| 1059 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1060 | c |
---|
| 1061 | do l = 1,lswitch-1 |
---|
| 1062 | rm(l,i_co2) = max(zycol(l, nqchem_min), 1.e-30) |
---|
| 1063 | rm(l,i_co) = max(zycol(l, nqchem_min + 1), 1.e-30) |
---|
| 1064 | rm(l,i_o) = max(zycol(l, nqchem_min + 2), 1.e-30) |
---|
| 1065 | rm(l,i_o1d) = max(zycol(l, nqchem_min + 3), 1.e-30) |
---|
| 1066 | rm(l,i_o2) = max(zycol(l, nqchem_min + 4), 1.e-30) |
---|
| 1067 | rm(l,i_o3) = max(zycol(l, nqchem_min + 5), 1.e-30) |
---|
| 1068 | rm(l,i_h) = max(zycol(l, nqchem_min + 6), 1.e-30) |
---|
| 1069 | rm(l,i_h2) = max(zycol(l, nqchem_min + 7), 1.e-30) |
---|
| 1070 | rm(l,i_oh) = max(zycol(l, nqchem_min + 8), 1.e-30) |
---|
| 1071 | rm(l,i_ho2) = max(zycol(l, nqchem_min + 9), 1.e-30) |
---|
| 1072 | rm(l,i_h2o2) = max(zycol(l, nqchem_min + 10), 1.e-30) |
---|
| 1073 | rm(l,i_n2) = max(zycol(l, nqchem_min + 11), 1.e-30) |
---|
| 1074 | rm(l,i_h2o) = max(zycol(l, nqmx), 1.e-30) |
---|
| 1075 | c |
---|
| 1076 | end do |
---|
| 1077 | c |
---|
| 1078 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1079 | c initialise chemical families c |
---|
| 1080 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1081 | c |
---|
| 1082 | do l = 1,lswitch-1 |
---|
| 1083 | rm(l,i_hox) = rm(l,i_h) |
---|
| 1084 | $ + rm(l,i_oh) |
---|
| 1085 | $ + rm(l,i_ho2) |
---|
| 1086 | rm(l,i_ox) = rm(l,i_o) |
---|
| 1087 | $ + rm(l,i_o3) |
---|
| 1088 | end do |
---|
| 1089 | c |
---|
| 1090 | return |
---|
| 1091 | end |
---|
| 1092 | c |
---|
| 1093 | c***************************************************************** |
---|
| 1094 | c |
---|
| 1095 | subroutine chimtogcm(zycol, lswitch, nesp, rm) |
---|
| 1096 | c |
---|
| 1097 | c***************************************************************** |
---|
| 1098 | c |
---|
| 1099 | implicit none |
---|
| 1100 | c |
---|
| 1101 | #include "dimensions.h" |
---|
| 1102 | #include "dimphys.h" |
---|
| 1103 | #include "callkeys.h" |
---|
| 1104 | c |
---|
| 1105 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1106 | c inputs: |
---|
| 1107 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1108 | c |
---|
| 1109 | integer nesp ! number of species in the chemistry |
---|
| 1110 | integer lswitch ! interface level between chemistries |
---|
| 1111 | c |
---|
| 1112 | real rm(nlayermx,nesp) ! species volume mixing ratio |
---|
| 1113 | c |
---|
| 1114 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1115 | c output: |
---|
| 1116 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1117 | c |
---|
| 1118 | real zycol(nlayermx,nqmx) ! species volume mixing ratio in the gcm |
---|
| 1119 | c |
---|
| 1120 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1121 | c local: |
---|
| 1122 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1123 | c |
---|
| 1124 | integer l |
---|
| 1125 | integer i_co2, i_co, i_o2, i_h2, i_h2o, i_h2o2, i_hox, i_ox, |
---|
| 1126 | $ i_o1d, i_o, i_o3, i_h, i_oh, i_ho2, i_n2 |
---|
| 1127 | c |
---|
| 1128 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1129 | c tracers numbering in the gcm |
---|
| 1130 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1131 | c |
---|
| 1132 | c co2 = nqchem_min |
---|
| 1133 | c co = nqchem_min + 1 |
---|
| 1134 | c o = nqchem_min + 2 |
---|
| 1135 | c o(1d) = nqchem_min + 3 |
---|
| 1136 | c o2 = nqchem_min + 4 |
---|
| 1137 | c o3 = nqchem_min + 5 |
---|
| 1138 | c h = nqchem_min + 6 |
---|
| 1139 | c h2 = nqchem_min + 7 |
---|
| 1140 | c oh = nqchem_min + 8 |
---|
| 1141 | c ho2 = nqchem_min + 9 |
---|
| 1142 | c h2o2 = nqchem_min + 10 |
---|
| 1143 | c n2 = nqchem_min + 11 |
---|
| 1144 | c ar = nqchem_min + 12 |
---|
| 1145 | c h2o = nqmx |
---|
| 1146 | c |
---|
| 1147 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1148 | c tracers numbering in the chemistry |
---|
| 1149 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1150 | c |
---|
| 1151 | i_co2 = 1 |
---|
| 1152 | i_co = 2 |
---|
| 1153 | i_o = 3 |
---|
| 1154 | i_o1d = 4 |
---|
| 1155 | i_o2 = 5 |
---|
| 1156 | i_o3 = 6 |
---|
| 1157 | i_h = 7 |
---|
| 1158 | i_h2 = 8 |
---|
| 1159 | i_oh = 9 |
---|
| 1160 | i_ho2 = 10 |
---|
| 1161 | i_h2o2 = 11 |
---|
| 1162 | i_h2o = 12 |
---|
| 1163 | i_n2 = 13 |
---|
| 1164 | i_hox = 14 |
---|
| 1165 | i_ox = 15 |
---|
| 1166 | c |
---|
| 1167 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1168 | c save mixing ratios for the gcm |
---|
| 1169 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1170 | c |
---|
| 1171 | do l = 1,lswitch-1 |
---|
| 1172 | zycol(l, nqchem_min) = rm(l,i_co2) |
---|
| 1173 | zycol(l, nqchem_min + 1) = rm(l,i_co) |
---|
| 1174 | zycol(l, nqchem_min + 2) = rm(l,i_o) |
---|
| 1175 | zycol(l, nqchem_min + 3) = rm(l,i_o1d) |
---|
| 1176 | zycol(l, nqchem_min + 4) = rm(l,i_o2) |
---|
| 1177 | zycol(l, nqchem_min + 5) = rm(l,i_o3) |
---|
| 1178 | zycol(l, nqchem_min + 6) = rm(l,i_h) |
---|
| 1179 | zycol(l, nqchem_min + 7) = rm(l,i_h2) |
---|
| 1180 | zycol(l, nqchem_min + 8) = rm(l,i_oh) |
---|
| 1181 | zycol(l, nqchem_min + 9) = rm(l,i_ho2) |
---|
| 1182 | zycol(l, nqchem_min + 10) = rm(l,i_h2o2) |
---|
| 1183 | zycol(l, nqchem_min + 11) = rm(l,i_n2) |
---|
| 1184 | zycol(l, nqmx) = rm(l,i_h2o) |
---|
| 1185 | end do |
---|
| 1186 | c |
---|
| 1187 | c water ice: zycol(nqmx-1) has not been changed if iceparty on |
---|
| 1188 | c |
---|
| 1189 | return |
---|
| 1190 | end |
---|
| 1191 | c |
---|
| 1192 | c***************************************************************** |
---|
| 1193 | c |
---|
| 1194 | subroutine chemrates(lswitch, dens, press, t, icesurf, |
---|
| 1195 | $ a001, a002, a003, |
---|
| 1196 | $ b001, b002, b003, b004, b005, b006, |
---|
| 1197 | $ c001, c002, c003, c004, c005, c006, |
---|
| 1198 | $ c007, c008, c009, c010, c011, c012, |
---|
| 1199 | $ c013, c014, c015, c016, c017, c018, |
---|
| 1200 | $ d001, d002, d003, |
---|
| 1201 | $ e001, e002, |
---|
| 1202 | $ h001, h002) |
---|
| 1203 | c |
---|
| 1204 | c***************************************************************** |
---|
| 1205 | c |
---|
| 1206 | implicit none |
---|
| 1207 | c |
---|
| 1208 | #include "dimensions.h" |
---|
| 1209 | #include "dimphys.h" |
---|
| 1210 | c |
---|
| 1211 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1212 | c inputs: c |
---|
| 1213 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1214 | c |
---|
| 1215 | integer lswitch ! interface level between chemistries |
---|
| 1216 | |
---|
| 1217 | real dens(nlayermx) ! density (cm-3) |
---|
| 1218 | real press(nlayermx) ! pressure (hpa) |
---|
| 1219 | real t(nlayermx) ! temperature (k) |
---|
| 1220 | real icesurf(nlayermx) ! ice surface area (cm^2/cm^3) |
---|
| 1221 | c |
---|
| 1222 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1223 | c outputs: c |
---|
| 1224 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1225 | c |
---|
| 1226 | real a001(nlayermx), a002(nlayermx), a003(nlayermx) |
---|
| 1227 | real b001(nlayermx), b002(nlayermx), b003(nlayermx), |
---|
| 1228 | $ b004(nlayermx), b005(nlayermx), b006(nlayermx) |
---|
| 1229 | real c001(nlayermx), c002(nlayermx), c003(nlayermx), |
---|
| 1230 | $ c004(nlayermx), c005(nlayermx), c006(nlayermx), |
---|
| 1231 | $ c007(nlayermx), c008(nlayermx), c009(nlayermx), |
---|
| 1232 | $ c010(nlayermx), c011(nlayermx), c012(nlayermx), |
---|
| 1233 | $ c013(nlayermx), c014(nlayermx), c015(nlayermx), |
---|
| 1234 | $ c016(nlayermx), c017(nlayermx), c018(nlayermx) |
---|
| 1235 | real d001(nlayermx), d002(nlayermx), d003(nlayermx) |
---|
| 1236 | real e001(nlayermx), e002(nlayermx) |
---|
| 1237 | real h001(nlayermx), h002(nlayermx) |
---|
| 1238 | c |
---|
| 1239 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1240 | c local: c |
---|
| 1241 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1242 | c |
---|
| 1243 | real ak0, ak1, rate, xpo |
---|
| 1244 | c |
---|
| 1245 | integer l |
---|
| 1246 | c |
---|
| 1247 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1248 | c compute reaction rates from JPL 1997, otherwise mentioned |
---|
| 1249 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1250 | c |
---|
| 1251 | do l = 1,lswitch-1 |
---|
| 1252 | c |
---|
| 1253 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1254 | c oxygen compounds |
---|
| 1255 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1256 | c |
---|
| 1257 | ccc a001: o + o2 + co2 -> o3 + co2 |
---|
| 1258 | c |
---|
| 1259 | c jpl 2003 |
---|
| 1260 | c |
---|
| 1261 | a001(l) = 2.5 |
---|
| 1262 | $ *6.0e-34*(t(l)/300.)**(-2.4)*dens(l) |
---|
| 1263 | c |
---|
| 1264 | c nair et al., 1994 |
---|
| 1265 | c |
---|
| 1266 | c a001(l) = 1.3e-34*exp(724./t(l))*dens(l) |
---|
| 1267 | c |
---|
| 1268 | ccc a002: o + o + co2 -> o2 + co2 |
---|
| 1269 | c |
---|
| 1270 | c Tsang and Hampson, J. Chem. Phys. Ref. Data, 15, 1087, 1986 |
---|
| 1271 | c |
---|
| 1272 | a002(l) = 2.5*5.2e-35*exp(900./t(l))*dens(l) |
---|
| 1273 | c |
---|
| 1274 | c Campbell and Gray, Chem. Phys. Lett., 18, 607, 1973 |
---|
| 1275 | c |
---|
| 1276 | c a002(l) = 1.2e-32*(300./t(l))**(2.0)*dens(l) |
---|
| 1277 | c |
---|
| 1278 | ccc a003: o + o3 -> o2 + o2 |
---|
| 1279 | c |
---|
| 1280 | c jpl 2003 |
---|
| 1281 | c |
---|
| 1282 | a003(l) = 8.0e-12*exp(-2060./t(l)) |
---|
| 1283 | c |
---|
| 1284 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1285 | c reactions with o(1d) |
---|
| 1286 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1287 | c |
---|
| 1288 | ccc b001: o(1d) + co2 -> o + co2 |
---|
| 1289 | c |
---|
| 1290 | c jpl 2003 |
---|
| 1291 | c |
---|
| 1292 | b001(l) = 7.4e-11*exp(120./t(l)) |
---|
| 1293 | c |
---|
| 1294 | ccc b002: o(1d) + h2o -> oh + oh |
---|
| 1295 | c |
---|
| 1296 | c jpl 2003 |
---|
| 1297 | c |
---|
| 1298 | b002(l) = 2.2e-10 |
---|
| 1299 | c |
---|
| 1300 | ccc b003: o(1d) + h2 -> oh + h |
---|
| 1301 | c |
---|
| 1302 | c jpl 2003 |
---|
| 1303 | c |
---|
| 1304 | b003(l) = 1.1e-10 |
---|
| 1305 | c |
---|
| 1306 | c nair et al., 1994 |
---|
| 1307 | c |
---|
| 1308 | c b003(l) = 1.0e-10 |
---|
| 1309 | c |
---|
| 1310 | ccc b004: o(1d) + o2 -> o + o2 |
---|
| 1311 | c |
---|
| 1312 | c jpl 2003 |
---|
| 1313 | c |
---|
| 1314 | b004(l) = 3.2e-11*exp(70./t(l)) |
---|
| 1315 | c |
---|
| 1316 | ccc b005: o(1d) + o3 -> o2 + o2 |
---|
| 1317 | c |
---|
| 1318 | c jpl 2003 |
---|
| 1319 | c |
---|
| 1320 | b005(l) = 1.2e-10 |
---|
| 1321 | c |
---|
| 1322 | ccc b006: o(1d) + o3 -> o2 + o + o |
---|
| 1323 | c |
---|
| 1324 | c jpl 2003 |
---|
| 1325 | c |
---|
| 1326 | b006(l) = 1.2e-10 |
---|
| 1327 | c |
---|
| 1328 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1329 | c hydrogen compounds |
---|
| 1330 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1331 | c |
---|
| 1332 | ccc c001: o + ho2 -> oh + o2 |
---|
| 1333 | c |
---|
| 1334 | c jpl 2003 |
---|
| 1335 | c |
---|
| 1336 | c001(l) = 3.0e-11*exp(200./t(l)) |
---|
| 1337 | c c001(l) = 0.75*3.0e-11*exp(200./t(l)) |
---|
| 1338 | c |
---|
| 1339 | ccc c002: o + oh -> o2 + h |
---|
| 1340 | c |
---|
| 1341 | c jpl 2003 |
---|
| 1342 | c |
---|
| 1343 | c002(l) = 2.2e-11*exp(120./t(l)) |
---|
| 1344 | c |
---|
| 1345 | ccc c003: h + o3 -> oh + o2 |
---|
| 1346 | c |
---|
| 1347 | c jpl 2003 |
---|
| 1348 | c |
---|
| 1349 | c003(l) = 1.4e-10*exp(-470./t(l)) |
---|
| 1350 | c |
---|
| 1351 | ccc c004: h + ho2 -> oh + oh |
---|
| 1352 | c |
---|
| 1353 | c jpl 2003 |
---|
| 1354 | c |
---|
| 1355 | c004(l) = 8.1e-11*0.90 |
---|
| 1356 | c |
---|
| 1357 | ccc c005: h + ho2 -> h2 + o2 |
---|
| 1358 | c |
---|
| 1359 | c jpl 2003 |
---|
| 1360 | c |
---|
| 1361 | c005(l) = 8.1e-11*0.08 |
---|
| 1362 | c |
---|
| 1363 | ccc c006: h + ho2 -> h2o + o |
---|
| 1364 | c |
---|
| 1365 | c jpl 2003 |
---|
| 1366 | c |
---|
| 1367 | c006(l) = 8.1e-11*0.02 |
---|
| 1368 | c |
---|
| 1369 | ccc c007: oh + ho2 -> h2o + o2 |
---|
| 1370 | c |
---|
| 1371 | c jpl 2003 |
---|
| 1372 | c |
---|
| 1373 | c007(l) = 4.8e-11*exp(250./t(l)) |
---|
| 1374 | c c007(l) = 0.75*4.8e-11*exp(250./t(l)) |
---|
| 1375 | c |
---|
| 1376 | ccc c008: ho2 + ho2 -> h2o2 + o2 |
---|
| 1377 | c |
---|
| 1378 | c jpl 2003 |
---|
| 1379 | c |
---|
| 1380 | c c008(l) = 2.3e-13*exp(600./t(l)) |
---|
| 1381 | c |
---|
| 1382 | c christensen et al., grl, 13, 2002 |
---|
| 1383 | c |
---|
| 1384 | c008(l) = 1.5e-12*exp(19./t(l)) |
---|
| 1385 | c |
---|
| 1386 | ccc c009: oh + h2o2 -> h2o + ho2 |
---|
| 1387 | c |
---|
| 1388 | c jpl 2003 |
---|
| 1389 | c |
---|
| 1390 | c009(l) = 2.9e-12*exp(-160./t(l)) |
---|
| 1391 | c |
---|
| 1392 | ccc c010: oh + h2 -> h2o + h |
---|
| 1393 | c |
---|
| 1394 | c jpl 2003 |
---|
| 1395 | c |
---|
| 1396 | c010(l) = 5.5e-12*exp(-2000./t(l)) |
---|
| 1397 | c |
---|
| 1398 | ccc c011: h + o2 + co2 -> ho2 + co2 |
---|
| 1399 | c |
---|
| 1400 | c jpl 2003 |
---|
| 1401 | c |
---|
| 1402 | ak0 = 2.5*5.7e-32*(t(l)/300.)**(-1.6) |
---|
| 1403 | ak1 = 7.5e-11*(t(l)/300.)**(0.0) |
---|
| 1404 | c |
---|
| 1405 | rate = (ak0*dens(l))/(1. + ak0*dens(l)/ak1) |
---|
| 1406 | xpo = 1./(1. + alog10((ak0*dens(l))/ak1)**2) |
---|
| 1407 | c011(l) = rate*0.6**xpo |
---|
| 1408 | c |
---|
| 1409 | ccc c012: o + h2o2 -> oh + ho2 |
---|
| 1410 | c |
---|
| 1411 | c jpl 2003 |
---|
| 1412 | c |
---|
| 1413 | c012(l) = 1.4e-12*exp(-2000./t(l)) |
---|
| 1414 | c |
---|
| 1415 | ccc c013: oh + oh -> h2o + o |
---|
| 1416 | c |
---|
| 1417 | c jpl 2003 |
---|
| 1418 | c |
---|
| 1419 | c013(l) = 4.2e-12*exp(-240./t(l)) |
---|
| 1420 | c |
---|
| 1421 | ccc c014: oh + o3 -> ho2 + o2 |
---|
| 1422 | c |
---|
| 1423 | c jpl 2003 |
---|
| 1424 | c |
---|
| 1425 | c014(l) = 1.7e-12*exp(-940./t(l)) |
---|
| 1426 | c |
---|
| 1427 | c jpl 2000 |
---|
| 1428 | c |
---|
| 1429 | c c014(l) = 1.5e-12*exp(-880./t(l)) |
---|
| 1430 | c |
---|
| 1431 | c nair et al., 1994 (jpl 1997) |
---|
| 1432 | c |
---|
| 1433 | c c014(l) = 1.6e-12*exp(-940./t(l)) |
---|
| 1434 | c |
---|
| 1435 | ccc c015: ho2 + o3 -> oh + o2 + o2 |
---|
| 1436 | c |
---|
| 1437 | c jpl 2003 |
---|
| 1438 | c |
---|
| 1439 | c015(l) = 1.0e-14*exp(-490./t(l)) |
---|
| 1440 | c |
---|
| 1441 | c jpl 2000 |
---|
| 1442 | c |
---|
| 1443 | c c015(l) = 2.0e-14*exp(-680./t(l)) |
---|
| 1444 | c |
---|
| 1445 | c nair et al., 1994 (jpl 1997) |
---|
| 1446 | c |
---|
| 1447 | c c015(l) = 1.1e-14*exp(-500./t(l)) |
---|
| 1448 | c |
---|
| 1449 | ccc c016: ho2 + ho2 + co2 -> h2o2 + o2 + co2 |
---|
| 1450 | c |
---|
| 1451 | c jpl 2003 |
---|
| 1452 | c |
---|
| 1453 | c016(l) = 2.5*1.7e-33 |
---|
| 1454 | $ *exp(1000./t(l))*dens(l) |
---|
| 1455 | c |
---|
| 1456 | ccc c017: oh + oh + co2 -> h2o2 + co2 |
---|
| 1457 | c |
---|
| 1458 | c jpl 2003 |
---|
| 1459 | c |
---|
| 1460 | ak0 = 2.5*6.9e-31*(t(l)/300.)**(-1.0) |
---|
| 1461 | ak1 = 2.6e-11*(t(l)/300.)**(0.0) |
---|
| 1462 | c |
---|
| 1463 | c jpl 1997 |
---|
| 1464 | c |
---|
| 1465 | c ak0 = 2.5*6.2e-31*(t(l)/300.)**(-1.0) |
---|
| 1466 | c ak1 = 2.6e-11*(t(l)/300.)**(0.0) |
---|
| 1467 | c |
---|
| 1468 | c nair et al., 1994 |
---|
| 1469 | c |
---|
| 1470 | c ak0 = 2.5*7.1e-31*(t(l)/300.)**(-0.8) |
---|
| 1471 | c ak1 = 1.5e-11*(t(l)/300.)**(0.0) |
---|
| 1472 | c |
---|
| 1473 | rate = (ak0*dens(l))/(1. + ak0*dens(l)/ak1) |
---|
| 1474 | xpo = 1./(1. + alog10((ak0*dens(l))/ak1)**2) |
---|
| 1475 | c017(l) = rate*0.6**xpo |
---|
| 1476 | c |
---|
| 1477 | ccc c018: h + h + co2 -> h2 + co2 |
---|
| 1478 | c |
---|
| 1479 | c baulch et al., 1992 |
---|
| 1480 | c |
---|
| 1481 | c018(l) = 2.5*8.85e-33*(t(l)/298.)**(-0.6)*dens(l) |
---|
| 1482 | c |
---|
| 1483 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1484 | c nitrogen compounds |
---|
| 1485 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1486 | c |
---|
| 1487 | ccc d001: no2 + o -> no + o2 |
---|
| 1488 | c |
---|
| 1489 | c jpl 2003 |
---|
| 1490 | c |
---|
| 1491 | d001(l) = 5.6e-12*exp(180./t(l)) |
---|
| 1492 | c |
---|
| 1493 | ccc d002: no + o3 -> no2 + o2 |
---|
| 1494 | c |
---|
| 1495 | c jpl 2003 |
---|
| 1496 | c |
---|
| 1497 | d002(l) = 3.0e-12*exp(-1500./t(l)) |
---|
| 1498 | c |
---|
| 1499 | ccc d003: no + ho2 -> no2 + oh |
---|
| 1500 | c |
---|
| 1501 | c jpl 2003 |
---|
| 1502 | c |
---|
| 1503 | d003(l) = 3.5e-12*exp(250./t(l)) |
---|
| 1504 | c |
---|
| 1505 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1506 | c carbon compounds |
---|
| 1507 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1508 | c |
---|
| 1509 | ccc e001: oh + co -> co2 + h |
---|
| 1510 | c |
---|
| 1511 | c jpl 2003 |
---|
| 1512 | c |
---|
| 1513 | c e001(l) = 1.5e-13*(1 + 0.6*press(l)/1013.) |
---|
| 1514 | c |
---|
| 1515 | c mccabe et al., grl, 28, 3135, 2001 |
---|
| 1516 | c |
---|
| 1517 | e001(l) = 1.57e-13 + 3.54e-33*dens(l) |
---|
| 1518 | c |
---|
| 1519 | ccc e002: o + co + m -> co2 + m |
---|
| 1520 | c |
---|
| 1521 | c tsang and hampson, 1986. |
---|
| 1522 | c |
---|
| 1523 | e002(l) = 2.5*6.5e-33*exp(-2184./t(l))*dens(l) |
---|
| 1524 | c |
---|
| 1525 | c baulch et al., butterworths, 1976. |
---|
| 1526 | c |
---|
| 1527 | c e002(l) = 1.6e-32*exp(-2184./t(l))*dens(l) |
---|
| 1528 | c |
---|
| 1529 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1530 | c heterogenous chemistry |
---|
| 1531 | cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 1532 | c |
---|
| 1533 | ccc h001: ho2 + ice surface -> (ho2)ads |
---|
| 1534 | ccc (ho2)ads + oh -> h2o + o2 |
---|
| 1535 | ccc total : oh + ho2 -> h2o + o2, first order |
---|
| 1536 | ccc K = h001*[ho2], |
---|
| 1537 | ccc h001 = (S*v*gamma1)/4. (s-1), v=100*sqrt(3kT/m_ho2) (cm s-1) |
---|
| 1538 | c |
---|
| 1539 | c cooper and abbatt, 1996 |
---|
| 1540 | c |
---|
| 1541 | h001(l) = icesurf(l) * 25.*sqrt(3*8.31*t(l)/33.) * 0.025 |
---|
| 1542 | c h001(l) = 0.0e0 |
---|
| 1543 | c |
---|
| 1544 | ccc h002: oh + ice surface -> (oh)ads |
---|
| 1545 | ccc (oh)ads + ho2 -> h2o + o2 |
---|
| 1546 | ccc total : oh + ho2 -> h2o + o2, first order |
---|
| 1547 | ccc K = h002*[oh], |
---|
| 1548 | ccc h002 = (S*v*gamma2)/4. (s-1), v=100*sqrt(3kT/m_oh) (cm s-1) |
---|
| 1549 | c |
---|
| 1550 | c cooper and abbatt, 1996 |
---|
| 1551 | c |
---|
| 1552 | h002(l) = icesurf(l) * 25.*sqrt(3*8.31*t(l)/17.) * 0.03 |
---|
| 1553 | c h002(l) = 0.0e0 |
---|
| 1554 | c |
---|
| 1555 | c write(*,*) "level ",l," / icesurf=",icesurf(l), |
---|
| 1556 | c $ " / 0.25vg1=",(25.*sqrt(3*8.31*t(l)/33.) * 0.025), |
---|
| 1557 | c $ " / 0.25vg2=",(25.*sqrt(3*8.31*t(l)/17.) * 0.03) |
---|
| 1558 | c |
---|
| 1559 | end do |
---|
| 1560 | c |
---|
| 1561 | return |
---|
| 1562 | end |
---|