PROGRAM testphys1d use comsoil_h, only: inertiedat, inertiesoil, nsoilmx, tsoil, nqsoil, qsoil use surfdat_h, only: albedodat, perennial_co2ice, watercap, tsurf, emis, qsurf use comslope_mod, only: def_slope, subslope_dist use phyredem, only: physdem0, physdem1 use watersat_mod, only: watersat use tracer_mod, only: igcm_h2o_vap, igcm_h2o_ice, noms use comcstfi_h, only: pi, rad, omeg, g, mugaz, rcp, r, cpp use time_phylmdz_mod, only: daysec, day_step use dimradmars_mod, only: tauvis, totcloudfrac, albedo use dust_param_mod, only: tauscaling use comvert_mod, only: ap, bp, aps, bps, pa, preff, sig use physiq_mod, only: physiq use turb_mod, only: q2 use write_output_mod, only: write_output use ioipsl_getincom, only: getin ! To use 'getin' use init_testphys1d_mod, only: init_testphys1d use writerestart1D_mod, only: writerestart1D ! Mostly for XIOS outputs: use mod_const_mpi, only: init_const_mpi use parallel_lmdz, only: init_parallel implicit none !======================================================================= ! subject: ! -------- ! PROGRAM useful to run physical part of the martian GCM in a 1D column ! ! Can be compiled with a command like (e.g. for 25 layers) ! "makegcm -p mars -d 25 testphys1d" ! It requires the files "testphys1d.def" "callphys.def" ! and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line) ! and a file describing the sigma layers (e.g. "z2sig.def") ! ! author: Frederic Hourdin, R. Fournier, F. Forget ! ------- ! ! update: 12/06/2003, including chemistry (S. Lebonnois) ! and water ice (F. Montmessin) ! 27/09/2023, upgrade to F90 (JB Clément) ! !======================================================================= include "dimensions.h" !#include "dimradmars.h" !#include "comgeomfi.h" !#include "surfdat.h" !#include "slope.h" !#include "comsoil.h" !#include "comdiurn.h" include "callkeys.h" !#include "comsaison.h" !#include "control.h" include "netcdf.inc" !#include "advtrac.h" !-------------------------------------------------------------- ! Declarations !-------------------------------------------------------------- integer, parameter :: ngrid = 1 ! (2+(jjm-1)*iim - 1/jjm) integer, parameter :: nlayer = llm real, parameter :: odpref = 610. ! DOD reference pressure (Pa) integer :: unitstart ! unite d'ecriture de "startfi" integer :: ndt, ilayer, ilevel, isoil, idt, iq logical :: firstcall, lastcall integer :: day0 ! initial (sol ; =0 at Ls=0) real :: day ! date during the run real :: time ! time (00, dry if =0 q(1,:,igcm_h2o_vap) = min(zqsat(:),atm_wat_profile*g/psurf) q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice else ! Relaxation towards the value atm_wat_profile with relaxation time atm_wat_tau 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) q(1,:,igcm_h2o_vap) = min(zqsat(:),q(1,:,igcm_h2o_vap)) q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice endif endif !EV profile ! IF(atm_wat_profile.eq.2) THEN ! DO ilayer=1,16 ! q(1,ilayer,igcm_h2o_vap)=min(zqsat(ilayer),6e-5) ! ENDDO! ilayer=1,16 ! DO ilayer=17,22 ! q(1,ilayer,igcm_h2o_vap)=min(zqsat(ilayer),3e-5) ! ENDDO! ilayer=17,22 ! DO ilayer=23,nlayer ! q(1,ilayer,igcm_h2o_vap)=0 ! ENDDO! ilayer=23,nlayer ! endif endif ! Call physics ! -------------------- call physiq(1,llm,nq,firstcall,lastcall,day,time,dttestphys,plev,play,phi,u,v,temp,q,w,du,dv,dtemp,dq,dpsurf) ! ^----- outputs -----^ ! Wind increment: specific for 1D ! -------------------------------- ! The physics compute the tendencies on u and v, ! here we just add Coriolos effect !do ilayer = 1,nlayer ! du(ilayer) = du(ilayer) + ptif*(v(ilayer) - grv) ! dv(ilayer) = dv(ilayer) + ptif*(-u(ilayer) + gru) !enddo ! For some tests: No coriolis force at equator !if (latitude(1) == 0.) then du(:) = du(:) + (gru - u(:))/1.e4 dv(:) = dv(:) + (grv - v(:))/1.e4 !endif ! Compute time for next time step ! ------------------------------- firstcall = .false. time = time + dttestphys/daysec if (time > 1.) then time = time - 1. day = day + 1 endif ! Compute winds and temperature for next time step ! ------------------------------------------------ u(:) = u(:) + dttestphys*du(:) v(:) = v(:) + dttestphys*dv(:) temp(:) = temp(:) + dttestphys*dtemp(:) ! Compute pressure for next time step ! ----------------------------------- psurf = psurf + dttestphys*dpsurf(1) ! surface pressure change plev(:) = ap(:) + psurf*bp(:) play(:) = aps(:) + psurf*bps(:) ! Increment tracers q(1,:,:) = q(1,:,:) + dttestphys*dq(1,:,:) enddo ! End of time stepping loop (idt=1,ndt) ! Writing the "restart1D.txt" file for the next run if (startfiles_1D) call writerestart1D('restart1D.txt',psurf,tsurf(1,:),nlayer,size(tsurf,2),temp,u,v,nq,noms,qsurf(1,:,:),q) write(*,*) "testphys1d: everything is cool." END PROGRAM testphys1d !*********************************************************************** !*********************************************************************** ! Dummy subroutine used only in 3D, but required to ! compile testphys1d (to cleanly use writediagfi) SUBROUTINE gr_fi_dyn(nfield,ngrid,im,jm,pfi,pdyn) implicit none integer :: im, jm, ngrid, nfield real, dimension(im,jm,nfield) :: pdyn real, dimension(ngrid,nfield) :: pfi if (ngrid /= 1) error stop 'gr_fi_dyn error: in 1D ngrid should be 1!!!' pdyn(1,1,1:nfield) = pfi(1,1:nfield) END SUBROUTINE gr_fi_dyn !*********************************************************************** !***********************************************************************