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