Ignore:
Timestamp:
Nov 17, 2009, 2:00:14 PM (15 years ago)
Author:
lguez
Message:

1) Reactivated ability to read ozone (that was deactivated because of
dependency on version of IOIPSL). Added ability to read a pressure
coordinate in Pa in "regr_lat_time_climoz".

2) Added the ability to read a second ozone climatology, corresponding to
daylight ozone:

-- "read_climoz" is now an integer variable, instead of a logical
variable.

-- Added argument "read_climoz" to "phys_state_var_init",
"phys_output_open" and "regr_lat_time_climoz".

-- Created new variable "ozone_daylight" for "hist*.nc" output files.

-- Added a third dimension to variable "wo" in module
"phys_state_var_mod" and variable "POZON" in "radlwsw": index 1 for
average day-night ozone, index 2 for daylight ozone.

-- Added a fourth dimension to variables "o3_in", "o3_regr_lat" and
"o3_out" in "regr_lat_time_climoz": index 1 for average day-night
ozone, index 2 for daylight ozone.

-- In "physiq", moved call to "conf_phys" before call to
"phys_state_var_init". Thus, "conf_phys" is now inside the block "if
(first)" instead of "IF (debut)". There were definitions of "bl95_b0"
and "bl95_b1" that were useless because the variables were overwritten
by "conf_phys". Removed those definitions.

-- In "radlwsw", we pass the average day-night ozone to "LW_LMDAR4"
and the daylight ozone, if we have it, to "SW_LMDAR4" or
"SW_AEROAR4". If we do not have a specific field for daylight ozone
then "SW_LMDAR4" or "SW_AEROAR4" just get the average day-night ozone.

-- "regr_lat_time_climoz" now manages latitudes where the input ozone
field is missing at all levels (polar night).

-- Encapsulated "radlwsw" in a module.

3) Modifications to make sequential and parallel versions of
"create_etat0_limit" almost identical:

-- In "dyn3dpar/create_etat0_limit.F". No need to call
"phys_state_var_init", removed "use phys_state_var_mod" statement. No
need for "clesphys.h", removed "include" statement.

-- In "dyn3dpar/etat0_netcdf.F". Added argument "tau_ratqs" in call to
"conf_phys" (this bug was already corrected in "dyn3d"). Moved call to
"inifilr" after call to "infotrac_init" (as in "dyn3d").

4) Other peripheral modifications:

-- Added procedures "nf95_get_att" and "nf95_def_var_scalar" in
NetCDF95 interface. Overloaded "nf95_put_var" with three more
procedures: "nf95_put_var_FourByteReal", "nf95_put_var_FourByteInt",
"nf95_put_var_1D_FourByteInt".

-- Overloaded "regr1_step_av" with one more procedure:
"regr14_step_av". Overloaded "regr3_lint" with one more procedure:
"regr34_lint".

-- Corrected call to "Init_Phys_lmdz" in "dyn3d/create_etat0_limit.F":
the last argument should be an array, not a scalar.

-- Encapsulated "conf_phys" in a module.

-- Splitted module "regr_pr" into "regr_pr_av_m" and "regr_pr_int_m".

5) Tests:

This revision was compared to revision 1259, with optimization options
"debug" and "dev", parallelization options "none", "mpi", "omp" and
"mpi_omp", 1 and 2 MPI processes, 1 and 2 OpenMP threads, with the
compiler "FORTRAN90/SX Version 2.0 for SX-8". Both programs
"create_etat0_limit" and "gcm" were tested. In all cases,
parallelization does not change the results. With "read_climoz = 0" in
the ".def" files, the results of revision 1259 and of this revision
are the same.

Location:
LMDZ4/branches/LMDZ4-dev/libf/bibio
Files:
1 added
6 edited

Legend:

