[38] | 1 | |
---|
| 2 | PROGRAM testphys1d |
---|
| 3 | ! to use 'getin' |
---|
[1036] | 4 | USE ioipsl_getincom, only: getin |
---|
[1543] | 5 | use dimphy, only : init_dimphy |
---|
| 6 | use mod_grid_phy_lmdz, only : regular_lonlat |
---|
[1130] | 7 | use infotrac, only: nqtot, tname |
---|
| 8 | use comsoil_h, only: volcapa, layer, mlayer, inertiedat, nsoilmx |
---|
[1541] | 9 | use comgeomfi_h, only: sinlat, ini_fillgeom |
---|
[1047] | 10 | use surfdat_h, only: albedodat, z0_default, emissiv, emisice, |
---|
| 11 | & albedice, iceradius, dtemisice, z0, |
---|
| 12 | & zmea, zstd, zsig, zgam, zthe, phisfi, |
---|
[2079] | 13 | & watercaptag, hmons, summit, base |
---|
[1047] | 14 | use slope_mod, only: theta_sl, psi_sl |
---|
[1130] | 15 | use phyredem, only: physdem0,physdem1 |
---|
[1543] | 16 | use geometry_mod, only: init_geometry |
---|
[1229] | 17 | use planete_h, only: year_day, periheli, aphelie, peri_day, |
---|
| 18 | & obliquit, emin_turb, lmixmin |
---|
[1524] | 19 | use comcstfi_h, only: pi, rad, omeg, g, mugaz, rcp, r, cpp |
---|
[1535] | 20 | use time_phylmdz_mod, only: daysec, dtphys, day_step, |
---|
| 21 | & ecritphy, iphysiq |
---|
[1919] | 22 | use dimradmars_mod, only: tauscaling,tauvis,totcloudfrac |
---|
[1944] | 23 | use co2cloud_mod, only: mem_Mccn_co2, mem_Mh2o_co2, |
---|
| 24 | & mem_Nccn_co2 |
---|
[1621] | 25 | USE comvert_mod, ONLY: ap,bp,aps,bps,pa,preff,sig, |
---|
| 26 | & presnivs,pseudoalt,scaleheight |
---|
| 27 | USE vertical_layers_mod, ONLY: init_vertical_layers |
---|
[1462] | 28 | USE logic_mod, ONLY: hybrid |
---|
[1543] | 29 | use physics_distribution_mod, only: init_physics_distribution |
---|
[1535] | 30 | use regular_lonlat_mod, only: init_regular_lonlat |
---|
[1543] | 31 | use mod_interface_dyn_phys, only: init_interface_dyn_phys |
---|
[1524] | 32 | USE phys_state_var_init_mod, ONLY: phys_state_var_init |
---|
[1549] | 33 | USE physiq_mod, ONLY: physiq |
---|
[38] | 34 | IMPLICIT NONE |
---|
| 35 | |
---|
| 36 | c======================================================================= |
---|
| 37 | c subject: |
---|
| 38 | c -------- |
---|
| 39 | c PROGRAM useful to run physical part of the martian GCM in a 1D column |
---|
| 40 | c |
---|
| 41 | c Can be compiled with a command like (e.g. for 25 layers) |
---|
| 42 | c "makegcm -p mars -d 25 testphys1d" |
---|
| 43 | c It requires the files "testphys1d.def" "callphys.def" |
---|
| 44 | c and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line) |
---|
| 45 | c and a file describing the sigma layers (e.g. "z2sig.def") |
---|
| 46 | c |
---|
| 47 | c author: Frederic Hourdin, R.Fournier,F.Forget |
---|
| 48 | c ------- |
---|
| 49 | c |
---|
| 50 | c update: 12/06/2003 including chemistry (S. Lebonnois) |
---|
| 51 | c and water ice (F. Montmessin) |
---|
| 52 | c |
---|
| 53 | c======================================================================= |
---|
| 54 | |
---|
[1621] | 55 | include "dimensions.h" |
---|
[1130] | 56 | integer, parameter :: ngrid = 1 !(2+(jjm-1)*iim - 1/jjm) |
---|
[1266] | 57 | integer, parameter :: nlayer = llm |
---|
[1047] | 58 | !#include "dimradmars.h" |
---|
| 59 | !#include "comgeomfi.h" |
---|
| 60 | !#include "surfdat.h" |
---|
| 61 | !#include "slope.h" |
---|
| 62 | !#include "comsoil.h" |
---|
| 63 | !#include "comdiurn.h" |
---|
[1621] | 64 | include "callkeys.h" |
---|
[1047] | 65 | !#include "comsaison.h" |
---|
[1130] | 66 | !#include "control.h" |
---|
[1621] | 67 | include "netcdf.inc" |
---|
| 68 | include "comg1d.h" |
---|
[1036] | 69 | !#include "advtrac.h" |
---|
[38] | 70 | |
---|
| 71 | c -------------------------------------------------------------- |
---|
| 72 | c Declarations |
---|
| 73 | c -------------------------------------------------------------- |
---|
| 74 | c |
---|
| 75 | INTEGER unitstart ! unite d'ecriture de "startfi" |
---|
[1273] | 76 | INTEGER nlevel,nsoil,ndt |
---|
[38] | 77 | INTEGER ilayer,ilevel,isoil,idt,iq |
---|
| 78 | LOGICAl firstcall,lastcall |
---|
| 79 | c |
---|
[627] | 80 | real,parameter :: odpref=610. ! DOD reference pressure (Pa) |
---|
| 81 | c |
---|
[38] | 82 | INTEGER day0 ! date initial (sol ; =0 a Ls=0) |
---|
| 83 | REAL day ! date durant le run |
---|
| 84 | REAL time ! time (0<time<1 ; time=0.5 a midi) |
---|
[1266] | 85 | REAL play(nlayer) ! Pressure at the middle of the layers (Pa) |
---|
| 86 | REAL plev(nlayer+1) ! intermediate pressure levels (pa) |
---|
[1130] | 87 | REAL psurf,tsurf(1) |
---|
[1266] | 88 | REAL u(nlayer),v(nlayer) ! zonal, meridional wind |
---|
[38] | 89 | REAL gru,grv ! prescribed "geostrophic" background wind |
---|
[1266] | 90 | REAL temp(nlayer) ! temperature at the middle of the layers |
---|
[1036] | 91 | REAL,ALLOCATABLE :: q(:,:) ! tracer mixing ratio (e.g. kg/kg) |
---|
| 92 | REAL,ALLOCATABLE :: qsurf(:) ! tracer surface budget (e.g. kg.m-2) |
---|
[38] | 93 | REAL tsoil(nsoilmx) ! subsurface soik temperature (K) |
---|
[1130] | 94 | REAL co2ice(1) ! co2ice layer (kg.m-2) |
---|
| 95 | REAL emis(1) ! surface layer |
---|
[1944] | 96 | REAL albedo(1,1) ! surface albedo |
---|
| 97 | REAL :: wstar(1)=0. ! Thermals vertical velocity |
---|
[1266] | 98 | REAL q2(nlayer+1) ! Turbulent Kinetic Energy |
---|
| 99 | REAL zlay(nlayer) ! altitude estimee dans les couches (km) |
---|
[38] | 100 | |
---|
| 101 | c Physical and dynamical tandencies (e.g. m.s-2, K/s, Pa/s) |
---|
[1266] | 102 | REAL du(nlayer),dv(nlayer),dtemp(nlayer) |
---|
| 103 | REAL dudyn(nlayer),dvdyn(nlayer),dtempdyn(nlayer) |
---|
[1549] | 104 | REAL dpsurf(1) |
---|
[1036] | 105 | REAL,ALLOCATABLE :: dq(:,:) |
---|
| 106 | REAL,ALLOCATABLE :: dqdyn(:,:) |
---|
[38] | 107 | |
---|
| 108 | c Various intermediate variables |
---|
| 109 | INTEGER thermo |
---|
| 110 | REAL zls |
---|
[1266] | 111 | REAL phi(nlayer),h(nlayer),s(nlayer) |
---|
| 112 | REAL pks, ptif, w(nlayer) |
---|
[1036] | 113 | REAL qtotinit,qtot |
---|
| 114 | real,allocatable :: mqtot(:) |
---|
[38] | 115 | INTEGER ierr, aslun |
---|
[1266] | 116 | REAL tmp1(0:nlayer),tmp2(0:nlayer) |
---|
[38] | 117 | integer :: nq=1 ! number of tracers |
---|
[1543] | 118 | real :: latitude(1), longitude(1), cell_area(1) |
---|
[38] | 119 | |
---|
| 120 | character*2 str2 |
---|
| 121 | character (len=7) :: str7 |
---|
| 122 | character(len=44) :: txt |
---|
| 123 | |
---|
[2228] | 124 | c New flag to compute paleo orbital configurations + few variables JN |
---|
| 125 | REAL halfaxe, excentric, Lsperi |
---|
| 126 | Logical paleomars |
---|
| 127 | |
---|
[38] | 128 | c======================================================================= |
---|
| 129 | |
---|
| 130 | c======================================================================= |
---|
| 131 | c INITIALISATION |
---|
| 132 | c======================================================================= |
---|
[1130] | 133 | ! initialize "serial/parallel" related stuff |
---|
| 134 | ! CALL init_phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/)) |
---|
[1543] | 135 | ! CALL init_phys_lmdz(1,1,llm,1,(/1/)) |
---|
| 136 | ! call initcomgeomphy |
---|
[38] | 137 | |
---|
| 138 | c ------------------------------------------------------ |
---|
[899] | 139 | c Prescribed constants to be set here |
---|
[38] | 140 | c ------------------------------------------------------ |
---|
| 141 | |
---|
| 142 | pi=2.E+0*asin(1.E+0) |
---|
| 143 | |
---|
[899] | 144 | c Mars planetary constants |
---|
[38] | 145 | c ---------------------------- |
---|
[899] | 146 | rad=3397200. ! mars radius (m) ~3397200 m |
---|
| 147 | daysec=88775. ! length of a sol (s) ~88775 s |
---|
| 148 | omeg=4.*asin(1.)/(daysec) ! rotation rate (rad.s-1) |
---|
| 149 | g=3.72 ! gravity (m.s-2) ~3.72 |
---|
| 150 | mugaz=43.49 ! atmosphere mola mass (g.mol-1) ~43.49 |
---|
[38] | 151 | rcp=.256793 ! = r/cp ~0.256793 |
---|
| 152 | r= 8.314511E+0 *1000.E+0/mugaz |
---|
| 153 | cpp= r/rcp |
---|
[899] | 154 | year_day = 669 ! lenght of year (sols) ~668.6 |
---|
| 155 | periheli = 206.66 ! minimum sun-mars distance (Mkm) ~206.66 |
---|
| 156 | aphelie = 249.22 ! maximum sun-mars distance (Mkm) ~249.22 |
---|
[2228] | 157 | halfaxe = 227.94 ! demi-grand axe de l'ellipse |
---|
[899] | 158 | peri_day = 485. ! perihelion date (sols since N. Spring) |
---|
| 159 | obliquit = 25.2 ! Obliquity (deg) ~25.2 |
---|
[2228] | 160 | excentric = 0.0934 ! Eccentricity (0.0934) |
---|
[38] | 161 | |
---|
[899] | 162 | c Planetary Boundary Layer and Turbulence parameters |
---|
| 163 | c -------------------------------------------------- |
---|
| 164 | z0_default = 1.e-2 ! surface roughness (m) ~0.01 |
---|
| 165 | emin_turb = 1.e-6 ! minimal turbulent energy ~1.e-8 |
---|
| 166 | lmixmin = 30 ! mixing length ~100 |
---|
[38] | 167 | |
---|
[899] | 168 | c cap properties and surface emissivities |
---|
[38] | 169 | c ---------------------------------------------------- |
---|
[899] | 170 | emissiv= 0.95 ! Bare ground emissivity ~.95 |
---|
| 171 | emisice(1)=0.95 ! Northern cap emissivity |
---|
| 172 | emisice(2)=0.95 ! Southern cap emisssivity |
---|
| 173 | albedice(1)=0.5 ! Northern cap albedo |
---|
| 174 | albedice(2)=0.5 ! Southern cap albedo |
---|
[38] | 175 | iceradius(1) = 100.e-6 ! mean scat radius of CO2 snow (north) |
---|
| 176 | iceradius(2) = 100.e-6 ! mean scat radius of CO2 snow (south) |
---|
| 177 | dtemisice(1) = 2. ! time scale for snow metamorphism (north) |
---|
| 178 | dtemisice(2) = 2. ! time scale for snow metamorphism (south |
---|
| 179 | |
---|
[1920] | 180 | c mesh surface (not a very usefull quantity in 1D) |
---|
| 181 | c ---------------------------------------------------- |
---|
| 182 | cell_area(1)=1.E+0 |
---|
| 183 | |
---|
[38] | 184 | c ------------------------------------------------------ |
---|
[899] | 185 | c Loading run parameters from "run.def" file |
---|
[38] | 186 | c ------------------------------------------------------ |
---|
| 187 | |
---|
| 188 | |
---|
| 189 | ! check if 'run.def' file is around (otherwise reading parameters |
---|
| 190 | ! from callphys.def via getin() routine won't work. |
---|
| 191 | open(99,file='run.def',status='old',form='formatted', |
---|
| 192 | & iostat=ierr) |
---|
| 193 | if (ierr.ne.0) then |
---|
| 194 | write(*,*) 'Cannot find required file "run.def"' |
---|
| 195 | write(*,*) ' (which should contain some input parameters' |
---|
| 196 | write(*,*) ' along with the following line:' |
---|
| 197 | write(*,*) ' INCLUDEDEF=callphys.def' |
---|
| 198 | write(*,*) ' )' |
---|
| 199 | write(*,*) ' ... might as well stop here ...' |
---|
| 200 | stop |
---|
| 201 | else |
---|
| 202 | close(99) |
---|
| 203 | endif |
---|
| 204 | |
---|
| 205 | ! check if we are going to run with or without tracers |
---|
| 206 | write(*,*) "Run with or without tracer transport ?" |
---|
| 207 | tracer=.false. ! default value |
---|
| 208 | call getin("tracer",tracer) |
---|
| 209 | write(*,*) " tracer = ",tracer |
---|
| 210 | |
---|
| 211 | ! while we're at it, check if there is a 'traceur.def' file |
---|
[1036] | 212 | ! and process it. |
---|
[38] | 213 | if (tracer) then |
---|
| 214 | ! load tracer names from file 'traceur.def' |
---|
| 215 | open(90,file='traceur.def',status='old',form='formatted', |
---|
| 216 | & iostat=ierr) |
---|
| 217 | if (ierr.ne.0) then |
---|
| 218 | write(*,*) 'Cannot find required file "traceur.def"' |
---|
| 219 | write(*,*) ' If you want to run with tracers, I need it' |
---|
| 220 | write(*,*) ' ... might as well stop here ...' |
---|
| 221 | stop |
---|
| 222 | else |
---|
| 223 | write(*,*) "testphys1d: Reading file traceur.def" |
---|
| 224 | ! read number of tracers: |
---|
| 225 | read(90,*,iostat=ierr) nq |
---|
[1036] | 226 | nqtot=nq ! set value of nqtot (in infotrac module) as nq |
---|
[38] | 227 | if (ierr.ne.0) then |
---|
| 228 | write(*,*) "testphys1d: error reading number of tracers" |
---|
| 229 | write(*,*) " (first line of traceur.def) " |
---|
| 230 | stop |
---|
| 231 | endif |
---|
| 232 | endif |
---|
[1036] | 233 | ! allocate arrays: |
---|
[1130] | 234 | allocate(tname(nq)) |
---|
[1266] | 235 | allocate(q(nlayer,nq)) |
---|
[1036] | 236 | allocate(qsurf(nq)) |
---|
[1266] | 237 | allocate(dq(nlayer,nq)) |
---|
| 238 | allocate(dqdyn(nlayer,nq)) |
---|
[1036] | 239 | allocate(mqtot(nq)) |
---|
| 240 | |
---|
[38] | 241 | ! read tracer names from file traceur.def |
---|
[1036] | 242 | do iq=1,nq |
---|
[1130] | 243 | read(90,*,iostat=ierr) tname(iq) |
---|
[38] | 244 | if (ierr.ne.0) then |
---|
| 245 | write(*,*) 'testphys1d: error reading tracer names...' |
---|
| 246 | stop |
---|
| 247 | endif |
---|
| 248 | enddo |
---|
| 249 | close(90) |
---|
| 250 | |
---|
| 251 | ! initialize tracers here: |
---|
| 252 | write(*,*) "testphys1d: initializing tracers" |
---|
| 253 | q(:,:)=0 ! default, set everything to zero |
---|
| 254 | qsurf(:)=0 |
---|
| 255 | ! "smarter" initialization of some tracers |
---|
| 256 | ! (get values from "profile_*" files, if these are available) |
---|
[1036] | 257 | do iq=1,nq |
---|
[38] | 258 | txt="" |
---|
[1130] | 259 | write(txt,"(a)") tname(iq) |
---|
[38] | 260 | write(*,*)" tracer:",trim(txt) |
---|
| 261 | ! CO2 |
---|
| 262 | if (txt.eq."co2") then |
---|
| 263 | q(:,iq)=0.95 ! kg /kg of atmosphere |
---|
| 264 | qsurf(iq)=0. ! kg/m2 (not used for CO2) |
---|
| 265 | ! even better, look for a "profile_co2" input file |
---|
| 266 | open(91,file='profile_co2',status='old', |
---|
| 267 | & form='formatted',iostat=ierr) |
---|
| 268 | if (ierr.eq.0) then |
---|
| 269 | read(91,*) qsurf(iq) |
---|
[1266] | 270 | do ilayer=1,nlayer |
---|
[38] | 271 | read(91,*) q(ilayer,iq) |
---|
| 272 | enddo |
---|
| 273 | endif |
---|
| 274 | close(91) |
---|
| 275 | endif ! of if (txt.eq."co2") |
---|
[161] | 276 | ! Allow for an initial profile of argon |
---|
| 277 | ! Can also be used to introduce a decaying tracer |
---|
| 278 | ! in the 1D (TBD) to study thermals |
---|
| 279 | if (txt.eq."ar") then |
---|
| 280 | !look for a "profile_ar" input file |
---|
| 281 | open(91,file='profile_ar',status='old', |
---|
| 282 | & form='formatted',iostat=ierr) |
---|
| 283 | if (ierr.eq.0) then |
---|
| 284 | read(91,*) qsurf(iq) |
---|
[1266] | 285 | do ilayer=1,nlayer |
---|
[161] | 286 | read(91,*) q(ilayer,iq) |
---|
| 287 | enddo |
---|
| 288 | else |
---|
| 289 | write(*,*) "No profile_ar file!" |
---|
| 290 | endif |
---|
| 291 | close(91) |
---|
| 292 | endif ! of if (txt.eq."ar") |
---|
| 293 | |
---|
[38] | 294 | ! WATER VAPOUR |
---|
| 295 | if (txt.eq."h2o_vap") then |
---|
| 296 | !look for a "profile_h2o_vap" input file |
---|
| 297 | open(91,file='profile_h2o_vap',status='old', |
---|
| 298 | & form='formatted',iostat=ierr) |
---|
| 299 | if (ierr.eq.0) then |
---|
| 300 | read(91,*) qsurf(iq) |
---|
[1266] | 301 | do ilayer=1,nlayer |
---|
[38] | 302 | read(91,*) q(ilayer,iq) |
---|
| 303 | enddo |
---|
| 304 | else |
---|
| 305 | write(*,*) "No profile_h2o_vap file!" |
---|
| 306 | endif |
---|
| 307 | close(91) |
---|
| 308 | endif ! of if (txt.eq."h2o_ice") |
---|
| 309 | ! WATER ICE |
---|
| 310 | if (txt.eq."h2o_ice") then |
---|
| 311 | !look for a "profile_h2o_vap" input file |
---|
| 312 | open(91,file='profile_h2o_ice',status='old', |
---|
| 313 | & form='formatted',iostat=ierr) |
---|
| 314 | if (ierr.eq.0) then |
---|
| 315 | read(91,*) qsurf(iq) |
---|
[1266] | 316 | do ilayer=1,nlayer |
---|
[38] | 317 | read(91,*) q(ilayer,iq) |
---|
| 318 | enddo |
---|
| 319 | else |
---|
| 320 | write(*,*) "No profile_h2o_ice file!" |
---|
| 321 | endif |
---|
| 322 | close(91) |
---|
| 323 | endif ! of if (txt.eq."h2o_ice") |
---|
| 324 | ! DUST |
---|
| 325 | !if (txt(1:4).eq."dust") then |
---|
| 326 | ! q(:,iq)=0.4 ! kg/kg of atmosphere |
---|
| 327 | ! qsurf(iq)=100 ! kg/m2 |
---|
| 328 | !endif |
---|
| 329 | ! DUST MMR |
---|
| 330 | if (txt.eq."dust_mass") then |
---|
| 331 | !look for a "profile_dust_mass" input file |
---|
| 332 | open(91,file='profile_dust_mass',status='old', |
---|
| 333 | & form='formatted',iostat=ierr) |
---|
| 334 | if (ierr.eq.0) then |
---|
| 335 | read(91,*) qsurf(iq) |
---|
[1266] | 336 | do ilayer=1,nlayer |
---|
[38] | 337 | read(91,*) q(ilayer,iq) |
---|
| 338 | ! write(*,*) "l=",ilayer," q(ilayer,iq)=",q(ilayer,iq) |
---|
| 339 | enddo |
---|
| 340 | else |
---|
| 341 | write(*,*) "No profile_dust_mass file!" |
---|
| 342 | endif |
---|
| 343 | close(91) |
---|
| 344 | endif ! of if (txt.eq."dust_mass") |
---|
| 345 | ! DUST NUMBER |
---|
| 346 | if (txt.eq."dust_number") then |
---|
| 347 | !look for a "profile_dust_number" input file |
---|
| 348 | open(91,file='profile_dust_number',status='old', |
---|
| 349 | & form='formatted',iostat=ierr) |
---|
| 350 | if (ierr.eq.0) then |
---|
| 351 | read(91,*) qsurf(iq) |
---|
[1266] | 352 | do ilayer=1,nlayer |
---|
[38] | 353 | read(91,*) q(ilayer,iq) |
---|
| 354 | enddo |
---|
| 355 | else |
---|
| 356 | write(*,*) "No profile_dust_number file!" |
---|
| 357 | endif |
---|
| 358 | close(91) |
---|
| 359 | endif ! of if (txt.eq."dust_number") |
---|
[358] | 360 | ! NB: some more initializations (chemistry) is done later |
---|
| 361 | ! CCN MASS |
---|
| 362 | if (txt.eq."ccn_mass") then |
---|
| 363 | !look for a "profile_ccn_mass" input file |
---|
| 364 | open(91,file='profile_ccn_mass',status='old', |
---|
| 365 | & form='formatted',iostat=ierr) |
---|
| 366 | if (ierr.eq.0) then |
---|
| 367 | read(91,*) qsurf(iq) |
---|
[1266] | 368 | do ilayer=1,nlayer |
---|
[358] | 369 | read(91,*) q(ilayer,iq) |
---|
| 370 | enddo |
---|
| 371 | else |
---|
| 372 | write(*,*) "No profile_ccn_mass file!" |
---|
| 373 | endif |
---|
| 374 | close(91) |
---|
| 375 | endif ! of if (txt.eq."ccn_mass") |
---|
| 376 | ! CCN NUMBER |
---|
| 377 | if (txt.eq."ccn_number") then |
---|
| 378 | !look for a "profile_ccn_number" input file |
---|
| 379 | open(91,file='profile_ccn_number',status='old', |
---|
| 380 | & form='formatted',iostat=ierr) |
---|
| 381 | if (ierr.eq.0) then |
---|
| 382 | read(91,*) qsurf(iq) |
---|
[1266] | 383 | do ilayer=1,nlayer |
---|
[358] | 384 | read(91,*) q(ilayer,iq) |
---|
| 385 | enddo |
---|
| 386 | else |
---|
| 387 | write(*,*) "No profile_ccn_number file!" |
---|
| 388 | endif |
---|
| 389 | close(91) |
---|
| 390 | endif ! of if (txt.eq."ccn_number") |
---|
[1036] | 391 | enddo ! of do iq=1,nq |
---|
[38] | 392 | |
---|
| 393 | else |
---|
[1036] | 394 | ! we still need to set (dummy) tracer number and names for physdem1 |
---|
| 395 | nq=1 |
---|
[1380] | 396 | nqtot=nq ! set value of nqtot (in infotrac module) as nq |
---|
[1036] | 397 | ! allocate arrays: |
---|
[1130] | 398 | allocate(tname(nq)) |
---|
[1266] | 399 | allocate(q(nlayer,nq)) |
---|
[1036] | 400 | allocate(qsurf(nq)) |
---|
[1266] | 401 | allocate(dq(nlayer,nq)) |
---|
| 402 | allocate(dqdyn(nlayer,nq)) |
---|
[1036] | 403 | allocate(mqtot(nq)) |
---|
[38] | 404 | do iq=1,nq |
---|
[1380] | 405 | write(str7,'(a1,i2.2)')'t',iq |
---|
[1130] | 406 | tname(iq)=str7 |
---|
[38] | 407 | enddo |
---|
[899] | 408 | ! and just to be clean, also initialize tracers to zero for physdem1 |
---|
| 409 | q(:,:)=0 |
---|
| 410 | qsurf(:)=0 |
---|
[38] | 411 | endif ! of if (tracer) |
---|
[358] | 412 | |
---|
| 413 | !write(*,*) "testphys1d q", q(1,:) |
---|
| 414 | !write(*,*) "testphys1d qsurf", qsurf |
---|
[38] | 415 | |
---|
[899] | 416 | c Date and local time at beginning of run |
---|
| 417 | c --------------------------------------- |
---|
| 418 | c Date (in sols since spring solstice) at beginning of run |
---|
[38] | 419 | day0 = 0 ! default value for day0 |
---|
| 420 | write(*,*) 'Initial date (in martian sols ; =0 at Ls=0)?' |
---|
| 421 | call getin("day0",day0) |
---|
| 422 | day=float(day0) |
---|
| 423 | write(*,*) " day0 = ",day0 |
---|
[899] | 424 | c Local time at beginning of run |
---|
[38] | 425 | time=0 ! default value for time |
---|
| 426 | write(*,*)'Initial local time (in hours, between 0 and 24)?' |
---|
| 427 | call getin("time",time) |
---|
| 428 | write(*,*)" time = ",time |
---|
| 429 | time=time/24.E+0 ! convert time (hours) to fraction of sol |
---|
| 430 | |
---|
[899] | 431 | c Discretization (Definition of grid and time steps) |
---|
[38] | 432 | c -------------- |
---|
| 433 | c |
---|
| 434 | nlevel=nlayer+1 |
---|
| 435 | nsoil=nsoilmx |
---|
| 436 | |
---|
| 437 | day_step=48 ! default value for day_step |
---|
| 438 | PRINT *,'Number of time steps per sol ?' |
---|
| 439 | call getin("day_step",day_step) |
---|
| 440 | write(*,*) " day_step = ",day_step |
---|
| 441 | |
---|
[1535] | 442 | ecritphy=day_step ! default value for ecritphy, output every time step |
---|
| 443 | |
---|
[38] | 444 | ndt=10 ! default value for ndt |
---|
| 445 | PRINT *,'Number of sols to run ?' |
---|
| 446 | call getin("ndt",ndt) |
---|
| 447 | write(*,*) " ndt = ",ndt |
---|
| 448 | |
---|
| 449 | ndt=ndt*day_step |
---|
| 450 | dtphys=daysec/day_step |
---|
[899] | 451 | |
---|
| 452 | c Imposed surface pressure |
---|
[38] | 453 | c ------------------------------------ |
---|
| 454 | c |
---|
| 455 | psurf=610. ! default value for psurf |
---|
| 456 | PRINT *,'Surface pressure (Pa) ?' |
---|
| 457 | call getin("psurf",psurf) |
---|
| 458 | write(*,*) " psurf = ",psurf |
---|
[899] | 459 | c Reference pressures |
---|
| 460 | pa=20. ! transition pressure (for hybrid coord.) |
---|
| 461 | preff=610. ! reference surface pressure |
---|
[38] | 462 | |
---|
[899] | 463 | c Aerosol properties |
---|
[38] | 464 | c -------------------------------- |
---|
[899] | 465 | tauvis=0.2 ! default value for tauvis (dust opacity) |
---|
[627] | 466 | write(*,'("Reference dust opacity at ",f4.0," Pa ?")')odpref |
---|
[38] | 467 | call getin("tauvis",tauvis) |
---|
| 468 | write(*,*) " tauvis = ",tauvis |
---|
[720] | 469 | |
---|
| 470 | c Orbital parameters |
---|
| 471 | c ------------------ |
---|
[2228] | 472 | paleomars=.false. ! Default: no water ice reservoir |
---|
| 473 | call getin("paleomars",paleomars) |
---|
| 474 | if (paleomars==.true.) then |
---|
| 475 | write(*,*) "paleomars=", paleomars |
---|
| 476 | write(*,*) "Orbital parameters from callphys.def" |
---|
| 477 | write(*,*) "Enter eccentricity & Lsperi" |
---|
| 478 | print *, 'Martian eccentricity (0<e<1) ?' |
---|
| 479 | call getin('excentric ',excentric) |
---|
| 480 | write(*,*)"excentric =",excentric |
---|
| 481 | print *, 'Solar longitude of perihelion (0<Ls<360) ?' |
---|
| 482 | call getin('Lsperi',Lsperi ) |
---|
| 483 | write(*,*)"Lsperi=",Lsperi |
---|
| 484 | Lsperi = Lsperi*pi/180.0 ! Put it in rad for peri_day |
---|
| 485 | periheli = halfaxe*(1-excentric) |
---|
| 486 | aphelie = halfaxe*(1+excentric) |
---|
| 487 | call call_dayperi(Lsperi,excentric,peri_day,year_day) |
---|
| 488 | write(*,*) "Corresponding orbital params for GCM" |
---|
| 489 | write(*,*) " periheli = ",periheli |
---|
| 490 | write(*,*) " aphelie = ",aphelie |
---|
| 491 | write(*,*) "date of perihelion (sol)",peri_day |
---|
| 492 | else |
---|
| 493 | write(*,*) "paleomars=", paleomars |
---|
| 494 | write(*,*) "Default present-day orbital parameters" |
---|
| 495 | write(*,*) "Unless specified otherwise" |
---|
| 496 | print *,'Min. distance Sun-Mars (Mkm)?' |
---|
| 497 | call getin("periheli",periheli) |
---|
| 498 | write(*,*) " periheli = ",periheli |
---|
[720] | 499 | |
---|
[2228] | 500 | print *,'Max. distance Sun-Mars (Mkm)?' |
---|
| 501 | call getin("aphelie",aphelie) |
---|
| 502 | write(*,*) " aphelie = ",aphelie |
---|
[720] | 503 | |
---|
[2228] | 504 | print *,'Day of perihelion?' |
---|
| 505 | call getin("periday",peri_day) |
---|
| 506 | write(*,*) " periday = ",peri_day |
---|
[720] | 507 | |
---|
[2228] | 508 | print *,'Obliquity?' |
---|
| 509 | call getin("obliquit",obliquit) |
---|
| 510 | write(*,*) " obliquit = ",obliquit |
---|
| 511 | end if |
---|
[38] | 512 | |
---|
| 513 | c latitude/longitude |
---|
| 514 | c ------------------ |
---|
[1311] | 515 | latitude(1)=0 ! default value for latitude |
---|
[38] | 516 | PRINT *,'latitude (in degrees) ?' |
---|
[1311] | 517 | call getin("latitude",latitude(1)) |
---|
[1047] | 518 | write(*,*) " latitude = ",latitude |
---|
| 519 | latitude=latitude*pi/180.E+0 |
---|
| 520 | longitude=0.E+0 |
---|
| 521 | longitude=longitude*pi/180.E+0 |
---|
[38] | 522 | |
---|
[1224] | 523 | ! some initializations (some of which have already been |
---|
[1047] | 524 | ! done above!) and loads parameters set in callphys.def |
---|
| 525 | ! and allocates some arrays |
---|
| 526 | !Mars possible matter with dtphys in input and include!!! |
---|
[1535] | 527 | ! Initializations below should mimick what is done in iniphysiq for 3D GCM |
---|
[1543] | 528 | call init_physics_distribution(regular_lonlat,4, |
---|
| 529 | & 1,1,1,nlayer,1) |
---|
| 530 | call init_interface_dyn_phys |
---|
[1535] | 531 | call init_regular_lonlat(1,1,longitude,latitude, |
---|
| 532 | & (/0.,0./),(/0.,0./)) |
---|
[1543] | 533 | call init_geometry(1,longitude,latitude, |
---|
| 534 | & (/0.,0.,0.,0./),(/0.,0.,0.,0./), |
---|
| 535 | & cell_area) |
---|
[1621] | 536 | ! Ehouarn: init_vertial_layers called later (because disvert not called yet) |
---|
| 537 | ! call init_vertical_layers(nlayer,preff,scaleheight, |
---|
| 538 | ! & ap,bp,aps,bps,presnivs,pseudoalt) |
---|
[1543] | 539 | call init_dimphy(1,nlayer) ! Initialize dimphy module |
---|
[1621] | 540 | call phys_state_var_init(1,llm,nq,tname, |
---|
[1524] | 541 | . day0,time,daysec,dtphys,rad,g,r,cpp) |
---|
[1311] | 542 | call ini_fillgeom(1,latitude,longitude,(/1.0/)) |
---|
[1246] | 543 | call conf_phys(1,llm,nq) |
---|
[1047] | 544 | |
---|
[1550] | 545 | ! in 1D model physics are called every time step |
---|
| 546 | ! ovverride iphysiq value that has been set by conf_phys |
---|
| 547 | if (iphysiq/=1) then |
---|
| 548 | write(*,*) "testphys1d: setting iphysiq=1" |
---|
| 549 | iphysiq=1 |
---|
| 550 | endif |
---|
[1047] | 551 | |
---|
[899] | 552 | c Initialize albedo / soil thermal inertia |
---|
| 553 | c ---------------------------------------- |
---|
[38] | 554 | c |
---|
| 555 | albedodat(1)=0.2 ! default value for albedodat |
---|
| 556 | PRINT *,'Albedo of bare ground ?' |
---|
| 557 | call getin("albedo",albedodat(1)) |
---|
| 558 | write(*,*) " albedo = ",albedodat(1) |
---|
[1944] | 559 | albedo(1,1)=albedodat(1) |
---|
[38] | 560 | |
---|
| 561 | inertiedat(1,1)=400 ! default value for inertiedat |
---|
| 562 | PRINT *,'Soil thermal inertia (SI) ?' |
---|
| 563 | call getin("inertia",inertiedat(1,1)) |
---|
| 564 | write(*,*) " inertia = ",inertiedat(1,1) |
---|
[224] | 565 | |
---|
| 566 | z0(1)=z0_default ! default value for roughness |
---|
| 567 | write(*,*) 'Surface roughness length z0 (m)?' |
---|
| 568 | call getin("z0",z0(1)) |
---|
| 569 | write(*,*) " z0 = ",z0(1) |
---|
[899] | 570 | |
---|
| 571 | ! Initialize local slope parameters (only matters if "callslope" |
---|
| 572 | ! is .true. in callphys.def) |
---|
| 573 | ! slope inclination angle (deg) 0: horizontal, 90: vertical |
---|
| 574 | theta_sl(1)=0.0 ! default: no inclination |
---|
| 575 | call getin("slope_inclination",theta_sl(1)) |
---|
| 576 | ! slope orientation (deg) |
---|
| 577 | ! 0 == Northward, 90 == Eastward, 180 == Southward, 270 == Westward |
---|
| 578 | psi_sl(1)=0.0 ! default value |
---|
| 579 | call getin("slope_orientation",psi_sl(1)) |
---|
| 580 | |
---|
[38] | 581 | c |
---|
[899] | 582 | c for the gravity wave scheme |
---|
[38] | 583 | c --------------------------------- |
---|
| 584 | c |
---|
| 585 | zmea(1)=0.E+0 |
---|
| 586 | zstd(1)=0.E+0 |
---|
| 587 | zsig(1)=0.E+0 |
---|
| 588 | zgam(1)=0.E+0 |
---|
[2199] | 589 | zthe(1)=0.E+0 |
---|
| 590 | c |
---|
| 591 | c for the slope wind scheme |
---|
| 592 | c --------------------------------- |
---|
| 593 | c |
---|
[1974] | 594 | hmons(1)=0.E+0 |
---|
[2199] | 595 | PRINT *,'hmons is initialized to ',hmons(1) |
---|
[2079] | 596 | summit(1)=0.E+0 |
---|
[2199] | 597 | PRINT *,'summit is initialized to ',summit(1) |
---|
[2079] | 598 | base(1)=0.E+0 |
---|
[2199] | 599 | c |
---|
| 600 | c Default values initializing the coefficients calculated later |
---|
| 601 | c --------------------------------- |
---|
| 602 | c |
---|
| 603 | tauscaling(1)=1. ! calculated in aeropacity_mod.F |
---|
| 604 | totcloudfrac(1)=1. ! calculated in watercloud_mod.F |
---|
| 605 | |
---|
[899] | 606 | c Specific initializations for "physiq" |
---|
| 607 | c ------------------------------------- |
---|
| 608 | c surface geopotential is not used (or useful) since in 1D |
---|
| 609 | c everything is controled by surface pressure |
---|
[38] | 610 | phisfi(1)=0.E+0 |
---|
| 611 | |
---|
[899] | 612 | c Initialization to take into account prescribed winds |
---|
[38] | 613 | c ------------------------------------------------------ |
---|
| 614 | ptif=2.E+0*omeg*sinlat(1) |
---|
| 615 | |
---|
[899] | 616 | c geostrophic wind |
---|
[38] | 617 | gru=10. ! default value for gru |
---|
| 618 | PRINT *,'zonal eastward component of the geostrophic wind (m/s) ?' |
---|
| 619 | call getin("u",gru) |
---|
| 620 | write(*,*) " u = ",gru |
---|
| 621 | grv=0. !default value for grv |
---|
| 622 | PRINT *,'meridional northward component of the geostrophic', |
---|
| 623 | &' wind (m/s) ?' |
---|
| 624 | call getin("v",grv) |
---|
| 625 | write(*,*) " v = ",grv |
---|
| 626 | |
---|
[899] | 627 | c Initialize winds for first time step |
---|
[38] | 628 | DO ilayer=1,nlayer |
---|
| 629 | u(ilayer)=gru |
---|
| 630 | v(ilayer)=grv |
---|
[1541] | 631 | w(ilayer)=0 ! default: no vertical wind |
---|
[38] | 632 | ENDDO |
---|
| 633 | |
---|
[899] | 634 | c Initialize turbulente kinetic energy |
---|
[38] | 635 | DO ilevel=1,nlevel |
---|
| 636 | q2(ilevel)=0.E+0 |
---|
| 637 | ENDDO |
---|
| 638 | |
---|
[899] | 639 | c CO2 ice on the surface |
---|
[38] | 640 | c ------------------- |
---|
[1130] | 641 | co2ice(1)=0.E+0 ! default value for co2ice |
---|
[38] | 642 | PRINT *,'Initial CO2 ice on the surface (kg.m-2)' |
---|
| 643 | call getin("co2ice",co2ice) |
---|
| 644 | write(*,*) " co2ice = ",co2ice |
---|
[1944] | 645 | ! Initialization for CO2 clouds (could be improved to read initial profiles) |
---|
| 646 | mem_Mccn_co2(:,:)=0 |
---|
| 647 | mem_Mh2o_co2(:,:)=0 |
---|
| 648 | mem_Nccn_co2(:,:)=0 |
---|
[38] | 649 | c |
---|
[899] | 650 | c emissivity |
---|
[38] | 651 | c ---------- |
---|
| 652 | emis=emissiv |
---|
[1130] | 653 | IF (co2ice(1).eq.1.E+0) THEN |
---|
[38] | 654 | emis=emisice(1) ! northern hemisphere |
---|
[1541] | 655 | IF(latitude(1).LT.0) emis=emisice(2) ! southern hemisphere |
---|
[38] | 656 | ENDIF |
---|
| 657 | |
---|
| 658 | |
---|
| 659 | |
---|
[899] | 660 | c Compute pressures and altitudes of atmospheric levels |
---|
[38] | 661 | c ---------------------------------------------------------------- |
---|
| 662 | |
---|
| 663 | c Vertical Coordinates |
---|
| 664 | c """""""""""""""""""" |
---|
| 665 | hybrid=.true. |
---|
| 666 | PRINT *,'Hybrid coordinates ?' |
---|
| 667 | call getin("hybrid",hybrid) |
---|
| 668 | write(*,*) " hybrid = ", hybrid |
---|
| 669 | |
---|
[2167] | 670 | CALL disvert_noterre |
---|
[1621] | 671 | ! now that disvert has been called, initialize module vertical_layers_mod |
---|
| 672 | call init_vertical_layers(nlayer,preff,scaleheight, |
---|
| 673 | & ap,bp,aps,bps,presnivs,pseudoalt) |
---|
[38] | 674 | |
---|
| 675 | DO ilevel=1,nlevel |
---|
| 676 | plev(ilevel)=ap(ilevel)+psurf*bp(ilevel) |
---|
| 677 | ENDDO |
---|
| 678 | |
---|
| 679 | DO ilayer=1,nlayer |
---|
| 680 | play(ilayer)=aps(ilayer)+psurf*bps(ilayer) |
---|
| 681 | ENDDO |
---|
| 682 | |
---|
| 683 | DO ilayer=1,nlayer |
---|
| 684 | zlay(ilayer)=-200.E+0 *r*log(play(ilayer)/plev(1)) |
---|
| 685 | & /g |
---|
| 686 | ENDDO |
---|
| 687 | |
---|
| 688 | |
---|
[899] | 689 | c Initialize temperature profile |
---|
[38] | 690 | c -------------------------------------- |
---|
| 691 | pks=psurf**rcp |
---|
| 692 | |
---|
[899] | 693 | c altitude in km in profile: divide zlay by 1000 |
---|
[38] | 694 | tmp1(0)=0.E+0 |
---|
| 695 | DO ilayer=1,nlayer |
---|
| 696 | tmp1(ilayer)=zlay(ilayer)/1000.E+0 |
---|
| 697 | ENDDO |
---|
| 698 | |
---|
| 699 | call profile(nlayer+1,tmp1,tmp2) |
---|
| 700 | |
---|
| 701 | tsurf=tmp2(0) |
---|
| 702 | DO ilayer=1,nlayer |
---|
| 703 | temp(ilayer)=tmp2(ilayer) |
---|
| 704 | ENDDO |
---|
| 705 | |
---|
| 706 | |
---|
| 707 | |
---|
| 708 | ! Initialize soil properties and temperature |
---|
| 709 | ! ------------------------------------------ |
---|
| 710 | volcapa=1.e6 ! volumetric heat capacity |
---|
| 711 | DO isoil=1,nsoil |
---|
| 712 | inertiedat(1,isoil)=inertiedat(1,1) ! soil thermal inertia |
---|
[1130] | 713 | tsoil(isoil)=tsurf(1) ! soil temperature |
---|
[38] | 714 | ENDDO |
---|
| 715 | |
---|
| 716 | ! Initialize depths |
---|
| 717 | ! ----------------- |
---|
| 718 | do isoil=0,nsoil-1 |
---|
| 719 | mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) ! mid-layer depth |
---|
| 720 | enddo |
---|
| 721 | do isoil=1,nsoil |
---|
| 722 | layer(isoil)=2.e-4*(2.**(isoil-1)) ! layer depth |
---|
| 723 | enddo |
---|
| 724 | |
---|
[899] | 725 | c Initialize traceurs |
---|
[38] | 726 | c --------------------------- |
---|
| 727 | |
---|
| 728 | if (photochem.or.callthermos) then |
---|
[899] | 729 | write(*,*) 'Initializing chemical species' |
---|
| 730 | ! thermo=0: initialize over all atmospheric layers |
---|
[38] | 731 | thermo=0 |
---|
[1541] | 732 | call inichim_newstart(q,psurf,sig,nq,latitude,longitude, |
---|
| 733 | $ cell_area,thermo,qsurf) |
---|
[38] | 734 | endif |
---|
[520] | 735 | |
---|
[899] | 736 | c Check if the surface is a water ice reservoir |
---|
[520] | 737 | c -------------------------------------------------- |
---|
[1047] | 738 | watercaptag(1)=.false. ! Default: no water ice reservoir |
---|
[520] | 739 | print *,'Water ice cap on ground ?' |
---|
| 740 | call getin("watercaptag",watercaptag) |
---|
| 741 | write(*,*) " watercaptag = ",watercaptag |
---|
[38] | 742 | |
---|
| 743 | |
---|
[899] | 744 | c Initialization for GRADS outputs in "g1d.dat" and "g1d.ctl" |
---|
[38] | 745 | c ---------------------------------------------------------------- |
---|
[899] | 746 | c (output done in "writeg1d", typically called by "physiq.F") |
---|
[38] | 747 | |
---|
| 748 | g1d_nlayer=nlayer |
---|
| 749 | g1d_nomfich='g1d.dat' |
---|
| 750 | g1d_unitfich=40 |
---|
| 751 | g1d_nomctl='g1d.ctl' |
---|
| 752 | g1d_unitctl=41 |
---|
| 753 | g1d_premier=.true. |
---|
| 754 | g2d_premier=.true. |
---|
| 755 | |
---|
[899] | 756 | c Write a "startfi" file |
---|
[38] | 757 | c -------------------- |
---|
[899] | 758 | c This file will be read during the first call to "physiq". |
---|
| 759 | c It is needed to transfert physics variables to "physiq"... |
---|
[38] | 760 | |
---|
[1541] | 761 | call physdem0("startfi.nc",longitude,latitude,nsoilmx,ngrid,llm, |
---|
[2214] | 762 | & nq,dtphys,float(day0),0.,cell_area, |
---|
[2079] | 763 | & albedodat,inertiedat,zmea,zstd,zsig,zgam,zthe, |
---|
| 764 | & hmons,summit,base) |
---|
[1112] | 765 | call physdem1("startfi.nc",nsoilmx,ngrid,llm,nq, |
---|
[1944] | 766 | & dtphys,time, |
---|
| 767 | & tsurf,tsoil,co2ice,albedo,emis,q2,qsurf,tauscaling, |
---|
| 768 | & totcloudfrac,wstar, |
---|
| 769 | & mem_Mccn_co2,mem_Nccn_co2, |
---|
| 770 | & mem_Mh2o_co2) |
---|
[899] | 771 | |
---|
[38] | 772 | c======================================================================= |
---|
[899] | 773 | c 1D MODEL TIME STEPPING LOOP |
---|
[38] | 774 | c======================================================================= |
---|
| 775 | c |
---|
| 776 | firstcall=.true. |
---|
| 777 | lastcall=.false. |
---|
| 778 | |
---|
| 779 | DO idt=1,ndt |
---|
| 780 | c IF (idt.eq.ndt) lastcall=.true. |
---|
| 781 | IF (idt.eq.ndt-day_step-1) then !test |
---|
| 782 | lastcall=.true. |
---|
| 783 | call solarlong(day*1.0,zls) |
---|
| 784 | write(103,*) 'Ls=',zls*180./pi |
---|
[1541] | 785 | write(103,*) 'Lat=', latitude(1)*180./pi |
---|
[627] | 786 | write(103,*) 'Tau=', tauvis/odpref*psurf |
---|
[38] | 787 | write(103,*) 'RunEnd - Atmos. Temp. File' |
---|
| 788 | write(103,*) 'RunEnd - Atmos. Temp. File' |
---|
| 789 | write(104,*) 'Ls=',zls*180./pi |
---|
[1541] | 790 | write(104,*) 'Lat=', latitude(1) |
---|
[627] | 791 | write(104,*) 'Tau=', tauvis/odpref*psurf |
---|
[38] | 792 | write(104,*) 'RunEnd - Atmos. Temp. File' |
---|
| 793 | ENDIF |
---|
| 794 | |
---|
[899] | 795 | c compute geopotential |
---|
[38] | 796 | c ~~~~~~~~~~~~~~~~~~~~~ |
---|
| 797 | DO ilayer=1,nlayer |
---|
| 798 | s(ilayer)=(aps(ilayer)/psurf+bps(ilayer))**rcp |
---|
| 799 | h(ilayer)=cpp*temp(ilayer)/(pks*s(ilayer)) |
---|
| 800 | ENDDO |
---|
| 801 | phi(1)=pks*h(1)*(1.E+0-s(1)) |
---|
| 802 | DO ilayer=2,nlayer |
---|
| 803 | phi(ilayer)=phi(ilayer-1)+ |
---|
| 804 | & pks*(h(ilayer-1)+h(ilayer))*.5E+0 |
---|
| 805 | & *(s(ilayer-1)-s(ilayer)) |
---|
| 806 | |
---|
| 807 | ENDDO |
---|
| 808 | |
---|
[899] | 809 | c call physics |
---|
[38] | 810 | c -------------------- |
---|
[358] | 811 | ! write(*,*) "testphys1d avant q", q(1,:) |
---|
[1036] | 812 | CALL physiq (1,llm,nq, |
---|
[38] | 813 | , firstcall,lastcall, |
---|
| 814 | , day,time,dtphys, |
---|
| 815 | , plev,play,phi, |
---|
| 816 | , u, v,temp, q, |
---|
| 817 | , w, |
---|
[899] | 818 | C - outputs |
---|
[1576] | 819 | s du, dv, dtemp, dq,dpsurf) |
---|
[358] | 820 | ! write(*,*) "testphys1d apres q", q(1,:) |
---|
[38] | 821 | |
---|
[358] | 822 | |
---|
[899] | 823 | c wind increment : specific for 1D |
---|
| 824 | c -------------------------------- |
---|
[38] | 825 | |
---|
[899] | 826 | c The physics compute the tendencies on u and v, |
---|
| 827 | c here we just add Coriolos effect |
---|
[38] | 828 | c |
---|
| 829 | c DO ilayer=1,nlayer |
---|
| 830 | c du(ilayer)=du(ilayer)+ptif*(v(ilayer)-grv) |
---|
| 831 | c dv(ilayer)=dv(ilayer)+ptif*(-u(ilayer)+gru) |
---|
| 832 | c ENDDO |
---|
| 833 | |
---|
[899] | 834 | c For some tests : No coriolis force at equator |
---|
[1541] | 835 | c if(latitude(1).eq.0.) then |
---|
[38] | 836 | DO ilayer=1,nlayer |
---|
| 837 | du(ilayer)=du(ilayer)+ (gru-u(ilayer))/1.e4 |
---|
| 838 | dv(ilayer)=dv(ilayer)+ (grv-v(ilayer))/1.e4 |
---|
| 839 | ENDDO |
---|
| 840 | c end if |
---|
| 841 | c |
---|
| 842 | c |
---|
[899] | 843 | c Compute time for next time step |
---|
[38] | 844 | c --------------------------------------- |
---|
| 845 | firstcall=.false. |
---|
| 846 | time=time+dtphys/daysec |
---|
| 847 | IF (time.gt.1.E+0) then |
---|
| 848 | time=time-1.E+0 |
---|
| 849 | day=day+1 |
---|
| 850 | ENDIF |
---|
| 851 | |
---|
[899] | 852 | c compute winds and temperature for next time step |
---|
[38] | 853 | c ---------------------------------------------------------- |
---|
| 854 | |
---|
| 855 | DO ilayer=1,nlayer |
---|
| 856 | u(ilayer)=u(ilayer)+dtphys*du(ilayer) |
---|
| 857 | v(ilayer)=v(ilayer)+dtphys*dv(ilayer) |
---|
| 858 | temp(ilayer)=temp(ilayer)+dtphys*dtemp(ilayer) |
---|
| 859 | ENDDO |
---|
| 860 | |
---|
[899] | 861 | c compute pressure for next time step |
---|
[38] | 862 | c ---------------------------------------------------------- |
---|
| 863 | |
---|
[1549] | 864 | psurf=psurf+dtphys*dpsurf(1) ! surface pressure change |
---|
[38] | 865 | DO ilevel=1,nlevel |
---|
| 866 | plev(ilevel)=ap(ilevel)+psurf*bp(ilevel) |
---|
| 867 | ENDDO |
---|
| 868 | DO ilayer=1,nlayer |
---|
| 869 | play(ilayer)=aps(ilayer)+psurf*bps(ilayer) |
---|
| 870 | ENDDO |
---|
| 871 | |
---|
[899] | 872 | ! increment tracers |
---|
[1036] | 873 | DO iq = 1, nq |
---|
[38] | 874 | DO ilayer=1,nlayer |
---|
| 875 | q(ilayer,iq)=q(ilayer,iq)+dtphys*dq(ilayer,iq) |
---|
| 876 | ENDDO |
---|
| 877 | ENDDO |
---|
| 878 | |
---|
[899] | 879 | ENDDO ! of idt=1,ndt ! end of time stepping loop |
---|
[38] | 880 | |
---|
| 881 | c ======================================================== |
---|
[899] | 882 | c OUTPUTS |
---|
[38] | 883 | c ======================================================== |
---|
| 884 | |
---|
[899] | 885 | c finalize and close grads files "g1d.dat" and "g1d.ctl" |
---|
[38] | 886 | |
---|
| 887 | c CALL endg1d(1,nlayer,zphi/(g*1000.),ndt) |
---|
| 888 | CALL endg1d(1,nlayer,zlay/1000.,ndt) |
---|
| 889 | |
---|
[1380] | 890 | write(*,*) "testphys1d: Everything is cool." |
---|
| 891 | |
---|
[38] | 892 | END |
---|
| 893 | |
---|
| 894 | c*********************************************************************** |
---|
| 895 | c*********************************************************************** |
---|
[899] | 896 | c Dummy subroutines used only in 3D, but required to |
---|
| 897 | c compile testphys1d (to cleanly use writediagfi) |
---|
[38] | 898 | |
---|
[899] | 899 | subroutine gr_fi_dyn(nfield,ngrid,im,jm,pfi,pdyn) |
---|
| 900 | |
---|
| 901 | IMPLICIT NONE |
---|
| 902 | |
---|
| 903 | INTEGER im,jm,ngrid,nfield |
---|
| 904 | REAL pdyn(im,jm,nfield) |
---|
| 905 | REAL pfi(ngrid,nfield) |
---|
| 906 | |
---|
| 907 | if (ngrid.ne.1) then |
---|
| 908 | write(*,*) "gr_fi_dyn error: in 1D ngrid should be 1!!!" |
---|
| 909 | stop |
---|
| 910 | endif |
---|
| 911 | |
---|
| 912 | pdyn(1,1,1:nfield)=pfi(1,1:nfield) |
---|
| 913 | |
---|
| 914 | end |
---|
[38] | 915 | |
---|
| 916 | c*********************************************************************** |
---|
| 917 | c*********************************************************************** |
---|
| 918 | |
---|