/[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 168 - (show annotations)
Wed Sep 9 10:41:47 2015 UTC (8 years, 9 months ago) by guez
Original Path: trunk/Sources/phylmd/Mobidic/regr_pr_int.f
File size: 3058 byte(s)
In order to be able to choose finer resolutions, set large memory
model in compiler options and use dynamic libraries.

Variables rlatd, rlond, cuphy and cvphy of module comgeomphy were
never used. (In LMDZ, they are used only for Orchid.)

There is a bug in PGI Fortran 13.10 that does not accept the
combination of forall, pack and spread in regr_pr_av and
regr_pr_int. In order to circumvent this bug, created the function
gr_dyn_phy.

In program test_inifilr, use a single latitude coordinate for north
and south.

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.
23 ! The input data does not depend on longitude, but the pressure
24 ! at LMDZ mid-layers does.
25 ! Therefore, the values on the LMDZ grid do depend on longitude.
26 ! Regridding is by linear interpolation.
27
28 use dimens_m, only: iim, jjm, llm
29 use dimphy, only: klon
30 use grid_change, only: gr_dyn_phy
31 use netcdf95, only: nf95_inq_varid, nf95_get_var
32 use nr_util, only: assert
33 use numer_rec_95, only: regr1_lint
34 use press_coefoz_m, only: plev
35
36 integer, intent(in):: ncid ! NetCDF ID of the file
37 character(len=*), intent(in):: name ! of the NetCDF variable
38 integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
39
40 real, intent(in):: pplay(:, :) ! (klon, llm)
41 ! (pression pour le mileu de chaque couche, en Pa)
42
43 real, intent(in):: top_value
44 ! (extra value of field at 0 pressure)
45
46 real, intent(out):: v3(:, :) ! (klon, llm)
47 ! Regridded field on the partial "physics" grid. "v3(i, k)" is at
48 ! longitude "xlon(i)", latitude "xlat(i)", middle of layer "k".
49
50 ! Variables local to the procedure:
51
52 integer varid ! for NetCDF
53
54 real v1(jjm + 1, 0:size(plev))
55 ! Input field at day "julien". "v1(j, k >=1)" is at latitude
56 ! "rlatu(j)" and pressure "plev(k)".
57
58 real v2(klon, 0:size(plev))
59 ! Field on the "physics" horizontal grid. "v2(i, k >= 1)" is at
60 ! longitude "xlon(i)", latitude "xlat(i)" and pressure "plev(k)".)
61
62 integer i
63
64 !--------------------------------------------
65
66 call assert(shape(v3) == (/klon, llm/), "regr_pr_int")
67
68 call nf95_inq_varid(ncid, name, varid)
69
70 ! Get data at the right day from the input file:
71 call nf95_get_var(ncid, varid, v1(:, 1:), start=(/1, 1, julien/))
72 ! Latitudes are in ascending order in the input file while
73 ! "rlatu" is in descending order so we need to invert order:
74 v1(:, 1:) = v1(jjm+1:1:-1, 1:)
75
76 ! Complete "v1" with the value at 0 pressure:
77 v1(:, 0) = top_value
78
79 v2 = gr_dyn_phy(spread(v1, dim = 1, ncopies = iim + 1))
80
81 ! Regrid in pressure at each horizontal position:
82 do i = 1, klon
83 v3(i, llm:1:-1) = regr1_lint(v2(i, :), (/0., plev/), pplay(i, llm:1:-1))
84 ! (invert order of indices because "pplay" is in descending order)
85 end do
86
87 end subroutine regr_pr_int
88
89 end module regr_pr_int_m

  ViewVC Help
Powered by ViewVC 1.1.21