source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phylmd/regr_pr_o3_m.F90 @ 5308

Last change on this file since 5308 was 2440, checked in by lguez, 9 years ago

For read_climoz = 1 or 2, replaced first order conservative regridding
of ozone by second order conservative regridding, with Van Leer
slope-limiting. The replacement is done for both latitude and pressure
regridding. The replacement is beneficial if the resolution of the
input data is coarser than the resolution of LMDZ. If the resolution
of the input data is finer, then the replacement is neutral, it does
not change much.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 3.8 KB
Line 
1! $Id$
2module regr_pr_o3_m
3
4  implicit none
5
6contains
7
8  subroutine regr_pr_o3(p3d, o3_mob_regr)
9
10    ! "regr_pr_o3" stands for "regrid pressure ozone".
11    ! This procedure reads Mobidic ozone mole fraction from
12    ! "coefoz_LMDZ.nc" at the initial day of the run and regrids it in
13    ! pressure.
14    ! Ozone mole fraction from "coefoz_LMDZ.nc" at the initial day is
15    ! a 2D latitude -- pressure variable.
16    ! The target horizontal LMDZ grid is the "scalar" grid: "rlonv", "rlatu".
17    ! The target vertical LMDZ grid is the grid of layer boundaries.
18    ! We assume that the input variable is already on the LMDZ "rlatu"
19    ! latitude grid.
20    ! The input variable does not depend on longitude, but the
21    ! pressure at LMDZ layers does.
22    ! Therefore, the values on the LMDZ grid do depend on longitude.
23    ! Regridding is by averaging, assuming a step function.
24    ! We assume that, in the input file, the pressure levels are in
25    ! hPa and strictly increasing.
26
27    use netcdf95, only: nf95_open, nf95_close, nf95_inq_varid, handle_err
28    use netcdf, only:  nf90_nowrite, nf90_get_var
29    use assert_m, only: assert
30    use regr1_conserv_m, only: regr1_conserv
31    use press_coefoz_m, only: press_in_edg
32    use time_phylmdz_mod, only: day_ref
33    use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev
34
35    REAL, intent(in):: p3d(:, :, :) ! pressure at layer interfaces, in Pa
36    ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
37    ! for interface "l")
38
39    real, intent(out):: o3_mob_regr(:, :, :) ! (iim + 1, jjm + 1, llm)
40    ! (ozone mole fraction from Mobidic adapted to the LMDZ grid)
41    ! ("o3_mob_regr(i, j, l)" is at longitude "rlonv(i)", latitude
42    ! "rlatu(j)" and pressure level "pls(i, j, l)")
43
44    ! Variables local to the procedure:
45
46    integer ncid, varid, ncerr ! for NetCDF
47    integer i, j
48
49    real r_mob(nbp_lat, size(press_in_edg) - 1)
50    ! (ozone mole fraction from Mobidic at day "day_ref")
51    ! (r_mob(j, k) is at latitude "rlatu(j)", in pressure interval
52    ! "[press_in_edg(k), press_in_edg(k+1)]".)
53
54    !------------------------------------------------------------
55
56    print *, "Call sequence information: regr_pr_o3"
57    call assert(shape(o3_mob_regr) == (/nbp_lon + 1, nbp_lat, nbp_lev/), &
58         "regr_pr_o3 o3_mob_regr")
59    call assert(shape(p3d) == (/nbp_lon + 1, nbp_lat, nbp_lev + 1/), "regr_pr_o3 p3d")
60
61    call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
62
63    call nf95_inq_varid(ncid, "r_Mob", varid)
64    ! Get data at the right day from the input file:
65    ncerr = nf90_get_var(ncid, varid, r_mob, start=(/1, 1, day_ref/))
66    call handle_err("nf90_get_var r_Mob", ncerr)
67    ! Latitudes are in ascending order in the input file while
68    ! "rlatu" is in descending order so we need to invert order:
69    r_mob = r_mob(nbp_lat:1:-1, :)
70
71    call nf95_close(ncid)
72
73    ! Regrid in pressure by averaging a step function of pressure:
74
75    ! Poles:
76    do j = 1, nbp_lat, nbp_lat-1
77       call regr1_conserv(r_mob(j, :), press_in_edg, &
78            p3d(1, j, nbp_lev + 1:1:-1), o3_mob_regr(1, j, nbp_lev:1:-1))
79       ! (invert order of indices because "p3d" is in descending order)
80    end do
81
82    ! Other latitudes:
83    do j = 2, nbp_lat-1
84       do i = 1, nbp_lon
85          call regr1_conserv(r_mob(j, :), press_in_edg, &
86               p3d(i, j, nbp_lev + 1:1:-1), o3_mob_regr(i, j, nbp_lev:1:-1))
87          ! (invert order of indices because "p3d" is in descending order)
88       end do
89    end do
90
91    ! Duplicate pole values on all longitudes:
92    o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies=nbp_lon)
93    o3_mob_regr(2:, nbp_lat, :) &
94         = spread(o3_mob_regr(1, nbp_lat, :), dim=1, ncopies=nbp_lon)
95
96    ! Duplicate first longitude to last longitude:
97    o3_mob_regr(nbp_lon + 1, 2:nbp_lat-1, :) = o3_mob_regr(1, 2:nbp_lat-1, :)
98
99  end subroutine regr_pr_o3
100
101end module regr_pr_o3_m
Note: See TracBrowser for help on using the repository browser.