PROGRAM testphys1d ! to use 'getin' USE ioipsl_getincom, only: getin use dimphy, only : init_dimphy use mod_grid_phy_lmdz, only : regular_lonlat use infotrac, only: nqtot, tname, nqperes,nqfils use comsoil_h, only: volcapa, layer, mlayer, inertiedat, & inertiesoil,nsoilmx,flux_geo use comgeomfi_h, only: sinlat, ini_fillgeom use surfdat_h, only: albedodat, z0_default, emissiv, emisice, & albedice, iceradius, dtemisice, z0, & zmea, zstd, zsig, zgam, zthe, phisfi, & watercaptag, watercap, hmons, summit, base use slope_mod, only: theta_sl, psi_sl use comslope_mod, only: def_slope,subslope_dist,def_slope_mean use phyredem, only: physdem0,physdem1 use geometry_mod, only: init_geometry use watersat_mod, only: watersat use tracer_mod, only: igcm_h2o_vap,igcm_h2o_ice,igcm_co2,noms use planete_h, only: year_day, periheli, aphelie, peri_day, & obliquit, emin_turb, lmixmin use comcstfi_h, only: pi, rad, omeg, g, mugaz, rcp, r, cpp use time_phylmdz_mod, only: daysec, dtphys, day_step, & ecritphy, iphysiq use dimradmars_mod, only: tauvis,totcloudfrac use dust_param_mod, only: tauscaling USE comvert_mod, ONLY: ap,bp,aps,bps,pa,preff,sig, & presnivs,pseudoalt,scaleheight USE vertical_layers_mod, ONLY: init_vertical_layers USE logic_mod, ONLY: hybrid use physics_distribution_mod, only: init_physics_distribution use regular_lonlat_mod, only: init_regular_lonlat use mod_interface_dyn_phys, only: init_interface_dyn_phys USE phys_state_var_init_mod, ONLY: phys_state_var_init USE physiq_mod, ONLY: physiq USE read_profile_mod, ONLY: read_profile use inichim_newstart_mod, only: inichim_newstart use phyetat0_mod, only: phyetat0 USE netcdf, only: NF90_OPEN, NF90_NOERR, NF90_NOWRITE, & nf90_strerror use iostart, only: open_startphy,get_var, close_startphy use write_output_mod, only: write_output IMPLICIT NONE c======================================================================= c subject: c -------- c PROGRAM useful to run physical part of the martian GCM in a 1D column c c Can be compiled with a command like (e.g. for 25 layers) c "makegcm -p mars -d 25 testphys1d" c It requires the files "testphys1d.def" "callphys.def" c and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line) c and a file describing the sigma layers (e.g. "z2sig.def") c c author: Frederic Hourdin, R.Fournier,F.Forget c ------- c c update: 12/06/2003 including chemistry (S. Lebonnois) c and water ice (F. Montmessin) c c======================================================================= include "dimensions.h" integer, parameter :: ngrid = 1 !(2+(jjm-1)*iim - 1/jjm) integer, parameter :: nlayer = llm !#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 "comg1d.h" !#include "advtrac.h" c -------------------------------------------------------------- c Declarations c -------------------------------------------------------------- c INTEGER unitstart ! unite d'ecriture de "startfi" INTEGER nlevel,nsoil,ndt INTEGER ilayer,ilevel,isoil,idt,iq LOGICAl firstcall,lastcall LOGICAL write_prof c real,parameter :: odpref=610. ! DOD reference pressure (Pa) c INTEGER day0,dayn ! date initial (sol ; =0 a Ls=0) and final REAL day ! date durant le run REAL time ! time (0=1!" stop endif endif ! allocate arrays: allocate(tname(nq)) allocate(q(nlayer,nq)) allocate(zqsat(nlayer)) allocate(qsurf(nq)) allocate(dq(nlayer,nq)) allocate(dqdyn(nlayer,nq)) allocate(mqtot(nq)) allocate(tnom_transp(nq)) ! read tracer names from file traceur.def do iq=1,nq read(90,'(80a)',iostat=ierr) line ! store the line from traceur.def if (ierr.ne.0) then write(*,*) 'testphys1d: error reading tracer names...' stop endif ! if format is tnom_0, tnom_transp (isotopes) read(line,*,iostat=ierr0) tname(iq),tnom_transp(iq) if (ierr0.ne.0) then read(line,*) tname(iq) tnom_transp(iq)='air' endif enddo close(90) ! Isotopes: as in the 3D case we have to determine father/son relations for isotopes and carrying fluid ALLOCATE(nqfils(nqtot)) nqperes=0 nqfils(:)=0 DO iq=1,nqtot if (tnom_transp(iq) == 'air') then ! ceci est un traceur père WRITE(*,*) 'Le traceur',iq,', appele ', & trim(tname(iq)),', est un pere' nqperes=nqperes+1 else !if (tnom_transp(iq) == 'air') then ! ceci est un fils. Qui est son père? WRITE(*,*) 'Le traceur',iq,', appele ', & trim(tname(iq)),', est un fils' continu=.true. ipere=1 do while (continu) if (tnom_transp(iq) .eq. tname(ipere)) then ! Son père est ipere WRITE(*,*) 'Le traceur',iq,'appele ', & trim(tname(iq)),' est le fils de ', & ipere,'appele ',trim(tname(ipere)) nqfils(ipere)=nqfils(ipere)+1 continu=.false. else !if (tnom_transp(iq) == tnom_0(ipere)) then ipere=ipere+1 if (ipere.gt.nqtot) then WRITE(*,*) 'Le traceur',iq,'appele ', & trim(tname(iq)),', est orpelin.' CALL abort_gcm('infotrac_init', & 'Un traceur est orphelin',1) endif !if (ipere.gt.nqtot) then endif !if (tnom_transp(iq) == tnom_0(ipere)) then enddo !do while (continu) endif !if (tnom_transp(iq) == 'air') then enddo !DO iq=1,nqtot WRITE(*,*) 'nqperes=',nqperes WRITE(*,*) 'nqfils=',nqfils ! initialize tracers here: write(*,*) "testphys1d: initializing tracers" call read_profile(nq, nlayer, qsurf, q) call init_physics_distribution(regular_lonlat,4, & 1,1,1,nlayer,1) if(.not. startfile_1D ) then c Date and local time at beginning of run c --------------------------------------- c Date (in sols since spring solstice) at beginning of run day0 = 0 ! default value for day0 write(*,*) 'Initial date (in martian sols ; =0 at Ls=0)?' call getin("day0",day0) day=float(day0) write(*,*) " day0 = ",day0 c Local time at beginning of run time=0 ! default value for time write(*,*)'Initial local time (in hours, between 0 and 24)?' call getin("time",time) write(*,*)" time = ",time time=time/24.E+0 ! convert time (hours) to fraction of sol else call open_startphy("startfi.nc") call get_var("controle",tab_cntrl,found) if (.not.found) then call abort_physic("open_startphy", & "tabfi: Failed reading array",1) else write(*,*)'tabfi: tab_cntrl',tab_cntrl endif day0 = tab_cntrl(3) day=float(day0) call get_var("Time",time,found) call close_startphy endif !startfile_1D c Discretization (Definition of grid and time steps) c -------------- c nlevel=nlayer+1 nsoil=nsoilmx day_step=48 ! default value for day_step PRINT *,'Number of time steps per sol ?' call getin("day_step",day_step) write(*,*) " day_step = ",day_step ecritphy=day_step ! default value for ecritphy, output every time step ndt=10 ! default value for ndt PRINT *,'Number of sols to run ?' call getin("ndt",ndt) write(*,*) " ndt = ",ndt dayn=day0+ndt ndt=ndt*day_step dtphys=daysec/day_step c Imposed surface pressure c ------------------------------------ c psurf=610. ! default value for psurf PRINT *,'Surface pressure (Pa) ?' call getin("psurf",psurf) write(*,*) " psurf = ",psurf c Reference pressures pa=20. ! transition pressure (for hybrid coord.) preff=610. ! reference surface pressure c Aerosol properties c -------------------------------- tauvis=0.2 ! default value for tauvis (dust opacity) write(*,'("Reference dust opacity at ",f4.0," Pa ?")')odpref call getin("tauvis",tauvis) write(*,*) " tauvis = ",tauvis c Orbital parameters c ------------------ if(.not. startfile_1D ) then paleomars=.false. ! Default: no water ice reservoir call getin("paleomars",paleomars) if (paleomars.eqv..true.) then write(*,*) "paleomars=", paleomars write(*,*) "Orbital parameters from callphys.def" write(*,*) "Enter eccentricity & Lsperi" print *, 'Martian eccentricity (0