| 1 | module init_phys_1d_mod |
|---|
| 2 | |
|---|
| 3 | implicit none |
|---|
| 4 | |
|---|
| 5 | contains |
|---|
| 6 | |
|---|
| 7 | #ifdef CPP_1D |
|---|
| 8 | |
|---|
| 9 | subroutine init_phys_1d(llm_in,nqtot_in,v,u,temp,q,psurf,time,ap_out,bp_out) |
|---|
| 10 | |
|---|
| 11 | !-------------- init_phys_1d ------------- |
|---|
| 12 | ! |
|---|
| 13 | ! This function is a copy of the test_phys_1d.F90 initialisation part. |
|---|
| 14 | ! It is meant to initialize all the variables in modules and for the ouput |
|---|
| 15 | ! that can't be initialzed via the startfi file. |
|---|
| 16 | ! |
|---|
| 17 | ! Basically it reads run.def and save variables values |
|---|
| 18 | ! |
|---|
| 19 | !-------------- ------------- |
|---|
| 20 | |
|---|
| 21 | USE ioipsl_getincom, only: getin |
|---|
| 22 | use comcstfi_h, only: pi, rad, omeg, g, mugaz, rcp, r, cpp |
|---|
| 23 | use time_phylmdz_mod, only: daysec, dtphys, day_step, & |
|---|
| 24 | ecritphy, iphysiq |
|---|
| 25 | use planete_h, only: year_day, periheli, aphelie, peri_day,& |
|---|
| 26 | obliquit, emin_turb, lmixmin |
|---|
| 27 | use surfdat_h, only: albedodat, z0_default, emissiv, emisice,& |
|---|
| 28 | albedice, iceradius, dtemisice, z0,& |
|---|
| 29 | zmea, zstd, zsig, zgam, zthe, phisfi,& |
|---|
| 30 | watercaptag, watercap, hmons, summit, base,& |
|---|
| 31 | tsurf, emis,qsurf |
|---|
| 32 | use infotrac, only: nqtot, tname, nqperes,nqfils |
|---|
| 33 | use regular_lonlat_mod, only: init_regular_lonlat |
|---|
| 34 | use mod_grid_phy_lmdz, only : regular_lonlat |
|---|
| 35 | USE comvert_mod, ONLY: ap,bp,aps,bps,pa,preff,sig,& |
|---|
| 36 | presnivs,pseudoalt,scaleheight |
|---|
| 37 | use dimradmars_mod, only: tauvis,totcloudfrac |
|---|
| 38 | use mod_interface_dyn_phys, only: init_interface_dyn_phys |
|---|
| 39 | use geometry_mod, only: init_geometry |
|---|
| 40 | use dimphy, only : init_dimphy |
|---|
| 41 | USE phys_state_var_init_mod, ONLY: phys_state_var_init |
|---|
| 42 | use comgeomfi_h, only: sinlat, ini_fillgeom |
|---|
| 43 | use slope_mod, only: theta_sl, psi_sl |
|---|
| 44 | use comslope_mod, only: def_slope,subslope_dist,def_slope_mean |
|---|
| 45 | use dimradmars_mod, only: tauvis,totcloudfrac |
|---|
| 46 | use dust_param_mod, only: tauscaling |
|---|
| 47 | use tracer_mod, only: igcm_h2o_vap,igcm_h2o_ice,igcm_co2,noms |
|---|
| 48 | USE logic_mod, ONLY: hybrid |
|---|
| 49 | USE vertical_layers_mod, ONLY: init_vertical_layers |
|---|
| 50 | use comsoil_h, only: volcapa, layer, mlayer, inertiedat, & |
|---|
| 51 | inertiesoil,nsoilmx,flux_geo |
|---|
| 52 | USE read_profile_mod, ONLY: read_profile |
|---|
| 53 | use inichim_newstart_mod, only: inichim_newstart |
|---|
| 54 | use physics_distribution_mod, only: init_physics_distribution |
|---|
| 55 | use iostart, only: open_startphy,get_var, close_startphy |
|---|
| 56 | #ifdef CPP_XIOS |
|---|
| 57 | use mod_const_mpi, only: COMM_LMDZ |
|---|
| 58 | #endif |
|---|
| 59 | include "dimensions.h" |
|---|
| 60 | integer, parameter :: ngrid = 1 !(2+(jjm-1)*iim - 1/jjm) |
|---|
| 61 | integer, parameter :: nlayer = llm |
|---|
| 62 | include "callkeys.h" |
|---|
| 63 | include "netcdf.inc" |
|---|
| 64 | |
|---|
| 65 | !-------------- INPUT VARIABLES ------------- |
|---|
| 66 | |
|---|
| 67 | integer, intent(in) :: llm_in,nqtot_in |
|---|
| 68 | |
|---|
| 69 | !-------------- OUTPUT VARIABLES ------------- |
|---|
| 70 | |
|---|
| 71 | REAL,INTENT(OUT):: u(nlayer),v(nlayer) ! zonal, meridional wind |
|---|
| 72 | REAL,INTENT(OUT):: temp(nlayer) ! temperature at the middle of the layers |
|---|
| 73 | REAL,INTENT(OUT) :: q(2,llm_in,nqtot_in) ! tracer mixing ratio (e.g. kg/kg) |
|---|
| 74 | REAL,INTENT(OUT):: psurf(2) |
|---|
| 75 | REAL,INTENT(OUT):: time ! time (0<time<1 ; time=0.5 a midi) |
|---|
| 76 | REAL,INTENT(OUT) :: ap_out(llm_in+1),bp_out(llm_in+1) |
|---|
| 77 | |
|---|
| 78 | |
|---|
| 79 | !-------------- LOCAL VARIABLES ------------- |
|---|
| 80 | |
|---|
| 81 | REAL halfaxe, excentric, Lsperi |
|---|
| 82 | integer :: nq=1 ! number of tracers |
|---|
| 83 | real :: latitude(1), longitude(1), cell_area(1) |
|---|
| 84 | |
|---|
| 85 | ! MVals: isotopes as in the dynamics (CRisi) |
|---|
| 86 | INTEGER :: ifils,ipere,generation |
|---|
| 87 | CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_transp ! transporting fluid short name |
|---|
| 88 | CHARACTER(len=80) :: line ! to store a line of text |
|---|
| 89 | INTEGER ierr0 |
|---|
| 90 | LOGICAL :: continu |
|---|
| 91 | |
|---|
| 92 | logical :: present,found |
|---|
| 93 | INTEGER nlevel,nsoil,ndt |
|---|
| 94 | REAL gru,grv ! prescribed "geostrophic" background wind |
|---|
| 95 | REAL pks, ptif, w(nlayer) |
|---|
| 96 | REAL q2(nlayer+1) ! Turbulent Kinetic Energy |
|---|
| 97 | REAL play(nlayer) ! Pressure at the middle of the layers (Pa) |
|---|
| 98 | REAL plev(nlayer+1) ! intermediate pressure levels (pa) |
|---|
| 99 | REAL zlay(nlayer) ! altitude estimee dans les couches (km) |
|---|
| 100 | REAL tmp1(0:nlayer),tmp2(0:nlayer) |
|---|
| 101 | INTEGER flagthermo,flagh2o |
|---|
| 102 | REAL atm_wat_profile |
|---|
| 103 | REAL,ALLOCATABLE :: zqsat(:) ! useful for (atm_wat_profile=2) |
|---|
| 104 | |
|---|
| 105 | INTEGER :: ierr,iq,ilayer,ilevel,isoil |
|---|
| 106 | INTEGER day0,dayn ! date initial (sol ; =0 a Ls=0) and final |
|---|
| 107 | REAL day ! date durant le run |
|---|
| 108 | REAL tab_cntrl(100) |
|---|
| 109 | real,parameter :: odpref=610. ! DOD reference pressure (Pa) |
|---|
| 110 | |
|---|
| 111 | ! LL: Subsurface geothermal flux |
|---|
| 112 | real :: flux_geo_tmp |
|---|
| 113 | |
|---|
| 114 | real,allocatable :: qdyn(:,:,:,:),psdyn(:,:) |
|---|
| 115 | |
|---|
| 116 | |
|---|
| 117 | ! check if 'run.def' file is around (otherwise reading parameters |
|---|
| 118 | ! from callphys.def via getin() routine won't work. |
|---|
| 119 | open(99,file='run.def',status='old',form='formatted',& |
|---|
| 120 | iostat=ierr) |
|---|
| 121 | if (ierr.ne.0) then |
|---|
| 122 | write(*,*) 'Cannot find required file "run.def"' |
|---|
| 123 | write(*,*) ' (which should contain some input parameters' |
|---|
| 124 | write(*,*) ' along with the following line:' |
|---|
| 125 | write(*,*) ' INCLUDEDEF=callphys.def' |
|---|
| 126 | write(*,*) ' )' |
|---|
| 127 | write(*,*) ' ... might as well stop here ...' |
|---|
| 128 | stop |
|---|
| 129 | else |
|---|
| 130 | close(99) |
|---|
| 131 | endif |
|---|
| 132 | |
|---|
| 133 | ! ------------------------------------------------------ |
|---|
| 134 | ! Prescribed constants to be set here |
|---|
| 135 | ! ------------------------------------------------------ |
|---|
| 136 | |
|---|
| 137 | ! if(.not. startfile_1D ) then |
|---|
| 138 | |
|---|
| 139 | pi=2.E+0*asin(1.E+0) |
|---|
| 140 | |
|---|
| 141 | ! Mars planetary constants |
|---|
| 142 | ! ---------------------------- |
|---|
| 143 | rad=3397200. ! mars radius (m) ~3397200 m |
|---|
| 144 | daysec=88775. ! length of a sol (s) ~88775 s |
|---|
| 145 | omeg=4.*asin(1.)/(daysec) ! rotation rate (rad.s-1) |
|---|
| 146 | g=3.72 ! gravity (m.s-2) ~3.72 |
|---|
| 147 | mugaz=43.49 ! atmosphere mola mass (g.mol-1) ~43.49 |
|---|
| 148 | rcp=.256793 ! = r/cp ~0.256793 |
|---|
| 149 | r= 8.314511E+0 *1000.E+0/mugaz |
|---|
| 150 | cpp= r/rcp |
|---|
| 151 | year_day = 669 ! lenght of year (sols) ~668.6 |
|---|
| 152 | periheli = 206.66 ! minimum sun-mars distance (Mkm) ~206.66 |
|---|
| 153 | aphelie = 249.22 ! maximum sun-mars distance (Mkm) ~249.22 |
|---|
| 154 | halfaxe = 227.94 ! demi-grand axe de l'ellipse |
|---|
| 155 | peri_day = 485. ! perihelion date (sols since N. Spring) |
|---|
| 156 | obliquit = 25.2 ! Obliquity (deg) ~25.2 |
|---|
| 157 | excentric = 0.0934 ! Eccentricity (0.0934) |
|---|
| 158 | |
|---|
| 159 | ! Planetary Boundary Layer and Turbulence parameters |
|---|
| 160 | ! -------------------------------------------------- |
|---|
| 161 | z0_default = 1.e-2 ! surface roughness (m) ~0.01 |
|---|
| 162 | emin_turb = 1.e-6 ! minimal turbulent energy ~1.e-8 |
|---|
| 163 | lmixmin = 30 ! mixing length ~100 |
|---|
| 164 | |
|---|
| 165 | ! cap properties and surface emissivities |
|---|
| 166 | ! ---------------------------------------------------- |
|---|
| 167 | emissiv= 0.95 ! Bare ground emissivity ~.95 |
|---|
| 168 | emisice(1)=0.95 ! Northern cap emissivity |
|---|
| 169 | emisice(2)=0.95 ! Southern cap emisssivity |
|---|
| 170 | albedice(1)=0.5 ! Northern cap albedo |
|---|
| 171 | albedice(2)=0.5 ! Southern cap albedo |
|---|
| 172 | iceradius(1) = 100.e-6 ! mean scat radius of CO2 snow (north) |
|---|
| 173 | iceradius(2) = 100.e-6 ! mean scat radius of CO2 snow (south) |
|---|
| 174 | dtemisice(1) = 2. ! time scale for snow metamorphism (north) |
|---|
| 175 | dtemisice(2) = 2. ! time scale for snow metamorphism (south) |
|---|
| 176 | |
|---|
| 177 | ! endif ! .not. startfile_1D |
|---|
| 178 | |
|---|
| 179 | ! mesh surface (not a very usefull quantity in 1D) |
|---|
| 180 | ! ---------------------------------------------------- |
|---|
| 181 | cell_area(1)=1.E+0 |
|---|
| 182 | |
|---|
| 183 | |
|---|
| 184 | ! check if there is a 'traceur.def' file |
|---|
| 185 | ! and process it. |
|---|
| 186 | ! load tracer names from file 'traceur.def' |
|---|
| 187 | open(90,file='traceur.def',status='old',form='formatted',& |
|---|
| 188 | iostat=ierr) |
|---|
| 189 | if (ierr.ne.0) then |
|---|
| 190 | write(*,*) 'Cannot find required file "traceur.def"' |
|---|
| 191 | write(*,*) ' If you want to run with tracers, I need it' |
|---|
| 192 | write(*,*) ' ... might as well stop here ...' |
|---|
| 193 | stop |
|---|
| 194 | else |
|---|
| 195 | nqtot=nqtot_in ! set value of nqtot (in infotrac module) as nq |
|---|
| 196 | nq=nqtot_in |
|---|
| 197 | endif |
|---|
| 198 | ! allocate arrays: |
|---|
| 199 | allocate(tname(nq)) |
|---|
| 200 | allocate(qsurf(1,1,nq)) |
|---|
| 201 | allocate(tnom_transp(nq)) |
|---|
| 202 | |
|---|
| 203 | ! read tracer names from file traceur.def |
|---|
| 204 | do iq=1,nq |
|---|
| 205 | read(90,'(80a)',iostat=ierr) line ! store the line from traceur.def |
|---|
| 206 | if (ierr.ne.0) then |
|---|
| 207 | write(*,*) 'testphys1d: error reading tracer names...' |
|---|
| 208 | stop |
|---|
| 209 | endif |
|---|
| 210 | ! if format is tnom_0, tnom_transp (isotopes) |
|---|
| 211 | read(line,*,iostat=ierr0) tname(iq),tnom_transp(iq) |
|---|
| 212 | if (ierr0.ne.0) then |
|---|
| 213 | read(line,*) tname(iq) |
|---|
| 214 | tnom_transp(iq)='air' |
|---|
| 215 | endif |
|---|
| 216 | |
|---|
| 217 | enddo |
|---|
| 218 | close(90) |
|---|
| 219 | |
|---|
| 220 | ! Isotopes: as in the 3D case we have to determine father/son relations for isotopes and carrying fluid |
|---|
| 221 | ALLOCATE(nqfils(nqtot)) |
|---|
| 222 | nqperes=0 |
|---|
| 223 | nqfils(:)=0 |
|---|
| 224 | DO iq=1,nqtot |
|---|
| 225 | if (tnom_transp(iq) == 'air') then |
|---|
| 226 | ! ceci est un traceur père |
|---|
| 227 | WRITE(*,*) 'Le traceur',iq,', appele ',& |
|---|
| 228 | trim(tname(iq)),', est un pere' |
|---|
| 229 | nqperes=nqperes+1 |
|---|
| 230 | else !if (tnom_transp(iq) == 'air') then |
|---|
| 231 | ! ceci est un fils. Qui est son père? |
|---|
| 232 | WRITE(*,*) 'Le traceur',iq,', appele ',& |
|---|
| 233 | trim(tname(iq)),', est un fils' |
|---|
| 234 | continu=.true. |
|---|
| 235 | ipere=1 |
|---|
| 236 | do while (continu) |
|---|
| 237 | if (tnom_transp(iq) .eq. tname(ipere)) then |
|---|
| 238 | ! Son père est ipere |
|---|
| 239 | WRITE(*,*) 'Le traceur',iq,'appele ',& |
|---|
| 240 | trim(tname(iq)),' est le fils de ',& |
|---|
| 241 | ipere,'appele ',trim(tname(ipere)) |
|---|
| 242 | nqfils(ipere)=nqfils(ipere)+1 |
|---|
| 243 | continu=.false. |
|---|
| 244 | else !if (tnom_transp(iq) == tnom_0(ipere)) then |
|---|
| 245 | ipere=ipere+1 |
|---|
| 246 | if (ipere.gt.nqtot) then |
|---|
| 247 | WRITE(*,*) 'Le traceur',iq,'appele ',& |
|---|
| 248 | trim(tname(iq)),', est orpelin.' |
|---|
| 249 | CALL abort_gcm('infotrac_init',& |
|---|
| 250 | 'Un traceur est orphelin',1) |
|---|
| 251 | endif !if (ipere.gt.nqtot) then |
|---|
| 252 | endif !if (tnom_transp(iq) == tnom_0(ipere)) then |
|---|
| 253 | enddo !do while (continu) |
|---|
| 254 | endif !if (tnom_transp(iq) == 'air') then |
|---|
| 255 | enddo !DO iq=1,nqtot |
|---|
| 256 | WRITE(*,*) 'nqperes=',nqperes |
|---|
| 257 | WRITE(*,*) 'nqfils=',nqfils |
|---|
| 258 | |
|---|
| 259 | ! initialize tracers here: |
|---|
| 260 | write(*,*) "testphys1d: initializing tracers" |
|---|
| 261 | call read_profile(nq, nlayer, qsurf, q(1,:,:)) |
|---|
| 262 | q(2,:,:)=q(1,:,:) |
|---|
| 263 | |
|---|
| 264 | #ifdef CPP_XIOS |
|---|
| 265 | call init_physics_distribution(regular_lonlat,4,& |
|---|
| 266 | 1,1,1,nlayer,COMM_LMDZ) |
|---|
| 267 | #else |
|---|
| 268 | call init_physics_distribution(regular_lonlat,4,& |
|---|
| 269 | 1,1,1,nlayer,1) |
|---|
| 270 | #endif |
|---|
| 271 | |
|---|
| 272 | call open_startphy("startfi_evol.nc") |
|---|
| 273 | call get_var("controle",tab_cntrl,found) |
|---|
| 274 | if (.not.found) then |
|---|
| 275 | call abort_physic("open_startphy", & |
|---|
| 276 | "tabfi: Failed reading <controle> array",1) |
|---|
| 277 | else |
|---|
| 278 | write(*,*)'tabfi: tab_cntrl',tab_cntrl |
|---|
| 279 | endif |
|---|
| 280 | day0 = tab_cntrl(3) |
|---|
| 281 | day=float(day0) |
|---|
| 282 | |
|---|
| 283 | call get_var("Time",time,found) |
|---|
| 284 | |
|---|
| 285 | call close_startphy |
|---|
| 286 | |
|---|
| 287 | |
|---|
| 288 | ! Discretization (Definition of grid and time steps) |
|---|
| 289 | ! -------------- |
|---|
| 290 | ! |
|---|
| 291 | nlevel=nlayer+1 |
|---|
| 292 | nsoil=nsoilmx |
|---|
| 293 | |
|---|
| 294 | day_step=48 ! default value for day_step |
|---|
| 295 | PRINT *,'Number of time steps per sol ?' |
|---|
| 296 | call getin("day_step",day_step) |
|---|
| 297 | write(*,*) " day_step = ",day_step |
|---|
| 298 | |
|---|
| 299 | ecritphy=day_step ! default value for ecritphy, output every time step |
|---|
| 300 | |
|---|
| 301 | ndt=10 ! default value for ndt |
|---|
| 302 | PRINT *,'Number of sols to run ?' |
|---|
| 303 | call getin("ndt",ndt) |
|---|
| 304 | write(*,*) " ndt = ",ndt |
|---|
| 305 | |
|---|
| 306 | dayn=day0+ndt |
|---|
| 307 | ndt=ndt*day_step |
|---|
| 308 | dtphys=daysec/day_step |
|---|
| 309 | |
|---|
| 310 | ! Imposed surface pressure |
|---|
| 311 | ! ------------------------------------ |
|---|
| 312 | ! |
|---|
| 313 | |
|---|
| 314 | psurf(1)=610. ! default value for psurf |
|---|
| 315 | PRINT *,'Surface pressure (Pa) ?' |
|---|
| 316 | call getin("psurf",psurf(1)) |
|---|
| 317 | write(*,*) " psurf = ",psurf(1) |
|---|
| 318 | psurf(2)=psurf(1) |
|---|
| 319 | |
|---|
| 320 | ! Reference pressures |
|---|
| 321 | pa=20. ! transition pressure (for hybrid coord.) |
|---|
| 322 | preff=610. ! reference surface pressure |
|---|
| 323 | |
|---|
| 324 | ! Aerosol properties |
|---|
| 325 | ! -------------------------------- |
|---|
| 326 | tauvis=0.2 ! default value for tauvis (dust opacity) |
|---|
| 327 | write(*,'("Reference dust opacity at ",f4.0," Pa ?")')odpref |
|---|
| 328 | call getin("tauvis",tauvis) |
|---|
| 329 | write(*,*) " tauvis = ",tauvis |
|---|
| 330 | |
|---|
| 331 | ! Orbital parameters |
|---|
| 332 | ! ------------------ |
|---|
| 333 | |
|---|
| 334 | ! latitude/longitude |
|---|
| 335 | ! ------------------ |
|---|
| 336 | latitude(1)=0 ! default value for latitude |
|---|
| 337 | PRINT *,'latitude (in degrees) ?' |
|---|
| 338 | call getin("latitude",latitude(1)) |
|---|
| 339 | write(*,*) " latitude = ",latitude |
|---|
| 340 | latitude=latitude*pi/180.E+0 |
|---|
| 341 | longitude=0.E+0 |
|---|
| 342 | longitude=longitude*pi/180.E+0 |
|---|
| 343 | |
|---|
| 344 | ! some initializations (some of which have already been |
|---|
| 345 | ! done above!) and loads parameters set in callphys.def |
|---|
| 346 | ! and allocates some arrays |
|---|
| 347 | !Mars possible matter with dtphys in input and include!!! |
|---|
| 348 | ! Initializations below should mimick what is done in iniphysiq for 3D GCM |
|---|
| 349 | call init_interface_dyn_phys |
|---|
| 350 | call init_regular_lonlat(1,1,longitude,latitude,& |
|---|
| 351 | (/0.,0./),(/0.,0./)) |
|---|
| 352 | call init_geometry(1,longitude,latitude,& |
|---|
| 353 | (/0.,0.,0.,0./),(/0.,0.,0.,0./),& |
|---|
| 354 | cell_area) |
|---|
| 355 | ! Ehouarn: init_vertial_layers called later (because disvert not called yet) |
|---|
| 356 | ! call init_vertical_layers(nlayer,preff,scaleheight, |
|---|
| 357 | ! & ap,bp,aps,bps,presnivs,pseudoalt) |
|---|
| 358 | call init_dimphy(1,nlayer) ! Initialize dimphy module |
|---|
| 359 | call phys_state_var_init(1,llm,nq,tname,& |
|---|
| 360 | day0,dayn,time,& |
|---|
| 361 | daysec,dtphys,& |
|---|
| 362 | rad,g,r,cpp,& |
|---|
| 363 | nqperes,nqfils)! MVals: variables isotopes |
|---|
| 364 | call ini_fillgeom(1,latitude,longitude,(/1.0/)) |
|---|
| 365 | call conf_phys(1,llm,nq) |
|---|
| 366 | |
|---|
| 367 | ! in 1D model physics are called every time step |
|---|
| 368 | ! ovverride iphysiq value that has been set by conf_phys |
|---|
| 369 | if (iphysiq/=1) then |
|---|
| 370 | write(*,*) "testphys1d: setting iphysiq=1" |
|---|
| 371 | iphysiq=1 |
|---|
| 372 | endif |
|---|
| 373 | |
|---|
| 374 | ! Initialize albedo / soil thermal inertia |
|---|
| 375 | ! ---------------------------------------- |
|---|
| 376 | ! |
|---|
| 377 | |
|---|
| 378 | |
|---|
| 379 | ! Initialize local slope parameters (only matters if "callslope" |
|---|
| 380 | ! is .true. in callphys.def) |
|---|
| 381 | ! slope inclination angle (deg) 0: horizontal, 90: vertical |
|---|
| 382 | theta_sl(1)=0.0 ! default: no inclination |
|---|
| 383 | call getin("slope_inclination",theta_sl(1)) |
|---|
| 384 | ! slope orientation (deg) |
|---|
| 385 | ! 0 == Northward, 90 == Eastward, 180 == Southward, 270 == Westward |
|---|
| 386 | psi_sl(1)=0.0 ! default value |
|---|
| 387 | call getin("slope_orientation",psi_sl(1)) |
|---|
| 388 | |
|---|
| 389 | ! sub-slopes parameters (assuming no sub-slopes distribution for now). |
|---|
| 390 | def_slope(1)=-90 ! minimum slope angle |
|---|
| 391 | def_slope(2)=90 ! maximum slope angle |
|---|
| 392 | subslope_dist(1,1)=1 ! fraction of subslopes in mesh |
|---|
| 393 | ! |
|---|
| 394 | ! for the gravity wave scheme |
|---|
| 395 | ! --------------------------------- |
|---|
| 396 | ! |
|---|
| 397 | zmea(1)=0.E+0 |
|---|
| 398 | zstd(1)=0.E+0 |
|---|
| 399 | zsig(1)=0.E+0 |
|---|
| 400 | zgam(1)=0.E+0 |
|---|
| 401 | zthe(1)=0.E+0 |
|---|
| 402 | ! |
|---|
| 403 | ! for the slope wind scheme |
|---|
| 404 | ! --------------------------------- |
|---|
| 405 | ! |
|---|
| 406 | hmons(1)=0.E+0 |
|---|
| 407 | PRINT *,'hmons is initialized to ',hmons(1) |
|---|
| 408 | summit(1)=0.E+0 |
|---|
| 409 | PRINT *,'summit is initialized to ',summit(1) |
|---|
| 410 | base(1)=0.E+0 |
|---|
| 411 | ! |
|---|
| 412 | ! Default values initializing the coefficients calculated later |
|---|
| 413 | ! --------------------------------- |
|---|
| 414 | ! |
|---|
| 415 | tauscaling(1)=1. ! calculated in aeropacity_mod.F |
|---|
| 416 | totcloudfrac(1)=1. ! calculated in watercloud_mod.F |
|---|
| 417 | |
|---|
| 418 | ! Specific initializations for "physiq" |
|---|
| 419 | ! ------------------------------------- |
|---|
| 420 | ! surface geopotential is not used (or useful) since in 1D |
|---|
| 421 | ! everything is controled by surface pressure |
|---|
| 422 | phisfi(1)=0.E+0 |
|---|
| 423 | |
|---|
| 424 | ! Initialization to take into account prescribed winds |
|---|
| 425 | ! ------------------------------------------------------ |
|---|
| 426 | ptif=2.E+0*omeg*sinlat(1) |
|---|
| 427 | |
|---|
| 428 | ! geostrophic wind |
|---|
| 429 | gru=10. ! default value for gru |
|---|
| 430 | PRINT *,'zonal eastward component of the geostrophic wind (m/s) ?' |
|---|
| 431 | call getin("u",gru) |
|---|
| 432 | write(*,*) " u = ",gru |
|---|
| 433 | grv=0. !default value for grv |
|---|
| 434 | PRINT *,'meridional northward component of the geostrophic',& |
|---|
| 435 | ' wind (m/s) ?' |
|---|
| 436 | call getin("v",grv) |
|---|
| 437 | write(*,*) " v = ",grv |
|---|
| 438 | |
|---|
| 439 | ! Initialize winds for first time step |
|---|
| 440 | DO ilayer=1,nlayer |
|---|
| 441 | u(ilayer)=gru |
|---|
| 442 | v(ilayer)=grv |
|---|
| 443 | w(ilayer)=0 ! default: no vertical wind |
|---|
| 444 | ENDDO |
|---|
| 445 | |
|---|
| 446 | ! Initialize turbulente kinetic energy |
|---|
| 447 | DO ilevel=1,nlevel |
|---|
| 448 | q2(ilevel)=0.E+0 |
|---|
| 449 | ENDDO |
|---|
| 450 | |
|---|
| 451 | ! CO2 ice on the surface |
|---|
| 452 | ! ------------------- |
|---|
| 453 | ! get the index of co2 tracer (not known at this stage) |
|---|
| 454 | igcm_co2=0 |
|---|
| 455 | do iq=1,nq |
|---|
| 456 | if (trim(tname(iq))=="co2") then |
|---|
| 457 | igcm_co2=iq |
|---|
| 458 | endif |
|---|
| 459 | enddo |
|---|
| 460 | if (igcm_co2==0) then |
|---|
| 461 | write(*,*) "testphys1d error, missing co2 tracer!" |
|---|
| 462 | stop |
|---|
| 463 | endif |
|---|
| 464 | |
|---|
| 465 | ! Compute pressures and altitudes of atmospheric levels |
|---|
| 466 | ! ---------------------------------------------------------------- |
|---|
| 467 | |
|---|
| 468 | ! Vertical Coordinates |
|---|
| 469 | ! """""""""""""""""""" |
|---|
| 470 | hybrid=.true. |
|---|
| 471 | PRINT *,'Hybrid coordinates ?' |
|---|
| 472 | call getin("hybrid",hybrid) |
|---|
| 473 | write(*,*) " hybrid = ", hybrid |
|---|
| 474 | |
|---|
| 475 | CALL disvert_noterre |
|---|
| 476 | ! now that disvert has been called, initialize module vertical_layers_mod |
|---|
| 477 | call init_vertical_layers(nlayer,preff,scaleheight,& |
|---|
| 478 | ap,bp,aps,bps,presnivs,pseudoalt) |
|---|
| 479 | |
|---|
| 480 | DO ilevel=1,nlevel |
|---|
| 481 | plev(ilevel)=ap(ilevel)+psurf(1)*bp(ilevel) |
|---|
| 482 | ENDDO |
|---|
| 483 | |
|---|
| 484 | DO ilayer=1,nlayer |
|---|
| 485 | play(ilayer)=aps(ilayer)+psurf(1)*bps(ilayer) |
|---|
| 486 | ENDDO |
|---|
| 487 | |
|---|
| 488 | DO ilayer=1,nlayer |
|---|
| 489 | zlay(ilayer)=-200.E+0 *r*log(play(ilayer)/plev(1))& |
|---|
| 490 | /g |
|---|
| 491 | ENDDO |
|---|
| 492 | |
|---|
| 493 | |
|---|
| 494 | ! Initialize temperature profile |
|---|
| 495 | ! -------------------------------------- |
|---|
| 496 | pks=psurf(1)**rcp |
|---|
| 497 | |
|---|
| 498 | ! altitude in km in profile: divide zlay by 1000 |
|---|
| 499 | tmp1(0)=0.E+0 |
|---|
| 500 | DO ilayer=1,nlayer |
|---|
| 501 | tmp1(ilayer)=zlay(ilayer)/1000.E+0 |
|---|
| 502 | ENDDO |
|---|
| 503 | |
|---|
| 504 | call profile(nlayer+1,tmp1,tmp2) |
|---|
| 505 | |
|---|
| 506 | tsurf=tmp2(0) |
|---|
| 507 | DO ilayer=1,nlayer |
|---|
| 508 | temp(ilayer)=tmp2(ilayer) |
|---|
| 509 | ENDDO |
|---|
| 510 | |
|---|
| 511 | |
|---|
| 512 | ! Initialize soil properties and temperature |
|---|
| 513 | ! ------------------------------------------ |
|---|
| 514 | volcapa=1.e6 ! volumetric heat capacity |
|---|
| 515 | |
|---|
| 516 | flux_geo_tmp=0. |
|---|
| 517 | call getin("flux_geo",flux_geo_tmp) |
|---|
| 518 | flux_geo(:,:) = flux_geo_tmp |
|---|
| 519 | |
|---|
| 520 | ! Initialize depths |
|---|
| 521 | ! ----------------- |
|---|
| 522 | do isoil=0,nsoil-1 |
|---|
| 523 | mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) ! mid-layer depth |
|---|
| 524 | enddo |
|---|
| 525 | do isoil=1,nsoil |
|---|
| 526 | layer(isoil)=2.e-4*(2.**(isoil-1)) ! layer depth |
|---|
| 527 | enddo |
|---|
| 528 | |
|---|
| 529 | ! Initialize traceurs |
|---|
| 530 | ! --------------------------- |
|---|
| 531 | |
|---|
| 532 | if (photochem.or.callthermos) then |
|---|
| 533 | write(*,*) 'Initializing chemical species' |
|---|
| 534 | ! flagthermo=0: initialize over all atmospheric layers |
|---|
| 535 | flagthermo=0 |
|---|
| 536 | ! check if "h2o_vap" has already been initialized |
|---|
| 537 | ! (it has been if there is a "profile_h2o_vap" file around) |
|---|
| 538 | inquire(file="profile_h2o_vap",exist=present) |
|---|
| 539 | if (present) then |
|---|
| 540 | flagh2o=0 ! 0: do not initialize h2o_vap |
|---|
| 541 | else |
|---|
| 542 | flagh2o=1 ! 1: initialize h2o_vap in inichim_newstart |
|---|
| 543 | endif |
|---|
| 544 | |
|---|
| 545 | ! hack to accomodate that inichim_newstart assumes that |
|---|
| 546 | ! q & psurf arrays are on the dynamics scalar grid |
|---|
| 547 | allocate(qdyn(2,1,llm,nq),psdyn(2,1)) |
|---|
| 548 | qdyn(1,1,1:llm,1:nq)=q(1,1:llm,1:nq) |
|---|
| 549 | psdyn(1:2,1)=psurf(1) |
|---|
| 550 | call inichim_newstart(ngrid,nq,qdyn,qsurf,psdyn,& |
|---|
| 551 | flagh2o,flagthermo) |
|---|
| 552 | q(1,1:llm,1:nq)=qdyn(1,1,1:llm,1:nq) |
|---|
| 553 | endif |
|---|
| 554 | |
|---|
| 555 | ! Check if the surface is a water ice reservoir |
|---|
| 556 | ! -------------------------------------------------- |
|---|
| 557 | watercaptag(1)=.false. ! Default: no water ice reservoir |
|---|
| 558 | print *,'Water ice cap on ground ?' |
|---|
| 559 | call getin("watercaptag",watercaptag) |
|---|
| 560 | write(*,*) " watercaptag = ",watercaptag |
|---|
| 561 | |
|---|
| 562 | ! Check if the atmospheric water profile is specified |
|---|
| 563 | ! --------------------------------------------------- |
|---|
| 564 | ! Adding an option to force atmospheric water values JN |
|---|
| 565 | atm_wat_profile=-1 ! Default: free atm wat profile |
|---|
| 566 | print *,'Force atmospheric water vapor profile ?' |
|---|
| 567 | call getin("atm_wat_profile",atm_wat_profile) |
|---|
| 568 | write(*,*) "atm_wat_profile = ", atm_wat_profile |
|---|
| 569 | if (atm_wat_profile.EQ.-1) then |
|---|
| 570 | write(*,*) "Free atmospheric water vapor profile" |
|---|
| 571 | write(*,*) "Total water is conserved on the column" |
|---|
| 572 | else if (atm_wat_profile.EQ.0) then |
|---|
| 573 | write(*,*) "Dry atmospheric water vapor profile" |
|---|
| 574 | else if ((atm_wat_profile.GT.0).and.(atm_wat_profile.LE.1)) then |
|---|
| 575 | write(*,*) "Prescribed atmospheric water vapor MMR profile" |
|---|
| 576 | write(*,*) "Unless it reaches saturation (maximal value)" |
|---|
| 577 | write(*,*) "MMR chosen : ", atm_wat_profile |
|---|
| 578 | endif |
|---|
| 579 | |
|---|
| 580 | ap_out=ap |
|---|
| 581 | bp_out=bp |
|---|
| 582 | |
|---|
| 583 | end subroutine init_phys_1d |
|---|
| 584 | |
|---|
| 585 | #endif |
|---|
| 586 | |
|---|
| 587 | end module init_phys_1d_mod |
|---|