source: LMDZ6/trunk/libf/phylmd/regr_pr_int_m.F90 @ 5249

Last change on this file since 5249 was 4489, checked in by lguez, 20 months ago

Merge LMDZ_ECRad branch back into trunk!

  • 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.6 KB
Line 
1! $Id$
2module regr_pr_int_m
3
4  ! Author: Lionel GUEZ
5
6  implicit none
7
8contains
9
10  subroutine regr_pr_int(ncid, name, julien, plev, pplay, top_value, v3)
11
12    ! "regr_pr_int" stands for "regrid pressure interpolation".
13    ! In this procedure:
14    ! -- the root process reads a 2D latitude-pressure field from a
15    !    NetCDF file, at a given day.
16    ! -- the field is packed to the LMDZ horizontal "physics"
17    !    grid and scattered to all threads of all processes;
18    ! -- in all the threads of all the processes, the field is regridded in
19    !    pressure to the LMDZ vertical grid.
20    ! We assume that, in the input file, the field has 3 dimensions:
21    ! latitude, pressure, julian day.
22    ! We assume that latitudes are in ascending order in the input file.
23    ! The target vertical LMDZ grid is the grid of mid-layers.
24    ! Regridding is by linear interpolation.
25
26    use dimphy, only: klon
27    use netcdf95, only: nf95_inq_varid, nf95_get_var
28    use assert_m, only: assert
29    use regr_lint_m, only: regr_lint
30    use mod_phys_lmdz_mpi_data, only: is_mpi_root
31    use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev
32    use mod_phys_lmdz_transfert_para, only: scatter2d
33    ! (pack to the LMDZ horizontal "physics" grid and scatter)
34
35    integer, intent(in):: ncid ! NetCDF ID of the file
36    character(len=*), intent(in):: name ! of the NetCDF variable
37    integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
38
39    real, intent(in):: plev(:)
40    ! (pressure level of input data, in Pa, in strictly ascending order)
41
42    real, intent(in):: pplay(:, :) ! (klon, nbp_lev)
43    ! (pression pour le mileu de chaque couche, en Pa)
44
45    real, intent(in):: top_value
46    ! (extra value of field at 0 pressure)
47
48    real, intent(out):: v3(:, :) ! (klon, nbp_lev)
49    ! (regridded field on the partial "physics" grid)
50    ! ("v3(i, k)" is at longitude "xlon(i)", latitude
51    ! "xlat(i)", middle of layer "k".)
52
53    ! Variables local to the procedure:
54
55    integer varid, ncerr ! for NetCDF
56
57    real  v1(nbp_lon, nbp_lat, 0:size(plev))
58    ! (input field at day "julien", on the global "dynamics" horizontal grid)
59    ! (First dimension is for longitude.
60    ! The value is the same for all longitudes.
61    ! "v1(:, j, k >=1)" is at latitude "rlatu(j)" and pressure "plev(k)".)
62
63    real v2(klon, 0:size(plev))
64    ! (field scattered to the partial "physics" horizontal grid)
65    ! "v2(i, k >= 1)" is at longitude "xlon(i)", latitude "xlat(i)"
66    ! and pressure "plev(k)".)
67
68    integer i
69
70    !--------------------------------------------
71
72    call assert(shape(v3) == (/klon, nbp_lev/), "regr_pr_int v3")
73    call assert(shape(pplay) == (/klon, nbp_lev/), "regr_pr_int pplay")
74
75    !$omp master
76    if (is_mpi_root) then
77       call nf95_inq_varid(ncid, name, varid)
78
79       ! Get data at the right day from the input file:
80       call nf95_get_var(ncid, varid, v1(1, :, 1:), start=(/1, 1, julien/))
81       ! Latitudes are in ascending order in the input file while
82       ! "rlatu" is in descending order so we need to invert order:
83       v1(1, :, 1:) = v1(1, nbp_lat:1:-1, 1:)
84
85       ! Complete "v1" with the value at 0 pressure:
86       v1(1, :, 0) = top_value
87
88       ! Duplicate on all longitudes:
89       v1(2:, :, :) = spread(v1(1, :, :), dim=1, ncopies=nbp_lon-1)
90    end if
91    !$omp end master
92
93    call scatter2d(v1, v2)
94
95    ! Regrid in pressure at each horizontal position:
96    do i = 1, klon
97       call regr_lint(1, v2(i, :), (/0., plev/), pplay(i, nbp_lev:1:-1), &
98                         v3(i, nbp_lev:1:-1))
99       ! (invert order of indices because "pplay" is in descending order)
100    end do
101
102  end subroutine regr_pr_int
103
104end module regr_pr_int_m
Note: See TracBrowser for help on using the repository browser.