Ignore:
Timestamp:
Oct 17, 2013, 4:19:26 PM (11 years ago)
Author:
tnavarro
Message:

Better handling of the first date of the file : Read and write the controle field, if possible, for zrecast, localtime and concatnc + cosmetic change in lslin. This is retrocompatible with previous versions.

File:
1 edited

Legend:

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

    r909 r1073  
    6060!              constants (radius, R, etc.) are now read from file
    6161! TN 01/2013 : Adapted for large output files with at least 2 variables > 2 GiB
     62! TN 10/2013 : Read and write controle field, if available
    6263!
    6364implicit none
     
    8081integer nbvar4dinfile ! # of 4D (lon,lat,alt,time) variables in input file
    8182integer outfid ! NetCDF output file ID
    82 integer lon_dimid,lat_dimid,alt_dimid,time_dimid ! NetCDF dimension IDs
    83 integer lon_varid,lat_varid,alt_varid,time_varid
     83integer lon_dimid,lat_dimid,alt_dimid,time_dimid,ctl_dimid ! NetCDF dimension IDs
     84integer lon_varid,lat_varid,alt_varid,time_varid,ctl_varid
    8485integer gcm_layers_dimid ! NetCDF dimension ID for # of layers in GCM
    8586integer sigma_varid,aps_varid,bps_varid
     
    103104integer latlength ! # of grid points along latitude
    104105integer altlength ! # of grid point along altitude (of input datasets)
     106real,dimension(:),allocatable :: ctl ! controle
     107integer ctllength ! # of grid points along controle
    105108real,dimension(:),allocatable :: time ! time
    106109integer timelength ! # of points along time
     
    439442  if (ierr.ne.NF_NOERR) then
    440443      stop "Error: Failed to get altitude length"
     444  endif
     445endif
     446
     447! controle
     448ierr=NF_INQ_DIMID(infid,"index",tmpdimid)
     449if (ierr.ne.NF_NOERR) then
     450  write(*,*) "Could not get controle dimension ID"
     451  ctllength = 0
     452else
     453  ierr=NF_INQ_VARID(infid,"controle",tmpvarid)
     454  if (ierr.ne.NF_NOERR) then
     455    stop "Error: Failed to get controle ID"
     456  else
     457    ierr=NF_INQ_DIMLEN(infid,tmpdimid,ctllength)
     458    if (ierr.ne.NF_NOERR) then
     459      stop "Error: Failed to get controle length"
     460    else
     461      allocate(ctl(ctllength),stat=ierr)
     462      if (ierr.ne.0) then
     463        write(*,*) "Error: Failed to allocate ctl(ctllength)"
     464        write(*,*) "     ctllength=",ctllength
     465        stop
     466      endif
     467      ierr=NF_GET_VAR_REAL(infid,tmpvarid,ctl)
     468      if (ierr.ne.NF_NOERR) then
     469        stop "Error: Failed reading controle"
     470      endif
     471    endif
    441472  endif
    442473endif
     
    10301061endif
    10311062
     1063! controle
     1064if (ctllength .ne. 0) then
     1065  ierr=NF_DEF_DIM(outfid,"index",ctllength,ctl_dimid)
     1066  if (ierr.ne.NF_NOERR) then
     1067    write(*,*) "Error: Could not define controle dimension"
     1068    write(*,*) NF_STRERROR(ierr)
     1069    stop
     1070  endif
     1071endif
     1072
    10321073! GCM layers (for sigma or aps and bps)
    10331074ierr=NF_DEF_DIM(outfid,"GCM_layers",altlength,gcm_layers_dimid)
     
    11621203  endif
    11631204endif ! of if (have_sigma)
     1205
     1206! controle
     1207if (ctllength .ne. 0) then
     1208  ierr=NF_DEF_VAR(outfid,"controle",NF_REAL,1,ctl_dimid,ctl_varid)
     1209  if (ierr.ne.NF_NOERR) then
     1210    write(*,*) "Error: Could not define controle variable"
     1211    write(*,*) NF_STRERROR(ierr)
     1212    stop
     1213  endif
     1214
     1215  ! controle attributes
     1216  text='Control parameters'
     1217  ierr=NF_PUT_ATT_TEXT(outfid,ctl_varid,'title',len_trim(text),text)
     1218  if (ierr.ne.NF_NOERR) then
     1219    stop "Error: Problem writing title for controle"
     1220  endif
     1221endif
    11641222
    11651223! GCM_layers
     
    15091567endif
    15101568
     1569! Write controle
     1570if (ctllength .ne. 0) then
     1571  ierr=NF_PUT_VAR_REAL(outfid,ctl_varid,ctl)
     1572  if (ierr.ne.NF_NOERR) then
     1573    write(*,*) "Error: Could not write controle data to output file"
     1574    write(*,*) NF_STRERROR(ierr)
     1575    stop
     1576  endif
     1577endif
     1578
    15111579! Write time
    15121580ierr=NF_PUT_VARA_REAL(outfid,time_varid,1,timelength,time)
     
    17711839subroutine init_planet_const(infid)
    17721840! initialize planetary constants using the "controle" array in the file
    1773 ! if "cointrole" array not found in file, look for it in "diagfi.nc"
     1841! if "controle" array not found in file, look for it in "diagfi.nc"
    17741842use planet_const
    17751843use netcdf
Note: See TracChangeset for help on using the changeset viewer.