| 1 | module physiq_mod |
|---|
| 2 | |
|---|
| 3 | implicit none |
|---|
| 4 | |
|---|
| 5 | contains |
|---|
| 6 | |
|---|
| 7 | subroutine physiq(ngrid,nlayer,nq, & |
|---|
| 8 | nametrac, & |
|---|
| 9 | firstcall,lastcall, & |
|---|
| 10 | pday,ptime,ptimestep, & |
|---|
| 11 | pplev,pplay,pphi, & |
|---|
| 12 | pu,pv,pt,pq, & |
|---|
| 13 | flxw, & |
|---|
| 14 | pdu,pdv,pdt,pdq,pdpsrf) |
|---|
| 15 | |
|---|
| 16 | use radinc_h, only : L_NSPECTI,L_NSPECTV |
|---|
| 17 | use radcommon_h, only: sigma, gzlat, gzlat_ig, Cmk, grav, BWNV |
|---|
| 18 | use surfdat_h, only: phisfi, zmea, zstd, zsig, zgam, zthe |
|---|
| 19 | use comchem_h, only: nkim, cnames, nlaykim_up, ykim_up, ykim_tot, botCH4 |
|---|
| 20 | use comdiurn_h, only: coslat, sinlat, coslon, sinlon |
|---|
| 21 | use comsaison_h, only: mu0, fract, dist_star, declin, right_ascen |
|---|
| 22 | use comsoil_h, only: nsoilmx, layer, mlayer, inertiedat |
|---|
| 23 | use datafile_mod, only: datadir, corrkdir, banddir |
|---|
| 24 | use geometry_mod, only: latitude, longitude, cell_area |
|---|
| 25 | USE comgeomfi_h, only: totarea, totarea_planet |
|---|
| 26 | USE tracer_h |
|---|
| 27 | use time_phylmdz_mod, only: ecritphy, iphysiq, nday |
|---|
| 28 | use phyetat0_mod, only: phyetat0 |
|---|
| 29 | use phyredem, only: physdem0, physdem1 |
|---|
| 30 | use planetwide_mod, only: planetwide_minval,planetwide_maxval,planetwide_sumval |
|---|
| 31 | use mod_phys_lmdz_para, only : is_master |
|---|
| 32 | use planete_mod, only: apoastr, periastr, year_day, peri_day, & |
|---|
| 33 | obliquit, nres, z0 |
|---|
| 34 | use comcstfi_mod, only: pi, g, rcp, r, rad, mugaz, cpp |
|---|
| 35 | use time_phylmdz_mod, only: daysec |
|---|
| 36 | use logic_mod, only: moyzon_ch |
|---|
| 37 | use moyzon_mod, only: zphibar, zphisbar, zplevbar, zplaybar, & |
|---|
| 38 | zzlevbar, zzlaybar, ztfibar, zqfibar |
|---|
| 39 | use callkeys_mod |
|---|
| 40 | use vertical_layers_mod, only: presnivs, pseudoalt |
|---|
| 41 | use ioipsl_getin_p_mod, only: getin_p |
|---|
| 42 | use mod_phys_lmdz_omp_data, ONLY: is_omp_master |
|---|
| 43 | #ifdef CPP_XIOS |
|---|
| 44 | use xios_output_mod, only: initialize_xios_output, & |
|---|
| 45 | update_xios_timestep, & |
|---|
| 46 | send_xios_field |
|---|
| 47 | use wxios, only: wxios_context_init, xios_context_finalize |
|---|
| 48 | #endif |
|---|
| 49 | use MMP_OPTICS |
|---|
| 50 | use muphy_diag |
|---|
| 51 | implicit none |
|---|
| 52 | |
|---|
| 53 | |
|---|
| 54 | !================================================================== |
|---|
| 55 | ! |
|---|
| 56 | ! Purpose |
|---|
| 57 | ! ------- |
|---|
| 58 | ! Central subroutine for all the physics parameterisations in the |
|---|
| 59 | ! universal model. Originally adapted from the Mars LMDZ model. |
|---|
| 60 | ! |
|---|
| 61 | ! The model can be run without or with tracer transport |
|---|
| 62 | ! depending on the value of "tracer" in file "callphys.def". |
|---|
| 63 | ! |
|---|
| 64 | ! |
|---|
| 65 | ! It includes: |
|---|
| 66 | ! |
|---|
| 67 | ! I. Initialization : |
|---|
| 68 | ! I.1 Firstcall initializations. |
|---|
| 69 | ! I.2 Initialization for every call to physiq. |
|---|
| 70 | ! |
|---|
| 71 | ! II. Compute radiative transfer tendencies (longwave and shortwave) : |
|---|
| 72 | ! II.a Option 1 : Call correlated-k radiative transfer scheme. |
|---|
| 73 | ! II.b Option 2 : Call Newtonian cooling scheme. |
|---|
| 74 | ! II.c Option 3 : Atmosphere has no radiative effect. |
|---|
| 75 | ! |
|---|
| 76 | ! III. Vertical diffusion (turbulent mixing) : |
|---|
| 77 | ! |
|---|
| 78 | ! IV. Dry Convective adjustment : |
|---|
| 79 | ! |
|---|
| 80 | ! V. Tracers |
|---|
| 81 | ! V.1. Microphysics |
|---|
| 82 | ! V.2. Chemistry |
|---|
| 83 | ! V.3. Updates (pressure variations, surface budget). |
|---|
| 84 | ! V.4. Surface Tracer Update. |
|---|
| 85 | ! |
|---|
| 86 | ! VI. Surface and sub-surface soil temperature. |
|---|
| 87 | ! |
|---|
| 88 | ! VII. Perform diagnostics and write output files. |
|---|
| 89 | ! |
|---|
| 90 | ! |
|---|
| 91 | ! arguments |
|---|
| 92 | ! --------- |
|---|
| 93 | ! |
|---|
| 94 | ! INPUT |
|---|
| 95 | ! ----- |
|---|
| 96 | ! |
|---|
| 97 | ! ngrid Size of the horizontal grid. |
|---|
| 98 | ! nlayer Number of vertical layers. |
|---|
| 99 | ! nq Number of advected fields. |
|---|
| 100 | ! nametrac Name of corresponding advected fields. |
|---|
| 101 | ! |
|---|
| 102 | ! firstcall True at the first call. |
|---|
| 103 | ! lastcall True at the last call. |
|---|
| 104 | ! |
|---|
| 105 | ! pday Number of days counted from the North. Spring equinoxe. |
|---|
| 106 | ! ptime Universal time (0<ptime<1): ptime=0.5 at 12:00 UT. |
|---|
| 107 | ! ptimestep timestep (s). |
|---|
| 108 | ! |
|---|
| 109 | ! pplay(ngrid,nlayer) Pressure at the middle of the layers (Pa). |
|---|
| 110 | ! pplev(ngrid,nlayer+1) Intermediate pressure levels (Pa). |
|---|
| 111 | ! pphi(ngrid,nlayer) Geopotential at the middle of the layers (m2.s-2). |
|---|
| 112 | ! |
|---|
| 113 | ! pu(ngrid,nlayer) u, zonal component of the wind (ms-1). |
|---|
| 114 | ! pv(ngrid,nlayer) v, meridional component of the wind (ms-1). |
|---|
| 115 | ! |
|---|
| 116 | ! pt(ngrid,nlayer) Temperature (K). |
|---|
| 117 | ! |
|---|
| 118 | ! pq(ngrid,nlayer,nq) Advected fields. |
|---|
| 119 | ! |
|---|
| 120 | ! pudyn(ngrid,nlayer) \ |
|---|
| 121 | ! pvdyn(ngrid,nlayer) \ Dynamical temporal derivative for the |
|---|
| 122 | ! ptdyn(ngrid,nlayer) / corresponding variables. |
|---|
| 123 | ! pqdyn(ngrid,nlayer,nq) / |
|---|
| 124 | ! flxw(ngrid,nlayer) vertical mass flux (kg/s) at layer lower boundary |
|---|
| 125 | ! |
|---|
| 126 | ! OUTPUT |
|---|
| 127 | ! ------ |
|---|
| 128 | ! |
|---|
| 129 | ! pdu(ngrid,nlayer) \ |
|---|
| 130 | ! pdv(ngrid,nlayer) \ Temporal derivative of the corresponding |
|---|
| 131 | ! pdt(ngrid,nlayer) / variables due to physical processes. |
|---|
| 132 | ! pdq(ngrid,nlayer) / |
|---|
| 133 | ! pdpsrf(ngrid) / |
|---|
| 134 | ! |
|---|
| 135 | ! |
|---|
| 136 | ! Authors |
|---|
| 137 | ! ------- |
|---|
| 138 | ! Frederic Hourdin 15/10/93 |
|---|
| 139 | ! Francois Forget 1994 |
|---|
| 140 | ! Christophe Hourdin 02/1997 |
|---|
| 141 | ! Subroutine completely rewritten by F. Forget (01/2000) |
|---|
| 142 | ! Water ice clouds: Franck Montmessin (update 06/2003) |
|---|
| 143 | ! Radiatively active tracers: J.-B. Madeleine (10/2008-06/2009) |
|---|
| 144 | ! New correlated-k radiative scheme: R. Wordsworth (2009) |
|---|
| 145 | ! Many specifically Martian subroutines removed: R. Wordsworth (2009) |
|---|
| 146 | ! Improved water cycle: R. Wordsworth / B. Charnay (2010) |
|---|
| 147 | ! To F90: R. Wordsworth (2010) |
|---|
| 148 | ! New turbulent diffusion scheme: J. Leconte (2012) |
|---|
| 149 | ! Loops converted to F90 matrix format: J. Leconte (2012) |
|---|
| 150 | ! No more ngridmx/nqmx, F90 commons and adaptation to parallel: A. Spiga (2012) |
|---|
| 151 | ! Purge of the code : M. Turbet (2015) |
|---|
| 152 | ! Fork for Titan : J. Vatant d'Ollone (2017) |
|---|
| 153 | ! + clean of all too-generic (ocean, water, co2 ...) routines |
|---|
| 154 | ! + Titan's chemistry |
|---|
| 155 | ! Microphysical moment model - J.Burgalat / J.Vatant d'Ollone (2017-2018) |
|---|
| 156 | !============================================================================================ |
|---|
| 157 | |
|---|
| 158 | |
|---|
| 159 | ! 0. Declarations : |
|---|
| 160 | ! ------------------ |
|---|
| 161 | |
|---|
| 162 | include "netcdf.inc" |
|---|
| 163 | |
|---|
| 164 | ! Arguments : |
|---|
| 165 | ! ----------- |
|---|
| 166 | |
|---|
| 167 | ! INPUTS: |
|---|
| 168 | ! ------- |
|---|
| 169 | |
|---|
| 170 | |
|---|
| 171 | integer,intent(in) :: ngrid ! Number of atmospheric columns. |
|---|
| 172 | integer,intent(in) :: nlayer ! Number of atmospheric layers. |
|---|
| 173 | integer,intent(in) :: nq ! Number of tracers. |
|---|
| 174 | character*30,intent(in) :: nametrac(nq) ! Names of the tracers taken from dynamics. |
|---|
| 175 | |
|---|
| 176 | logical,intent(in) :: firstcall ! Signals first call to physics. |
|---|
| 177 | logical,intent(in) :: lastcall ! Signals last call to physics. |
|---|
| 178 | |
|---|
| 179 | real,intent(in) :: pday ! Number of elapsed sols since reference Ls=0. |
|---|
| 180 | real,intent(in) :: ptime ! "Universal time", given as fraction of sol (e.g.: 0.5 for noon). |
|---|
| 181 | real,intent(in) :: ptimestep ! Physics timestep (s). |
|---|
| 182 | real,intent(in) :: pplev(ngrid,nlayer+1) ! Inter-layer pressure (Pa). |
|---|
| 183 | real,intent(in) :: pplay(ngrid,nlayer) ! Mid-layer pressure (Pa). |
|---|
| 184 | real,intent(in) :: pphi(ngrid,nlayer) ! Geopotential at mid-layer (m2s-2). |
|---|
| 185 | real,intent(in) :: pu(ngrid,nlayer) ! Zonal wind component (m/s). |
|---|
| 186 | real,intent(in) :: pv(ngrid,nlayer) ! Meridional wind component (m/s). |
|---|
| 187 | real,intent(in) :: pt(ngrid,nlayer) ! Temperature (K). |
|---|
| 188 | real,intent(in) :: pq(ngrid,nlayer,nq) ! Tracers (kg/kg_of_air). |
|---|
| 189 | real,intent(in) :: flxw(ngrid,nlayer) ! Vertical mass flux (ks/s) at lower boundary of layer |
|---|
| 190 | |
|---|
| 191 | ! OUTPUTS: |
|---|
| 192 | ! -------- |
|---|
| 193 | |
|---|
| 194 | ! Physical tendencies : |
|---|
| 195 | |
|---|
| 196 | real,intent(out) :: pdu(ngrid,nlayer) ! Zonal wind tendencies (m/s/s). |
|---|
| 197 | real,intent(out) :: pdv(ngrid,nlayer) ! Meridional wind tendencies (m/s/s). |
|---|
| 198 | real,intent(out) :: pdt(ngrid,nlayer) ! Temperature tendencies (K/s). |
|---|
| 199 | real,intent(out) :: pdq(ngrid,nlayer,nq) ! Tracer tendencies (kg/kg_of_air/s). |
|---|
| 200 | real,intent(out) :: pdpsrf(ngrid) ! Surface pressure tendency (Pa/s). |
|---|
| 201 | |
|---|
| 202 | ! Local saved variables: |
|---|
| 203 | ! ---------------------- |
|---|
| 204 | |
|---|
| 205 | integer,save :: day_ini ! Initial date of the run (sol since Ls=0). |
|---|
| 206 | integer,save :: icount ! Counter of calls to physiq during the run. |
|---|
| 207 | !$OMP THREADPRIVATE(day_ini,icount) |
|---|
| 208 | |
|---|
| 209 | real, dimension(:),allocatable,save :: tsurf ! Surface temperature (K). |
|---|
| 210 | real, dimension(:,:),allocatable,save :: tsoil ! Sub-surface temperatures (K). |
|---|
| 211 | real, dimension(:,:),allocatable,save :: albedo ! Surface Spectral albedo. By MT2015. |
|---|
| 212 | real, dimension(:),allocatable,save :: albedo_equivalent ! Spectral Mean albedo. |
|---|
| 213 | |
|---|
| 214 | !$OMP THREADPRIVATE(tsurf,tsoil,albedo,albedo_equivalent) |
|---|
| 215 | |
|---|
| 216 | real,dimension(:),allocatable,save :: albedo_bareground ! Bare Ground Albedo. By MT 2015. |
|---|
| 217 | |
|---|
| 218 | !$OMP THREADPRIVATE(albedo_bareground) |
|---|
| 219 | |
|---|
| 220 | real,dimension(:),allocatable,save :: emis ! Thermal IR surface emissivity. |
|---|
| 221 | real,dimension(:,:),allocatable,save :: dtrad ! Net atmospheric radiative heating rate (K.s-1). |
|---|
| 222 | real,dimension(:),allocatable,save :: fluxrad_sky ! Radiative flux from sky absorbed by surface (W.m-2). |
|---|
| 223 | real,dimension(:),allocatable,save :: fluxrad ! Net radiative surface flux (W.m-2). |
|---|
| 224 | real,dimension(:),allocatable,save :: capcal ! Surface heat capacity (J m-2 K-1). |
|---|
| 225 | real,dimension(:),allocatable,save :: fluxgrd ! Surface conduction flux (W.m-2). |
|---|
| 226 | real,dimension(:,:),allocatable,save :: qsurf ! Tracer on surface (e.g. kg.m-2). |
|---|
| 227 | real,dimension(:,:),allocatable,save :: q2 ! Turbulent Kinetic Energy. |
|---|
| 228 | |
|---|
| 229 | !$OMP THREADPRIVATE(emis,dtrad,fluxrad_sky,fluxrad,capcal,fluxgrd,qsurf,q2) |
|---|
| 230 | |
|---|
| 231 | |
|---|
| 232 | ! Local variables : |
|---|
| 233 | ! ----------------- |
|---|
| 234 | |
|---|
| 235 | real zh(ngrid,nlayer) ! Potential temperature (K). |
|---|
| 236 | real pw(ngrid,nlayer) ! Vertical velocity (m/s). (NOTE : >0 WHEN DOWNWARDS !!) |
|---|
| 237 | |
|---|
| 238 | integer l,ig,ierr,iq,nw,isoil |
|---|
| 239 | |
|---|
| 240 | ! FOR DIAGNOSTIC : |
|---|
| 241 | |
|---|
| 242 | real,dimension(:),allocatable,save :: fluxsurf_lw ! Incident Long Wave (IR) surface flux (W.m-2). |
|---|
| 243 | real,dimension(:),allocatable,save :: fluxsurf_sw ! Incident Short Wave (stellar) surface flux (W.m-2). |
|---|
| 244 | real,dimension(:),allocatable,save :: fluxsurfabs_sw ! Absorbed Short Wave (stellar) flux by the surface (W.m-2). |
|---|
| 245 | real,dimension(:),allocatable,save :: fluxtop_lw ! Outgoing LW (IR) flux to space (W.m-2). |
|---|
| 246 | real,dimension(:),allocatable,save :: fluxabs_sw ! Absorbed SW (stellar) flux (W.m-2). |
|---|
| 247 | real,dimension(:),allocatable,save :: fluxtop_dn ! Incoming SW (stellar) radiation at the top of the atmosphere (W.m-2). |
|---|
| 248 | real,dimension(:),allocatable,save :: fluxdyn ! Horizontal heat transport by dynamics (W.m-2). |
|---|
| 249 | real,dimension(:,:),allocatable,save :: OLR_nu ! Outgoing LW radiation in each band (Normalized to the band width (W/m2/cm-1)). |
|---|
| 250 | real,dimension(:,:),allocatable,save :: OSR_nu ! Outgoing SW radiation in each band (Normalized to the band width (W/m2/cm-1)). |
|---|
| 251 | real,dimension(:,:),allocatable,save :: zdtlw ! LW heating tendencies (K/s). |
|---|
| 252 | real,dimension(:,:),allocatable,save :: zdtsw ! SW heating tendencies (K/s). |
|---|
| 253 | real,dimension(:),allocatable,save :: sensibFlux ! Turbulent flux given by the atmosphere to the surface (W.m-2). |
|---|
| 254 | |
|---|
| 255 | !$OMP THREADPRIVATE(fluxsurf_lw,fluxsurf_sw,fluxsurfabs_sw,fluxtop_lw,fluxabs_sw,fluxtop_dn,fluxdyn,OLR_nu,OSR_nu,& |
|---|
| 256 | |
|---|
| 257 | !$OMP zdtlw,zdtsw,sensibFlux) |
|---|
| 258 | |
|---|
| 259 | real zls ! Solar longitude (radians). |
|---|
| 260 | real zlss ! Sub solar point longitude (radians). |
|---|
| 261 | real zday ! Date (time since Ls=0, calculated in sols). |
|---|
| 262 | real zzlay(ngrid,nlayer) ! Altitude at the middle of the atmospheric layers. |
|---|
| 263 | real zzlev(ngrid,nlayer+1) ! Altitude at the atmospheric layer boundaries. |
|---|
| 264 | |
|---|
| 265 | ! TENDENCIES due to various processes : |
|---|
| 266 | |
|---|
| 267 | ! For Surface Temperature : (K/s) |
|---|
| 268 | real zdtsurf(ngrid) ! Cumulated tendencies. |
|---|
| 269 | real zdtsurfmr(ngrid) ! Mass_redistribution routine. |
|---|
| 270 | real zdtsdif(ngrid) ! Turbdiff/vdifc routines. |
|---|
| 271 | |
|---|
| 272 | ! For Atmospheric Temperatures : (K/s) |
|---|
| 273 | real zdtdif(ngrid,nlayer) ! Turbdiff/vdifc routines. |
|---|
| 274 | real zdtmr(ngrid,nlayer) ! Mass_redistribution routine. |
|---|
| 275 | real zdtsw1(ngrid,nlayer), zdtlw1(ngrid,nlayer) ! Callcorrk routine. |
|---|
| 276 | |
|---|
| 277 | ! For Surface Tracers : (kg/m2/s) |
|---|
| 278 | real dqsurf(ngrid,nq) ! Cumulated tendencies. |
|---|
| 279 | real zdqsdif(ngrid,nq) ! Turbdiff/vdifc routines. |
|---|
| 280 | real zdqsurfmr(ngrid,nq) ! Mass_redistribution routine. |
|---|
| 281 | |
|---|
| 282 | ! For Tracers : (kg/kg_of_air/s) |
|---|
| 283 | real zdqadj(ngrid,nlayer,nq) ! Convadj routine. |
|---|
| 284 | real zdqdif(ngrid,nlayer,nq) ! Turbdiff/vdifc routines. |
|---|
| 285 | real zdqevap(ngrid,nlayer) ! Turbdiff routine. |
|---|
| 286 | real zdqmr(ngrid,nlayer,nq) ! Mass_redistribution routine. |
|---|
| 287 | |
|---|
| 288 | real zdqchi(ngrid,nlayer,nq) ! Chemical tendency ( chemistry routine ). |
|---|
| 289 | |
|---|
| 290 | real zdqmufi(ngrid,nlayer,nq) ! Microphysical tendency. |
|---|
| 291 | |
|---|
| 292 | real zdqfibar(ngrid,nlayer,nq) ! For 2D chemistry |
|---|
| 293 | real zdqmufibar(ngrid,nlayer,nq) ! For 2D chemistry |
|---|
| 294 | |
|---|
| 295 | ! For Winds : (m/s/s) |
|---|
| 296 | real zdvadj(ngrid,nlayer),zduadj(ngrid,nlayer) ! Convadj routine. |
|---|
| 297 | real zdumr(ngrid,nlayer),zdvmr(ngrid,nlayer) ! Mass_redistribution routine. |
|---|
| 298 | real zdvdif(ngrid,nlayer),zdudif(ngrid,nlayer) ! Turbdiff/vdifc routines. |
|---|
| 299 | real zdhdif(ngrid,nlayer) ! Turbdiff/vdifc routines. |
|---|
| 300 | real zdhadj(ngrid,nlayer) ! Convadj routine. |
|---|
| 301 | |
|---|
| 302 | ! For Pressure and Mass : |
|---|
| 303 | real zdmassmr(ngrid,nlayer) ! Atmospheric Mass tendency for mass_redistribution (kg_of_air/m2/s). |
|---|
| 304 | real zdmassmr_col(ngrid) ! Atmospheric Column Mass tendency for mass_redistribution (kg_of_air/m2/s). |
|---|
| 305 | real zdpsrfmr(ngrid) ! Pressure tendency for mass_redistribution routine (Pa/s). |
|---|
| 306 | |
|---|
| 307 | |
|---|
| 308 | |
|---|
| 309 | ! Local variables for LOCAL CALCULATIONS: |
|---|
| 310 | ! --------------------------------------- |
|---|
| 311 | real zflubid(ngrid) |
|---|
| 312 | real zplanck(ngrid),zpopsk(ngrid,nlayer) |
|---|
| 313 | real ztim1,ztim2,ztim3, z1,z2 |
|---|
| 314 | real ztime_fin |
|---|
| 315 | real zdh(ngrid,nlayer) |
|---|
| 316 | real gmplanet |
|---|
| 317 | real taux(ngrid),tauy(ngrid) |
|---|
| 318 | |
|---|
| 319 | |
|---|
| 320 | |
|---|
| 321 | ! local variables for DIAGNOSTICS : (diagfi & stat) |
|---|
| 322 | ! ------------------------------------------------- |
|---|
| 323 | real ps(ngrid) ! Surface Pressure. |
|---|
| 324 | real zt(ngrid,nlayer) ! Atmospheric Temperature. |
|---|
| 325 | real zu(ngrid,nlayer),zv(ngrid,nlayer) ! Zonal and Meridional Winds. |
|---|
| 326 | real zq(ngrid,nlayer,nq) ! Atmospheric Tracers. |
|---|
| 327 | real zdtadj(ngrid,nlayer) ! Convadj Diagnostic. |
|---|
| 328 | real zdtdyn(ngrid,nlayer) ! Dynamical Heating (K/s). |
|---|
| 329 | real zdudyn(ngrid,nlayer) ! Dynamical Zonal Wind tendency (m.s-2). |
|---|
| 330 | real,allocatable,dimension(:,:),save :: ztprevious ! Previous loop Atmospheric Temperature (K) ! Useful for Dynamical Heating calculation. |
|---|
| 331 | real,allocatable,dimension(:,:),save :: zuprevious ! Previous loop Zonal Wind (m.s-1) ! Useful for Zonal Wind tendency calculation. |
|---|
| 332 | !$OMP THREADPRIVATE(ztprevious,zuprevious) |
|---|
| 333 | |
|---|
| 334 | real zhorizwind(ngrid,nlayer) ! Horizontal Wind ( sqrt(u**+v*v)) |
|---|
| 335 | |
|---|
| 336 | real vmr(ngrid,nlayer) ! volume mixing ratio |
|---|
| 337 | real time_phys |
|---|
| 338 | |
|---|
| 339 | real ISR,ASR,OLR,GND,DYN,GSR,Ts1,Ts2,Ts3,TsS ! for Diagnostic. |
|---|
| 340 | |
|---|
| 341 | ! to test energy conservation (RW+JL) |
|---|
| 342 | real mass(ngrid,nlayer),massarea(ngrid,nlayer) |
|---|
| 343 | real dEtot, dEtots, AtmToSurf_TurbFlux |
|---|
| 344 | real,save :: dEtotSW, dEtotsSW, dEtotLW, dEtotsLW |
|---|
| 345 | !$OMP THREADPRIVATE(dEtotSW, dEtotsSW, dEtotLW, dEtotsLW) |
|---|
| 346 | real dEzRadsw(ngrid,nlayer),dEzRadlw(ngrid,nlayer),dEzdiff(ngrid,nlayer) |
|---|
| 347 | real dEdiffs(ngrid),dEdiff(ngrid) |
|---|
| 348 | |
|---|
| 349 | !JL12 conservation test for mean flow kinetic energy has been disabled temporarily |
|---|
| 350 | |
|---|
| 351 | real dItot, dItot_tmp, dVtot, dVtot_tmp |
|---|
| 352 | |
|---|
| 353 | real dWtot, dWtot_tmp, dWtots, dWtots_tmp |
|---|
| 354 | |
|---|
| 355 | |
|---|
| 356 | ! For Clear Sky Case. |
|---|
| 357 | real fluxsurf_lw1(ngrid), fluxsurf_sw1(ngrid), fluxsurfabs_sw1(ngrid) ! For SW/LW flux. |
|---|
| 358 | real fluxtop_lw1(ngrid), fluxabs_sw1(ngrid) ! For SW/LW flux. |
|---|
| 359 | real albedo_equivalent1(ngrid) ! For Equivalent albedo calculation. |
|---|
| 360 | real tf, ntf |
|---|
| 361 | |
|---|
| 362 | real,allocatable,dimension(:,:),save :: qsurf_hist |
|---|
| 363 | !$OMP THREADPRIVATE(qsurf_hist) |
|---|
| 364 | |
|---|
| 365 | ! Miscellaneous : |
|---|
| 366 | character(len=10) :: tmp1 |
|---|
| 367 | character(len=10) :: tmp2 |
|---|
| 368 | |
|---|
| 369 | ! Local variables for Titan chemistry and microphysics |
|---|
| 370 | ! ---------------------------------------------------- |
|---|
| 371 | |
|---|
| 372 | real,save :: ctimestep ! Chemistry timestep (s) |
|---|
| 373 | !$OMP THREADPRIVATE(ctimestep) |
|---|
| 374 | |
|---|
| 375 | ! Chemical tracers in molar fraction |
|---|
| 376 | real, dimension(ngrid,nlayer,nkim) :: ychim ! (mol/mol) |
|---|
| 377 | real, dimension(ngrid,nlayer,nkim) :: ychimbar ! For 2D chemistry |
|---|
| 378 | |
|---|
| 379 | ! Molar fraction tendencies ( chemistry, condensation and evaporation ) for tracers (mol/mol/s) |
|---|
| 380 | real, dimension(:,:,:), allocatable, save :: dycchi ! NB : Only for chem tracers. Saved since chemistry is not called every step. |
|---|
| 381 | !$OMP THREADPRIVATE(dycchi) |
|---|
| 382 | real, dimension(ngrid,nlayer,nq) :: dyccond ! Condensation rate. NB : for all tracers, as we want to use indx on it. |
|---|
| 383 | real, dimension(ngrid,nlayer,nq) :: dyccondbar ! For 2D chemistry |
|---|
| 384 | real, dimension(ngrid) :: dycevapCH4 ! Surface "pseudo-evaporation" rate (forcing constant surface humidity). |
|---|
| 385 | |
|---|
| 386 | ! Saturation profiles |
|---|
| 387 | real, dimension(ngrid,nlayer,nkim) :: ysat ! (mol/mol) |
|---|
| 388 | |
|---|
| 389 | ! Surface methane |
|---|
| 390 | real, dimension(:), allocatable, save :: tankCH4 ! Depth of surface methane tank (m) |
|---|
| 391 | !$OMP THREADPRIVATE(tankCH4) |
|---|
| 392 | |
|---|
| 393 | real :: i2e(ngrid,nlayer) ! int 2 ext factor ( X.kg-1 -> X.m-3 for diags ) |
|---|
| 394 | |
|---|
| 395 | real,save,dimension(:,:,:), allocatable :: tpq ! Tracers for decoupled microphysical tests ( temporary in 01/18 ) |
|---|
| 396 | !$OMP THREADPRIVATE(tpq) |
|---|
| 397 | real,dimension(ngrid,nlayer,nq) :: dtpq ! (temporary in 01/18) |
|---|
| 398 | |
|---|
| 399 | |
|---|
| 400 | !----------------------------------------------------------------------------- |
|---|
| 401 | ! Interface to calmufi |
|---|
| 402 | ! --> needed in order to pass assumed-shape arrays. Otherwise we must put calmufi in a module |
|---|
| 403 | ! (to have an explicit interface generated by the compiler). |
|---|
| 404 | ! Or one can put calmufi in MMP_GCM module (in muphytitan). |
|---|
| 405 | INTERFACE |
|---|
| 406 | SUBROUTINE calmufi(dt, plev, zlev, play, zlay, g3d, temp, pq, zdqfi, zdq) |
|---|
| 407 | REAL(kind=8), INTENT(IN) :: dt !! Physics timestep (s). |
|---|
| 408 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: plev !! Pressure levels (Pa). |
|---|
| 409 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: zlev !! Altitude levels (m). |
|---|
| 410 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: play !! Pressure layers (Pa). |
|---|
| 411 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: zlay !! Altitude at the center of each layer (m). |
|---|
| 412 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: g3d !! Latitude-Altitude depending gravitational acceleration (m.s-2). |
|---|
| 413 | REAL(kind=8), DIMENSION(:,:), INTENT(IN) :: temp !! Temperature at the center of each layer (K). |
|---|
| 414 | REAL(kind=8), DIMENSION(:,:,:), INTENT(IN) :: pq !! Tracers (\(kg.kg^{-1}}\)). |
|---|
| 415 | REAL(kind=8), DIMENSION(:,:,:), INTENT(IN) :: zdqfi !! Tendency from former processes for tracers (\(kg.kg^{-1}}\)). |
|---|
| 416 | REAL(kind=8), DIMENSION(:,:,:), INTENT(OUT) :: zdq !! Microphysical tendency for tracers (\(kg.kg^{-1}}\)). |
|---|
| 417 | END SUBROUTINE calmufi |
|---|
| 418 | END INTERFACE |
|---|
| 419 | |
|---|
| 420 | !================================================================================================== |
|---|
| 421 | |
|---|
| 422 | ! ----------------- |
|---|
| 423 | ! I. INITIALISATION |
|---|
| 424 | ! ----------------- |
|---|
| 425 | |
|---|
| 426 | ! -------------------------------- |
|---|
| 427 | ! I.1 First Call Initialisation. |
|---|
| 428 | ! -------------------------------- |
|---|
| 429 | if (firstcall) then |
|---|
| 430 | allocate(tpq(ngrid,nlayer,nq)) |
|---|
| 431 | tpq(:,:,:) = pq(:,:,:) |
|---|
| 432 | |
|---|
| 433 | ! Initialisation of nmicro as well as tracers names, indexes ... |
|---|
| 434 | if (ngrid.ne.1) then ! Already done in rcm1d |
|---|
| 435 | call initracer2(nq,nametrac) ! WARNING JB (27/03/2018): should be wrapped in an OMP SINGLE statement (see module notes) |
|---|
| 436 | endif |
|---|
| 437 | |
|---|
| 438 | ! Allocate saved arrays. |
|---|
| 439 | ALLOCATE(tsurf(ngrid)) |
|---|
| 440 | ALLOCATE(tsoil(ngrid,nsoilmx)) |
|---|
| 441 | ALLOCATE(albedo(ngrid,L_NSPECTV)) |
|---|
| 442 | ALLOCATE(albedo_equivalent(ngrid)) |
|---|
| 443 | ALLOCATE(albedo_bareground(ngrid)) |
|---|
| 444 | ALLOCATE(emis(ngrid)) |
|---|
| 445 | ALLOCATE(dtrad(ngrid,nlayer)) |
|---|
| 446 | ALLOCATE(fluxrad_sky(ngrid)) |
|---|
| 447 | ALLOCATE(fluxrad(ngrid)) |
|---|
| 448 | ALLOCATE(capcal(ngrid)) |
|---|
| 449 | ALLOCATE(fluxgrd(ngrid)) |
|---|
| 450 | ALLOCATE(qsurf(ngrid,nq)) |
|---|
| 451 | ALLOCATE(q2(ngrid,nlayer+1)) |
|---|
| 452 | ALLOCATE(tankCH4(ngrid)) |
|---|
| 453 | ALLOCATE(ztprevious(ngrid,nlayer)) |
|---|
| 454 | ALLOCATE(zuprevious(ngrid,nlayer)) |
|---|
| 455 | ALLOCATE(qsurf_hist(ngrid,nq)) |
|---|
| 456 | ALLOCATE(fluxsurf_lw(ngrid)) |
|---|
| 457 | ALLOCATE(fluxsurf_sw(ngrid)) |
|---|
| 458 | ALLOCATE(fluxsurfabs_sw(ngrid)) |
|---|
| 459 | ALLOCATE(fluxtop_lw(ngrid)) |
|---|
| 460 | ALLOCATE(fluxabs_sw(ngrid)) |
|---|
| 461 | ALLOCATE(fluxtop_dn(ngrid)) |
|---|
| 462 | ALLOCATE(fluxdyn(ngrid)) |
|---|
| 463 | ALLOCATE(OLR_nu(ngrid,L_NSPECTI)) |
|---|
| 464 | ALLOCATE(OSR_nu(ngrid,L_NSPECTV)) |
|---|
| 465 | ALLOCATE(sensibFlux(ngrid)) |
|---|
| 466 | ALLOCATE(zdtlw(ngrid,nlayer)) |
|---|
| 467 | ALLOCATE(zdtsw(ngrid,nlayer)) |
|---|
| 468 | |
|---|
| 469 | ! This is defined in comsaison_h |
|---|
| 470 | ALLOCATE(mu0(ngrid)) |
|---|
| 471 | ALLOCATE(fract(ngrid)) |
|---|
| 472 | ! This is defined in radcommon_h |
|---|
| 473 | ALLOCATE(gzlat(ngrid,nlayer)) |
|---|
| 474 | ALLOCATE(gzlat_ig(nlayer)) |
|---|
| 475 | ALLOCATE(Cmk(nlayer)) |
|---|
| 476 | |
|---|
| 477 | ! Variables set to 0 |
|---|
| 478 | ! ~~~~~~~~~~~~~~~~~~ |
|---|
| 479 | dtrad(:,:) = 0.0 |
|---|
| 480 | fluxrad(:) = 0.0 |
|---|
| 481 | zdtsw(:,:) = 0.0 |
|---|
| 482 | zdtlw(:,:) = 0.0 |
|---|
| 483 | |
|---|
| 484 | ! Initialize setup for correlated-k radiative transfer |
|---|
| 485 | ! JVO 17 : Was in callcorrk firstcall, but we need spectral intervals for microphysics |
|---|
| 486 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 487 | |
|---|
| 488 | if (corrk) then |
|---|
| 489 | |
|---|
| 490 | call system('rm -f surf_vals_long.out') |
|---|
| 491 | |
|---|
| 492 | write( tmp1, '(i3)' ) L_NSPECTI |
|---|
| 493 | write( tmp2, '(i3)' ) L_NSPECTV |
|---|
| 494 | banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2)) |
|---|
| 495 | banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir)) |
|---|
| 496 | |
|---|
| 497 | call setspi ! Basic infrared properties. |
|---|
| 498 | call setspv ! Basic visible properties. |
|---|
| 499 | call sugas_corrk ! Set up gaseous absorption properties. |
|---|
| 500 | |
|---|
| 501 | OLR_nu(:,:) = 0. |
|---|
| 502 | OSR_nu(:,:) = 0. |
|---|
| 503 | |
|---|
| 504 | endif |
|---|
| 505 | |
|---|
| 506 | ! Initialize names and timestep for chemistry |
|---|
| 507 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 508 | |
|---|
| 509 | if ( callchim ) then |
|---|
| 510 | |
|---|
| 511 | if ( moyzon_ch .and. ngrid.eq.1 ) then |
|---|
| 512 | print *, "moyzon_ch=",moyzon_ch," and ngrid=1" |
|---|
| 513 | print *, "Please desactivate zonal mean for 1D !" |
|---|
| 514 | stop |
|---|
| 515 | endif |
|---|
| 516 | |
|---|
| 517 | allocate(dycchi(ngrid,nlayer,nkim)) ! only for chemical tracers |
|---|
| 518 | |
|---|
| 519 | ! Chemistry timestep |
|---|
| 520 | ctimestep = ptimestep*REAL(ichim) |
|---|
| 521 | |
|---|
| 522 | zdqchi(:,:,:) = 0.0 |
|---|
| 523 | |
|---|
| 524 | endif |
|---|
| 525 | |
|---|
| 526 | ! Initialize microphysics. |
|---|
| 527 | ! ~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 528 | |
|---|
| 529 | IF ( callmufi ) THEN |
|---|
| 530 | ! WARNING JB (27/03/2018): inimufi AND mmp_initialize_optics should be wrapped in an OMP SINGLE statement. |
|---|
| 531 | call inimufi(ptimestep) |
|---|
| 532 | ! Optical coupling of YAMMS is plugged but inactivated for now |
|---|
| 533 | ! as long as the microphysics only isn't fully debugged -- JVO 01/18 |
|---|
| 534 | ! NOTE JB (03/18): mmp_optic_file is initialized in inimufi => mmp_initialize_optics must be called AFTER inimufi |
|---|
| 535 | IF (.NOT.uncoupl_optic_haze) call mmp_initialize_optics(mmp_optic_file) |
|---|
| 536 | |
|---|
| 537 | ! initialize microphysics diagnostics arrays. |
|---|
| 538 | call ini_diag_arrays(ngrid,nlayer,nice) |
|---|
| 539 | |
|---|
| 540 | ENDIF |
|---|
| 541 | |
|---|
| 542 | |
|---|
| 543 | #ifdef CPP_XIOS |
|---|
| 544 | ! Initialize XIOS context |
|---|
| 545 | write(*,*) "physiq: call wxios_context_init" |
|---|
| 546 | CALL wxios_context_init |
|---|
| 547 | #endif |
|---|
| 548 | |
|---|
| 549 | ! Read 'startfi.nc' file. |
|---|
| 550 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 551 | call phyetat0(startphy_file,ngrid,nlayer,"startfi.nc",0,0,nsoilmx,nq, & |
|---|
| 552 | day_ini,time_phys,tsurf,tsoil,emis,q2,qsurf,tankCH4) |
|---|
| 553 | if (.not.startphy_file) then |
|---|
| 554 | ! additionnal "academic" initialization of physics |
|---|
| 555 | if (is_master) write(*,*) "Physiq: initializing tsurf(:) to pt(:,1) !!" |
|---|
| 556 | tsurf(:)=pt(:,1) |
|---|
| 557 | if (is_master) write(*,*) "Physiq: initializing tsoil(:) to pt(:,1) !!" |
|---|
| 558 | do isoil=1,nsoilmx |
|---|
| 559 | tsoil(:,isoil)=tsurf(:) |
|---|
| 560 | enddo |
|---|
| 561 | if (is_master) write(*,*) "Physiq: initializing day_ini to pdat !" |
|---|
| 562 | day_ini=pday |
|---|
| 563 | endif |
|---|
| 564 | |
|---|
| 565 | if (pday.ne.day_ini) then |
|---|
| 566 | write(*,*) "ERROR in physiq.F90:" |
|---|
| 567 | write(*,*) "bad synchronization between physics and dynamics" |
|---|
| 568 | write(*,*) "dynamics day: ",pday |
|---|
| 569 | write(*,*) "physics day: ",day_ini |
|---|
| 570 | stop |
|---|
| 571 | endif |
|---|
| 572 | write (*,*) 'In physiq day_ini =', day_ini |
|---|
| 573 | |
|---|
| 574 | ! Initialize albedo calculation. |
|---|
| 575 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 576 | albedo(:,:)=0.0 |
|---|
| 577 | albedo_bareground(:)=0.0 |
|---|
| 578 | call surfini(ngrid,nq,qsurf,albedo,albedo_bareground) |
|---|
| 579 | |
|---|
| 580 | ! Initialize orbital calculation. |
|---|
| 581 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 582 | call iniorbit(apoastr,periastr,year_day,peri_day,obliquit) |
|---|
| 583 | |
|---|
| 584 | if(tlocked)then |
|---|
| 585 | print*,'Planet is tidally locked at resonance n=',nres |
|---|
| 586 | print*,'Make sure you have the right rotation rate!!!' |
|---|
| 587 | endif |
|---|
| 588 | |
|---|
| 589 | |
|---|
| 590 | ! Initialize soil. |
|---|
| 591 | ! ~~~~~~~~~~~~~~~~ |
|---|
| 592 | if (callsoil) then |
|---|
| 593 | |
|---|
| 594 | call soil(ngrid,nsoilmx,firstcall,lastcall,inertiedat, & |
|---|
| 595 | ptimestep,tsurf,tsoil,capcal,fluxgrd) |
|---|
| 596 | |
|---|
| 597 | else ! else of 'callsoil'. |
|---|
| 598 | |
|---|
| 599 | print*,'WARNING! Thermal conduction in the soil turned off' |
|---|
| 600 | capcal(:)=1.e6 |
|---|
| 601 | fluxgrd(:)=intheat |
|---|
| 602 | print*,'Flux from ground = ',intheat,' W m^-2' |
|---|
| 603 | |
|---|
| 604 | endif ! end of 'callsoil'. |
|---|
| 605 | |
|---|
| 606 | icount=1 |
|---|
| 607 | |
|---|
| 608 | |
|---|
| 609 | ! Initialize surface history variable. |
|---|
| 610 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 611 | qsurf_hist(:,:)=qsurf(:,:) |
|---|
| 612 | |
|---|
| 613 | ! Initialize variable for dynamical heating and zonal wind tendency diagnostic |
|---|
| 614 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 615 | ztprevious(:,:)=pt(:,:) |
|---|
| 616 | zuprevious(:,:)=pu(:,:) |
|---|
| 617 | |
|---|
| 618 | |
|---|
| 619 | if(meanOLR)then |
|---|
| 620 | call system('rm -f rad_bal.out') ! to record global radiative balance. |
|---|
| 621 | call system('rm -f tem_bal.out') ! to record global mean/max/min temperatures. |
|---|
| 622 | call system('rm -f h2o_bal.out') ! to record global hydrological balance. |
|---|
| 623 | endif |
|---|
| 624 | |
|---|
| 625 | if (ngrid.ne.1) then |
|---|
| 626 | ! Note : no need to create a restart file in 1d. |
|---|
| 627 | call physdem0("restartfi.nc",longitude,latitude,nsoilmx,ngrid,nlayer,nq, & |
|---|
| 628 | ptimestep,pday+nday,time_phys,cell_area, & |
|---|
| 629 | albedo_bareground,inertiedat,zmea,zstd,zsig,zgam,zthe) |
|---|
| 630 | endif |
|---|
| 631 | |
|---|
| 632 | ! XIOS outputs |
|---|
| 633 | #ifdef CPP_XIOS |
|---|
| 634 | |
|---|
| 635 | write(*,*) "physiq: call initialize_xios_output" |
|---|
| 636 | call initialize_xios_output(pday,ptime,ptimestep,daysec, & |
|---|
| 637 | presnivs,pseudoalt) |
|---|
| 638 | #endif |
|---|
| 639 | write(*,*) "physiq: end of firstcall" |
|---|
| 640 | endif ! end of 'firstcall' |
|---|
| 641 | |
|---|
| 642 | ! ------------------------------------------------------ |
|---|
| 643 | ! I.2 Initializations done at every physical timestep: |
|---|
| 644 | ! ------------------------------------------------------ |
|---|
| 645 | |
|---|
| 646 | #ifdef CPP_XIOS |
|---|
| 647 | ! update XIOS time/calendar |
|---|
| 648 | call update_xios_timestep |
|---|
| 649 | #endif |
|---|
| 650 | |
|---|
| 651 | ! Initialize various variables |
|---|
| 652 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 653 | |
|---|
| 654 | pdt(:,:) = 0.0 |
|---|
| 655 | zdtsurf(:) = 0.0 |
|---|
| 656 | pdq(:,:,:) = 0.0 |
|---|
| 657 | dqsurf(:,:) = 0.0 |
|---|
| 658 | pdu(:,:) = 0.0 |
|---|
| 659 | pdv(:,:) = 0.0 |
|---|
| 660 | pdpsrf(:) = 0.0 |
|---|
| 661 | zflubid(:) = 0.0 |
|---|
| 662 | taux(:) = 0.0 |
|---|
| 663 | tauy(:) = 0.0 |
|---|
| 664 | |
|---|
| 665 | zday=pday+ptime ! Compute time, in sols (and fraction thereof). |
|---|
| 666 | |
|---|
| 667 | ! Compute Stellar Longitude (Ls), and orbital parameters. |
|---|
| 668 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 669 | if (season) then |
|---|
| 670 | call stellarlong(zday,zls) |
|---|
| 671 | else |
|---|
| 672 | call stellarlong(float(day_ini),zls) |
|---|
| 673 | end if |
|---|
| 674 | |
|---|
| 675 | call orbite(zls,dist_star,declin,right_ascen) |
|---|
| 676 | |
|---|
| 677 | if (tlocked) then |
|---|
| 678 | zlss=Mod(-(2.*pi*(zday/year_day)*nres - right_ascen),2.*pi) |
|---|
| 679 | elseif (diurnal) then |
|---|
| 680 | zlss=-2.*pi*(zday-.5) |
|---|
| 681 | else if(diurnal .eqv. .false.) then |
|---|
| 682 | zlss=9999. |
|---|
| 683 | endif |
|---|
| 684 | |
|---|
| 685 | |
|---|
| 686 | ! Compute variations of g with latitude (oblate case). |
|---|
| 687 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 688 | if (oblate .eqv. .false.) then |
|---|
| 689 | gzlat(:,:) = g |
|---|
| 690 | else if (flatten .eq. 0.0 .or. J2 .eq. 0.0 .or. Rmean .eq. 0.0 .or. MassPlanet .eq. 0.0) then |
|---|
| 691 | print*,'I need values for flatten, J2, Rmean and MassPlanet to compute gzlat (else set oblate=.false.)' |
|---|
| 692 | call abort |
|---|
| 693 | else |
|---|
| 694 | gmplanet = MassPlanet*grav*1e24 |
|---|
| 695 | do ig=1,ngrid |
|---|
| 696 | gzlat(ig,:)= gmplanet/(Rmean**2) * (1.D0 + 0.75 *J2 - 2.0*flatten/3. + (2.*flatten - 15./4.* J2) * cos(2. * (pi/2. - latitude(ig)))) |
|---|
| 697 | end do |
|---|
| 698 | endif |
|---|
| 699 | |
|---|
| 700 | ! Compute altitudes with the geopotential coming from the dynamics. |
|---|
| 701 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 702 | |
|---|
| 703 | if (eff_gz .eqv. .false.) then |
|---|
| 704 | |
|---|
| 705 | do l=1,nlayer |
|---|
| 706 | zzlay(:,l) = pphi(:,l) / gzlat(:,l) |
|---|
| 707 | enddo |
|---|
| 708 | |
|---|
| 709 | else ! In this case we consider variations of g with altitude |
|---|
| 710 | |
|---|
| 711 | do l=1,nlayer |
|---|
| 712 | zzlay(:,l) = g*rad*rad / ( g*rad - pphi(:,l) ) - rad |
|---|
| 713 | gzlat(:,l) = g*rad*rad / ( rad + zzlay(:,l) )**2 |
|---|
| 714 | end do |
|---|
| 715 | |
|---|
| 716 | endif ! if eff_gz |
|---|
| 717 | |
|---|
| 718 | zzlev(:,1)=0. |
|---|
| 719 | zzlev(:,nlayer+1)=1.e7 ! Dummy top of last layer above 10000 km... |
|---|
| 720 | |
|---|
| 721 | do l=2,nlayer |
|---|
| 722 | do ig=1,ngrid |
|---|
| 723 | z1=(pplay(ig,l-1)+pplev(ig,l))/(pplay(ig,l-1)-pplev(ig,l)) |
|---|
| 724 | z2=(pplev(ig,l)+pplay(ig,l))/(pplev(ig,l)-pplay(ig,l)) |
|---|
| 725 | zzlev(ig,l)=(z1*zzlay(ig,l-1)+z2*zzlay(ig,l))/(z1+z2) |
|---|
| 726 | enddo |
|---|
| 727 | enddo |
|---|
| 728 | |
|---|
| 729 | ! Zonal averages needed for chemistry |
|---|
| 730 | ! ~~~~~~~ Taken from old Titan ~~~~~~ |
|---|
| 731 | if (moyzon_ch) then |
|---|
| 732 | |
|---|
| 733 | if (eff_gz .eqv. .false.) then |
|---|
| 734 | zzlaybar(1,:)=(zphibar(1,:)+zphisbar(1))/g |
|---|
| 735 | else ! if we take into account effective altitude-dependent gravity field |
|---|
| 736 | zzlaybar(1,:)=g*rad*rad/(g*rad-(zphibar(1,:)+zphisbar(1)))-rad |
|---|
| 737 | endif |
|---|
| 738 | zzlevbar(1,1)=zphisbar(1)/g |
|---|
| 739 | DO l=2,nlayer |
|---|
| 740 | z1=(zplaybar(1,l-1)+zplevbar(1,l))/(zplaybar(1,l-1)-zplevbar(1,l)) |
|---|
| 741 | z2=(zplevbar(1,l) +zplaybar(1,l))/(zplevbar(1,l) -zplaybar(1,l)) |
|---|
| 742 | zzlevbar(1,l)=(z1*zzlaybar(1,l-1)+z2*zzlaybar(1,l))/(z1+z2) |
|---|
| 743 | ENDDO |
|---|
| 744 | zzlevbar(1,nlayer+1)=zzlaybar(1,nlayer)+(zzlaybar(1,nlayer)-zzlevbar(1,nlayer)) |
|---|
| 745 | |
|---|
| 746 | DO ig=2,ngrid |
|---|
| 747 | if (latitude(ig).ne.latitude(ig-1)) then |
|---|
| 748 | DO l=1,nlayer |
|---|
| 749 | if (eff_gz .eqv. .false.) then |
|---|
| 750 | zzlaybar(ig,l)=(zphibar(ig,l)+zphisbar(ig))/g |
|---|
| 751 | else ! if we take into account effective altitude-dependent gravity field |
|---|
| 752 | zzlaybar(ig,l)=g*rad*rad/(g*rad-(zphibar(ig,l)+zphisbar(ig)))-rad |
|---|
| 753 | endif |
|---|
| 754 | ENDDO |
|---|
| 755 | zzlevbar(ig,1)=zphisbar(ig)/g |
|---|
| 756 | DO l=2,nlayer |
|---|
| 757 | z1=(zplaybar(ig,l-1)+zplevbar(ig,l))/ (zplaybar(ig,l-1)-zplevbar(ig,l)) |
|---|
| 758 | z2=(zplevbar(ig,l) +zplaybar(ig,l))/(zplevbar(ig,l) -zplaybar(ig,l)) |
|---|
| 759 | zzlevbar(ig,l)=(z1*zzlaybar(ig,l-1)+z2*zzlaybar(ig,l))/(z1+z2) |
|---|
| 760 | ENDDO |
|---|
| 761 | zzlevbar(ig,nlayer+1)=zzlaybar(ig,nlayer)+(zzlaybar(ig,nlayer)-zzlevbar(ig,nlayer)) |
|---|
| 762 | else |
|---|
| 763 | zzlaybar(ig,:)=zzlaybar(ig-1,:) |
|---|
| 764 | zzlevbar(ig,:)=zzlevbar(ig-1,:) |
|---|
| 765 | endif |
|---|
| 766 | ENDDO |
|---|
| 767 | |
|---|
| 768 | endif ! moyzon |
|---|
| 769 | |
|---|
| 770 | ! ------------------------------------------------------------------------------------- |
|---|
| 771 | ! Compute potential temperature |
|---|
| 772 | ! Note : Potential temperature calculation may not be the same in physiq and dynamic... |
|---|
| 773 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 774 | do l=1,nlayer |
|---|
| 775 | do ig=1,ngrid |
|---|
| 776 | zpopsk(ig,l)=(pplay(ig,l)/pplev(ig,1))**rcp |
|---|
| 777 | zh(ig,l)=pt(ig,l)/zpopsk(ig,l) |
|---|
| 778 | mass(ig,l) = (pplev(ig,l) - pplev(ig,l+1))/gzlat(ig,l) |
|---|
| 779 | massarea(ig,l)=mass(ig,l)*cell_area(ig) |
|---|
| 780 | enddo |
|---|
| 781 | enddo |
|---|
| 782 | |
|---|
| 783 | ! Compute vertical velocity (m/s) from vertical mass flux |
|---|
| 784 | ! w = F / (rho*area) and rho = P/(r*T) |
|---|
| 785 | ! But first linearly interpolate mass flux to mid-layers |
|---|
| 786 | do l=1,nlayer-1 |
|---|
| 787 | pw(:,l)=0.5*(flxw(:,l)+flxw(:,l+1)) |
|---|
| 788 | enddo |
|---|
| 789 | pw(:,nlayer)=0.5*flxw(:,nlayer) ! since flxw(nlayer+1)=0 |
|---|
| 790 | do l=1,nlayer |
|---|
| 791 | pw(:,l)=(pw(:,l)*r*pt(:,l)) / (pplay(:,l)*cell_area(:)) |
|---|
| 792 | enddo |
|---|
| 793 | |
|---|
| 794 | !--------------------------------- |
|---|
| 795 | ! II. Compute radiative tendencies |
|---|
| 796 | !--------------------------------- |
|---|
| 797 | |
|---|
| 798 | if (callrad) then |
|---|
| 799 | if( mod(icount-1,iradia).eq.0.or.lastcall) then |
|---|
| 800 | |
|---|
| 801 | ! Compute local stellar zenith angles |
|---|
| 802 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 803 | if (tlocked) then |
|---|
| 804 | ! JL14 corrects tidally resonant (and inclined) cases. nres=omega_rot/omega_orb |
|---|
| 805 | ztim1=SIN(declin) |
|---|
| 806 | ztim2=COS(declin)*COS(zlss) |
|---|
| 807 | ztim3=COS(declin)*SIN(zlss) |
|---|
| 808 | |
|---|
| 809 | call stelang(ngrid,sinlon,coslon,sinlat,coslat, & |
|---|
| 810 | ztim1,ztim2,ztim3,mu0,fract, flatten) |
|---|
| 811 | |
|---|
| 812 | elseif (diurnal) then |
|---|
| 813 | ztim1=SIN(declin) |
|---|
| 814 | ztim2=COS(declin)*COS(2.*pi*(zday-.5)) |
|---|
| 815 | ztim3=-COS(declin)*SIN(2.*pi*(zday-.5)) |
|---|
| 816 | |
|---|
| 817 | call stelang(ngrid,sinlon,coslon,sinlat,coslat, & |
|---|
| 818 | ztim1,ztim2,ztim3,mu0,fract, flatten) |
|---|
| 819 | else if(diurnal .eqv. .false.) then |
|---|
| 820 | |
|---|
| 821 | call mucorr(ngrid,declin,latitude,mu0,fract,10000.,rad,flatten) |
|---|
| 822 | ! WARNING: this function appears not to work in 1D |
|---|
| 823 | |
|---|
| 824 | endif |
|---|
| 825 | |
|---|
| 826 | ! Eclipse incoming sunlight (e.g. Saturn ring shadowing). |
|---|
| 827 | if(rings_shadow) then |
|---|
| 828 | call call_rings(ngrid, ptime, pday, diurnal) |
|---|
| 829 | endif |
|---|
| 830 | |
|---|
| 831 | |
|---|
| 832 | if (corrk) then |
|---|
| 833 | |
|---|
| 834 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 835 | ! II.a Call correlated-k radiative transfer scheme |
|---|
| 836 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 837 | |
|---|
| 838 | call call_profilgases(nlayer) |
|---|
| 839 | |
|---|
| 840 | ! standard callcorrk |
|---|
| 841 | call callcorrk(ngrid,nlayer,pq,nq,qsurf,zday, & |
|---|
| 842 | albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt, & |
|---|
| 843 | tsurf,fract,dist_star, & |
|---|
| 844 | zdtlw,zdtsw,fluxsurf_lw,fluxsurf_sw, & |
|---|
| 845 | fluxsurfabs_sw,fluxtop_lw, & |
|---|
| 846 | fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu, & |
|---|
| 847 | lastcall) |
|---|
| 848 | |
|---|
| 849 | ! Radiative flux from the sky absorbed by the surface (W.m-2). |
|---|
| 850 | GSR=0.0 |
|---|
| 851 | fluxrad_sky(:)=emis(:)*fluxsurf_lw(:)+fluxsurfabs_sw(:) |
|---|
| 852 | |
|---|
| 853 | !if(noradsurf)then ! no lower surface; SW flux just disappears |
|---|
| 854 | ! GSR = SUM(fluxsurf_sw(:)*cell_area(:))/totarea |
|---|
| 855 | ! fluxrad_sky(:)=emis(:)*fluxsurf_lw(:) |
|---|
| 856 | ! print*,'SW lost in deep atmosphere = ',GSR,' W m^-2' |
|---|
| 857 | !endif |
|---|
| 858 | |
|---|
| 859 | ! Net atmospheric radiative heating rate (K.s-1) |
|---|
| 860 | dtrad(:,:)=zdtsw(:,:)+zdtlw(:,:) |
|---|
| 861 | |
|---|
| 862 | elseif(newtonian)then |
|---|
| 863 | |
|---|
| 864 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 865 | ! II.b Call Newtonian cooling scheme |
|---|
| 866 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 867 | call newtrelax(ngrid,nlayer,mu0,sinlat,zpopsk,pt,pplay,pplev,dtrad,firstcall) |
|---|
| 868 | |
|---|
| 869 | zdtsurf(:) = +(pt(:,1)-tsurf(:))/ptimestep |
|---|
| 870 | ! e.g. surface becomes proxy for 1st atmospheric layer ? |
|---|
| 871 | |
|---|
| 872 | else |
|---|
| 873 | |
|---|
| 874 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 875 | ! II.c Atmosphere has no radiative effect |
|---|
| 876 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 877 | fluxtop_dn(:) = fract(:)*mu0(:)*Fat1AU/dist_star**2 |
|---|
| 878 | if(ngrid.eq.1)then ! / by 4 globally in 1D case... |
|---|
| 879 | fluxtop_dn(1) = fract(1)*Fat1AU/dist_star**2/2.0 |
|---|
| 880 | endif |
|---|
| 881 | fluxsurf_sw(:) = fluxtop_dn(:) |
|---|
| 882 | print*,'------------WARNING---WARNING------------' ! by MT2015. |
|---|
| 883 | print*,'You are in corrk=false mode, ' |
|---|
| 884 | print*,'and the surface albedo is taken equal to the first visible spectral value' |
|---|
| 885 | |
|---|
| 886 | fluxsurfabs_sw(:) = fluxtop_dn(:)*(1.-albedo(:,1)) |
|---|
| 887 | fluxrad_sky(:) = fluxsurfabs_sw(:) |
|---|
| 888 | fluxtop_lw(:) = emis(:)*sigma*tsurf(:)**4 |
|---|
| 889 | |
|---|
| 890 | dtrad(:,:)=0.0 ! no atmospheric radiative heating |
|---|
| 891 | |
|---|
| 892 | endif ! end of corrk |
|---|
| 893 | |
|---|
| 894 | endif ! of if(mod(icount-1,iradia).eq.0) |
|---|
| 895 | |
|---|
| 896 | |
|---|
| 897 | ! Transformation of the radiative tendencies |
|---|
| 898 | ! ------------------------------------------ |
|---|
| 899 | zplanck(:)=tsurf(:)*tsurf(:) |
|---|
| 900 | zplanck(:)=emis(:)*sigma*zplanck(:)*zplanck(:) |
|---|
| 901 | fluxrad(:)=fluxrad_sky(:)-zplanck(:) |
|---|
| 902 | pdt(:,:)=pdt(:,:)+dtrad(:,:) |
|---|
| 903 | |
|---|
| 904 | ! Test of energy conservation |
|---|
| 905 | !---------------------------- |
|---|
| 906 | if(enertest)then |
|---|
| 907 | call planetwide_sumval(cpp*massarea(:,:)*zdtsw(:,:)/totarea_planet,dEtotSW) |
|---|
| 908 | call planetwide_sumval(cpp*massarea(:,:)*zdtlw(:,:)/totarea_planet,dEtotLW) |
|---|
| 909 | !call planetwide_sumval(fluxsurf_sw(:)*(1.-albedo_equivalent(:))*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk |
|---|
| 910 | call planetwide_sumval(fluxsurfabs_sw(:)*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk |
|---|
| 911 | call planetwide_sumval((fluxsurf_lw(:)*emis(:)-zplanck(:))*cell_area(:)/totarea_planet,dEtotsLW) |
|---|
| 912 | dEzRadsw(:,:)=cpp*mass(:,:)*zdtsw(:,:) |
|---|
| 913 | dEzRadlw(:,:)=cpp*mass(:,:)*zdtlw(:,:) |
|---|
| 914 | if (is_master) then |
|---|
| 915 | print*,'---------------------------------------------------------------' |
|---|
| 916 | print*,'In corrk SW atmospheric heating =',dEtotSW,' W m-2' |
|---|
| 917 | print*,'In corrk LW atmospheric heating =',dEtotLW,' W m-2' |
|---|
| 918 | print*,'atmospheric net rad heating (SW+LW) =',dEtotLW+dEtotSW,' W m-2' |
|---|
| 919 | print*,'In corrk SW surface heating =',dEtotsSW,' W m-2' |
|---|
| 920 | print*,'In corrk LW surface heating =',dEtotsLW,' W m-2' |
|---|
| 921 | print*,'surface net rad heating (SW+LW) =',dEtotsLW+dEtotsSW,' W m-2' |
|---|
| 922 | endif |
|---|
| 923 | endif ! end of 'enertest' |
|---|
| 924 | |
|---|
| 925 | endif ! of if (callrad) |
|---|
| 926 | |
|---|
| 927 | |
|---|
| 928 | |
|---|
| 929 | ! -------------------------------------------- |
|---|
| 930 | ! III. Vertical diffusion (turbulent mixing) : |
|---|
| 931 | ! -------------------------------------------- |
|---|
| 932 | |
|---|
| 933 | if (calldifv) then |
|---|
| 934 | |
|---|
| 935 | zflubid(:)=fluxrad(:)+fluxgrd(:) |
|---|
| 936 | |
|---|
| 937 | ! JL12 the following if test is temporarily there to allow us to compare the old vdifc with turbdiff. |
|---|
| 938 | if (UseTurbDiff) then |
|---|
| 939 | |
|---|
| 940 | call turbdiff(ngrid,nlayer,nq, & |
|---|
| 941 | ptimestep,capcal,lwrite, & |
|---|
| 942 | pplay,pplev,zzlay,zzlev,z0, & |
|---|
| 943 | pu,pv,pt,zpopsk,pq,tsurf,emis,qsurf, & |
|---|
| 944 | pdt,pdq,zflubid, & |
|---|
| 945 | zdudif,zdvdif,zdtdif,zdtsdif, & |
|---|
| 946 | sensibFlux,q2,zdqdif,zdqsdif, & |
|---|
| 947 | taux,tauy,lastcall) |
|---|
| 948 | |
|---|
| 949 | else |
|---|
| 950 | |
|---|
| 951 | zdh(:,:)=pdt(:,:)/zpopsk(:,:) |
|---|
| 952 | |
|---|
| 953 | call vdifc(ngrid,nlayer,nq,zpopsk, & |
|---|
| 954 | ptimestep,capcal,lwrite, & |
|---|
| 955 | pplay,pplev,zzlay,zzlev,z0, & |
|---|
| 956 | pu,pv,zh,pq,tsurf,emis,qsurf, & |
|---|
| 957 | zdh,pdq,zflubid, & |
|---|
| 958 | zdudif,zdvdif,zdhdif,zdtsdif, & |
|---|
| 959 | sensibFlux,q2,zdqdif,zdqsdif, & |
|---|
| 960 | taux,tauy,lastcall) |
|---|
| 961 | |
|---|
| 962 | zdtdif(:,:)=zdhdif(:,:)*zpopsk(:,:) ! for diagnostic only |
|---|
| 963 | zdqevap(:,:)=0. |
|---|
| 964 | |
|---|
| 965 | end if !end of 'UseTurbDiff' |
|---|
| 966 | |
|---|
| 967 | |
|---|
| 968 | pdv(:,:)=pdv(:,:)+zdvdif(:,:) |
|---|
| 969 | pdu(:,:)=pdu(:,:)+zdudif(:,:) |
|---|
| 970 | pdt(:,:)=pdt(:,:)+zdtdif(:,:) |
|---|
| 971 | zdtsurf(:)=zdtsurf(:)+zdtsdif(:) |
|---|
| 972 | |
|---|
| 973 | if (tracer) then |
|---|
| 974 | pdq(:,:,:)=pdq(:,:,:)+ zdqdif(:,:,:) |
|---|
| 975 | dqsurf(:,:)=dqsurf(:,:) + zdqsdif(:,:) |
|---|
| 976 | end if ! of if (tracer) |
|---|
| 977 | |
|---|
| 978 | |
|---|
| 979 | ! test energy conservation |
|---|
| 980 | !------------------------- |
|---|
| 981 | if(enertest)then |
|---|
| 982 | |
|---|
| 983 | dEzdiff(:,:)=cpp*mass(:,:)*zdtdif(:,:) |
|---|
| 984 | do ig = 1, ngrid |
|---|
| 985 | dEdiff(ig)=SUM(dEzdiff (ig,:))+ sensibFlux(ig)! subtract flux to the ground |
|---|
| 986 | dEzdiff(ig,1)= dEzdiff(ig,1)+ sensibFlux(ig)! subtract flux to the ground |
|---|
| 987 | enddo |
|---|
| 988 | |
|---|
| 989 | call planetwide_sumval(dEdiff(:)*cell_area(:)/totarea_planet,dEtot) |
|---|
| 990 | dEdiffs(:)=capcal(:)*zdtsdif(:)-zflubid(:)-sensibFlux(:) |
|---|
| 991 | call planetwide_sumval(dEdiffs(:)*cell_area(:)/totarea_planet,dEtots) |
|---|
| 992 | call planetwide_sumval(sensibFlux(:)*cell_area(:)/totarea_planet,AtmToSurf_TurbFlux) |
|---|
| 993 | |
|---|
| 994 | if (is_master) then |
|---|
| 995 | |
|---|
| 996 | if (UseTurbDiff) then |
|---|
| 997 | print*,'In TurbDiff sensible flux (atm=>surf) =',AtmToSurf_TurbFlux,' W m-2' |
|---|
| 998 | print*,'In TurbDiff non-cons atm nrj change =',dEtot,' W m-2' |
|---|
| 999 | print*,'In TurbDiff (correc rad+latent heat) surf nrj change =',dEtots,' W m-2' |
|---|
| 1000 | else |
|---|
| 1001 | print*,'In vdifc sensible flux (atm=>surf) =',AtmToSurf_TurbFlux,' W m-2' |
|---|
| 1002 | print*,'In vdifc non-cons atm nrj change =',dEtot,' W m-2' |
|---|
| 1003 | print*,'In vdifc (correc rad+latent heat) surf nrj change =',dEtots,' W m-2' |
|---|
| 1004 | end if |
|---|
| 1005 | endif ! end of 'is_master' |
|---|
| 1006 | |
|---|
| 1007 | ! JL12 : note that the black body radiative flux emitted by the surface has been updated by the implicit scheme but not given back elsewhere. |
|---|
| 1008 | endif ! end of 'enertest' |
|---|
| 1009 | |
|---|
| 1010 | else ! calldifv |
|---|
| 1011 | |
|---|
| 1012 | if(.not.newtonian)then |
|---|
| 1013 | |
|---|
| 1014 | zdtsurf(:) = zdtsurf(:) + (fluxrad(:) + fluxgrd(:))/capcal(:) |
|---|
| 1015 | |
|---|
| 1016 | endif |
|---|
| 1017 | |
|---|
| 1018 | endif ! end of 'calldifv' |
|---|
| 1019 | |
|---|
| 1020 | |
|---|
| 1021 | !---------------------------------- |
|---|
| 1022 | ! IV. Dry convective adjustment : |
|---|
| 1023 | !---------------------------------- |
|---|
| 1024 | |
|---|
| 1025 | if(calladj) then |
|---|
| 1026 | |
|---|
| 1027 | zdh(:,:) = pdt(:,:)/zpopsk(:,:) |
|---|
| 1028 | zduadj(:,:)=0.0 |
|---|
| 1029 | zdvadj(:,:)=0.0 |
|---|
| 1030 | zdhadj(:,:)=0.0 |
|---|
| 1031 | |
|---|
| 1032 | |
|---|
| 1033 | call convadj(ngrid,nlayer,nq,ptimestep, & |
|---|
| 1034 | pplay,pplev,zpopsk, & |
|---|
| 1035 | pu,pv,zh,pq, & |
|---|
| 1036 | pdu,pdv,zdh,pdq, & |
|---|
| 1037 | zduadj,zdvadj,zdhadj, & |
|---|
| 1038 | zdqadj) |
|---|
| 1039 | |
|---|
| 1040 | pdu(:,:) = pdu(:,:) + zduadj(:,:) |
|---|
| 1041 | pdv(:,:) = pdv(:,:) + zdvadj(:,:) |
|---|
| 1042 | pdt(:,:) = pdt(:,:) + zdhadj(:,:)*zpopsk(:,:) |
|---|
| 1043 | zdtadj(:,:) = zdhadj(:,:)*zpopsk(:,:) ! for diagnostic only |
|---|
| 1044 | |
|---|
| 1045 | if(tracer) then |
|---|
| 1046 | pdq(:,:,:) = pdq(:,:,:) + zdqadj(:,:,:) |
|---|
| 1047 | end if |
|---|
| 1048 | |
|---|
| 1049 | ! Test energy conservation |
|---|
| 1050 | if(enertest)then |
|---|
| 1051 | call planetwide_sumval(cpp*massarea(:,:)*zdtadj(:,:)/totarea_planet,dEtot) |
|---|
| 1052 | if (is_master) print*,'In convadj atmospheric energy change =',dEtot,' W m-2' |
|---|
| 1053 | endif |
|---|
| 1054 | |
|---|
| 1055 | |
|---|
| 1056 | endif ! end of 'calladj' |
|---|
| 1057 | |
|---|
| 1058 | |
|---|
| 1059 | !--------------------------------------------- |
|---|
| 1060 | ! V. Specific parameterizations for tracers |
|---|
| 1061 | !--------------------------------------------- |
|---|
| 1062 | |
|---|
| 1063 | if (tracer) then |
|---|
| 1064 | |
|---|
| 1065 | ! ------------------- |
|---|
| 1066 | ! V.1 Microphysics |
|---|
| 1067 | ! ------------------- |
|---|
| 1068 | |
|---|
| 1069 | ! JVO 05/18 : We must call microphysics before chemistry, for condensation ! |
|---|
| 1070 | |
|---|
| 1071 | if (callmufi) then |
|---|
| 1072 | #ifdef USE_QTEST |
|---|
| 1073 | dtpq(:,:,:) = 0.0 ! we want tpq to go only through mufi |
|---|
| 1074 | call calmufi(ptimestep,pplev,zzlev,pplay,zzlay,gzlat,pt,tpq,dtpq,zdqmufi) |
|---|
| 1075 | tpq(:,:,:) = tpq(:,:,:) + zdqmufi(:,:,:)*ptimestep ! only manipulation of tpq->*ptimestep here |
|---|
| 1076 | #else |
|---|
| 1077 | call calmufi(ptimestep,pplev,zzlev,pplay,zzlay,gzlat,pt,pq,pdq,zdqmufi) |
|---|
| 1078 | pdq(:,:,:) = pdq(:,:,:) + zdqmufi(:,:,:) |
|---|
| 1079 | #endif |
|---|
| 1080 | |
|---|
| 1081 | ! Microphysics condensation for 2D fields to sent non-saturated fields to photochem |
|---|
| 1082 | if ( callclouds .and. moyzon_ch .and. mod(icount-1,ichim).eq.0 ) then |
|---|
| 1083 | zdqfibar(:,:,:) = 0.0 ! We work in zonal average -> forget processes other than condensation |
|---|
| 1084 | call calmufi(ptimestep,zplevbar,zzlevbar,zplaybar,zzlaybar, & |
|---|
| 1085 | gzlat,ztfibar,zqfibar,zdqfibar,zdqmufibar) |
|---|
| 1086 | endif |
|---|
| 1087 | |
|---|
| 1088 | endif |
|---|
| 1089 | |
|---|
| 1090 | ! ----------------- |
|---|
| 1091 | ! V.2. Chemistry |
|---|
| 1092 | ! ----------------- |
|---|
| 1093 | ! NB : Must be call last ( brings fields back to an equilibrium ) |
|---|
| 1094 | |
|---|
| 1095 | ! Known bug ? ( JVO 18 ) : If you try to use chimi_indx instead of iq+nmicro |
|---|
| 1096 | ! it leads to weird results / crash on dev mod ( ok in debug ) ... Why ? Idk ... |
|---|
| 1097 | |
|---|
| 1098 | if (callchim) then |
|---|
| 1099 | |
|---|
| 1100 | ! o. Convert updated tracers to molar fraction |
|---|
| 1101 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 1102 | do iq = 1,nkim |
|---|
| 1103 | ychim(:,:,iq) = ( pq(:,:,iq+nmicro) + pdq(:,:,iq+nmicro) ) / rat_mmol(iq+nmicro) |
|---|
| 1104 | enddo |
|---|
| 1105 | |
|---|
| 1106 | ! JVO 05/18 : We update zonally averaged fields with condensation |
|---|
| 1107 | ! as it is compulsory to have correct photochem production. But for other |
|---|
| 1108 | ! processes ( convadj ... ) we miss them in any case as we work in zonally/diurnal |
|---|
| 1109 | ! mean -> no fine diurnal/short time evolution, only seasonal evolution only. |
|---|
| 1110 | if ( moyzon_ch .and. mod(icount-1,ichim).eq. 0 ) then |
|---|
| 1111 | do iq = 1,nkim |
|---|
| 1112 | ychimbar(:,:,iq) = zqfibar(:,:,iq+nmicro) / rat_mmol(iq+nmicro) |
|---|
| 1113 | if ( callclouds ) then |
|---|
| 1114 | ychimbar(:,:,iq) = ychimbar(:,:,iq) + ( zdqmufibar(:,:,iq+nmicro) / rat_mmol(iq+nmicro) ) |
|---|
| 1115 | endif |
|---|
| 1116 | enddo |
|---|
| 1117 | endif |
|---|
| 1118 | |
|---|
| 1119 | ! i. Condensation of the 3D tracers after the transport |
|---|
| 1120 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 1121 | |
|---|
| 1122 | call calc_ysat(ngrid,nlayer,pplay/100.0,pt,ysat) ! Compute saturation profiles for every grid point |
|---|
| 1123 | |
|---|
| 1124 | dyccond(:,:,:) = 0.0 ! Default value -> no condensation |
|---|
| 1125 | |
|---|
| 1126 | do iq=1,nkim |
|---|
| 1127 | where ( ychim(:,:,iq).gt.ysat(:,:,iq) ) & |
|---|
| 1128 | dyccond(:,:,iq+nmicro) = ( -ychim(:,:,iq)+ysat(:,:,iq) ) / ptimestep |
|---|
| 1129 | enddo |
|---|
| 1130 | |
|---|
| 1131 | if ( callclouds ) dyccond(:,:,ices_indx) = 0.0 ! Condensation have been calculated in the cloud microphysics |
|---|
| 1132 | |
|---|
| 1133 | do iq=1,nkim |
|---|
| 1134 | ychim(:,:,iq) = ychim(:,:,iq) + dyccond(:,:,iq+nmicro) ! update molar ychim for following calchim |
|---|
| 1135 | |
|---|
| 1136 | pdq(:,:,iq+nmicro) = pdq(:,:,iq+nmicro) + dyccond(:,:,iq+nmicro)*rat_mmol(iq+nmicro) ! convert tendencies to mass mixing ratio |
|---|
| 1137 | enddo |
|---|
| 1138 | |
|---|
| 1139 | |
|---|
| 1140 | ! ii. 2D zonally averaged fields need to condense and evap before photochem |
|---|
| 1141 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 1142 | if ( moyzon_ch .and. mod(icount-1,ichim).eq. 0 ) then |
|---|
| 1143 | |
|---|
| 1144 | call calc_ysat(ngrid,nlayer,zplaybar/100.0,ztfibar,ysat) ! Compute saturation profiles for every grid point for the zon-ave fields |
|---|
| 1145 | |
|---|
| 1146 | dyccondbar(:,:,:) = 0.0 ! Default value -> no condensation |
|---|
| 1147 | |
|---|
| 1148 | do iq = 1,nkim |
|---|
| 1149 | where ( ychimbar(:,:,iq).gt.ysat(:,:,iq) ) & |
|---|
| 1150 | dyccondbar(:,:,iq+nmicro) = ( -ychimbar(:,:,iq)+ysat(:,:,iq) ) / ptimestep |
|---|
| 1151 | enddo |
|---|
| 1152 | |
|---|
| 1153 | if ( callclouds ) dyccondbar(:,:,ices_indx) = 0.0 ! Condensation have been calculated in the cloud microphysics |
|---|
| 1154 | |
|---|
| 1155 | do iq=1,nkim |
|---|
| 1156 | ychimbar(:,:,iq) = ychimbar(:,:,iq) + dyccondbar(:,:,iq+nmicro) |
|---|
| 1157 | enddo |
|---|
| 1158 | |
|---|
| 1159 | ! Pseudo-evap ( forcing constant surface humidity ) |
|---|
| 1160 | do ig=1,ngrid |
|---|
| 1161 | if ( ychimbar(ig,1,7+nmicro) .lt. botCH4 ) ychimbar(ig,1,7+nmicro) = botCH4 |
|---|
| 1162 | enddo |
|---|
| 1163 | |
|---|
| 1164 | endif ! if ( moyzon_ch .and. mod(icount-1,ichim).eq. 0 ) |
|---|
| 1165 | |
|---|
| 1166 | ! iii. Photochemistry ( must be call after condensation (and evap of 2D) ) |
|---|
| 1167 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 1168 | if( mod(icount-1,ichim).eq.0. ) then |
|---|
| 1169 | |
|---|
| 1170 | print *, "We enter in the chemistry ..." |
|---|
| 1171 | |
|---|
| 1172 | if (moyzon_ch) then ! 2D zonally averaged chemistry |
|---|
| 1173 | |
|---|
| 1174 | ! Here we send zonal average fields ( corrected with cond ) from dynamics to chem. module |
|---|
| 1175 | call calchim(ngrid,ychimbar,declin,ctimestep,ztfibar,zphibar, & |
|---|
| 1176 | zplaybar,zplevbar,zzlaybar,zzlevbar,dycchi) |
|---|
| 1177 | else ! 3D chemistry (or 1D run) |
|---|
| 1178 | call calchim(ngrid,ychim,declin,ctimestep,pt,pphi, & |
|---|
| 1179 | pplay,pplev,zzlay,zzlev,dycchi) |
|---|
| 1180 | endif ! if moyzon |
|---|
| 1181 | |
|---|
| 1182 | endif |
|---|
| 1183 | |
|---|
| 1184 | ! iv. Surface pseudo-evaporation |
|---|
| 1185 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 1186 | do ig=1,ngrid |
|---|
| 1187 | if ( (ychim(ig,1,7+nmicro)+dycchi(ig,1,7+nmicro)*ptimestep) .lt. botCH4 ) then ! +dycchi because ychim not yet updated |
|---|
| 1188 | dycevapCH4(ig) = ( -ychim(ig,1,7+nmicro)+botCH4 ) / ptimestep - dycchi(ig,1,7+nmicro) |
|---|
| 1189 | else |
|---|
| 1190 | dycevapCH4(ig) = 0.0 |
|---|
| 1191 | endif |
|---|
| 1192 | enddo |
|---|
| 1193 | |
|---|
| 1194 | pdq(:,1,7+nmicro) = pdq(:,1,7+nmicro) + dycevapCH4(:)*rat_mmol(7+nmicro) |
|---|
| 1195 | |
|---|
| 1196 | ! v. Updates and positivity check |
|---|
| 1197 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 1198 | do iq=1,nkim |
|---|
| 1199 | zdqchi(:,:,iq+nmicro) = dycchi(:,:,iq)*rat_mmol(iq+nmicro) ! convert tendencies to mass mixing ratio |
|---|
| 1200 | |
|---|
| 1201 | where( (pq(:,:,iq+nmicro) + ( pdq(:,:,iq+nmicro)+zdqchi(:,:,iq+nmicro) )*ptimestep ) .LT. 0.) & ! When using zonal means we set the same tendency |
|---|
| 1202 | zdqchi(:,:,iq+nmicro) = 1.e-30 - pdq(:,:,iq+nmicro) - pq(:,:,iq+nmicro)/ptimestep ! everywhere in longitude -> could lead to negs ! |
|---|
| 1203 | enddo |
|---|
| 1204 | |
|---|
| 1205 | pdq(:,:,:) = pdq(:,:,:) + zdqchi(:,:,:) |
|---|
| 1206 | |
|---|
| 1207 | endif ! end of 'callchim' |
|---|
| 1208 | |
|---|
| 1209 | ! --------------- |
|---|
| 1210 | ! V.3 Updates |
|---|
| 1211 | ! --------------- |
|---|
| 1212 | |
|---|
| 1213 | ! Updating Atmospheric Mass and Tracers budgets. |
|---|
| 1214 | if(mass_redistrib) then |
|---|
| 1215 | |
|---|
| 1216 | zdmassmr(:,:) = mass(:,:) * zdqevap(:,:) |
|---|
| 1217 | |
|---|
| 1218 | do ig = 1, ngrid |
|---|
| 1219 | zdmassmr_col(ig)=SUM(zdmassmr(ig,:)) |
|---|
| 1220 | enddo |
|---|
| 1221 | |
|---|
| 1222 | call writediagfi(ngrid,"mass_evap","mass gain"," ",3,zdmassmr) |
|---|
| 1223 | call writediagfi(ngrid,"mass_evap_col","mass gain col"," ",2,zdmassmr_col) |
|---|
| 1224 | call writediagfi(ngrid,"mass","mass","kg/m2",3,mass) |
|---|
| 1225 | |
|---|
| 1226 | call mass_redistribution(ngrid,nlayer,nq,ptimestep, & |
|---|
| 1227 | capcal,pplay,pplev,pt,tsurf,pq,qsurf, & |
|---|
| 1228 | pu,pv,pdt,zdtsurf,pdq,pdu,pdv,zdmassmr, & |
|---|
| 1229 | zdtmr,zdtsurfmr,zdpsrfmr,zdumr,zdvmr,zdqmr,zdqsurfmr) |
|---|
| 1230 | |
|---|
| 1231 | pdq(:,:,:) = pdq(:,:,:) + zdqmr(:,:,:) |
|---|
| 1232 | dqsurf(:,:) = dqsurf(:,:) + zdqsurfmr(:,:) |
|---|
| 1233 | pdt(:,:) = pdt(:,:) + zdtmr(:,:) |
|---|
| 1234 | pdu(:,:) = pdu(:,:) + zdumr(:,:) |
|---|
| 1235 | pdv(:,:) = pdv(:,:) + zdvmr(:,:) |
|---|
| 1236 | pdpsrf(:) = pdpsrf(:) + zdpsrfmr(:) |
|---|
| 1237 | zdtsurf(:) = zdtsurf(:) + zdtsurfmr(:) |
|---|
| 1238 | |
|---|
| 1239 | endif |
|---|
| 1240 | |
|---|
| 1241 | ! ----------------------------- |
|---|
| 1242 | ! V.4. Surface Tracer Update |
|---|
| 1243 | ! ----------------------------- |
|---|
| 1244 | |
|---|
| 1245 | qsurf(:,:) = qsurf(:,:) + ptimestep*dqsurf(:,:) |
|---|
| 1246 | |
|---|
| 1247 | ! Add qsurf to qsurf_hist, which is what we save in diagfi.nc. At the same time, we set the water |
|---|
| 1248 | ! content of ocean gridpoints back to zero, in order to avoid rounding errors in vdifc, rain. |
|---|
| 1249 | qsurf_hist(:,:) = qsurf(:,:) |
|---|
| 1250 | |
|---|
| 1251 | endif! end of if 'tracer' |
|---|
| 1252 | |
|---|
| 1253 | |
|---|
| 1254 | !------------------------------------------------ |
|---|
| 1255 | ! VI. Surface and sub-surface soil temperature |
|---|
| 1256 | !------------------------------------------------ |
|---|
| 1257 | |
|---|
| 1258 | |
|---|
| 1259 | ! Increment surface temperature |
|---|
| 1260 | |
|---|
| 1261 | tsurf(:)=tsurf(:)+ptimestep*zdtsurf(:) |
|---|
| 1262 | |
|---|
| 1263 | ! Compute soil temperatures and subsurface heat flux. |
|---|
| 1264 | if (callsoil) then |
|---|
| 1265 | call soil(ngrid,nsoilmx,.false.,lastcall,inertiedat, & |
|---|
| 1266 | ptimestep,tsurf,tsoil,capcal,fluxgrd) |
|---|
| 1267 | endif |
|---|
| 1268 | |
|---|
| 1269 | |
|---|
| 1270 | ! Test energy conservation |
|---|
| 1271 | if(enertest)then |
|---|
| 1272 | call planetwide_sumval(cell_area(:)*capcal(:)*zdtsurf(:)/totarea_planet,dEtots) |
|---|
| 1273 | if (is_master) print*,'Surface energy change =',dEtots,' W m-2' |
|---|
| 1274 | endif |
|---|
| 1275 | |
|---|
| 1276 | |
|---|
| 1277 | !--------------------------------------------------- |
|---|
| 1278 | ! VII. Perform diagnostics and write output files |
|---|
| 1279 | !--------------------------------------------------- |
|---|
| 1280 | |
|---|
| 1281 | ! Note : For output only: the actual model integration is performed in the dynamics. |
|---|
| 1282 | |
|---|
| 1283 | ! Temperature, zonal and meridional winds. |
|---|
| 1284 | zt(:,:) = pt(:,:) + pdt(:,:)*ptimestep |
|---|
| 1285 | zu(:,:) = pu(:,:) + pdu(:,:)*ptimestep |
|---|
| 1286 | zv(:,:) = pv(:,:) + pdv(:,:)*ptimestep |
|---|
| 1287 | |
|---|
| 1288 | ! Diagnostic. |
|---|
| 1289 | zdtdyn(:,:) = (pt(:,:)-ztprevious(:,:)) / ptimestep |
|---|
| 1290 | ztprevious(:,:) = zt(:,:) |
|---|
| 1291 | |
|---|
| 1292 | zdudyn(:,:) = (pu(:,:)-zuprevious(:,:)) / ptimestep |
|---|
| 1293 | zuprevious(:,:) = zu(:,:) |
|---|
| 1294 | |
|---|
| 1295 | if(firstcall)then |
|---|
| 1296 | zdtdyn(:,:)=0.0 |
|---|
| 1297 | zdudyn(:,:)=0.0 |
|---|
| 1298 | endif |
|---|
| 1299 | |
|---|
| 1300 | ! Horizotal wind |
|---|
| 1301 | zhorizwind(:,:) = sqrt( zu(:,:)*zu(:,:) + zv(:,:)*zv(:,:) ) |
|---|
| 1302 | |
|---|
| 1303 | ! Dynamical heating diagnostic. |
|---|
| 1304 | do ig=1,ngrid |
|---|
| 1305 | fluxdyn(ig)= SUM(zdtdyn(ig,:) *mass(ig,:))*cpp |
|---|
| 1306 | enddo |
|---|
| 1307 | |
|---|
| 1308 | ! Tracers. |
|---|
| 1309 | zq(:,:,:) = pq(:,:,:) + pdq(:,:,:)*ptimestep |
|---|
| 1310 | |
|---|
| 1311 | ! Surface pressure. |
|---|
| 1312 | ps(:) = pplev(:,1) + pdpsrf(:)*ptimestep |
|---|
| 1313 | |
|---|
| 1314 | |
|---|
| 1315 | |
|---|
| 1316 | ! Surface and soil temperature information |
|---|
| 1317 | call planetwide_sumval(cell_area(:)*tsurf(:)/totarea_planet,Ts1) |
|---|
| 1318 | call planetwide_minval(tsurf(:),Ts2) |
|---|
| 1319 | call planetwide_maxval(tsurf(:),Ts3) |
|---|
| 1320 | if(callsoil)then |
|---|
| 1321 | TsS = SUM(cell_area(:)*tsoil(:,nsoilmx))/totarea ! mean temperature at bottom soil layer |
|---|
| 1322 | if (is_master) then |
|---|
| 1323 | print*,' ave[Tsurf] min[Tsurf] max[Tsurf] ave[Tdeep]' |
|---|
| 1324 | print*,Ts1,Ts2,Ts3,TsS |
|---|
| 1325 | end if |
|---|
| 1326 | else |
|---|
| 1327 | if (is_master) then |
|---|
| 1328 | print*,' ave[Tsurf] min[Tsurf] max[Tsurf]' |
|---|
| 1329 | print*,Ts1,Ts2,Ts3 |
|---|
| 1330 | endif |
|---|
| 1331 | end if |
|---|
| 1332 | |
|---|
| 1333 | |
|---|
| 1334 | ! Check the energy balance of the simulation during the run |
|---|
| 1335 | if(corrk)then |
|---|
| 1336 | |
|---|
| 1337 | call planetwide_sumval(cell_area(:)*fluxtop_dn(:)/totarea_planet,ISR) |
|---|
| 1338 | call planetwide_sumval(cell_area(:)*fluxabs_sw(:)/totarea_planet,ASR) |
|---|
| 1339 | call planetwide_sumval(cell_area(:)*fluxtop_lw(:)/totarea_planet,OLR) |
|---|
| 1340 | call planetwide_sumval(cell_area(:)*fluxgrd(:)/totarea_planet,GND) |
|---|
| 1341 | call planetwide_sumval(cell_area(:)*fluxdyn(:)/totarea_planet,DYN) |
|---|
| 1342 | do ig=1,ngrid |
|---|
| 1343 | if(fluxtop_dn(ig).lt.0.0)then |
|---|
| 1344 | print*,'fluxtop_dn has gone crazy' |
|---|
| 1345 | print*,'fluxtop_dn=',fluxtop_dn(ig) |
|---|
| 1346 | print*,'temp= ',pt(ig,:) |
|---|
| 1347 | print*,'pplay= ',pplay(ig,:) |
|---|
| 1348 | call abort |
|---|
| 1349 | endif |
|---|
| 1350 | end do |
|---|
| 1351 | |
|---|
| 1352 | if(ngrid.eq.1)then |
|---|
| 1353 | DYN=0.0 |
|---|
| 1354 | endif |
|---|
| 1355 | |
|---|
| 1356 | if (is_master) then |
|---|
| 1357 | print*,' ISR ASR OLR GND DYN [W m^-2]' |
|---|
| 1358 | print*, ISR,ASR,OLR,GND,DYN |
|---|
| 1359 | endif |
|---|
| 1360 | |
|---|
| 1361 | if(enertest .and. is_master)then |
|---|
| 1362 | print*,'SW flux/heating difference SW++ - ASR = ',dEtotSW+dEtotsSW-ASR,' W m-2' |
|---|
| 1363 | print*,'LW flux/heating difference LW++ - OLR = ',dEtotLW+dEtotsLW+OLR,' W m-2' |
|---|
| 1364 | print*,'LW energy balance LW++ + ASR = ',dEtotLW+dEtotsLW+ASR,' W m-2' |
|---|
| 1365 | endif |
|---|
| 1366 | |
|---|
| 1367 | if(meanOLR .and. is_master)then |
|---|
| 1368 | if((ngrid.gt.1) .or. (mod(icount-1,ecritphy).eq.0))then |
|---|
| 1369 | ! to record global radiative balance |
|---|
| 1370 | open(92,file="rad_bal.out",form='formatted',position='append') |
|---|
| 1371 | write(92,*) zday,ISR,ASR,OLR |
|---|
| 1372 | close(92) |
|---|
| 1373 | open(93,file="tem_bal.out",form='formatted',position='append') |
|---|
| 1374 | if(callsoil)then |
|---|
| 1375 | write(93,*) zday,Ts1,Ts2,Ts3,TsS |
|---|
| 1376 | else |
|---|
| 1377 | write(93,*) zday,Ts1,Ts2,Ts3 |
|---|
| 1378 | endif |
|---|
| 1379 | close(93) |
|---|
| 1380 | endif |
|---|
| 1381 | endif |
|---|
| 1382 | |
|---|
| 1383 | endif ! end of 'corrk' |
|---|
| 1384 | |
|---|
| 1385 | |
|---|
| 1386 | ! Diagnostic to test radiative-convective timescales in code. |
|---|
| 1387 | if(testradtimes)then |
|---|
| 1388 | open(38,file="tau_phys.out",form='formatted',position='append') |
|---|
| 1389 | ig=1 |
|---|
| 1390 | do l=1,nlayer |
|---|
| 1391 | write(38,*) -1./pdt(ig,l),pt(ig,l),pplay(ig,l) |
|---|
| 1392 | enddo |
|---|
| 1393 | close(38) |
|---|
| 1394 | print*,'As testradtimes enabled,' |
|---|
| 1395 | print*,'exiting physics on first call' |
|---|
| 1396 | call abort |
|---|
| 1397 | endif |
|---|
| 1398 | |
|---|
| 1399 | |
|---|
| 1400 | if (is_master) print*,'--> Ls =',zls*180./pi |
|---|
| 1401 | |
|---|
| 1402 | |
|---|
| 1403 | !---------------------------------------------------------------------- |
|---|
| 1404 | ! Writing NetCDF file "RESTARTFI" at the end of the run |
|---|
| 1405 | !---------------------------------------------------------------------- |
|---|
| 1406 | |
|---|
| 1407 | ! Note: 'restartfi' is stored just before dynamics are stored |
|---|
| 1408 | ! in 'restart'. Between now and the writting of 'restart', |
|---|
| 1409 | ! there will have been the itau=itau+1 instruction and |
|---|
| 1410 | ! a reset of 'time' (lastacll = .true. when itau+1= itaufin) |
|---|
| 1411 | ! thus we store for time=time+dtvr |
|---|
| 1412 | |
|---|
| 1413 | |
|---|
| 1414 | |
|---|
| 1415 | if(lastcall) then |
|---|
| 1416 | ztime_fin = ptime + ptimestep/(float(iphysiq)*daysec) |
|---|
| 1417 | |
|---|
| 1418 | if (ngrid.ne.1) then |
|---|
| 1419 | write(*,*)'PHYSIQ: for physdem ztime_fin =',ztime_fin |
|---|
| 1420 | |
|---|
| 1421 | call physdem1("restartfi.nc",nsoilmx,ngrid,nlayer,nq, & |
|---|
| 1422 | ptimestep,ztime_fin, & |
|---|
| 1423 | tsurf,tsoil,emis,q2,qsurf_hist,tankCH4) |
|---|
| 1424 | endif |
|---|
| 1425 | endif ! end of 'lastcall' |
|---|
| 1426 | |
|---|
| 1427 | |
|---|
| 1428 | !----------------------------------------------------------------------------------------------------- |
|---|
| 1429 | ! OUTPUT in netcdf file "DIAGFI.NC", containing any variable for diagnostic |
|---|
| 1430 | ! |
|---|
| 1431 | ! Note 1 : output with period "ecritphy", set in "run.def" |
|---|
| 1432 | ! |
|---|
| 1433 | ! Note 2 : writediagfi can also be called from any other subroutine for any variable, |
|---|
| 1434 | ! but its preferable to keep all the calls in one place ... |
|---|
| 1435 | !----------------------------------------------------------------------------------------------------- |
|---|
| 1436 | |
|---|
| 1437 | |
|---|
| 1438 | call writediagfi(ngrid,"Ls","solar longitude","deg",0,zls*180./pi) |
|---|
| 1439 | call writediagfi(ngrid,"Lss","sub solar longitude","deg",0,zlss*180./pi) |
|---|
| 1440 | call writediagfi(ngrid,"RA","right ascension","deg",0,right_ascen*180./pi) |
|---|
| 1441 | call writediagfi(ngrid,"Declin","solar declination","deg",0,declin*180./pi) |
|---|
| 1442 | call writediagfi(ngrid,"tsurf","Surface temperature","K",2,tsurf) |
|---|
| 1443 | call writediagfi(ngrid,"ps","Surface pressure","Pa",2,ps) |
|---|
| 1444 | call writediagfi(ngrid,"temp","temperature","K",3,zt) |
|---|
| 1445 | call writediagfi(ngrid,"teta","potential temperature","K",3,zh) |
|---|
| 1446 | call writediagfi(ngrid,"u","Zonal wind","m.s-1",3,zu) |
|---|
| 1447 | call writediagfi(ngrid,"v","Meridional wind","m.s-1",3,zv) |
|---|
| 1448 | call writediagfi(ngrid,"w","Vertical wind","m.s-1",3,pw) |
|---|
| 1449 | call writediagfi(ngrid,"p","Pressure","Pa",3,pplay) |
|---|
| 1450 | |
|---|
| 1451 | ! Subsurface temperatures |
|---|
| 1452 | ! call writediagsoil(ngrid,"tsurf","Surface temperature","K",2,tsurf) |
|---|
| 1453 | ! call writediagsoil(ngrid,"temp","temperature","K",3,tsoil) |
|---|
| 1454 | |
|---|
| 1455 | |
|---|
| 1456 | |
|---|
| 1457 | ! Total energy balance diagnostics |
|---|
| 1458 | if(callrad.and.(.not.newtonian))then |
|---|
| 1459 | |
|---|
| 1460 | !call writediagfi(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent) |
|---|
| 1461 | call writediagfi(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn) |
|---|
| 1462 | call writediagfi(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw) |
|---|
| 1463 | call writediagfi(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw) |
|---|
| 1464 | |
|---|
| 1465 | ! call writediagfi(ngrid,"ASRcs","absorbed stellar rad (cs).","W m-2",2,fluxabs_sw1) |
|---|
| 1466 | ! call writediagfi(ngrid,"OLRcs","outgoing longwave rad (cs).","W m-2",2,fluxtop_lw1) |
|---|
| 1467 | ! call writediagfi(ngrid,"fluxsurfsw","sw surface flux.","W m-2",2,fluxsurf_sw) |
|---|
| 1468 | ! call writediagfi(ngrid,"fluxsurflw","lw back radiation.","W m-2",2,fluxsurf_lw) |
|---|
| 1469 | ! call writediagfi(ngrid,"fluxsurfswcs","sw surface flux (cs).","W m-2",2,fluxsurf_sw1) |
|---|
| 1470 | ! call writediagfi(ngrid,"fluxsurflwcs","lw back radiation (cs).","W m-2",2,fluxsurf_lw1) |
|---|
| 1471 | |
|---|
| 1472 | |
|---|
| 1473 | call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrd) |
|---|
| 1474 | |
|---|
| 1475 | call writediagfi(ngrid,"DYN","dynamical heat input","W m-2",2,fluxdyn) |
|---|
| 1476 | |
|---|
| 1477 | endif ! end of 'callrad' |
|---|
| 1478 | |
|---|
| 1479 | if(enertest) then |
|---|
| 1480 | |
|---|
| 1481 | if (calldifv) then |
|---|
| 1482 | |
|---|
| 1483 | call writediagfi(ngrid,"q2","turbulent kinetic energy","J.kg^-1",3,q2) |
|---|
| 1484 | call writediagfi(ngrid,"sensibFlux","sensible heat flux","w.m^-2",2,sensibFlux) |
|---|
| 1485 | |
|---|
| 1486 | ! call writediagfi(ngrid,"dEzdiff","turbulent diffusion heating (-sensible flux)","w.m^-2",3,dEzdiff) |
|---|
| 1487 | ! call writediagfi(ngrid,"dEdiff","integrated turbulent diffusion heating (-sensible flux)","w.m^-2",2,dEdiff) |
|---|
| 1488 | ! call writediagfi(ngrid,"dEdiffs","In TurbDiff (correc rad+latent heat) surf nrj change","w.m^-2",2,dEdiffs) |
|---|
| 1489 | |
|---|
| 1490 | endif |
|---|
| 1491 | |
|---|
| 1492 | if (corrk) then |
|---|
| 1493 | call writediagfi(ngrid,"dEzradsw","radiative heating","w.m^-2",3,dEzradsw) |
|---|
| 1494 | call writediagfi(ngrid,"dEzradlw","radiative heating","w.m^-2",3,dEzradlw) |
|---|
| 1495 | endif |
|---|
| 1496 | |
|---|
| 1497 | endif ! end of 'enertest' |
|---|
| 1498 | |
|---|
| 1499 | ! Temporary inclusions for winds diagnostics. |
|---|
| 1500 | call writediagfi(ngrid,"zdudif","Turbdiff tend. zon. wind","m s-2",3,zdudif) |
|---|
| 1501 | call writediagfi(ngrid,"zdudyn","Dyn. tend. zon. wind","m s-2",3,zdudyn) |
|---|
| 1502 | |
|---|
| 1503 | ! Temporary inclusions for heating diagnostics. |
|---|
| 1504 | call writediagfi(ngrid,"zdtsw","SW heating","T s-1",3,zdtsw) |
|---|
| 1505 | call writediagfi(ngrid,"zdtlw","LW heating","T s-1",3,zdtlw) |
|---|
| 1506 | call writediagfi(ngrid,"dtrad","radiative heating","K s-1",3,dtrad) |
|---|
| 1507 | call writediagfi(ngrid,"zdtdyn","Dyn. heating","T s-1",3,zdtdyn) |
|---|
| 1508 | |
|---|
| 1509 | ! For Debugging. |
|---|
| 1510 | call writediagfi(ngrid,'pphi','Geopotential',' ',3,pphi) |
|---|
| 1511 | |
|---|
| 1512 | |
|---|
| 1513 | ! Output tracers. |
|---|
| 1514 | if (tracer) then |
|---|
| 1515 | |
|---|
| 1516 | if (callmufi) then |
|---|
| 1517 | ! Microphysical tracers are expressed in unit/m3. |
|---|
| 1518 | ! convert X.kg-1 --> X.m-3 (whereas for optics was -> X.m-2) |
|---|
| 1519 | i2e(:,:) = ( pplev(:,1:nlayer)-pplev(:,2:nlayer+1) ) / gzlat(:,1:nlayer) /(zzlev(:,2:nlayer+1)-zzlev(:,1:nlayer)) |
|---|
| 1520 | |
|---|
| 1521 | #ifdef USE_QTEST |
|---|
| 1522 | ! Microphysical tracers passed through dyn+phys(except mufi) |
|---|
| 1523 | call writediagfi(ngrid,"mu_m0as_dp","Dynphys only spherical mode 0th order moment",'m-3',3,zq(:,:,micro_indx(1))*i2e) |
|---|
| 1524 | call writediagfi(ngrid,"mu_m3as_dp","Dynphys only spherical mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(2))*i2e) |
|---|
| 1525 | call writediagfi(ngrid,"mu_m0af_dp","Dynphys only fractal mode 0th order moment",'m-3',3,zq(:,:,micro_indx(3))*i2e) |
|---|
| 1526 | call writediagfi(ngrid,"mu_m3af_dp","Dynphys only fractal mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(4))*i2e) |
|---|
| 1527 | ! Microphysical tracers passed through mufi only |
|---|
| 1528 | call writediagfi(ngrid,"mu_m0as_mo","Mufi only spherical mode 0th order moment",'m-3',3,tpq(:,:,micro_indx(1))*i2e) |
|---|
| 1529 | call writediagfi(ngrid,"mu_m3as_mo","Mufi only spherical mode 3rd order moment",'m3/m3',3,tpq(:,:,micro_indx(2))*i2e) |
|---|
| 1530 | call writediagfi(ngrid,"mu_m0af_mo","Mufi only fractal mode 0th order moment",'m-3',3,tpq(:,:,micro_indx(3))*i2e) |
|---|
| 1531 | call writediagfi(ngrid,"mu_m3af_mo","Mufi only fractal mode 3rd order moment",'m3/m3',3,tpq(:,:,micro_indx(4))*i2e) |
|---|
| 1532 | #else |
|---|
| 1533 | call writediagfi(ngrid,"mu_m0as","Spherical mode 0th order moment",'m-3',3,zq(:,:,micro_indx(1))*i2e) |
|---|
| 1534 | call writediagfi(ngrid,"mu_m3as","Spherical mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(2))*i2e) |
|---|
| 1535 | call writediagfi(ngrid,"mu_m0af","Fractal mode 0th order moment",'m-3',3,zq(:,:,micro_indx(3))*i2e) |
|---|
| 1536 | call writediagfi(ngrid,"mu_m3af","Fractal mode 3rd order moment",'m3/m3',3,zq(:,:,micro_indx(4))*i2e) |
|---|
| 1537 | #endif |
|---|
| 1538 | |
|---|
| 1539 | ! Microphysical diagnostics |
|---|
| 1540 | call writediagfi(ngrid,"mmd_aer_prec","Total aerosols precipitations",'m',2,mmd_aer_prec) |
|---|
| 1541 | call writediagfi(ngrid,"mmd_aer_s_flux","Spherical aerosols sedimentation flux",'kg.m-2.s-1',3,mmd_aer_s_flux) |
|---|
| 1542 | call writediagfi(ngrid,"mmd_aer_f_flux","Fractal aerosols sedimentation flux",'kg.m-2.s-1',3,mmd_aer_f_flux) |
|---|
| 1543 | call writediagfi(ngrid,"mmd_rc_sph","Spherical mode caracteristic radius",'m',3,mmd_rc_sph) |
|---|
| 1544 | call writediagfi(ngrid,"mmd_rc_fra","Fractal mode caracteristic radius",'m',3,mmd_rc_fra) |
|---|
| 1545 | |
|---|
| 1546 | endif ! end of 'callmufi' |
|---|
| 1547 | |
|---|
| 1548 | ! Chemical tracers |
|---|
| 1549 | if (callchim) then |
|---|
| 1550 | do iq=1,nkim |
|---|
| 1551 | call writediagfi(ngrid,cnames(iq),cnames(iq),'mol/mol',3,zq(:,:,iq+nmicro)/rat_mmol(iq+nmicro)) |
|---|
| 1552 | enddo |
|---|
| 1553 | call writediagfi(ngrid,"evapCH4","Surface CH4 pseudo-evaporation rate",'mol/mol/s',2,dycevapCH4) |
|---|
| 1554 | endif |
|---|
| 1555 | |
|---|
| 1556 | endif ! end of 'tracer' |
|---|
| 1557 | |
|---|
| 1558 | ! XIOS outputs |
|---|
| 1559 | #ifdef CPP_XIOS |
|---|
| 1560 | ! Send fields to XIOS: (NB these fields must also be defined as |
|---|
| 1561 | ! <field id="..." /> in context_lmdz_physics.xml to be correctly used) |
|---|
| 1562 | CALL send_xios_field("ls",zls*180./pi) |
|---|
| 1563 | CALL send_xios_field("lss",zlss*180./pi) |
|---|
| 1564 | CALL send_xios_field("RA",right_ascen*180./pi) |
|---|
| 1565 | CALL send_xios_field("Declin",declin*180./pi) |
|---|
| 1566 | |
|---|
| 1567 | ! Total energy balance diagnostics |
|---|
| 1568 | if (callrad.and.(.not.newtonian)) then |
|---|
| 1569 | CALL send_xios_field("ISR_TOA",fluxtop_dn) |
|---|
| 1570 | CALL send_xios_field("OLR_TOA",fluxtop_lw) |
|---|
| 1571 | endif |
|---|
| 1572 | |
|---|
| 1573 | CALL send_xios_field("area",cell_area) |
|---|
| 1574 | CALL send_xios_field("pphi",pphi) |
|---|
| 1575 | CALL send_xios_field("pphis",phisfi) |
|---|
| 1576 | |
|---|
| 1577 | CALL send_xios_field("ps",ps) |
|---|
| 1578 | CALL send_xios_field("tsurf",tsurf) |
|---|
| 1579 | |
|---|
| 1580 | if(enertest) then |
|---|
| 1581 | if (calldifv) then |
|---|
| 1582 | CALL send_xios_field("sensibFlux",sensibFlux) |
|---|
| 1583 | endif |
|---|
| 1584 | endif |
|---|
| 1585 | |
|---|
| 1586 | CALL send_xios_field("temp",zt) |
|---|
| 1587 | CALL send_xios_field("teta",zh) |
|---|
| 1588 | CALL send_xios_field("u",zu) |
|---|
| 1589 | CALL send_xios_field("v",zv) |
|---|
| 1590 | CALL send_xios_field("w",pw) |
|---|
| 1591 | CALL send_xios_field("p",pplay) |
|---|
| 1592 | |
|---|
| 1593 | ! Winds diagnostics. |
|---|
| 1594 | CALL send_xios_field("dudif",zdudif) |
|---|
| 1595 | CALL send_xios_field("dudyn",zdudyn) |
|---|
| 1596 | |
|---|
| 1597 | CALL send_xios_field("horizwind",zhorizwind) |
|---|
| 1598 | |
|---|
| 1599 | ! Heating diagnostics. |
|---|
| 1600 | CALL send_xios_field("dtsw",zdtsw) |
|---|
| 1601 | CALL send_xios_field("dtlw",zdtlw) |
|---|
| 1602 | CALL send_xios_field("dtrad",dtrad) |
|---|
| 1603 | CALL send_xios_field("dtdyn",zdtdyn) |
|---|
| 1604 | |
|---|
| 1605 | ! Chemical tracers |
|---|
| 1606 | if (callchim) then |
|---|
| 1607 | |
|---|
| 1608 | ! Advected fields |
|---|
| 1609 | do iq=1,nkim |
|---|
| 1610 | CALL send_xios_field(trim(cnames(iq)),zq(:,:,iq+nmicro)/rat_mmol(iq+nmicro)) ! kg/kg -> mol/mol |
|---|
| 1611 | enddo |
|---|
| 1612 | |
|---|
| 1613 | ! Upper chemistry fields |
|---|
| 1614 | do iq=1,nkim |
|---|
| 1615 | CALL send_xios_field(trim(cnames(iq))//"_up",ykim_up(iq,:,:)) ! mol/mol |
|---|
| 1616 | enddo |
|---|
| 1617 | |
|---|
| 1618 | ! Append fields in ykim_tot for output on the total vertical grid (0->1300km) |
|---|
| 1619 | do iq=1,nkim |
|---|
| 1620 | |
|---|
| 1621 | ! GCM levels |
|---|
| 1622 | do l=1,nlayer |
|---|
| 1623 | ykim_tot(iq,:,l) = zq(:,l,iq+nmicro)/rat_mmol(iq+nmicro) |
|---|
| 1624 | enddo |
|---|
| 1625 | ! Upper levels |
|---|
| 1626 | do l=1,nlaykim_up |
|---|
| 1627 | ykim_tot(iq,:,nlayer+l) = ykim_up(iq,:,l) |
|---|
| 1628 | enddo |
|---|
| 1629 | |
|---|
| 1630 | CALL send_xios_field(trim(cnames(iq))//"_tot",ykim_tot(iq,:,:)) ! mol/mol |
|---|
| 1631 | |
|---|
| 1632 | enddo |
|---|
| 1633 | |
|---|
| 1634 | ! Condensation tendencies ( mol/mol/s ) |
|---|
| 1635 | CALL send_xios_field("dqcond_CH4",dyccond(:,:,7+nmicro)) |
|---|
| 1636 | CALL send_xios_field("dqcond_C2H2",dyccond(:,:,10+nmicro)) |
|---|
| 1637 | CALL send_xios_field("dqcond_C2H4",dyccond(:,:,12+nmicro)) |
|---|
| 1638 | CALL send_xios_field("dqcond_C2H6",dyccond(:,:,14+nmicro)) |
|---|
| 1639 | CALL send_xios_field("dqcond_C3H6",dyccond(:,:,17+nmicro)) |
|---|
| 1640 | CALL send_xios_field("dqcond_C4H4",dyccond(:,:,21+nmicro)) |
|---|
| 1641 | CALL send_xios_field("dqcond_CH3CCH",dyccond(:,:,24+nmicro)) |
|---|
| 1642 | CALL send_xios_field("dqcond_C3H8",dyccond(:,:,25+nmicro)) |
|---|
| 1643 | CALL send_xios_field("dqcond_C4H2",dyccond(:,:,26+nmicro)) |
|---|
| 1644 | CALL send_xios_field("dqcond_C4H6",dyccond(:,:,27+nmicro)) |
|---|
| 1645 | CALL send_xios_field("dqcond_C4H10",dyccond(:,:,28+nmicro)) |
|---|
| 1646 | CALL send_xios_field("dqcond_AC6H6",dyccond(:,:,29+nmicro)) |
|---|
| 1647 | CALL send_xios_field("dqcond_HCN",dyccond(:,:,36+nmicro)) |
|---|
| 1648 | CALL send_xios_field("dqcond_CH3CN",dyccond(:,:,40+nmicro)) |
|---|
| 1649 | CALL send_xios_field("dqcond_HC3N",dyccond(:,:,42+nmicro)) |
|---|
| 1650 | CALL send_xios_field("dqcond_NCCN",dyccond(:,:,43+nmicro)) |
|---|
| 1651 | CALL send_xios_field("dqcond_C4N2",dyccond(:,:,44+nmicro)) |
|---|
| 1652 | |
|---|
| 1653 | ! Pseudo-evaporation flux (mol/mol/s) |
|---|
| 1654 | CALL send_xios_field("evapCH4",dycevapCH4(:)) |
|---|
| 1655 | |
|---|
| 1656 | endif ! of 'if callchim' |
|---|
| 1657 | |
|---|
| 1658 | ! Microphysical tracers |
|---|
| 1659 | if (callmufi) then |
|---|
| 1660 | CALL send_xios_field("mu_m0as",zq(:,:,micro_indx(1))*i2e) |
|---|
| 1661 | CALL send_xios_field("mu_m3as",zq(:,:,micro_indx(2))*i2e) |
|---|
| 1662 | CALL send_xios_field("mu_m0af",zq(:,:,micro_indx(3))*i2e) |
|---|
| 1663 | CALL send_xios_field("mu_m3af",zq(:,:,micro_indx(4))*i2e) |
|---|
| 1664 | endif |
|---|
| 1665 | |
|---|
| 1666 | if (lastcall.and.is_omp_master) then |
|---|
| 1667 | write(*,*) "physiq: call xios_context_finalize" |
|---|
| 1668 | call xios_context_finalize |
|---|
| 1669 | endif |
|---|
| 1670 | #endif |
|---|
| 1671 | |
|---|
| 1672 | icount=icount+1 |
|---|
| 1673 | |
|---|
| 1674 | end subroutine physiq |
|---|
| 1675 | |
|---|
| 1676 | end module physiq_mod |
|---|