1 | ! $Id$ |
---|
2 | module regr_pr_int_m |
---|
3 | |
---|
4 | ! Author: Lionel GUEZ |
---|
5 | |
---|
6 | IMPLICIT NONE |
---|
7 | |
---|
8 | CONTAINS |
---|
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 lmdz_assert, ONLY: assert |
---|
29 | USE lmdz_regr_lint, ONLY: regr_lint |
---|
30 | USE lmdz_phys_mpi_data, ONLY: is_mpi_root |
---|
31 | USE lmdz_grid_phy, ONLY: nbp_lon, nbp_lat, nbp_lev |
---|
32 | USE lmdz_phys_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 | |
---|
104 | END MODULE regr_pr_int_m |
---|