/[lmdze]/trunk/phylmd/Mobidic/regr_pr_int.f
ViewVC logotype

Contents of /trunk/phylmd/Mobidic/regr_pr_int.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 265 - (show annotations)
Tue Mar 20 09:35:59 2018 UTC (6 years, 2 months ago) by guez
File size: 3054 byte(s)
Rename module dimens_m to dimensions.
1 module regr_pr_int_m
2
3 ! Author: Lionel GUEZ
4
5 implicit none
6
7 contains
8
9 subroutine regr_pr_int(ncid, name, julien, pplay, top_value, v3)
10
11 ! "regr_pr_int" stands for "regrid pressure interpolation".
12
13 ! This procedure reads a 2D latitude-pressure field from a NetCDF
14 ! file, at a given day, regrids the field in pressure to the LMDZ
15 ! vertical grid and packs it to the LMDZ horizontal "physics"
16 ! grid. We assume that, in the input file, the field has 3
17 ! dimensions: latitude, pressure, julian day. We assume that
18 ! latitudes are in ascending order in the input file. We assume
19 ! that the input data is already on the LMDZ "rlatu" latitude
20 ! grid.
21
22 ! The target vertical LMDZ grid is the grid of mid-layers. The
23 ! input data does not depend on longitude, but the pressure at
24 ! LMDZ mid-layers does. Therefore, the values on the LMDZ grid do
25 ! depend on longitude. Regridding is by linear interpolation.
26
27 use dimensions, only: iim, jjm, llm
28 use dimphy, only: klon
29 use grid_change, only: gr_dyn_phy
30 use netcdf95, only: nf95_inq_varid, nf95_get_var
31 use nr_util, only: assert
32 use numer_rec_95, only: regr1_lint
33 use press_coefoz_m, only: plev
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):: pplay(:, :) ! (klon, llm)
40 ! (pression pour le mileu de chaque couche, en Pa)
41
42 real, intent(in):: top_value
43 ! (extra value of field at 0 pressure)
44
45 real, intent(out):: v3(:, :) ! (klon, llm)
46 ! Regridded field on the partial "physics" grid. "v3(i, k)" is at
47 ! longitude "xlon(i)", latitude "xlat(i)", middle of layer "k".
48
49 ! Variables local to the procedure:
50
51 integer varid ! for NetCDF
52
53 real v1(jjm + 1, 0:size(plev))
54 ! Input field at day "julien". "v1(j, k >=1)" is at latitude
55 ! "rlatu(j)" and pressure "plev(k)".
56
57 real v2(klon, 0:size(plev))
58 ! Field on the "physics" horizontal grid. "v2(i, k >= 1)" is at
59 ! longitude "xlon(i)", latitude "xlat(i)" and pressure "plev(k)".)
60
61 integer i
62
63 !--------------------------------------------
64
65 call assert(shape(v3) == (/klon, llm/), "regr_pr_int")
66
67 call nf95_inq_varid(ncid, name, varid)
68
69 ! Get data at the right day from the input file:
70 call nf95_get_var(ncid, varid, v1(:, 1:), start=(/1, 1, julien/))
71 ! Latitudes are in ascending order in the input file while
72 ! "rlatu" is in descending order so we need to invert order:
73 v1(:, 1:) = v1(jjm+1:1:-1, 1:)
74
75 ! Complete "v1" with the value at 0 pressure:
76 v1(:, 0) = top_value
77
78 v2 = gr_dyn_phy(spread(v1, dim = 1, ncopies = iim + 1))
79
80 ! Regrid in pressure at each horizontal position:
81 do i = 1, klon
82 v3(i, llm:1:-1) = regr1_lint(v2(i, :), (/0., plev/), pplay(i, llm:1:-1))
83 ! (invert order of indices because "pplay" is in descending order)
84 end do
85
86 end subroutine regr_pr_int
87
88 end module regr_pr_int_m

  ViewVC Help
Powered by ViewVC 1.1.21