Unmodified
Added
Removed
  • LMDZ4/branches/LMDZ4-dev/libf/bibio/netcdf95.F90

    r1157 r1263  
    2525  ! criticisms for some (not all) procedures.
    2626
     27  ! "nf95_get_att" is more secure than "nf90_get_att" because it
     28  ! checks that the "values" argument is long enough and removes the
     29  ! null terminator, if any.
     30
    2731  ! This module replaces some of the official NetCDF procedures.
    2832  ! This module also provides the procedures "handle_err" and "nf95_gw_var".
     
    3539  use nf95_gw_var_m
    3640  use nf95_put_att_m
     41  use nf95_get_att_m
    3742  use simple
    3843  use handle_err_m
  • LMDZ4/branches/LMDZ4-dev/libf/bibio/nf95_def_var_m.F90

    r1157 r1263  
    11! $Id$
    22module nf95_def_var_m
     3
     4  ! The generic procedure name "nf90_def_var" applies to
     5  ! "nf90_def_var_Scalar" but we cannot apply the generic procedure name
     6  ! "nf95_def_var" to "nf95_def_var_scalar" because of the additional
     7  ! optional argument.
     8  ! "nf95_def_var_scalar" cannot be distinguished from "nf95_def_var_oneDim".
    39
    410  implicit none
     
    915
    1016  private
    11   public nf95_def_var
     17  public nf95_def_var, nf95_def_var_scalar
    1218
    1319contains
     20
     21  subroutine nf95_def_var_scalar(ncid, name, xtype, varid, ncerr)
     22
     23    use netcdf, only: nf90_def_var
     24    use handle_err_m, only: handle_err
     25
     26    integer,               intent( in) :: ncid
     27    character (len = *),   intent( in) :: name
     28    integer,               intent( in) :: xtype
     29    integer,               intent(out) :: varid
     30    integer, intent(out), optional:: ncerr
     31
     32    ! Variable local to the procedure:
     33    integer ncerr_not_opt
     34
     35    !-------------------
     36
     37    ncerr_not_opt = nf90_def_var(ncid, name, xtype, varid)
     38    if (present(ncerr)) then
     39       ncerr = ncerr_not_opt
     40    else
     41       call handle_err("nf95_def_var_scalar " // name, ncerr_not_opt, ncid)
     42    end if
     43
     44  end subroutine nf95_def_var_scalar
     45
     46  !***********************
    1447
    1548  subroutine nf95_def_var_oneDim(ncid, name, xtype, dimids, varid, ncerr)
  • LMDZ4/branches/LMDZ4-dev/libf/bibio/nf95_put_var_m.F90

    r1157 r1263  
    55
    66  interface nf95_put_var
    7      module procedure nf95_put_var_1D_FourByteReal, &
     7     module procedure nf95_put_var_FourByteReal, nf95_put_var_FourByteInt, &
     8          nf95_put_var_1D_FourByteReal, nf95_put_var_1D_FourByteInt, &
    89          nf95_put_var_2D_FourByteReal, nf95_put_var_3D_FourByteReal, &
    910          nf95_put_var_4D_FourByteReal
     
    1920contains
    2021
     22  subroutine nf95_put_var_FourByteReal(ncid, varid, values, start, ncerr)
     23
     24    use netcdf, only: nf90_put_var
     25    use handle_err_m, only: handle_err
     26
     27    integer, intent( in) :: ncid, varid
     28    real, intent( in) :: values
     29    integer, dimension(:), optional, intent( in) :: start
     30    integer, intent(out), optional:: ncerr
     31
     32    ! Variable local to the procedure:
     33    integer ncerr_not_opt
     34
     35    !-------------------
     36
     37    ncerr_not_opt = nf90_put_var(ncid, varid, values, start)
     38    if (present(ncerr)) then
     39       ncerr = ncerr_not_opt
     40    else
     41       call handle_err("nf95_put_var_FourByteReal", ncerr_not_opt, ncid, &
     42            varid)
     43    end if
     44
     45  end subroutine nf95_put_var_FourByteReal
     46
     47  !***********************
     48
     49  subroutine nf95_put_var_FourByteInt(ncid, varid, values, start, ncerr)
     50
     51    use netcdf, only: nf90_put_var
     52    use handle_err_m, only: handle_err
     53
     54    integer, intent( in) :: ncid, varid
     55    integer, intent( in) :: values
     56    integer, dimension(:), optional, intent( in) :: start
     57    integer, intent(out), optional:: ncerr
     58
     59    ! Variable local to the procedure:
     60    integer ncerr_not_opt
     61
     62    !-------------------
     63
     64    ncerr_not_opt = nf90_put_var(ncid, varid, values, start)
     65    if (present(ncerr)) then
     66       ncerr = ncerr_not_opt
     67    else
     68       call handle_err("nf95_put_var_FourByteInt", ncerr_not_opt, ncid, &
     69            varid)
     70    end if
     71
     72  end subroutine nf95_put_var_FourByteInt
     73
     74  !***********************
     75
    2176  subroutine nf95_put_var_1D_FourByteReal(ncid, varid, values, start, count, &
    2277       stride, map, ncerr)
     
    45100
    46101  end subroutine nf95_put_var_1D_FourByteReal
     102
     103  !***********************
     104
     105  subroutine nf95_put_var_1D_FourByteInt(ncid, varid, values, start, count, &
     106       stride, map, ncerr)
     107
     108    use netcdf, only: nf90_put_var
     109    use handle_err_m, only: handle_err
     110
     111    integer,                         intent(in) :: ncid, varid
     112    integer, intent(in) :: values(:)
     113    integer, dimension(:), optional, intent(in) :: start, count, stride, map
     114    integer, intent(out), optional:: ncerr
     115
     116    ! Variable local to the procedure:
     117    integer ncerr_not_opt
     118
     119    !-------------------
     120
     121    ncerr_not_opt = nf90_put_var(ncid, varid, values, start, count, stride, &
     122         map)
     123    if (present(ncerr)) then
     124       ncerr = ncerr_not_opt
     125    else
     126       call handle_err("nf95_put_var_1D_FourByteInt", ncerr_not_opt, ncid, &
     127            varid)
     128    end if
     129
     130  end subroutine nf95_put_var_1D_FourByteInt
    47131
    48132  !***********************
  • LMDZ4/branches/LMDZ4-dev/libf/bibio/regr1_step_av_m.F90

    r1157 r1263  
    1717     ! The difference between the procedures is the rank of the first argument.
    1818
    19      module procedure regr11_step_av, regr12_step_av, regr13_step_av
     19     module procedure regr11_step_av, regr12_step_av, regr13_step_av, &
     20          regr14_step_av
    2021  end interface
    2122
     
    203204  end function regr13_step_av
    204205
     206  !********************************************
     207
     208  function regr14_step_av(vs, xs, xt) result(vt)
     209
     210    ! "vs" has rank 4.
     211
     212    use assert_eq_m, only: assert_eq
     213    use assert_m, only: assert
     214    use interpolation, only: locate
     215
     216    real, intent(in):: vs(:, :, :, :) ! values of steps on the source grid
     217    ! (Step "is" is between "xs(is)" and "xs(is + 1)".)
     218
     219    real, intent(in):: xs(:)
     220    ! (edges of steps on the source grid, in strictly increasing order)
     221
     222    real, intent(in):: xt(:)
     223    ! (edges of cells of the target grid, in strictly increasing order)
     224
     225    real vt(size(xt) - 1, size(vs, 2), size(vs, 3), size(vs, 4))
     226    ! (average values on the target grid)
     227    ! (Cell "it" is between "xt(it)" and "xt(it + 1)".)
     228
     229    ! Variables local to the procedure:
     230    integer is, it, ns, nt
     231    real left_edge
     232
     233    !---------------------------------------------
     234
     235    ns = assert_eq(size(vs, 1), size(xs) - 1, "regr14_step_av ns")
     236    nt = size(xt) - 1
     237
     238    ! Quick check on sort order:
     239    call assert(xs(1) < xs(2), "regr14_step_av xs bad order")
     240    call assert(xt(1) < xt(2), "regr14_step_av xt bad order")
     241
     242    call assert(xs(1) <= xt(1) .and. xt(nt + 1) <= xs(ns + 1), &
     243         "regr14_step_av extrapolation")
     244
     245    is = locate(xs, xt(1)) ! 1 <= is <= ns, because we forbid extrapolation
     246    do it = 1, nt
     247       ! 1 <= is <= ns
     248       ! xs(is) <= xt(it) < xs(is + 1)
     249       ! Compute "vt(it, :, :, :)":
     250       left_edge = xt(it)
     251       vt(it, :, :, :) = 0.
     252       do while (xs(is + 1) < xt(it + 1))
     253          ! 1 <= is <= ns - 1
     254          vt(it, :, :, :) = vt(it, :, :, :) + (xs(is + 1) - left_edge) &
     255               * vs(is, :, :, :)
     256          is = is + 1
     257          left_edge = xs(is)
     258       end do
     259       ! 1 <= is <= ns
     260       vt(it, :, :, :) = (vt(it, :, :, :) + (xt(it + 1) - left_edge) &
     261            * vs(is, :, :, :)) / (xt(it + 1) - xt(it))
     262       if (xs(is + 1) == xt(it + 1)) is = is + 1
     263       ! 1 <= is <= ns .or. it == nt
     264    end do
     265
     266  end function regr14_step_av
     267
    205268end module regr1_step_av_m
  • LMDZ4/branches/LMDZ4-dev/libf/bibio/regr3_lint_m.F90

    r1157 r1263  
    1111     ! input array.
    1212     ! The difference betwwen the procedures is the rank of the first argument.
    13      module procedure regr33_lint
     13     module procedure regr33_lint, regr34_lint
    1414  end interface
    1515
     
    5757  end function regr33_lint
    5858
     59  !*********************************************************
     60
     61  function regr34_lint(vs, xs, xt) result(vt)
     62
     63    ! "vs" has rank 4.
     64
     65    use assert_eq_m, only: assert_eq
     66    use interpolation, only: hunt
     67
     68    real, intent(in):: vs(:, :, :, :)
     69    ! (values of the function at source points "xs")
     70
     71    real, intent(in):: xs(:)
     72    ! (abscissas of points in source grid, in strictly monotonic order)
     73
     74    real, intent(in):: xt(:)
     75    ! (abscissas of points in target grid)
     76
     77    real vt(size(vs, 1), size(vs, 2), size(xt), size(vs, 4))
     78    ! (values of the function on the target grid)
     79
     80    ! Variables local to the procedure:
     81    integer is, it, ns
     82    integer is_b ! "is" bound between 1 and "ns - 1"
     83
     84    !--------------------------------------
     85
     86    ns = assert_eq(size(vs, 3), size(xs), "regr34_lint ns")
     87
     88    is = -1 ! go immediately to bisection on first call to "hunt"
     89
     90    do it = 1, size(xt)
     91       call hunt(xs, xt(it), is)
     92       is_b = min(max(is, 1), ns - 1)
     93       vt(:, :, it, :) = ((xs(is_b+1) - xt(it)) * vs(:, :, is_b, :) &
     94            + (xt(it) - xs(is_b)) * vs(:, :, is_b+1, :)) &
     95            / (xs(is_b+1) - xs(is_b))
     96    end do
     97
     98  end function regr34_lint
     99
    59100end module regr3_lint_m
  • LMDZ4/branches/LMDZ4-dev/libf/bibio/simple.F90

    r1157 r1263  
    118118    ! unknown in the calling procedure, before the call.
    119119    ! Here we use a better solution: a pointer argument array.
    120     ! This procedure associates and defines '"dimids" if it is present.
     120    ! This procedure associates and defines "dimids" if it is present.
    121121
    122122    use netcdf, only: nf90_inquire_variable, nf90_max_var_dims
Note: See TracChangeset for help on using the changeset viewer.