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

  ViewVC Help
Powered by ViewVC 1.1.21