source: LMDZ6/trunk/libf/phylmd/regr_pr_o3_m.F90 @ 5075

Last change on this file since 5075 was 5075, checked in by abarral, 4 months ago

[continued & end] replace netcdf by lmdz_netcdf.F90 wrapper
"use netcdf" is now only used in lmdz_netcdf.F90 (except ecrad and obsolete/)
<include "netcdf.inc"> is now likewise only used in lmdz_netcdf.F90.

systematically specify explicitely <USE lmdz_netcdf, ONLY:> (probably left some missing, to correct later on)

Further replacement of nf_put_* by nf90_put_* (same for _get_)

[minor] replace deprecated boolean operators along the way

  • 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
RevLine 
[1379]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
[4489]27    use netcdf95, only: nf95_open, nf95_close, nf95_inq_varid, nf95_get_var
[5075]28    use lmdz_netcdf, only:  nf90_nowrite
[1379]29    use assert_m, only: assert
[2788]30    use regr_conserv_m, only: regr_conserv
[1379]31    use press_coefoz_m, only: press_in_edg
[2345]32    use time_phylmdz_mod, only: day_ref
[2346]33    use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev
[1379]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
[2346]49    real r_mob(nbp_lat, size(press_in_edg) - 1)
[2345]50    ! (ozone mole fraction from Mobidic at day "day_ref")
[1379]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"
[2346]57    call assert(shape(o3_mob_regr) == (/nbp_lon + 1, nbp_lat, nbp_lev/), &
[1379]58         "regr_pr_o3 o3_mob_regr")
[2346]59    call assert(shape(p3d) == (/nbp_lon + 1, nbp_lat, nbp_lev + 1/), "regr_pr_o3 p3d")
[1379]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:
[4489]65    call nf95_get_var(ncid, varid, r_mob, start=(/1, 1, day_ref/))
[1379]66    ! Latitudes are in ascending order in the input file while
67    ! "rlatu" is in descending order so we need to invert order:
[2346]68    r_mob = r_mob(nbp_lat:1:-1, :)
[1379]69
70    call nf95_close(ncid)
71
72    ! Regrid in pressure by averaging a step function of pressure:
73
74    ! Poles:
[2346]75    do j = 1, nbp_lat, nbp_lat-1
[2788]76       call regr_conserv(1, r_mob(j, :), press_in_edg, &
[2440]77            p3d(1, j, nbp_lev + 1:1:-1), o3_mob_regr(1, j, nbp_lev:1:-1))
[1379]78       ! (invert order of indices because "p3d" is in descending order)
79    end do
80
81    ! Other latitudes:
[2346]82    do j = 2, nbp_lat-1
83       do i = 1, nbp_lon
[2788]84          call regr_conserv(1, r_mob(j, :), press_in_edg, &
[2440]85               p3d(i, j, nbp_lev + 1:1:-1), o3_mob_regr(i, j, nbp_lev:1:-1))
86          ! (invert order of indices because "p3d" is in descending order)
[1379]87       end do
88    end do
89
90    ! Duplicate pole values on all longitudes:
[2346]91    o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies=nbp_lon)
92    o3_mob_regr(2:, nbp_lat, :) &
93         = spread(o3_mob_regr(1, nbp_lat, :), dim=1, ncopies=nbp_lon)
[1379]94
95    ! Duplicate first longitude to last longitude:
[2346]96    o3_mob_regr(nbp_lon + 1, 2:nbp_lat-1, :) = o3_mob_regr(1, 2:nbp_lat-1, :)
[1379]97
98  end subroutine regr_pr_o3
99
100end module regr_pr_o3_m
Note: See TracBrowser for help on using the repository browser.