Ignore:
Timestamp:
Aug 5, 2009, 4:38:34 PM (15 years ago)
Author:
lguez
Message:

-- Replaced "integer*4" declarations by "integer", "real*8" by

"real(kind=8)" and "real*4" by "real". Note that these are the only
modifications in the files "radiation_AR4.F" and "sw_aeroAR4.F90".

-- Corrected the kind of arguments to "max" and "min".

-- Replaced "nH" edit descriptors, which is a deleted feature in

Fortran 95, by character strings.

-- "regr_lat_time_climoz" now allows the pressure coordinate in the

input file to be in descending order.

-- Replaced call to not standard function "float" by call to intrinsic

function "real".

-- Included file "radepsi.h" in "physiq" was not used. Removed it.

The following set of modifications is related to the management of time.

-- In "gcm", "leapfrog" and "sortvarc0", "day_ini" was defined as 1

plus number of days between the reference date "(annee_ref,
day_ref)" and the first day of the current simulation. Changed
definition: "(annee_ref, day_ini)" is the first day of the current
simulation. There is an accompanying modification for "day_end".

-- Corrected bug in call to "ioconf_startdate" in "gcm".

-- Added call to "ioconf_calendar" in "create_etat0_limit".

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ4/branches/LMDZ4-dev/libf/phylmd/regr_lat_time_climoz_m.F90

    r1158 r1220  
    4444
    4545    ! We assume that in the input file:
    46     ! -- the latitude is in degrees and strictly monotonic (as all
    47     ! NetCDF coordinate variables should be);
     46    ! -- latitude is in degrees;
     47    ! -- pressure is in hPa (even though we do not use pressure values
     48    ! here, we write the unit of pressure in the NetCDF header, and we
     49    ! will use the assumption later, when we regrid in pressure).
     50    ! -- latitude and pressure are strictly monotonic (as all NetCDF
     51    ! coordinate variables should be);
    4852    ! -- time increases (even though we do not use values of the input
    4953    ! time coordinate);
    50     ! -- pressure is in hPa and in strictly ascending order (even
    51     ! though we do not use pressure values here, we write the unit of
    52     ! pressure in the NetCDF header, and we will use the assumptions later,
    53     ! when we regrid in pressure).
    5454
    5555    use regr1_step_av_m, only: regr1_step_av
     
    8181    ! ascending order)
    8282
    83     real, pointer:: plev(:) ! pressure level of input data
     83    real, pointer:: plev(:)
     84    ! pressure levels of input data, sorted in strictly ascending order
     85
    8486    logical desc_lat ! latitude in descending order in the input file
     87    logical desc_plev ! pressure levels in descending order in the input file
    8588
    8689    real, pointer:: o3_in(:, :, :) ! (n_lat, n_plev, 12 or 0:13)
     
    110113
    111114    integer j
    112 
    113     integer varid_in, varid_out, varid_plev, varid_time
    114     integer varid
    115     ! (for NetCDF)
     115    integer varid_in, varid_out, varid_plev, varid_time, varid ! (for NetCDF)
    116116
    117117    real, parameter:: tmidmonth(0:13) = (/(-15. + 30. * j, j = 0, 13)/)
     
    153153    call nf95_gw_var(ncid_in, varid, plev)
    154154    n_plev = size(plev)
    155     ! (We only need the pressure coordinate to copy it to the output file.)
     155    ! We only need the pressure coordinate to copy it to the output file.
     156    ! The program "gcm" will assume that pressure levels are in
     157    ! ascending order in the regridded climatology so invert order if
     158    ! necessary:
     159    desc_plev = plev(1) > plev(n_plev)
     160    if (desc_plev) plev = plev(n_plev:1:-1)
    156161
    157162    ! Create the output file and get the variable IDs:
     
    168173    call nf95_gw_var(ncid_in, varid_in, o3_in)
    169174    if (desc_lat) o3_in = o3_in(n_lat:1:-1, :, :)
     175    if (desc_plev) o3_in = o3_in(:, n_plev:1:-1, :)
    170176
    171177    call nf95_close(ncid_in)
     
    219225    use netcdf95, only: nf95_create, nf95_def_dim, nf95_def_var, &
    220226         nf95_put_att, nf95_enddef, nf95_copy_att, nf95_put_var
    221     use netcdf, only: nf90_clobber, nf90_float, nf90_copy_att, nf90_global
     227    use netcdf, only: nf90_clobber, nf90_float, nf90_global
    222228
    223229    integer, intent(in):: ncid_in, n_plev
     
    276282    ! Global attributes:
    277283
    278     call nf95_put_att(ncid_out, nf90_global, "comment", "Regridded for LMDZ")
    279 
    280284    ! The following commands, copying attributes, may fail.
    281285    ! That is OK.
     
    290294    call handle_err_copy_att("title")
    291295
     296    call nf95_copy_att(ncid_in, nf90_global, "institution", ncid_out, &
     297         nf90_global, ncerr)
     298    call handle_err_copy_att("institution")
     299
    292300    call nf95_copy_att(ncid_in, nf90_global, "source", ncid_out, nf90_global, &
    293301         ncerr)
    294302    call handle_err_copy_att("source")
    295303
     304    call nf95_put_att(ncid_out, nf90_global, "comment", "Regridded for LMDZ")
     305
    296306    call nf95_enddef(ncid_out)
    297307
     
    311321
    312322      if (ncerr /= nf90_noerr) then
    313          print *, "regr_lat_time_climoz_m prepare_out nf90_copy_att " &
     323         print *, "regr_lat_time_climoz_m prepare_out nf95_copy_att " &
    314324              // att_name // " -- " // trim(nf90_strerror(ncerr))
    315325      end if
Note: See TracChangeset for help on using the changeset viewer.