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