[3053] | 1 | PROGRAM testphys1d |
---|
[38] | 2 | |
---|
[3113] | 3 | use comsoil_h, only: inertiedat, inertiesoil, nsoilmx, tsoil, nqsoil, qsoil |
---|
[3130] | 4 | use surfdat_h, only: albedodat, perennial_co2ice, watercap, tsurf, emis, qsurf |
---|
[3067] | 5 | use comslope_mod, only: def_slope, subslope_dist |
---|
| 6 | use phyredem, only: physdem0, physdem1 |
---|
| 7 | use watersat_mod, only: watersat |
---|
| 8 | use tracer_mod, only: igcm_h2o_vap, igcm_h2o_ice, noms |
---|
| 9 | use comcstfi_h, only: pi, rad, omeg, g, mugaz, rcp, r, cpp |
---|
| 10 | use time_phylmdz_mod, only: daysec, day_step |
---|
| 11 | use dimradmars_mod, only: tauvis, totcloudfrac, albedo |
---|
| 12 | use dust_param_mod, only: tauscaling |
---|
| 13 | use comvert_mod, only: ap, bp, aps, bps, pa, preff, sig |
---|
| 14 | use physiq_mod, only: physiq |
---|
| 15 | use turb_mod, only: q2 |
---|
| 16 | use write_output_mod, only: write_output |
---|
[3129] | 17 | use ioipsl_getincom, only: getin ! To use 'getin' |
---|
[3067] | 18 | use init_testphys1d_mod, only: init_testphys1d |
---|
[3074] | 19 | use writerestart1D_mod, only: writerestart1D |
---|
[3054] | 20 | ! Mostly for XIOS outputs: |
---|
[3067] | 21 | use mod_const_mpi, only: init_const_mpi |
---|
| 22 | use parallel_lmdz, only: init_parallel |
---|
[3003] | 23 | |
---|
[3053] | 24 | implicit none |
---|
[38] | 25 | |
---|
[3053] | 26 | !======================================================================= |
---|
| 27 | ! subject: |
---|
| 28 | ! -------- |
---|
| 29 | ! PROGRAM useful to run physical part of the martian GCM in a 1D column |
---|
| 30 | ! |
---|
| 31 | ! Can be compiled with a command like (e.g. for 25 layers) |
---|
| 32 | ! "makegcm -p mars -d 25 testphys1d" |
---|
| 33 | ! It requires the files "testphys1d.def" "callphys.def" |
---|
| 34 | ! and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line) |
---|
| 35 | ! and a file describing the sigma layers (e.g. "z2sig.def") |
---|
| 36 | ! |
---|
[3092] | 37 | ! author: Frederic Hourdin, R. Fournier, F. Forget |
---|
[3053] | 38 | ! ------- |
---|
| 39 | ! |
---|
| 40 | ! update: 12/06/2003, including chemistry (S. Lebonnois) |
---|
[3056] | 41 | ! and water ice (F. Montmessin) |
---|
| 42 | ! 27/09/2023, upgrade to F90 (JB Clément) |
---|
[3053] | 43 | ! |
---|
| 44 | !======================================================================= |
---|
[38] | 45 | |
---|
[3053] | 46 | include "dimensions.h" |
---|
[1047] | 47 | !#include "dimradmars.h" |
---|
| 48 | !#include "comgeomfi.h" |
---|
| 49 | !#include "surfdat.h" |
---|
| 50 | !#include "slope.h" |
---|
| 51 | !#include "comsoil.h" |
---|
| 52 | !#include "comdiurn.h" |
---|
[3053] | 53 | include "callkeys.h" |
---|
[1047] | 54 | !#include "comsaison.h" |
---|
[1130] | 55 | !#include "control.h" |
---|
[3053] | 56 | include "netcdf.inc" |
---|
[1036] | 57 | !#include "advtrac.h" |
---|
[38] | 58 | |
---|
[3053] | 59 | !-------------------------------------------------------------- |
---|
| 60 | ! Declarations |
---|
| 61 | !-------------------------------------------------------------- |
---|
[3067] | 62 | integer, parameter :: ngrid = 1 ! (2+(jjm-1)*iim - 1/jjm) |
---|
| 63 | integer, parameter :: nlayer = llm |
---|
| 64 | real, parameter :: odpref = 610. ! DOD reference pressure (Pa) |
---|
| 65 | integer :: unitstart ! unite d'ecriture de "startfi" |
---|
| 66 | integer :: ndt, ilayer, ilevel, isoil, idt, iq |
---|
| 67 | logical :: firstcall, lastcall |
---|
| 68 | integer :: day0 ! initial (sol ; =0 at Ls=0) |
---|
| 69 | real :: day ! date during the run |
---|
| 70 | real :: time ! time (0<time<1 ; time=0.5 a midi) |
---|
| 71 | real :: dttestphys ! testphys1d timestep |
---|
| 72 | real, dimension(nlayer) :: play ! Pressure at the middle of the layers (Pa) |
---|
| 73 | real, dimension(nlayer + 1) :: plev ! intermediate pressure levels (pa) |
---|
| 74 | real :: psurf ! Surface pressure |
---|
| 75 | real, dimension(nlayer) :: u, v ! zonal, meridional wind |
---|
| 76 | real :: gru, grv ! prescribed "geostrophic" background wind |
---|
| 77 | real, dimension(nlayer) :: temp ! temperature at the middle of the layers |
---|
| 78 | real, dimension(:,:,:), allocatable :: q ! tracer mixing ratio (e.g. kg/kg) |
---|
| 79 | real, dimension(1) :: wstar = 0. ! Thermals vertical velocity |
---|
[38] | 80 | |
---|
[3056] | 81 | ! Physical and dynamical tandencies (e.g. m.s-2, K/s, Pa/s) |
---|
[3067] | 82 | real, dimension(nlayer) :: du, dv, dtemp, dudyn, dvdyn, dtempdyn |
---|
| 83 | real, dimension(1) :: dpsurf |
---|
| 84 | real, dimension(:,:,:), allocatable :: dq, dqdyn |
---|
[38] | 85 | |
---|
[3053] | 86 | ! Various intermediate variables |
---|
[3056] | 87 | integer :: aslun |
---|
| 88 | real :: zls, pks, ptif, qtotinit, qtot |
---|
| 89 | real, dimension(nlayer) :: phi, h, s, w |
---|
| 90 | integer :: nq = 1 ! number of tracers |
---|
| 91 | real, dimension(1) :: latitude, longitude, cell_area |
---|
[3092] | 92 | character(2) :: str2 |
---|
| 93 | character(7) :: str7 |
---|
| 94 | character(44) :: txt |
---|
[38] | 95 | |
---|
[3065] | 96 | ! RV & JBC: Use of starting files for 1D |
---|
[3129] | 97 | logical :: startfiles_1D, therestart1D, therestartfi, there |
---|
[2789] | 98 | |
---|
[3053] | 99 | ! JN & JBC: Force atmospheric water profiles |
---|
| 100 | real :: atm_wat_profile, atm_wat_tau |
---|
| 101 | real, dimension(:), allocatable :: zqsat ! useful for (atm_wat_profile=2) |
---|
| 102 | !======================================================================= |
---|
| 103 | |
---|
| 104 | !======================================================================= |
---|
| 105 | ! INITIALISATION |
---|
| 106 | !======================================================================= |
---|
[2997] | 107 | #ifdef CPP_XIOS |
---|
[3053] | 108 | call init_const_mpi |
---|
| 109 | call init_parallel |
---|
[2997] | 110 | #endif |
---|
| 111 | |
---|
[3053] | 112 | ! Initialize "serial/parallel" related stuff |
---|
| 113 | !call init_phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/)) |
---|
| 114 | !call init_phys_lmdz(1,1,llm,1,(/1/)) |
---|
| 115 | !call initcomgeomphy |
---|
[38] | 116 | |
---|
[3129] | 117 | startfiles_1D = .false. |
---|
| 118 | !------------------------------------------------------ |
---|
| 119 | ! Loading run parameters from "run.def" file |
---|
| 120 | !------------------------------------------------------ |
---|
| 121 | ! check if 'run.def' file is around. Otherwise reading parameters |
---|
| 122 | ! from callphys.def via getin() routine won't work. |
---|
| 123 | inquire(file = 'run.def',exist = there) |
---|
| 124 | if (.not. there) then |
---|
| 125 | write(*,*) 'Cannot find required file "run.def"' |
---|
| 126 | write(*,*) ' (which should contain some input parameters along with the following line: INCLUDEDEF=callphys.def)' |
---|
| 127 | write(*,*) ' ... might as well stop here ...' |
---|
| 128 | error stop |
---|
| 129 | endif |
---|
| 130 | |
---|
| 131 | startfiles_1D = .false. ! Default value |
---|
| 132 | write(*,*) 'Do you want to use starting files and/or to write restarting files?' |
---|
| 133 | call getin('startfiles_1D',startfiles_1D) |
---|
| 134 | write(*,*) 'startfiles_1D =', startfiles_1D |
---|
| 135 | |
---|
| 136 | therestart1D = .false. ! Default value |
---|
| 137 | inquire(file = 'start1D.txt',exist = therestart1D) |
---|
| 138 | if (startfiles_1D .and. .not. therestart1D) then |
---|
| 139 | write(*,*) 'There is no "start1D.txt" file!' |
---|
| 140 | write(*,*) 'Initialization is done with default values.' |
---|
| 141 | endif |
---|
| 142 | therestartfi = .false. ! Default value |
---|
| 143 | inquire(file = 'startfi.nc',exist = therestartfi) |
---|
| 144 | if (.not. therestartfi) then |
---|
| 145 | write(*,*) 'There is no "startfi.nc" file!' |
---|
| 146 | write(*,*) 'Initialization is done with default values.' |
---|
| 147 | endif |
---|
| 148 | |
---|
| 149 | call init_testphys1d('start1D.txt','startfi.nc',startfiles_1D,therestart1D,therestartfi,ngrid,nlayer,odpref, & |
---|
| 150 | nq,q,time,psurf,u,v,temp,ndt,ptif,pks,dttestphys,zqsat,dq,dqdyn,day0,day,gru,grv,w, & |
---|
[3067] | 151 | play,plev,latitude,longitude,cell_area,atm_wat_profile,atm_wat_tau) |
---|
[3060] | 152 | |
---|
[3053] | 153 | ! Write a "startfi" file |
---|
| 154 | ! ---------------------- |
---|
| 155 | ! This file will be read during the first call to "physiq". |
---|
| 156 | ! It is needed to transfert physics variables to "physiq"... |
---|
[3060] | 157 | if (.not. therestartfi) then |
---|
[3053] | 158 | call physdem0("startfi.nc",longitude,latitude,nsoilmx,ngrid, & |
---|
| 159 | llm,nq,dttestphys,float(day0),0.,cell_area, & |
---|
| 160 | albedodat,inertiedat,def_slope,subslope_dist) |
---|
[3117] | 161 | call physdem1("startfi.nc",nsoilmx,ngrid,llm,nq,nqsoil,dttestphys,time, & |
---|
[3113] | 162 | tsurf,tsoil,inertiesoil,albedo,emis,q2,qsurf,qsoil,tauscaling, & |
---|
[3130] | 163 | totcloudfrac,wstar,watercap,perennial_co2ice) |
---|
[3060] | 164 | endif !(.not. therestartfi) |
---|
[38] | 165 | |
---|
[3053] | 166 | !======================================================================= |
---|
| 167 | ! 1D MODEL TIME STEPPING LOOP |
---|
| 168 | !======================================================================= |
---|
| 169 | firstcall = .true. |
---|
| 170 | lastcall = .false. |
---|
| 171 | do idt = 1,ndt |
---|
| 172 | if (idt == ndt) lastcall = .true. |
---|
| 173 | ! if (idt == ndt - day_step - 1) then ! test |
---|
| 174 | ! lastcall = .true. |
---|
| 175 | ! call solarlong(day*1.0,zls) |
---|
| 176 | ! write(103,*) 'Ls=',zls*180./pi |
---|
| 177 | ! write(103,*) 'Lat=', latitude(1)*180./pi |
---|
| 178 | ! write(103,*) 'Tau=', tauvis/odpref*psurf |
---|
| 179 | ! write(103,*) 'RunEnd - Atmos. Temp. File' |
---|
| 180 | ! write(103,*) 'RunEnd - Atmos. Temp. File' |
---|
| 181 | ! write(104,*) 'Ls=',zls*180./pi |
---|
| 182 | ! write(104,*) 'Lat=', latitude(1) |
---|
| 183 | ! write(104,*) 'Tau=', tauvis/odpref*psurf |
---|
| 184 | ! write(104,*) 'RunEnd - Atmos. Temp. File' |
---|
| 185 | ! endif |
---|
[899] | 186 | |
---|
[3053] | 187 | ! Compute geopotential |
---|
| 188 | ! ~~~~~~~~~~~~~~~~~~~~ |
---|
| 189 | s(:) = (aps(:)/psurf + bps(:))**rcp |
---|
| 190 | h(:) = cpp*temp(:)/(pks*s(:)) |
---|
[38] | 191 | |
---|
[3053] | 192 | phi(1) = pks*h(1)*(1. - s(1)) |
---|
| 193 | do ilayer = 2,nlayer |
---|
| 194 | phi(ilayer) = phi(ilayer - 1) + pks*(h(ilayer - 1) + h(ilayer))*.5*(s(ilayer - 1) - s(ilayer)) |
---|
| 195 | enddo |
---|
[38] | 196 | |
---|
[3053] | 197 | ! Modify atmospheric water if asked |
---|
| 198 | ! Added "atm_wat_profile" flag (JN + JBC) |
---|
| 199 | ! --------------------------------------- |
---|
| 200 | if (water) then |
---|
| 201 | call watersat(nlayer,temp,play,zqsat) |
---|
| 202 | if (0. <= atm_wat_profile .and. atm_wat_profile <= 1.) then |
---|
| 203 | ! If atmospheric water is monitored |
---|
| 204 | if (atm_wat_tau < 0.) then ! Prescribed atm_wat_profile: wet if >0, dry if =0 |
---|
[3067] | 205 | q(1,:,igcm_h2o_vap) = min(zqsat(:),atm_wat_profile*g/psurf) |
---|
| 206 | q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice |
---|
[3053] | 207 | else ! Relaxation towards the value atm_wat_profile with relaxation time atm_wat_tau |
---|
[3067] | 208 | q(1,:,igcm_h2o_vap) = atm_wat_profile*g/psurf + (q(1,:,igcm_h2o_vap) - atm_wat_profile*g/psurf)*dexp(-dttestphys/atm_wat_tau) |
---|
| 209 | q(1,:,igcm_h2o_vap) = min(zqsat(:),q(1,:,igcm_h2o_vap)) |
---|
| 210 | q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice |
---|
[3053] | 211 | endif |
---|
| 212 | endif |
---|
[3098] | 213 | !EV profile |
---|
| 214 | ! IF(atm_wat_profile.eq.2) THEN |
---|
| 215 | ! DO ilayer=1,16 |
---|
| 216 | ! q(1,ilayer,igcm_h2o_vap)=min(zqsat(ilayer),6e-5) |
---|
| 217 | ! ENDDO! ilayer=1,16 |
---|
| 218 | ! DO ilayer=17,22 |
---|
| 219 | ! q(1,ilayer,igcm_h2o_vap)=min(zqsat(ilayer),3e-5) |
---|
| 220 | ! ENDDO! ilayer=17,22 |
---|
| 221 | ! DO ilayer=23,nlayer |
---|
| 222 | ! q(1,ilayer,igcm_h2o_vap)=0 |
---|
| 223 | ! ENDDO! ilayer=23,nlayer |
---|
| 224 | ! endif |
---|
| 225 | endif |
---|
[38] | 226 | |
---|
[3053] | 227 | ! Call physics |
---|
| 228 | ! -------------------- |
---|
| 229 | call physiq(1,llm,nq,firstcall,lastcall,day,time,dttestphys,plev,play,phi,u,v,temp,q,w,du,dv,dtemp,dq,dpsurf) |
---|
| 230 | ! ^----- outputs -----^ |
---|
[2789] | 231 | |
---|
[3053] | 232 | ! Wind increment: specific for 1D |
---|
| 233 | ! -------------------------------- |
---|
| 234 | ! The physics compute the tendencies on u and v, |
---|
| 235 | ! here we just add Coriolos effect |
---|
| 236 | !do ilayer = 1,nlayer |
---|
| 237 | ! du(ilayer) = du(ilayer) + ptif*(v(ilayer) - grv) |
---|
| 238 | ! dv(ilayer) = dv(ilayer) + ptif*(-u(ilayer) + gru) |
---|
| 239 | !enddo |
---|
| 240 | ! For some tests: No coriolis force at equator |
---|
[3067] | 241 | !if (latitude(1) == 0.) then |
---|
[3053] | 242 | du(:) = du(:) + (gru - u(:))/1.e4 |
---|
| 243 | dv(:) = dv(:) + (grv - v(:))/1.e4 |
---|
| 244 | !endif |
---|
[38] | 245 | |
---|
[3053] | 246 | ! Compute time for next time step |
---|
| 247 | ! ------------------------------- |
---|
| 248 | firstcall = .false. |
---|
| 249 | time = time + dttestphys/daysec |
---|
| 250 | if (time > 1.) then |
---|
| 251 | time = time - 1. |
---|
| 252 | day = day + 1 |
---|
| 253 | endif |
---|
[38] | 254 | |
---|
[3053] | 255 | ! Compute winds and temperature for next time step |
---|
| 256 | ! ------------------------------------------------ |
---|
| 257 | u(:) = u(:) + dttestphys*du(:) |
---|
| 258 | v(:) = v(:) + dttestphys*dv(:) |
---|
| 259 | temp(:) = temp(:) + dttestphys*dtemp(:) |
---|
[3021] | 260 | |
---|
[3053] | 261 | ! Compute pressure for next time step |
---|
| 262 | ! ----------------------------------- |
---|
| 263 | psurf = psurf + dttestphys*dpsurf(1) ! surface pressure change |
---|
| 264 | plev(:) = ap(:) + psurf*bp(:) |
---|
| 265 | play(:) = aps(:) + psurf*bps(:) |
---|
[38] | 266 | |
---|
[3053] | 267 | ! Increment tracers |
---|
[3067] | 268 | q(1,:,:) = q(1,:,:) + dttestphys*dq(1,:,:) |
---|
[3053] | 269 | enddo ! End of time stepping loop (idt=1,ndt) |
---|
[38] | 270 | |
---|
[3053] | 271 | ! Writing the "restart1D.txt" file for the next run |
---|
[3069] | 272 | if (startfiles_1D) call writerestart1D('restart1D.txt',psurf,tsurf(1,:),nlayer,size(tsurf,2),temp,u,v,nq,noms,qsurf(1,:,:),q) |
---|
[38] | 273 | |
---|
[3053] | 274 | write(*,*) "testphys1d: everything is cool." |
---|
[38] | 275 | |
---|
[3053] | 276 | END PROGRAM testphys1d |
---|
[3054] | 277 | |
---|
[3053] | 278 | !*********************************************************************** |
---|
| 279 | !*********************************************************************** |
---|
| 280 | ! Dummy subroutine used only in 3D, but required to |
---|
| 281 | ! compile testphys1d (to cleanly use writediagfi) |
---|
| 282 | SUBROUTINE gr_fi_dyn(nfield,ngrid,im,jm,pfi,pdyn) |
---|
[38] | 283 | |
---|
[3053] | 284 | implicit none |
---|
[1380] | 285 | |
---|
[3053] | 286 | integer :: im, jm, ngrid, nfield |
---|
| 287 | real, dimension(im,jm,nfield) :: pdyn |
---|
| 288 | real, dimension(ngrid,nfield) :: pfi |
---|
[38] | 289 | |
---|
[3069] | 290 | if (ngrid /= 1) error stop 'gr_fi_dyn error: in 1D ngrid should be 1!!!' |
---|
[899] | 291 | |
---|
[3053] | 292 | pdyn(1,1,1:nfield) = pfi(1,1:nfield) |
---|
[899] | 293 | |
---|
[3053] | 294 | END SUBROUTINE gr_fi_dyn |
---|
[3054] | 295 | |
---|
[3053] | 296 | !*********************************************************************** |
---|
| 297 | !*********************************************************************** |
---|