subroutine readtesassim(ngrid,nlayer,zday,pplev,tauref) ! Reading of the dust assimilation file use netcdf implicit none #include "dimensions.h" #include "dimphys.h" #include "comgeomfi.h" #include "datafile.h" #include "callkeys.h" integer, intent(in) :: ngrid,nlayer real, intent(in) :: zday ! date in martian days real, dimension(ngrid,nlayer+1), intent(in) :: pplev real, dimension(ngrid), intent(out) :: tauref ! Local variables real :: realday integer nid,nvarid,ierr integer ltloop,lsloop,iloop,jloop,varloop,ig real, dimension(2) :: taubuf real tau1(4),tau real alt(4) integer latp(4),lonp(4) integer yinf,ysup,xinf,xsup,tinf,tsup real latinf,latsup,loninf,lonsup real latintmp,lonintmp,latdeg,londeg real colat,dlat,dlon,colattmp logical, save :: firstcall=.true. logical :: timeflag real,save :: radeg,pi integer :: timedim,londim,latdim real, dimension(:), allocatable, save :: lat,lon,time real, dimension(:,:,:), allocatable, save :: tautes integer, save :: timelen,lonlen,latlen character(len=33) :: filename realday=mod(zday,669.) if (firstcall) then firstcall=.false. pi=acos(-1.) radeg=180/pi ! assimilated dust file: (NB: iaervar is a common in "callkeys.h") ! iaervar=4 means read dust_tes.nc file ! iaervar=24 means read dust_tes_MY24.nc file ! iaervar=25 means read dust_tes_MY25.nc file ! iaervar=26 means read dust_tes_MY26.nc file if (iaervar.eq.4) then filename="dust_tes.nc" else if (iaervar.eq.24) then filename="dust_tes_MY24.nc" ! which currently is identical to dust_tes.nc else if (iaervar.eq.25) then filename="dust_tes_MY25.nc" else if (iaervar.eq.26) then filename="dust_tes_MY26.nc" endif ! Note: datafile() is defined in "datafile.h" !ierr=NF_OPEN(trim(datafile)//"/"//trim(filename),NF_NOWRITE,nid) ierr=nf90_open(trim(datafile)//"/"//trim(filename),NF90_NOWRITE,nid) IF (ierr.NE.nf90_noerr) THEN write(*,*)'Problem opening ',trim(filename),' (in phymars/readtesassim.F90)' write(*,*)'It should be in :',trim(datafile),'/' write(*,*)'1) You can change this directory address in ' write(*,*)' file phymars/datafile.h' write(*,*)'2) If necessary, dust.nc (and other datafiles)' write(*,*)' can be obtained online on:' write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile' CALL ABORT ENDIF ierr=nf90_inq_dimid(nid,"Time",timedim) ierr=nf90_inquire_dimension(nid,timedim,len=timelen) ierr=nf90_inq_dimid(nid,"latitude",latdim) ierr=nf90_inquire_dimension(nid,latdim,len=latlen) ierr=nf90_inq_dimid(nid,"longitude",londim) ierr=nf90_inquire_dimension(nid,londim,len=lonlen) allocate(tautes(lonlen,latlen,timelen)) allocate(lat(latlen), lon(lonlen), time(timelen)) ierr=nf90_inq_varid(nid,"dustop",nvarid) ierr=nf90_get_var(nid,nvarid,tautes) IF (ierr .NE. nf90_noerr) THEN PRINT*, "Error: Readtesassim not found" write(*,*)trim(nf90_strerror(ierr)) stop ENDIF ierr=nf90_inq_varid(nid,"Time",nvarid) ierr=nf90_get_var(nid,nvarid,time) IF (ierr .NE. nf90_noerr) THEN PRINT*, "Error: Readtesassim