Changeset 632 for trunk/LMDZ.MARS/util


Ignore:
Timestamp:
Apr 25, 2012, 10:23:52 AM (13 years ago)
Author:
emillour
Message:

Mars GCM utilities: minor improvements:

  • zrecast also looks for a 'diagfi1.nc' file when looking for 'phisinit' variable
  • concatnc now also includes 'aires' in its output file

EM

Location:
trunk/LMDZ.MARS/util
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/util/concatnc.F90

    r410 r632  
    713713!==============================================================================
    714714! Purpose:
    715 ! Copy ap() , bp(), aps(), bps() and phisinit() from input file to outpout file
     715! Copy ap() , bp(), aps(), bps(), aire() and phisinit()
     716! from input file to outpout file
    716717!==============================================================================
    717718! Remarks:
     
    743744real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates
    744745real,dimension(:),allocatable :: sigma ! sigma levels
     746real,dimension(:,:),allocatable :: aire ! mesh areas
    745747real,dimension(:,:),allocatable :: phisinit ! Ground geopotential
    746748integer :: apsid,bpsid,sigmaid,phisinitid
     
    748750integer :: ierr
    749751integer :: tmpvarid ! temporary variable ID
     752logical :: area ! is "aire" available ?
    750753logical :: phis ! is "phisinit" available ?
    751754logical :: hybrid ! are "aps" and "bps" available ?
     
    843846endif ! of if (.not.hybrid)
    844847
     848! mesh area
     849allocate(aire(lonlen,latlen),stat=ierr)
     850if (ierr.ne.0) then
     851  write(*,*) "init2: failed to allocate aire!"
     852  stop
     853endif
     854ierr=NF_INQ_VARID(infid,"aire",tmpvarid)
     855if (ierr.ne.NF_NOERR) then
     856  write(*,*)"init2 warning: Failed to get aire ID."
     857  area = .false.
     858else
     859  ierr=NF_GET_VAR_REAL(infid,tmpvarid,aire)
     860  if (ierr.ne.NF_NOERR) then
     861    stop "init2 Error: Failed reading aire"
     862  endif
     863  area = .true.
     864endif
     865
    845866! ground geopotential phisinit
    846867allocate(phisinit(lonlen,latlen),stat=ierr)
     
    958979
    959980!==============================================================================
    960 ! 2.2. phisinit()
    961 !==============================================================================
     981! 2.2. aire() and phisinit()
     982!==============================================================================
     983
     984if (area) then
     985  ! define aire
     986  call def_var(nout,"aire","Mesh area","m2",2,&
     987           (/londimout,latdimout/),tmpvarid,ierr)
     988  if (ierr.ne.NF_NOERR) then
     989     stop "init2 Error: Failed to def_var aire"
     990  endif
     991 
     992  ! write aire
     993#ifdef NC_DOUBLE
     994  ierr=NF_PUT_VAR_DOUBLE(outfid,tmpvarid,aire)
     995#else
     996  ierr=NF_PUT_VAR_REAL(outfid,tmpvarid,aire)
     997#endif
     998  if (ierr.ne.NF_NOERR) then
     999    stop "init2 Error: Failed to write aire"
     1000  endif
     1001endif ! of if (area)
    9621002
    9631003IF (phis) THEN
     
    9901030if (allocated(sigma)) deallocate(sigma)
    9911031if (allocated(phisinit)) deallocate(phisinit)
     1032if (allocated(aire)) deallocate(aire)
    9921033
    9931034end subroutine init2
  • trunk/LMDZ.MARS/util/zrecast.F90

    r410 r632  
    499499  if (ierr.ne.NF_NOERR) then
    500500    write(*,*) "Problem: Could not find/open that file"
    501     infile2="phisinit.nc"
     501    infile2="diagfi1.nc"
    502502    write(*,*) "         Trying file ",trim(infile2)
    503503    ierr=NF_OPEN(infile2,NF_NOWRITE,infid2)
    504504    if (ierr.ne.NF_NOERR) then
    505       write(*,*) "Error: Could not open that file either"
    506       stop "Might as well stop here"
     505      write(*,*) "Problem: Could not find/open that file"
     506      infile2="phisinit.nc"
     507      write(*,*) "         Trying file ",trim(infile2)
     508      ierr=NF_OPEN(infile2,NF_NOWRITE,infid2)
     509      if (ierr.ne.NF_NOERR) then
     510        write(*,*) "Error: Could not open that file either"
     511        stop "Might as well stop here"
     512      endif
    507513    endif
    508514  endif
Note: See TracChangeset for help on using the changeset viewer.