[2162] | 1 | MODULE calchim_mod |
---|
| 2 | |
---|
| 3 | IMPLICIT NONE |
---|
| 4 | |
---|
[2164] | 5 | INTEGER, SAVE :: ichemistry ! compute chemistry every ichemistry physics step |
---|
[2162] | 6 | REAL,SAVE,ALLOCATABLE :: zdqchim(:,:,:) ! Tendancy on pq due to photochemistry |
---|
| 7 | REAL,SAVE,ALLOCATABLE :: zdqschim(:,:) ! Tendancy on qsurf due to photochemistry |
---|
| 8 | |
---|
| 9 | CONTAINS |
---|
| 10 | |
---|
[2158] | 11 | subroutine calchim(ngrid,nlayer,nq, & |
---|
[1495] | 12 | ptimestep,pplay,pplev,pt,pdt,dist_sol,mu0, & |
---|
| 13 | zzlev,zzlay,zday,pq,pdq,dqchim,dqschim,dqcloud, & |
---|
[2031] | 14 | dqscloud,tau,co2ice, & |
---|
[1495] | 15 | pu,pdu,pv,pdv,surfdust,surfice) |
---|
| 16 | |
---|
| 17 | use tracer_mod, only: igcm_co2, igcm_co, igcm_o, igcm_o1d, igcm_o2, & |
---|
| 18 | igcm_o3, igcm_h, igcm_h2, igcm_oh, igcm_ho2, & |
---|
| 19 | igcm_h2o2, igcm_ch4, igcm_n2, igcm_h2o_vap, & |
---|
| 20 | igcm_no, igcm_n, igcm_no2, igcm_n2d, & |
---|
| 21 | igcm_o2plus, igcm_co2plus, igcm_oplus, & |
---|
| 22 | igcm_coplus, igcm_cplus, igcm_nplus, & |
---|
| 23 | igcm_noplus, igcm_n2plus, igcm_hplus, & |
---|
| 24 | igcm_hco2plus, igcm_elec, mmol |
---|
| 25 | |
---|
| 26 | use conc_mod, only: mmean ! mean molecular mass of the atmosphere |
---|
[2162] | 27 | use comcstfi_h, only: pi |
---|
| 28 | use photolysis_mod, only: jonline, init_photolysis |
---|
[2158] | 29 | use iono_h, only: temp_elect |
---|
[1495] | 30 | |
---|
| 31 | implicit none |
---|
| 32 | |
---|
| 33 | !======================================================================= |
---|
| 34 | ! |
---|
| 35 | ! subject: |
---|
| 36 | ! -------- |
---|
| 37 | ! |
---|
| 38 | ! Prepare the call for the photochemical module, and send back the |
---|
| 39 | ! tendencies from photochemistry in the chemical species mass mixing ratios |
---|
| 40 | ! |
---|
| 41 | ! Author: Sebastien Lebonnois (08/11/2002) |
---|
| 42 | ! ------- |
---|
| 43 | ! update 12/06/2003 for water ice clouds and compatibility with dust |
---|
| 44 | ! update 07/2003 for coupling with thermosphere (Monica Angelats-i-Coll) |
---|
| 45 | ! update 03/05/2005 cosmetic changes (Franck Lefevre) |
---|
| 46 | ! update sept. 2008 identify tracers by their names (Ehouarn Millour) |
---|
| 47 | ! update 05/12/2011 synchronize with latest version of chemistry (Franck Lefevre) |
---|
| 48 | ! update 16/03/2012 optimization (Franck Lefevre) |
---|
| 49 | ! update 11/12/2013 optimization (Franck Lefevre) |
---|
| 50 | ! |
---|
| 51 | ! Arguments: |
---|
| 52 | ! ---------- |
---|
| 53 | ! |
---|
| 54 | ! Input: |
---|
| 55 | ! |
---|
[2031] | 56 | ! ptimestep timestep (s) |
---|
[1495] | 57 | ! pplay(ngrid,nlayer) Pressure at the middle of the layers (Pa) |
---|
| 58 | ! pplev(ngrid,nlayer+1) Intermediate pressure levels (Pa) |
---|
| 59 | ! pt(ngrid,nlayer) Temperature (K) |
---|
| 60 | ! pdt(ngrid,nlayer) Temperature tendency (K) |
---|
| 61 | ! pu(ngrid,nlayer) u component of the wind (ms-1) |
---|
| 62 | ! pdu(ngrid,nlayer) u component tendency (K) |
---|
| 63 | ! pv(ngrid,nlayer) v component of the wind (ms-1) |
---|
| 64 | ! pdv(ngrid,nlayer) v component tendency (K) |
---|
[2031] | 65 | ! dist_sol distance of the sun (AU) |
---|
| 66 | ! mu0(ngrid) cos of solar zenith angle (=1 when sun at zenith) |
---|
| 67 | ! pq(ngrid,nlayer,nq) advected fields, ie chemical species here |
---|
| 68 | ! pdq(ngrid,nlayer,nq) previous tendencies on pq |
---|
| 69 | ! tau(ngrid) dust optical depth |
---|
| 70 | ! co2ice(ngrid) co2 ice surface layer (kg.m-2) |
---|
[1495] | 71 | ! surfdust(ngrid,nlayer) dust surface area (m2/m3) |
---|
| 72 | ! surfice(ngrid,nlayer) ice surface area (m2/m3) |
---|
| 73 | ! |
---|
| 74 | ! Output: |
---|
| 75 | ! |
---|
[2031] | 76 | ! dqchim(ngrid,nlayer,nq) tendencies on pq due to chemistry |
---|
| 77 | ! dqschim(ngrid,nq) tendencies on qsurf |
---|
[1495] | 78 | ! |
---|
| 79 | !======================================================================= |
---|
| 80 | |
---|
[2162] | 81 | include "callkeys.h" |
---|
[1495] | 82 | |
---|
| 83 | ! input: |
---|
| 84 | |
---|
| 85 | integer,intent(in) :: ngrid ! number of atmospheric columns |
---|
| 86 | integer,intent(in) :: nlayer ! number of atmospheric layers |
---|
| 87 | integer,intent(in) :: nq ! number of tracers |
---|
| 88 | real :: ptimestep |
---|
| 89 | real :: pplay(ngrid,nlayer) ! pressure at the middle of the layers |
---|
| 90 | real :: zzlay(ngrid,nlayer) ! pressure at the middle of the layers |
---|
| 91 | real :: pplev(ngrid,nlayer+1) ! intermediate pressure levels |
---|
| 92 | real :: zzlev(ngrid,nlayer+1) ! altitude at layer boundaries |
---|
| 93 | real :: pt(ngrid,nlayer) ! temperature |
---|
| 94 | real :: pdt(ngrid,nlayer) ! temperature tendency |
---|
| 95 | real :: pu(ngrid,nlayer) ! u component of the wind (m.s-1) |
---|
| 96 | real :: pdu(ngrid,nlayer) ! u component tendency |
---|
| 97 | real :: pv(ngrid,nlayer) ! v component of the wind (m.s-1) |
---|
| 98 | real :: pdv(ngrid,nlayer) ! v component tendency |
---|
| 99 | real :: dist_sol ! distance of the sun (AU) |
---|
| 100 | real :: mu0(ngrid) ! cos of solar zenith angle (=1 when sun at zenith) |
---|
| 101 | real :: pq(ngrid,nlayer,nq) ! tracers mass mixing ratio |
---|
| 102 | real :: pdq(ngrid,nlayer,nq) ! previous tendencies |
---|
| 103 | real :: zday ! date (time since Ls=0, in martian days) |
---|
[2031] | 104 | real :: tau(ngrid) ! dust optical depth |
---|
[1495] | 105 | real :: co2ice(ngrid) ! co2 ice surface layer (kg.m-2) |
---|
| 106 | real :: surfdust(ngrid,nlayer) ! dust surface area (m2/m3) |
---|
[2031] | 107 | real :: surfice(ngrid,nlayer) ! ice surface area (m2/m3) |
---|
[1495] | 108 | |
---|
| 109 | ! output: |
---|
| 110 | |
---|
| 111 | real :: dqchim(ngrid,nlayer,nq) ! tendencies on pq due to chemistry |
---|
| 112 | real :: dqschim(ngrid,nq) ! tendencies on qsurf |
---|
| 113 | real :: dqcloud(ngrid,nlayer,nq) ! tendencies on pq due to condensation |
---|
| 114 | real :: dqscloud(ngrid,nq) ! tendencies on qsurf |
---|
| 115 | |
---|
| 116 | ! local variables: |
---|
| 117 | |
---|
| 118 | integer,save :: nbq ! number of tracers used in the chemistry |
---|
| 119 | integer,allocatable,save :: niq(:) ! array storing the indexes of the tracers |
---|
| 120 | integer :: iloc(1) ! index of major species |
---|
| 121 | integer :: ig,l,i,iq,iqmax |
---|
| 122 | integer :: foundswitch, lswitch |
---|
| 123 | integer,save :: chemthermod |
---|
| 124 | |
---|
| 125 | integer,save :: i_co2 = 0 |
---|
| 126 | integer,save :: i_co = 0 |
---|
| 127 | integer,save :: i_o = 0 |
---|
| 128 | integer,save :: i_o1d = 0 |
---|
| 129 | integer,save :: i_o2 = 0 |
---|
| 130 | integer,save :: i_o3 = 0 |
---|
| 131 | integer,save :: i_h = 0 |
---|
| 132 | integer,save :: i_h2 = 0 |
---|
| 133 | integer,save :: i_oh = 0 |
---|
| 134 | integer,save :: i_ho2 = 0 |
---|
| 135 | integer,save :: i_h2o2 = 0 |
---|
| 136 | integer,save :: i_ch4 = 0 |
---|
| 137 | integer,save :: i_n2 = 0 |
---|
| 138 | integer,save :: i_h2o = 0 |
---|
| 139 | integer,save :: i_n = 0 |
---|
| 140 | integer,save :: i_no = 0 |
---|
| 141 | integer,save :: i_no2 = 0 |
---|
| 142 | integer,save :: i_n2d = 0 |
---|
| 143 | integer,save :: i_co2plus=0 |
---|
| 144 | integer,save :: i_oplus=0 |
---|
| 145 | integer,save :: i_o2plus=0 |
---|
| 146 | integer,save :: i_coplus=0 |
---|
| 147 | integer,save :: i_cplus=0 |
---|
| 148 | integer,save :: i_nplus=0 |
---|
| 149 | integer,save :: i_noplus=0 |
---|
| 150 | integer,save :: i_n2plus=0 |
---|
| 151 | integer,save :: i_hplus=0 |
---|
| 152 | integer,save :: i_hco2plus=0 |
---|
| 153 | integer,save :: i_elec=0 |
---|
| 154 | |
---|
| 155 | integer :: ig_vl1 |
---|
| 156 | |
---|
| 157 | real :: latvl1, lonvl1 |
---|
| 158 | real :: zq(ngrid,nlayer,nq) ! pq+pdq*ptimestep before chemistry |
---|
| 159 | ! new mole fraction after |
---|
| 160 | real :: zt(ngrid,nlayer) ! temperature |
---|
| 161 | real :: zu(ngrid,nlayer) ! u component of the wind |
---|
| 162 | real :: zv(ngrid,nlayer) ! v component of the wind |
---|
[2027] | 163 | real :: taucol ! dust optical depth at the surface |
---|
[2030] | 164 | real :: kb ! boltzmann constant |
---|
[1495] | 165 | |
---|
| 166 | logical,save :: firstcall = .true. |
---|
| 167 | logical,save :: depos = .false. ! switch for dry deposition |
---|
| 168 | |
---|
| 169 | ! for each column of atmosphere: |
---|
| 170 | |
---|
| 171 | real :: zpress(nlayer) ! Pressure (mbar) |
---|
| 172 | real :: zdens(nlayer) ! Density (cm-3) |
---|
| 173 | real :: ztemp(nlayer) ! Temperature (K) |
---|
[2158] | 174 | real :: ztelec(nlayer) ! Electronic temperature (K) |
---|
[1495] | 175 | real :: zlocal(nlayer) ! Altitude (km) |
---|
| 176 | real :: zycol(nlayer,nq) ! Composition (mole fractions) |
---|
| 177 | real :: zmmean(nlayer) ! Mean molecular mass (g.mole-1) |
---|
| 178 | real :: szacol ! Solar zenith angle |
---|
| 179 | real :: surfice1d(nlayer) ! Ice surface area (cm2/cm3) |
---|
| 180 | real :: surfdust1d(nlayer) ! Dust surface area (cm2/cm3) |
---|
| 181 | real :: jo3(nlayer) ! Photodissociation rate O3->O1D (s-1) |
---|
[2030] | 182 | real :: jh2o(nlayer) ! Photodissociation rate H2O->H+OH (s-1) |
---|
[1888] | 183 | real :: em_no(nlayer) ! NO nightglow emission rate |
---|
| 184 | real :: em_o2(nlayer) ! O2 nightglow emission rate |
---|
[1495] | 185 | |
---|
| 186 | integer :: iter(nlayer) ! Number of chemical iterations |
---|
| 187 | ! within one physical timestep |
---|
| 188 | |
---|
| 189 | ! for output: |
---|
| 190 | |
---|
| 191 | logical :: output ! to issue calls to writediagfi and stats |
---|
| 192 | parameter (output = .true.) |
---|
| 193 | real :: jo3_3d(ngrid,nlayer) ! Photodissociation rate O3->O1D (s-1) |
---|
[2030] | 194 | real :: jh2o_3d(ngrid,nlayer) ! Photodissociation rate H2O->H+OH (s-1) |
---|
[1888] | 195 | real :: emission_no(ngrid,nlayer) !NO emission rate |
---|
| 196 | real :: emission_o2(ngrid,nlayer) !O2 emission rate |
---|
[1495] | 197 | real :: iter_3d(ngrid,nlayer) ! Number of chemical iterations |
---|
| 198 | ! within one physical timestep |
---|
| 199 | |
---|
[2158] | 200 | logical :: unichim ! only one, unified chemical scheme at all |
---|
| 201 | ! layers (default), or upper atmospheric |
---|
| 202 | ! scheme in the thermosphere? |
---|
| 203 | |
---|
[1495] | 204 | !======================================================================= |
---|
| 205 | ! initialization of the chemistry (first call only) |
---|
| 206 | !======================================================================= |
---|
| 207 | |
---|
[2158] | 208 | !call only unified chemistry (default), or also upper atmospheric model |
---|
| 209 | !unichim = .true. unified chemistry |
---|
| 210 | !unichim = .false. 2 different models |
---|
| 211 | unichim = .true. |
---|
| 212 | |
---|
[1495] | 213 | if (firstcall) then |
---|
| 214 | |
---|
| 215 | if (photochem) then |
---|
[2044] | 216 | if (jonline) then |
---|
| 217 | print*,'calchim: Read UV absorption cross-sections' |
---|
| 218 | call init_photolysis ! on-line photolysis |
---|
| 219 | else |
---|
| 220 | print*,'calchim: Read photolysis lookup table' |
---|
| 221 | call read_phototable ! off-line photolysis |
---|
| 222 | end if |
---|
[1495] | 223 | end if |
---|
[2158] | 224 | |
---|
| 225 | if(.not.unichim) then |
---|
| 226 | !Read reaction rates from external file if the upper atmospheric |
---|
| 227 | !chemistry is called |
---|
| 228 | call chemthermos_readini |
---|
| 229 | endif |
---|
| 230 | |
---|
[1495] | 231 | ! find index of chemical tracers to use |
---|
| 232 | allocate(niq(nq)) |
---|
| 233 | ! Listed here are all tracers that can go into photochemistry |
---|
| 234 | nbq = 0 ! to count number of tracers |
---|
[2158] | 235 | ! Species ALWAYS present if photochem=.T. |
---|
[1495] | 236 | i_co2 = igcm_co2 |
---|
| 237 | if (i_co2 == 0) then |
---|
| 238 | write(*,*) "calchim: Error; no CO2 tracer !!!" |
---|
| 239 | stop |
---|
| 240 | else |
---|
| 241 | nbq = nbq + 1 |
---|
| 242 | niq(nbq) = i_co2 |
---|
| 243 | end if |
---|
| 244 | i_co = igcm_co |
---|
| 245 | if (i_co == 0) then |
---|
| 246 | write(*,*) "calchim: Error; no CO tracer !!!" |
---|
| 247 | stop |
---|
| 248 | else |
---|
| 249 | nbq = nbq + 1 |
---|
| 250 | niq(nbq) = i_co |
---|
| 251 | end if |
---|
| 252 | i_o = igcm_o |
---|
| 253 | if (i_o == 0) then |
---|
| 254 | write(*,*) "calchim: Error; no O tracer !!!" |
---|
| 255 | stop |
---|
| 256 | else |
---|
| 257 | nbq = nbq + 1 |
---|
| 258 | niq(nbq) = i_o |
---|
| 259 | end if |
---|
| 260 | i_o1d = igcm_o1d |
---|
| 261 | if (i_o1d == 0) then |
---|
| 262 | write(*,*) "calchim: Error; no O1D tracer !!!" |
---|
| 263 | stop |
---|
| 264 | else |
---|
| 265 | nbq = nbq + 1 |
---|
| 266 | niq(nbq) = i_o1d |
---|
| 267 | end if |
---|
| 268 | i_o2 = igcm_o2 |
---|
| 269 | if (i_o2 == 0) then |
---|
| 270 | write(*,*) "calchim: Error; no O2 tracer !!!" |
---|
| 271 | stop |
---|
| 272 | else |
---|
| 273 | nbq = nbq + 1 |
---|
| 274 | niq(nbq) = i_o2 |
---|
| 275 | end if |
---|
| 276 | i_o3 = igcm_o3 |
---|
| 277 | if (i_o3 == 0) then |
---|
| 278 | write(*,*) "calchim: Error; no O3 tracer !!!" |
---|
| 279 | stop |
---|
| 280 | else |
---|
| 281 | nbq = nbq + 1 |
---|
| 282 | niq(nbq) = i_o3 |
---|
| 283 | end if |
---|
| 284 | i_h = igcm_h |
---|
| 285 | if (i_h == 0) then |
---|
| 286 | write(*,*) "calchim: Error; no H tracer !!!" |
---|
| 287 | stop |
---|
| 288 | else |
---|
| 289 | nbq = nbq + 1 |
---|
| 290 | niq(nbq) = i_h |
---|
| 291 | end if |
---|
| 292 | i_h2 = igcm_h2 |
---|
| 293 | if (i_h2 == 0) then |
---|
| 294 | write(*,*) "calchim: Error; no H2 tracer !!!" |
---|
| 295 | stop |
---|
| 296 | else |
---|
| 297 | nbq = nbq + 1 |
---|
| 298 | niq(nbq) = i_h2 |
---|
| 299 | end if |
---|
| 300 | i_oh = igcm_oh |
---|
| 301 | if (i_oh == 0) then |
---|
| 302 | write(*,*) "calchim: Error; no OH tracer !!!" |
---|
| 303 | stop |
---|
| 304 | else |
---|
| 305 | nbq = nbq + 1 |
---|
| 306 | niq(nbq) = i_oh |
---|
| 307 | end if |
---|
| 308 | i_ho2 = igcm_ho2 |
---|
| 309 | if (i_ho2 == 0) then |
---|
| 310 | write(*,*) "calchim: Error; no HO2 tracer !!!" |
---|
| 311 | stop |
---|
| 312 | else |
---|
| 313 | nbq = nbq + 1 |
---|
| 314 | niq(nbq) = i_ho2 |
---|
| 315 | end if |
---|
| 316 | i_h2o2 = igcm_h2o2 |
---|
| 317 | if (i_h2o2 == 0) then |
---|
| 318 | write(*,*) "calchim: Error; no H2O2 tracer !!!" |
---|
| 319 | stop |
---|
| 320 | else |
---|
| 321 | nbq = nbq + 1 |
---|
| 322 | niq(nbq) = i_h2o2 |
---|
| 323 | end if |
---|
| 324 | i_ch4 = igcm_ch4 |
---|
| 325 | if (i_ch4 == 0) then |
---|
| 326 | write(*,*) "calchim: Error; no CH4 tracer !!!" |
---|
| 327 | write(*,*) "CH4 will be ignored in the chemistry" |
---|
| 328 | else |
---|
| 329 | nbq = nbq + 1 |
---|
| 330 | niq(nbq) = i_ch4 |
---|
| 331 | end if |
---|
| 332 | i_n2 = igcm_n2 |
---|
| 333 | if (i_n2 == 0) then |
---|
| 334 | write(*,*) "calchim: Error; no N2 tracer !!!" |
---|
| 335 | stop |
---|
| 336 | else |
---|
| 337 | nbq = nbq + 1 |
---|
| 338 | niq(nbq) = i_n2 |
---|
| 339 | end if |
---|
| 340 | i_n = igcm_n |
---|
| 341 | if (i_n == 0) then |
---|
| 342 | if (photochem) then |
---|
| 343 | write(*,*) "calchim: Error; no N tracer !!!" |
---|
| 344 | stop |
---|
| 345 | end if |
---|
| 346 | else |
---|
| 347 | nbq = nbq + 1 |
---|
| 348 | niq(nbq) = i_n |
---|
| 349 | end if |
---|
| 350 | i_n2d = igcm_n2d |
---|
| 351 | if (i_n2d == 0) then |
---|
| 352 | if (photochem) then |
---|
| 353 | write(*,*) "calchim: Error; no N2D tracer !!!" |
---|
| 354 | stop |
---|
| 355 | end if |
---|
| 356 | else |
---|
| 357 | nbq = nbq + 1 |
---|
| 358 | niq(nbq) = i_n2d |
---|
| 359 | end if |
---|
| 360 | i_no = igcm_no |
---|
| 361 | if (i_no == 0) then |
---|
| 362 | if (photochem) then |
---|
| 363 | write(*,*) "calchim: Error; no NO tracer !!!" |
---|
| 364 | stop |
---|
| 365 | end if |
---|
| 366 | else |
---|
| 367 | nbq = nbq + 1 |
---|
| 368 | niq(nbq) = i_no |
---|
| 369 | end if |
---|
| 370 | i_no2 = igcm_no2 |
---|
| 371 | if (i_no2 == 0) then |
---|
| 372 | if (photochem) then |
---|
| 373 | write(*,*) "calchim: Error; no NO2 tracer !!!" |
---|
| 374 | stop |
---|
| 375 | end if |
---|
| 376 | else |
---|
| 377 | nbq = nbq + 1 |
---|
| 378 | niq(nbq) = i_no2 |
---|
| 379 | end if |
---|
| 380 | i_h2o = igcm_h2o_vap |
---|
| 381 | if (i_h2o == 0) then |
---|
| 382 | write(*,*) "calchim: Error; no water vapor tracer !!!" |
---|
| 383 | stop |
---|
| 384 | else |
---|
| 385 | nbq = nbq + 1 |
---|
| 386 | niq(nbq) = i_h2o |
---|
| 387 | end if |
---|
[2158] | 388 | i_o2plus = igcm_o2plus |
---|
| 389 | if (i_o2plus == 0) then |
---|
| 390 | write(*,*) "calchim: Error, no O2+ tracer !!!" |
---|
| 391 | stop |
---|
| 392 | else |
---|
| 393 | nbq = nbq + 1 |
---|
| 394 | niq(nbq) = i_o2plus |
---|
| 395 | end if |
---|
| 396 | i_co2plus = igcm_co2plus |
---|
| 397 | if (i_co2plus == 0) then |
---|
| 398 | write(*,*) "calchim: Error, no CO2+ tracer !!!" |
---|
| 399 | stop |
---|
| 400 | else |
---|
| 401 | nbq = nbq + 1 |
---|
| 402 | niq(nbq) = i_co2plus |
---|
| 403 | end if |
---|
| 404 | i_oplus=igcm_oplus |
---|
| 405 | if (i_oplus == 0) then |
---|
| 406 | write(*,*) "calchim: Error, no O+ tracer !!!" |
---|
| 407 | stop |
---|
| 408 | else |
---|
| 409 | nbq = nbq + 1 |
---|
| 410 | niq(nbq) = i_oplus |
---|
| 411 | end if |
---|
| 412 | i_noplus=igcm_noplus |
---|
| 413 | if (i_noplus == 0) then |
---|
| 414 | write(*,*) "calchim: Error, no NO+ tracer !!!" |
---|
| 415 | stop |
---|
| 416 | else |
---|
| 417 | nbq = nbq + 1 |
---|
| 418 | niq(nbq) = i_noplus |
---|
| 419 | end if |
---|
| 420 | i_coplus=igcm_coplus |
---|
| 421 | if (i_coplus == 0) then |
---|
| 422 | write(*,*) "calchim: Error, no CO+ tracer !!!" |
---|
| 423 | stop |
---|
| 424 | else |
---|
| 425 | nbq = nbq + 1 |
---|
| 426 | niq(nbq) = i_coplus |
---|
| 427 | end if |
---|
| 428 | i_cplus=igcm_cplus |
---|
| 429 | if (i_cplus == 0) then |
---|
| 430 | write(*,*) "calchim: Error, no C+ tracer !!!" |
---|
| 431 | stop |
---|
| 432 | else |
---|
| 433 | nbq = nbq + 1 |
---|
| 434 | niq(nbq) = i_cplus |
---|
| 435 | end if |
---|
| 436 | i_n2plus=igcm_n2plus |
---|
| 437 | if (i_n2plus == 0) then |
---|
| 438 | write(*,*) "calchim: Error, no N2+ tracer !!!" |
---|
| 439 | stop |
---|
| 440 | else |
---|
| 441 | nbq = nbq + 1 |
---|
| 442 | niq(nbq) = i_n2plus |
---|
| 443 | end if |
---|
| 444 | i_nplus=igcm_nplus |
---|
| 445 | if (i_nplus == 0) then |
---|
| 446 | write(*,*) "calchim: Error, no N+ tracer !!!" |
---|
| 447 | stop |
---|
| 448 | else |
---|
| 449 | nbq = nbq + 1 |
---|
| 450 | niq(nbq) = i_nplus |
---|
| 451 | end if |
---|
| 452 | i_hplus=igcm_hplus |
---|
| 453 | if (i_hplus == 0) then |
---|
| 454 | write(*,*) "calchim: Error, no H+ tracer !!!" |
---|
| 455 | stop |
---|
| 456 | else |
---|
| 457 | nbq = nbq + 1 |
---|
| 458 | niq(nbq) = i_hplus |
---|
| 459 | end if |
---|
| 460 | i_hco2plus=igcm_hco2plus |
---|
| 461 | if (i_hco2plus == 0) then |
---|
| 462 | write(*,*) "calchim: Error, no HCO2+ tracer !!!" |
---|
| 463 | stop |
---|
| 464 | else |
---|
| 465 | nbq = nbq + 1 |
---|
| 466 | niq(nbq) = i_hco2plus |
---|
| 467 | end if |
---|
| 468 | i_elec = igcm_elec |
---|
| 469 | if (i_elec == 0) then |
---|
| 470 | write(*,*) "calchim: Error, no e- tracer !!!" |
---|
| 471 | stop |
---|
| 472 | else |
---|
| 473 | nbq = nbq + 1 |
---|
| 474 | niq(nbq) = i_elec |
---|
| 475 | end if |
---|
| 476 | |
---|
[1495] | 477 | write(*,*) 'calchim: found nbq = ',nbq,' tracers' |
---|
| 478 | |
---|
| 479 | firstcall = .false. |
---|
| 480 | end if ! if (firstcall) |
---|
| 481 | |
---|
| 482 | ! Initializations |
---|
| 483 | |
---|
| 484 | zycol(:,:) = 0. |
---|
| 485 | dqchim(:,:,:) = 0. |
---|
| 486 | dqschim(:,:) = 0. |
---|
| 487 | |
---|
[2030] | 488 | kb = 1.3806e-23 |
---|
| 489 | |
---|
[1495] | 490 | ! latvl1= 22.27 |
---|
| 491 | ! lonvl1= -47.94 |
---|
[2029] | 492 | ! ig_vl1= 1+ int( (1.5-(latvl1-90.)*48./180.) -2 )*64. + & |
---|
| 493 | ! int(1.5+(lonvl1+180)*64./360.) |
---|
[1495] | 494 | |
---|
| 495 | !======================================================================= |
---|
| 496 | ! loop over grid |
---|
| 497 | !======================================================================= |
---|
| 498 | |
---|
| 499 | do ig = 1,ngrid |
---|
| 500 | |
---|
| 501 | foundswitch = 0 |
---|
| 502 | do l = 1,nlayer |
---|
| 503 | do i = 1,nbq |
---|
| 504 | iq = niq(i) ! get tracer index |
---|
| 505 | zq(ig,l,iq) = pq(ig,l,iq) + pdq(ig,l,iq)*ptimestep |
---|
| 506 | zycol(l,iq) = zq(ig,l,iq)*mmean(ig,l)/mmol(iq) |
---|
| 507 | end do |
---|
| 508 | zt(ig,l) = pt(ig,l) + pdt(ig,l)*ptimestep |
---|
| 509 | zu(ig,l) = pu(ig,l) + pdu(ig,l)*ptimestep |
---|
| 510 | zv(ig,l) = pv(ig,l) + pdv(ig,l)*ptimestep |
---|
| 511 | zpress(l) = pplay(ig,l)/100. |
---|
| 512 | ztemp(l) = zt(ig,l) |
---|
| 513 | zdens(l) = zpress(l)/(kb*1.e4*ztemp(l)) |
---|
| 514 | zlocal(l) = zzlay(ig,l)/1000. |
---|
| 515 | zmmean(l) = mmean(ig,l) |
---|
[2158] | 516 | ztelec(l) = temp_elect(zlocal(l),ztemp(l),1) |
---|
| 517 | !Electronic temperature. Index 1 -> Viking; Index 2-> MAVEN |
---|
[1495] | 518 | |
---|
| 519 | ! surfdust1d and surfice1d: conversion from m2/m3 to cm2/cm3 |
---|
| 520 | |
---|
| 521 | surfdust1d(l) = surfdust(ig,l)*1.e-2 |
---|
| 522 | surfice1d(l) = surfice(ig,l)*1.e-2 |
---|
| 523 | |
---|
[2158] | 524 | ! search for switch index between regions |
---|
| 525 | if(unichim) then |
---|
| 526 | lswitch=nlayer+1 |
---|
| 527 | else |
---|
[1495] | 528 | if (foundswitch == 0 .and. pplay(ig,l) < 1.e-1) then |
---|
| 529 | lswitch = l |
---|
| 530 | foundswitch = 1 |
---|
| 531 | end if |
---|
[2158] | 532 | endif |
---|
[1495] | 533 | end do ! of do l=1,nlayer |
---|
| 534 | |
---|
| 535 | szacol = acos(mu0(ig))*180./pi |
---|
[2031] | 536 | taucol = tau(ig) |
---|
[1495] | 537 | |
---|
| 538 | !======================================================================= |
---|
| 539 | ! call chemical subroutines |
---|
| 540 | !======================================================================= |
---|
| 541 | |
---|
[2158] | 542 | ! chemistry: only one scheme at all layers |
---|
| 543 | |
---|
[1495] | 544 | if (photochem) then |
---|
[2007] | 545 | call photochemistry(nlayer,nq, & |
---|
[2027] | 546 | ig,lswitch,zycol,szacol,ptimestep, & |
---|
[2158] | 547 | zpress,zlocal,ztemp,ztelec,zdens,zmmean, & |
---|
| 548 | dist_sol,zday,surfdust1d,surfice1d, & |
---|
| 549 | jo3,jh2o,taucol,iter) |
---|
[1495] | 550 | |
---|
| 551 | ! ozone photolysis, for output |
---|
| 552 | |
---|
| 553 | do l = 1,nlayer |
---|
| 554 | jo3_3d(ig,l) = jo3(l) |
---|
[2030] | 555 | jh2o_3d(ig,l) = jh2o(l) |
---|
[1495] | 556 | iter_3d(ig,l) = iter(l) |
---|
| 557 | end do |
---|
| 558 | ! condensation of h2o2 |
---|
| 559 | |
---|
| 560 | call perosat(ngrid, nlayer, nq, & |
---|
| 561 | ig,ptimestep,pplev,pplay, & |
---|
| 562 | ztemp,zycol,dqcloud,dqscloud) |
---|
| 563 | |
---|
[2158] | 564 | if(.not.unichim) then |
---|
| 565 | chemthermod=3 !C/O/H/N/ions |
---|
| 566 | call chemthermos(ig,nlayer,lswitch,chemthermod,zycol,ztemp,& |
---|
| 567 | zdens,zpress,zlocal,szacol,ptimestep,zday,& |
---|
| 568 | em_no,em_o2) |
---|
| 569 | do l=1,nlayer |
---|
| 570 | emission_no(ig,l)=em_no(l) |
---|
| 571 | emission_o2(ig,l)=em_o2(l) |
---|
| 572 | enddo |
---|
| 573 | end if |
---|
[1495] | 574 | end if |
---|
| 575 | |
---|
| 576 | ! dry deposition |
---|
| 577 | |
---|
| 578 | if (depos) then |
---|
| 579 | call deposition(ngrid, nlayer, nq, & |
---|
| 580 | ig, ig_vl1, pplay, pplev, zzlay, zzlev, & |
---|
| 581 | zu, zv, zt, zycol, ptimestep, co2ice) |
---|
| 582 | end if |
---|
| 583 | |
---|
| 584 | !======================================================================= |
---|
| 585 | ! tendencies |
---|
| 586 | !======================================================================= |
---|
| 587 | |
---|
| 588 | ! index of the most abundant species at each level |
---|
| 589 | |
---|
| 590 | ! major(:) = maxloc(zycol, dim = 2) |
---|
| 591 | |
---|
| 592 | ! tendency for the most abundant species = - sum of others |
---|
| 593 | |
---|
| 594 | do l = 1,nlayer |
---|
| 595 | iloc=maxloc(zycol(l,:)) |
---|
| 596 | iqmax=iloc(1) |
---|
| 597 | do i = 1,nbq |
---|
| 598 | iq = niq(i) ! get tracer index |
---|
| 599 | if (iq /= iqmax) then |
---|
| 600 | dqchim(ig,l,iq) = (zycol(l,iq)*mmol(iq)/mmean(ig,l) & |
---|
| 601 | - zq(ig,l,iq))/ptimestep |
---|
| 602 | dqchim(ig,l,iqmax) = dqchim(ig,l,iqmax) & |
---|
| 603 | - dqchim(ig,l,iq) |
---|
| 604 | end if |
---|
| 605 | end do |
---|
| 606 | end do ! of do l = 1,nlayer |
---|
| 607 | |
---|
| 608 | !======================================================================= |
---|
| 609 | ! end of loop over grid |
---|
| 610 | !======================================================================= |
---|
| 611 | |
---|
| 612 | end do ! of do ig=1,ngrid |
---|
| 613 | |
---|
| 614 | !======================================================================= |
---|
| 615 | ! write outputs |
---|
| 616 | !======================================================================= |
---|
| 617 | |
---|
| 618 | ! value of parameter 'output' to trigger writting of outputs |
---|
| 619 | ! is set above at the declaration of the variable. |
---|
| 620 | |
---|
| 621 | if (photochem .and. output) then |
---|
| 622 | if (ngrid > 1) then |
---|
| 623 | call writediagfi(ngrid,'jo3','j o3->o1d', & |
---|
| 624 | 's-1',3,jo3_3d(1,1)) |
---|
[2030] | 625 | call writediagfi(ngrid,'jh2o','jh2o', & |
---|
| 626 | 's-1',3,jh2o_3d(1,1)) |
---|
[1495] | 627 | call writediagfi(ngrid,'iter','iterations', & |
---|
| 628 | ' ',3,iter_3d(1,1)) |
---|
[2158] | 629 | if(.not.unichim) then |
---|
| 630 | call writediagfi(ngrid,'emission_no', & |
---|
| 631 | 'NO nightglow emission rate','cm-3 s-1',3,emission_no) |
---|
| 632 | call writediagfi(ngrid,'emission_o2', & |
---|
| 633 | 'O2 nightglow emission rate','cm-3 s-1',3,emission_o2) |
---|
| 634 | endif |
---|
[2042] | 635 | |
---|
[1495] | 636 | if (callstats) then |
---|
| 637 | call wstats(ngrid,'jo3','j o3->o1d', & |
---|
| 638 | 's-1',3,jo3_3d(1,1)) |
---|
[1888] | 639 | call wstats(ngrid,'emission_no', & |
---|
| 640 | 'NO nightglow emission rate','cm-3 s-1',3,emission_no) |
---|
| 641 | call wstats(ngrid,'emission_o2', & |
---|
| 642 | 'O2 nightglow emission rate','cm-3 s-1',3,emission_o2) |
---|
[1495] | 643 | call wstats(ngrid,'mmean','mean molecular mass', & |
---|
| 644 | 'g.mole-1',3,mmean(1,1)) |
---|
| 645 | endif |
---|
| 646 | end if ! of if (ngrid.gt.1) |
---|
| 647 | end if ! of if (output) |
---|
| 648 | |
---|
[2162] | 649 | end subroutine calchim |
---|
| 650 | |
---|
| 651 | |
---|
| 652 | subroutine ini_calchim_mod(ngrid,nlayer,nq) |
---|
| 653 | |
---|
| 654 | implicit none |
---|
| 655 | |
---|
| 656 | integer,intent(in) :: ngrid ! number of atmospheric columns |
---|
| 657 | integer,intent(in) :: nlayer ! number of atmospheric layers |
---|
| 658 | integer,intent(in) :: nq ! number of tracers |
---|
| 659 | |
---|
| 660 | allocate(zdqchim(ngrid,nlayer,nq)) |
---|
| 661 | zdqchim(:,:,:)=0 |
---|
| 662 | allocate(zdqschim(ngrid,nq)) |
---|
| 663 | zdqschim(:,:)=0 |
---|
| 664 | |
---|
| 665 | end subroutine ini_calchim_mod |
---|
| 666 | |
---|
| 667 | |
---|
| 668 | subroutine end_calchim_mod |
---|
| 669 | |
---|
| 670 | implicit none |
---|
| 671 | |
---|
| 672 | if (allocated(zdqchim)) deallocate(zdqchim) |
---|
| 673 | if (allocated(zdqschim)) deallocate(zdqschim) |
---|
| 674 | |
---|
| 675 | end subroutine end_calchim_mod |
---|
| 676 | |
---|
| 677 | END MODULE calchim_mod |
---|
| 678 | |
---|