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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 162 - (show annotations)
Fri Jul 24 16:54:30 2015 UTC (8 years, 10 months ago) by guez
File size: 3085 byte(s)
Variable pls of module pressure_var was only modified by calfis but I
could not move it to module calfis_m because it was used by a
procedure called by calfis (it would have been a cyclic
dependency). In the same way, variable p3d of module pressure_var was
only modified by leapfrog but I could not move it to module
leapfrog_m. So removed module pressure_var. p3d becomes a local
variable of leapfrog and an argument of calfis. Use paprs and play in
regr_pr_int and regr_pr_av (following LMDZ). The idea in regr_pr_int
and regr_pr_av is to spread and pack before the regridding instead of
packing afterward. The cost in memory should only be a two-dimensional
temporary array created by spread. The cost in clarity is only the
transiting of paprs and pplay through regr_pr_comb_coefoz, but this is
more than compensated by removing the side effect on module variables.

1 module regr_pr_int_m
2
3 implicit none
4
5 contains
6
7 subroutine regr_pr_int(ncid, name, julien, pplay, 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
34 integer, intent(in):: ncid ! NetCDF ID of the file
35 character(len=*), intent(in):: name ! of the NetCDF variable
36 integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
37
38 real, intent(in):: pplay(:, :) ! (klon, llm)
39 ! (pression pour le mileu de chaque couche, en Pa)
40
41 real, intent(in):: top_value
42 ! (extra value of field at 0 pressure)
43
44 real, intent(out):: v3(:, :) ! (klon, llm)
45 ! Regridded field on the partial "physics" grid. "v3(i, k)" is at
46 ! longitude "xlon(i)", latitude "xlat(i)", middle of layer "k".
47
48 ! Variables local to the procedure:
49
50 integer varid ! for NetCDF
51
52 real v1(jjm + 1, 0:size(plev))
53 ! Input field at day "julien". "v1(j, k >=1)" is at latitude
54 ! "rlatu(j)" and pressure "plev(k)".
55
56 real v2(klon, 0:size(plev))
57 ! Field on the "physics" horizontal grid. "v2(i, k >= 1)" is at
58 ! longitude "xlon(i)", latitude "xlat(i)" and pressure "plev(k)".)
59
60 integer i, k
61
62 !--------------------------------------------
63
64 call assert(shape(v3) == (/klon, llm/), "regr_pr_int")
65
66 call nf95_inq_varid(ncid, name, varid)
67
68 ! Get data at the right day from the input file:
69 call nf95_get_var(ncid, varid, v1(:, 1:), start=(/1, 1, julien/))
70 ! Latitudes are in ascending order in the input file while
71 ! "rlatu" is in descending order so we need to invert order:
72 v1(:, 1:) = v1(jjm+1:1:-1, 1:)
73
74 ! Complete "v1" with the value at 0 pressure:
75 v1(:, 0) = top_value
76
77 forall (k = 0:size(plev)) v2(:, k) = pack(spread(v1(:, k), dim = 1, &
78 ncopies = iim + 1), dyn_phy)
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