subroutine physdem0(filename,lonfi,latfi,nsoil,ngrid,nlay,nq, . phystep,day_ini,time,airefi, . alb,ith,pzmea,pzstd,pzsig,pzgam,pzthe) use infotrac, only: nqtot, tnom use comsoil_h, only: inertiedat, volcapa, mlayer use comgeomfi_h, only: area use surfdat_h, only: albedodat, zmea, zstd, zsig, zgam, zthe, & z0_default, albedice, emisice, emissiv, & iceradius, dtemisice, phisfi, z0 use yomaer_h, only: tauvis implicit none c c Variables qui NE possèdent PAS un axe des temps : c airefi c alb c ith c pzmea c pzstd c pzsig c pzgam c pzthe c c------------------------------------------------------------- c c create physics (re-)start data file "restartfi.nc" c c c #include "dimensions.h" #include "paramet.h" c----------------------------------------------------------------------- #include "comvert.h" #include "comgeom2.h" #include "control.h" #include "comdissnew.h" #include "logic.h" #include "ener.h" #include "netcdf.inc" !#include "dimphys.h" !#include "advtrac.h" #include "callkeys.h" c INTEGER nid,iq INTEGER, parameter :: ivap=1 REAL, parameter :: qsolmax= 150.0 character (len=*) :: filename character (len=7) :: str7 REAL day_ini INTEGER,INTENT(IN) :: nsoil,ngrid,nlay,nq integer ierr,idim1,idim2,idim3,idim4,idim5,idim6,nvarid c REAL phystep,time REAL latfi(ngrid), lonfi(ngrid) ! REAL champhys(ngrid) INTEGER length PARAMETER (length=100) REAL tab_cntrl(length) c ! EXTERNAL defrun_new,iniconst,geopot,inigeom,massdair,pression ! EXTERNAL exner_hyb , SSUM c #include "serre.h" #include "clesph0.h" #include "fxyprim.h" !#include "comgeomfi.h" !#include "surfdat.h" !#include "comsoil.h" #include "planete.h" !#include "dimradmars.h" !#include "yomaer.h" #include "comcstfi.h" real airefi(ngrid) real alb(ngrid),ith(ngrid,nsoil) real pzmea(ngrid),pzstd(ngrid) real pzsig(ngrid),pzgam(ngrid),pzthe(ngrid) integer ig ! flag which identifies if we are using old tracer names (qsurf01,...) logical :: oldtracernames=.false. integer :: count character(len=30) :: txt ! to store some text ! indexes of water vapour & water ice tracers (if any): integer :: i_h2o_vap=0 integer :: i_h2o_ice=0 c----------------------------------------------------------------------- ! copy airefi(:) to area(:) CALL SCOPY(ngrid,airefi,1,area,1) ! note: area() is defined in comgeomfi_h DO ig=1,ngrid albedodat(ig)=alb(ig) ! note: albedodat() is defined in surfdat_h zmea(ig)=pzmea(ig) ! note: zmea() is defined in surfdat_h zstd(ig)=pzstd(ig) ! note: zstd() is defined in surfdat_h zsig(ig)=pzsig(ig) ! note: zsig() is defined in surfdat_h zgam(ig)=pzgam(ig) ! note: zgam() is defined in surfdat_h zthe(ig)=pzthe(ig) ! note: zthe() is defined in surfdat_h ENDDO inertiedat(:,:)=ith(:,:) ! note inertiedat() is defined in comsoil_h c c things to store in the physics start file: c ierr = NF_CREATE(adjustl(filename), & IOR(NF_CLOBBER,NF_64BIT_OFFSET), nid) IF (ierr.NE.NF_NOERR) THEN WRITE(6,*)'physdem0: Problem creating file ',adjustl(filename) write(6,*) NF_STRERROR(ierr) CALL ABORT ENDIF c ierr = NF_PUT_ATT_TEXT (nid, NF_GLOBAL, "title", 18, . "Physics start file") c ierr = NF_DEF_DIM (nid,"index",length,idim1) if (ierr.ne.NF_NOERR) then WRITE(6,*)'physdem0: Problem defining index dimension' write(6,*) NF_STRERROR(ierr) call abort endif c ierr = NF_DEF_DIM (nid,"physical_points",ngrid,idim2) if (ierr.ne.NF_NOERR) then WRITE(6,*)'physdem0: Problem defining physical_points dimension' write(6,*) NF_STRERROR(ierr) call abort endif c ierr = NF_DEF_DIM (nid,"subsurface_layers",nsoil,idim3) if (ierr.ne.NF_NOERR) then WRITE(6,*)'physdem0: Problem defining subsurface_layers dimension' write(6,*) NF_STRERROR(ierr) call abort endif c ! ierr = NF_DEF_DIM (nid,"nlayer+1",llm+1,idim4) ierr = NF_DEF_DIM (nid,"nlayer_plus_1",llm+1,idim4) if (ierr.ne.NF_NOERR) then WRITE(6,*)'physdem0: Problem defining nlayer+1 dimension' write(6,*) NF_STRERROR(ierr) call abort endif c ierr = NF_DEF_DIM (nid,"number_of_advected_fields",nq,idim5) if (ierr.ne.NF_NOERR) then WRITE(6,*)'physdem0: Problem defining advected fields dimension' WRITE(6,*)' nq = ',nq,' and ierr = ', ierr write(6,*) NF_STRERROR(ierr) endif c ierr = NF_DEF_DIM (nid,"Time",NF_UNLIMITED,idim6) if (ierr.ne.NF_NOERR) then WRITE(6,*)'physdem0: Problem defining time dimension' WRITE(6,*)' ierr = ', ierr write(6,*) NF_STRERROR(ierr) endif ierr = NF_ENDDEF(nid) ! exit NetCDF define mode c clear tab_cntrl(:) array DO ierr = 1, length tab_cntrl(ierr) = 0.0 ENDDO write(*,*) "physdem0: ngrid: ",ngrid ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc c Fill control array tab_cntrl(:) with parameters for this run ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc c Informations on the physics grid tab_cntrl(1) = float(ngrid) ! number of nodes on physics grid tab_cntrl(2) = float(nlay) ! number of atmospheric layers tab_cntrl(3) = day_ini + int(time) ! initial day tab_cntrl(4) = time -int(time) ! initial time of day c Informations about Mars, used by dynamics and physics tab_cntrl(5) = rad ! radius of Mars (m) ~3397200 tab_cntrl(6) = omeg ! rotation rate (rad.s-1) tab_cntrl(7) = g ! gravity (m.s-2) ~3.72 tab_cntrl(8) = mugaz ! Molar mass of the atmosphere (g.mol-1) ~43.49 tab_cntrl(9) = rcp ! = r/cp ~0.256793 (=kappa dans dynamique) tab_cntrl(10) = daysec ! length of a sol (s) ~88775 tab_cntrl(11) = phystep ! time step in the physics tab_cntrl(12) = 0. tab_cntrl(13) = 0. c Informations about Mars, only for physics tab_cntrl(14) = year_day ! length of year (sols) ~668.6 tab_cntrl(15) = periheli ! min. Sun-Mars distance (Mkm) ~206.66 tab_cntrl(16) = aphelie ! max. SUn-Mars distance (Mkm) ~249.22 tab_cntrl(17) = peri_day ! date of perihelion (sols since N. spring) tab_cntrl(18) = obliquit ! Obliquity of the planet (deg) ~23.98 c Boundary layer and turbulence tab_cntrl(19) = z0_default ! default surface roughness (m) ~0.01 tab_cntrl(20) = lmixmin ! mixing length ~100 tab_cntrl(21) = emin_turb ! minimal energy ~1.e-8 c Optical properties of polar caps and ground emissivity tab_cntrl(22) = albedice(1) ! Albedo of northern cap ~0.5 tab_cntrl(23) = albedice(2) ! Albedo of southern cap ~0.5 tab_cntrl(24) = emisice(1) ! Emissivity of northern cap ~0.95 tab_cntrl(25) = emisice(2) ! Emissivity of southern cap ~0.95 tab_cntrl(26) = emissiv ! Emissivity of martian soil ~.95 tab_cntrl(31) = iceradius(1) ! mean scat radius of CO2 snow (north) tab_cntrl(32) = iceradius(2) ! mean scat radius of CO2 snow (south) tab_cntrl(33) = dtemisice(1) ! time scale for snow metamorphism (north) tab_cntrl(34) = dtemisice(2) ! time scale for snow metamorphism (south) c dust aerosol properties tab_cntrl(27) = tauvis ! mean visible optical depth tab_cntrl(28) = 0. tab_cntrl(29) = 0. tab_cntrl(30) = 0. ! Soil properties: tab_cntrl(35) = volcapa ! soil volumetric heat capacity c ierr = NF_REDEF (nid) ! Enter NetCDF (re-)define mode IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to swich to NetCDF define mode' CALL abort ENDIF ! define variable #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "controle", NF_DOUBLE, 1, idim1,nvarid) #else ierr = NF_DEF_VAR (nid, "controle", NF_FLOAT, 1, idim1,nvarid) #endif IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to define controle' CALL abort ENDIF ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 18, . "Control parameters") IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to define controle title attribute' CALL abort ENDIF ierr = NF_ENDDEF(nid) ! Leave NetCDF define mode IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to swich out of NetCDF define mode' CALL abort ENDIF ! write variable #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,tab_cntrl) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,tab_cntrl) #endif IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to write controle data' CALL abort ENDIF ! write mid-layer depths mlayer() !known from comsoil_h ierr = NF_REDEF (nid) ! Enter NetCDF (re-)define mode ! define variable #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid,"soildepth",NF_DOUBLE,1,idim3,nvarid) #else ierr = NF_DEF_VAR (nid,"soildepth",NF_FLOAT,1,idim3,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 20, . "Soil mid-layer depth") ierr = NF_ENDDEF(nid) ! Leave NetCDF define mode ! write variable #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,mlayer) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,mlayer) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "longitude", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "longitude", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 26, . "Longitudes of physics grid") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,lonfi) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,lonfi) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "latitude", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "latitude", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 25, . "Latitudes of physics grid") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,latfi) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,latfi) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "area", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "area", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 9, . "Mesh area") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,area) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,area) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "phisfi", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "phisfi", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 27, . "Geopotential at the surface") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,phisfi) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,phisfi) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "albedodat", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "albedodat", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 21, . "Albedo of bare ground") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,albedodat) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,albedodat) #endif c c some data for Francois Lott's programs c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "ZMEA", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "ZMEA", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 19, . "Relief: mean relief") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,zmea) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,zmea) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "ZSTD", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "ZSTD", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 26, . "Relief: standard deviation") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,zstd) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,zstd) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "ZSIG", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "ZSIG", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 23, . "Relief: sigma parameter") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,zsig) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,zsig) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "ZGAM", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "ZGAM", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 23, . "Relief: gamma parameter") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,zgam) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,zgam) #endif c ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "ZTHE", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "ZTHE", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 23, . "Relief: theta parameter") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,zthe) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,zthe) #endif ! Soil Thermal inertia ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid,"inertiedat",NF_DOUBLE, & 2,(/idim2,idim3/),nvarid) #else ierr = NF_DEF_VAR (nid,"inertiedat",NF_FLOAT, & 2,(/idim2,idim3/),nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 20, . "Soil thermal inertia") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,inertiedat) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,inertiedat) #endif ! surface roughness length (z0 is a common in surfdat_h) ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "z0", NF_DOUBLE, 1, idim2,nvarid) #else ierr = NF_DEF_VAR (nid, "z0", NF_FLOAT, 1, idim2,nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 24, . "Surface roughness length") ierr = NF_ENDDEF(nid) #ifdef NC_DOUBLE ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,z0) #else ierr = NF_PUT_VAR_REAL (nid,nvarid,z0) #endif c--------------------- Time dependent variables --------------------- c------------------------------------------------------------------- c------------------------------------------------------------------- ! Time ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "Time", NF_DOUBLE, 1, . (/idim6/), nvarid) #else ierr = NF_DEF_VAR (nid, "Time", NF_FLOAT, 1, . (/idim6/), nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 19, . "Temps de simulation") ierr = NF_ENDDEF(nid) ! CO2 Ice Cover ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "co2ice", NF_DOUBLE, 2, . (/idim2,idim6/), nvarid) #else ierr = NF_DEF_VAR (nid, "co2ice", NF_FLOAT, 2, . (/idim2,idim6/), nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 13, . "CO2 ice cover") ierr = NF_ENDDEF(nid) ! Surface temperature ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "tsurf", NF_DOUBLE, 2, . (/idim2,idim6/), nvarid) #else ierr = NF_DEF_VAR (nid, "tsurf", NF_FLOAT, 2, . (/idim2,idim6/), nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 19, . "Surface temperature") ierr = NF_ENDDEF(nid) ! Soil temperature ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "tsoil", NF_DOUBLE, 3, . (/idim2,idim3,idim6/), nvarid) #else ierr = NF_DEF_VAR (nid, "tsoil", NF_FLOAT, 3, . (/idim2,idim3,idim6/), nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 16, . "Soil temperature") ierr = NF_ENDDEF(nid) c emissivity ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid,"emis", NF_DOUBLE, 2, . (/idim2,idim6/), nvarid) #else ierr = NF_DEF_VAR (nid, "emis", NF_FLOAT, 2, . (/idim2,idim6/), nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 18, . "Surface emissivity") ierr = NF_ENDDEF(nid) c planetary boundary layer ierr = NF_REDEF (nid) #ifdef NC_DOUBLE ierr = NF_DEF_VAR (nid, "q2", NF_DOUBLE, 3, . (/idim2,idim4,idim6/), nvarid) #else ierr = NF_DEF_VAR (nid, "q2", NF_FLOAT, 3, . (/idim2,idim4,idim6/), nvarid) #endif ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 17, . "pbl wind variance") ierr = NF_ENDDEF(nid) c tracers ! ! Preliminary stuff: check if tracers follow old naming convention (qsurf01, ! qsurf02, ...) count=0 do iq=1,nqtot txt= " " write(txt,'(a1,i2.2)')'q',iq if (txt.ne.tnom(iq)) then ! use tracer names stored in dynamics ! did not find old tracer name exit ! might as well stop here else ! found old tracer name count=count+1 endif enddo if (count.eq.nqtot) then write(*,*) "physdem0:tracers seem to follow old naming ", & "convention (qsurf01,qsurf02,...)" write(*,*) " => will work for now ... " write(*,*) " but you should run newstart to rename them" oldtracernames=.true. endif ! of if (count.eq.nqtot) IF(nq.GE.1) THEN ! preliminary stuff: look for water vapour & water ice tracers (if any) if (.not.oldtracernames) then do iq=1,nq if (tnom(iq).eq."h2o_vap") then i_h2o_vap=iq endif if (tnom(iq).eq."h2o_ice") then i_h2o_ice=iq endif enddo ! of iq=1,nq ! handle special case of only water vapour tracer (no ice) if ((i_h2o_vap.ne.0).and.(i_h2o_ice.eq.0)) then ! then the index of (surface) ice is i_h2o_vap i_h2o_ice=i_h2o_vap endif endif ! of if (.not.oldtracernames) DO iq=1,nq IF (oldtracernames) THEN txt=" " write(txt,'(a5,i2.2)')'qsurf',iq ELSE txt=tnom(iq) ! Exception: there is no water vapour surface tracer if (txt.eq."h2o_vap") then write(*,*)"physdem0: skipping water vapour tracer" if (i_h2o_ice.eq.i_h2o_vap) then ! then there is no "water ice" tracer; but still ! there is some water ice on the surface write(*,*)" writing water ice instead" txt="h2o_ice" else ! there is a "water ice" tracer which has been / will be ! delt with in due time cycle endif ! of if (igcm_h2o_ice.eq.igcm_h2o_vap) endif ! of if (txt.eq."h2o_vap") ENDIF ! of IF (oldtracernames) ierr=NF_REDEF(nid) IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to swich to NetCDF define mode' CALL abort ENDIF #ifdef NC_DOUBLE ierr=NF_DEF_VAR(nid,txt,NF_DOUBLE,2,(/idim2,idim6/),nvarid) #else ierr=NF_DEF_VAR(nid,txt,NF_FLOAT,2,(/idim2,idim6/),nvarid) #endif IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to define ',trim(txt) CALL abort ENDIF ierr=NF_PUT_ATT_TEXT (nid, nvarid, "title", 17, & "tracer on surface") IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to define ',trim(txt), & ' title attribute' CALL abort ENDIF ierr=NF_ENDDEF(nid) IF (ierr.NE.NF_NOERR) THEN PRINT*, 'physdem0: Failed to swich out of define mode' CALL abort ENDIF ENDDO ENDIF c close file ierr = NF_CLOSE(nid) END cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc subroutine physdem1(filename,nsoil,ngrid,nlay,nq, . phystep,time, . tsurf,tsoil,co2ice,emis,q2,qsurf) use infotrac, only: nqtot, tnom implicit none c c Variables qui possèdent un axe des temps : c tsurf c tsoil c co2ice c q2 c qsurf c emis (modifié dans newcondens --> co2snow) c c------------------------------------------------------------- c c write physics (re-)start data file "restartfi.nc" c c c #include "dimensions.h" #include "paramet.h" c----------------------------------------------------------------------- #include "comvert.h" #include "comgeom2.h" #include "control.h" #include "comdissnew.h" #include "logic.h" #include "ener.h" #include "netcdf.inc" !#include "dimphys.h" !#include "advtrac.h" #include "callkeys.h" c INTEGER nid,iq INTEGER, parameter :: ivap=1 REAL, parameter :: qsolmax= 150.0 character (len=*) :: filename character (len=7) :: str7 REAL day_ini INTEGER,INTENT(IN) :: nsoil,ngrid,nlay,nq integer ierr,nvarid c REAL phystep,time ! REAL champhys(ngrid) REAL tsurf(ngrid) INTEGER nb INTEGER edges(3),corner(3) c ! EXTERNAL defrun_new,iniconst,geopot,inigeom,massdair,pression ! EXTERNAL exner_hyb , SSUM c #include "serre.h" #include "clesph0.h" #include "fxyprim.h" !#include "comgeomfi.h" !#include "surfdat.h" !#include "comsoil.h" #include "planete.h" !#include "dimradmars.h" !#include "yomaer.h" #include "comcstfi.h" real co2ice(ngrid),tsoil(ngrid,nsoil),emis(ngrid) real q2(ngrid,nlay+1),qsurf(ngrid,nq) integer ig ! flag which identifies if we are using old tracer names (qsurf01,...) logical,save :: oldtracernames=.false. logical,save :: firstcall = .true. integer :: count character(len=30) :: txt ! to store some text ! indexes of water vapour & water ice tracers (if any): integer,save :: i_h2o_vap=0 integer,save :: i_h2o_ice=0 c----------------------------------------------------------------------- ierr = NF_OPEN(filename, NF_WRITE, nid) IF (ierr .NE. NF_NOERR) THEN PRINT*, "Pb. d ouverture "//filename CALL abort ENDIF c On a single run, different files can be written with physdem1. c Therefore, get the last time index from the file itself: ierr = NF_INQ_DIMID(nid,"Time",nvarid) ierr = NF_INQ_DIMLEN(nid,nvarid,nb) c Ecriture/extension de la coordonnee temps nb = nb + 1 ierr = NF_INQ_VARID(nid, "Time", nvarid) IF (ierr .NE. NF_NOERR) THEN print *, NF_STRERROR(ierr) print*,'Variable Time n est pas definie' CALL abort ENDIF #ifdef NC_DOUBLE ierr = NF_PUT_VAR1_DOUBLE (nid,nvarid,nb,time) #else ierr = NF_PUT_VAR1_REAL (nid,nvarid,nb,time) #endif IF (ierr .NE. NF_NOERR) THEN print*, "Erreur ecriture temps!!" print*, NF_STRERROR(ierr) ENDIF !PRINT*, "Enregistrement pour ", nb, time c Write the physical fields ! Soil temperature corner(1)=1 corner(2)=1 corner(3)=nb edges(1)=ngrid edges(2)=nsoil edges(3)=1 ierr = NF_INQ_VARID(nid, "tsoil", nvarid) IF (ierr .NE. NF_NOERR) THEN PRINT*, "Variable tsoil n est pas definie" CALL abort ENDIF #ifdef NC_DOUBLE ierr = NF_PUT_VARA_DOUBLE (nid,nvarid,corner,edges,tsoil) #else ierr = NF_PUT_VARA_REAL (nid,nvarid,corner,edges,tsoil) #endif IF (ierr .NE. NF_NOERR) THEN print*, "Erreur ecriture tsoil!!" print*, NF_STRERROR(ierr) ENDIF c planetary boundary layer corner(1)=1 corner(2)=1 corner(3)=nb edges(1)=ngrid edges(2)=llm+1 edges(3)=1 ierr = NF_INQ_VARID(nid, "q2", nvarid) IF (ierr .NE. NF_NOERR) THEN PRINT*, "Variable q2 n est pas definie" CALL abort ENDIF #ifdef NC_DOUBLE ierr = NF_PUT_VARA_DOUBLE (nid,nvarid,corner,edges,q2) #else ierr = NF_PUT_VARA_REAL (nid,nvarid,corner,edges,q2) #endif IF (ierr .NE. NF_NOERR) THEN print*, "Erreur ecriture q2!!" print*, NF_STRERROR(ierr) ENDIF c Following corner and egdes are the same for co2ice,tsurf,qusrf,emis and tracers corner(1)=1 corner(2)=nb edges(1)=ngrid edges(2)=1 ! CO2 Ice Cover ierr = NF_INQ_VARID(nid, "co2ice", nvarid) IF (ierr .NE. NF_NOERR) THEN PRINT*, "Variable co2ice n est pas definie" CALL abort ENDIF #ifdef NC_DOUBLE ierr = NF_PUT_VARA_DOUBLE (nid,nvarid,corner,edges,co2ice) #else ierr = NF_PUT_VARA_REAL (nid,nvarid,corner,edges,co2ice) #endif IF (ierr .NE. NF_NOERR) THEN print*, "Erreur ecriture co2ice!!" print*, NF_STRERROR(ierr) ENDIF ! Surface temperature ierr = NF_INQ_VARID(nid, "tsurf", nvarid) IF (ierr .NE. NF_NOERR) THEN PRINT*, "Variable tsurf n est pas definie" CALL abort ENDIF #ifdef NC_DOUBLE ierr = NF_PUT_VARA_DOUBLE (nid,nvarid,corner,edges,tsurf) #else ierr = NF_PUT_VARA_REAL (nid,nvarid,corner,edges,tsurf) #endif IF (ierr .NE. NF_NOERR) THEN print*, "Erreur ecriture tsurf!!" print*, NF_STRERROR(ierr) ENDIF c emissivity ierr = NF_INQ_VARID(nid, "emis", nvarid) IF (ierr .NE. NF_NOERR) THEN PRINT*, "Variable emis n est pas definie" CALL abort ENDIF #ifdef NC_DOUBLE ierr = NF_PUT_VARA_DOUBLE (nid,nvarid,corner,edges,emis) #else ierr = NF_PUT_VARA_REAL (nid,nvarid,corner,edges,emis) #endif IF (ierr .NE. NF_NOERR) THEN print*, "Erreur ecriture emis!!" print*, NF_STRERROR(ierr) ENDIF c tracers ! Preliminary stuff: check if tracers follow old naming convention (qsurf01, ! qsurf02, ...) IF (firstcall) THEN count=0 do iq=1,nqtot txt= " " write(txt,'(a1,i2.2)')'q',iq if (txt.ne.tnom(iq)) then ! use tracer names stored in dynamics ! did not find old tracer name exit ! might as well stop here else ! found old tracer name count=count+1 endif enddo if (count.eq.nqtot) then write(*,*) "physdem1:tracers seem to follow old naming ", & "convention (qsurf01,qsurf02,...)" write(*,*) " => will work for now ... " write(*,*) " but you should run newstart to rename them" oldtracernames=.true. endif ! of if (count.eq.nqtot) ENDIF ! of if(firstcall) ! If computing water cycle with ice, move surface ice ! back to qsurf(nqtot) IF (oldtracernames .and. water) THEN !"loop" to avoid potential out-of-bounds on arrays write(*,*)'physdem1: moving surface water ice to index ',nqtot do iq=nqtot,nqtot qsurf(1:ngrid,iq)=qsurf(1:ngrid,iq-1) qsurf(1:ngrid,iq-1)=0 enddo ENDIF IF(nq.GE.1) THEN IF (firstcall) THEN ! preliminary stuff: look for water vapour & water ice tracers (if any) if (.not.oldtracernames) then do iq=1,nq if (tnom(iq).eq."h2o_vap") then i_h2o_vap=iq endif if (tnom(iq).eq."h2o_ice") then i_h2o_ice=iq endif enddo ! of iq=1,nq ! handle special case of only water vapour tracer (no ice) if ((i_h2o_vap.ne.0).and.(i_h2o_ice.eq.0)) then ! then the index of (surface) ice is i_h2o_vap i_h2o_ice=i_h2o_vap endif endif ! of if (.not.oldtracernames) ENDIF ! of if(firstcall) firstcall = .false. DO iq=1,nq IF (oldtracernames) THEN txt=" " write(txt,'(a5,i2.2)')'qsurf',iq ELSE txt=tnom(iq) ! Exception: there is no water vapour surface tracer if (txt.eq."h2o_vap") then if (i_h2o_ice.eq.i_h2o_vap) then ! then there is no "water ice" tracer; but still ! there is some water ice on the surface txt="h2o_ice" else ! there is a "water ice" tracer which has been / will be ! delt with in due time cycle endif ! of if (igcm_h2o_ice.eq.igcm_h2o_vap) endif ! of if (txt.eq."h2o_vap") ENDIF ! of IF (oldtracernames) ierr = NF_INQ_VARID(nid, trim(txt), nvarid) IF (ierr .NE. NF_NOERR) THEN PRINT*, "Variable "//trim(txt)//" n est pas definie" CALL abort ENDIF #ifdef NC_DOUBLE ierr = NF_PUT_VARA_DOUBLE (nid,nvarid,corner,edges, . qsurf(:,iq)) #else ierr = NF_PUT_VARA_REAL (nid,nvarid,corner,edges, . qsurf(:,iq)) #endif IF (ierr .NE. NF_NOERR) THEN print*, "Erreur ecriture "//trim(txt)//"!!" print*, NF_STRERROR(ierr) ENDIF ENDDO ! of DO iq=1,nq ENDIF ! of IF(nq.GE.1) c close file ierr = NF_CLOSE(nid) RETURN END