source: LMDZ4/branches/LMDZ4-dev/libf/phylmd/regr_pr.F90 @ 1235

Last change on this file since 1235 was 1225, checked in by lguez, 15 years ago

Added some "intent" attributes in declarations.

In "phyredem", "dtime" is not declared. It is not in any included
file. Probably accepted by compilers as an intrinsic non-standard
function. Removed this element of "tab_cntrl".

Added some "only" clauses in "use" statements.

If the ozone field is read from a file, it is now updated every
360th of the length of the current year, regardless of that length.

In "physiq", "omega" was output before it was defined. Moved the
output instruction after the definition.

File size: 6.7 KB
RevLine 
[1158]1! $Id$
2module regr_pr
3
[1220]4  ! Author: Lionel GUEZ
5
[1158]6  ! In both procedures of this module:
7  ! -- the root process reads a 2D latitude-pressure field from a
8  !    NetCDF file, at a given day.
9  ! -- the field is packed to the LMDZ horizontal "physics"
10  !    grid and scattered to all threads of all processes;
11  ! -- in all the threads of all the processes, the field is regridded in
12  !    pressure to the LMDZ vertical grid.
13  ! We assume that, in the input file, the field has 3 dimensions:
14  ! latitude, pressure, julian day.
15  ! We assume that latitudes are in ascending order in the input file.
16
17
18  implicit none
19
20contains
21
22  subroutine regr_pr_av(ncid, name, julien, press_in_edg, paprs, v3)
23
24    ! "regr_pr_av" stands for "regrid pressure averaging".
25    ! The target vertical LMDZ grid is the grid of layer boundaries.
26    ! Regridding in pressure is done by averaging a step function of pressure.
27
28    use dimphy, only: klon
29    use netcdf95, only: nf95_inq_varid, handle_err
30    use netcdf, only: nf90_get_var
31    use assert_m, only: assert
32    use regr1_step_av_m, only: regr1_step_av
33    use mod_phys_lmdz_mpi_data, only: is_mpi_root
34
35    use mod_phys_lmdz_transfert_para, only: scatter2d
36    ! (pack to the LMDZ horizontal "physics" grid and scatter)
37
38    integer, intent(in):: ncid ! NetCDF ID of the file
39    character(len=*), intent(in):: name ! of the NetCDF variable
40    integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
41
42    real, intent(in):: press_in_edg(:)
[1225]43    ! edges of pressure intervals for input data, in Pa, in strictly
44    ! ascending order
[1158]45
46    real, intent(in):: paprs(:, :) ! (klon, llm + 1)
47    ! (pression pour chaque inter-couche, en Pa)
48
49    real, intent(out):: v3(:, :) ! (klon, llm)
50    ! (regridded field on the partial "physics" grid)
51    ! ("v3(i, k)" is at longitude "xlon(i)", latitude
52    ! "xlat(i)", in pressure interval "[paprs(i, l+1), paprs(i, l)]".)
53
54    ! Variables local to the procedure:
55
56    include "dimensions.h"
57    integer varid, ncerr ! for NetCDF
58
59    real  v1(iim, jjm + 1, size(press_in_edg) - 1)
60    ! (input field at day "julien", on the global "dynamics" horizontal grid)
61    ! (First dimension is for longitude.
62    ! The value is the same for all longitudes.
63    ! "v1(:, j, k)" is at latitude "rlatu(j)" and for
64    ! pressure interval "[press_in_edg(k), press_in_edg(k+1)]".)
65
66    real v2(klon, size(press_in_edg) - 1)
67    ! (field scattered to the partial "physics" horizontal grid)
68    ! ("v2(i, k)" is at longitude "xlon(i)", latitude "xlat(i)" and
69    ! for pressure interval "[press_in_edg(k), press_in_edg(k+1)]".)
70
71    integer i
72
73    !--------------------------------------------
74
75    call assert(shape(v3) == (/klon, llm/), "regr_pr_av v3")
76    call assert(shape(paprs) == (/klon, llm+1/), "regr_pr_av paprs")
77
78    !$omp master
79    if (is_mpi_root) then
80       call nf95_inq_varid(ncid, name, varid)
81
82       ! Get data at the right day from the input file:
83       ncerr = nf90_get_var(ncid, varid, v1(1, :, :), start=(/1, 1, julien/))
84       call handle_err("regr_pr_av nf90_get_var " // name, ncerr, ncid)
85       ! Latitudes are in ascending order in the input file while
86       ! "rlatu" is in descending order so we need to invert order:
87       v1(1, :, :) = v1(1, jjm+1:1:-1, :)
88
89       ! Duplicate on all longitudes:
90       v1(2:, :, :) = spread(v1(1, :, :), dim=1, ncopies=iim-1)
91    end if
92    !$omp end master
93
94    call scatter2d(v1, v2)
95
96    ! Regrid in pressure at each horizontal position:
97    do i = 1, klon
98       v3(i, llm:1:-1) = regr1_step_av(v2(i, :), press_in_edg, &
99            paprs(i, llm+1:1:-1))
100       ! (invert order of indices because "paprs" is in descending order)
101    end do
102
103  end subroutine regr_pr_av
104
105  !***************************************************************
106
107  subroutine regr_pr_int(ncid, name, julien, plev, pplay, top_value, v3)
108
109    ! "regr_pr_int" stands for "regrid pressure interpolation".
110    ! The target vertical LMDZ grid is the grid of mid-layers.
111    ! Regridding is by linear interpolation.
112
113    use dimphy, only: klon
114    use netcdf95, only: nf95_inq_varid, handle_err
115    use netcdf, only: nf90_get_var
116    use assert_m, only: assert
117    use regr1_lint_m, only: regr1_lint
118    use mod_phys_lmdz_mpi_data, only: is_mpi_root
119
120    use mod_phys_lmdz_transfert_para, only: scatter2d
121    ! (pack to the LMDZ horizontal "physics" grid and scatter)
122
123    integer, intent(in):: ncid ! NetCDF ID of the file
124    character(len=*), intent(in):: name ! of the NetCDF variable
125    integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
126
127    real, intent(in):: plev(:)
128    ! (pressure level of input data, in Pa, in strictly ascending order)
129
130    real, intent(in):: pplay(:, :) ! (klon, llm)
131    ! (pression pour le mileu de chaque couche, en Pa)
132
133    real, intent(in):: top_value
134    ! (extra value of field at 0 pressure)
135
136    real, intent(out):: v3(:, :) ! (klon, llm)
137    ! (regridded field on the partial "physics" grid)
138    ! ("v3(i, k)" is at longitude "xlon(i)", latitude
139    ! "xlat(i)", middle of layer "k".)
140
141    ! Variables local to the procedure:
142
143    include "dimensions.h"
144    integer varid, ncerr ! for NetCDF
145
146    real  v1(iim, jjm + 1, 0:size(plev))
147    ! (input field at day "julien", on the global "dynamics" horizontal grid)
148    ! (First dimension is for longitude.
149    ! The value is the same for all longitudes.
150    ! "v1(:, j, k >=1)" is at latitude "rlatu(j)" and pressure "plev(k)".)
151
152    real v2(klon, 0:size(plev))
153    ! (field scattered to the partial "physics" horizontal grid)
154    ! "v2(i, k >= 1)" is at longitude "xlon(i)", latitude "xlat(i)"
155    ! and pressure "plev(k)".)
156
157    integer i
158
159    !--------------------------------------------
160
161    call assert(shape(v3) == (/klon, llm/), "regr_pr_int v3")
162    call assert(shape(pplay) == (/klon, llm/), "regr_pr_int pplay")
163
164    !$omp master
165    if (is_mpi_root) then
166       call nf95_inq_varid(ncid, name, varid)
167
168       ! Get data at the right day from the input file:
169       ncerr = nf90_get_var(ncid, varid, v1(1, :, 1:), start=(/1, 1, julien/))
170       call handle_err("regr_pr_int nf90_get_var " // name, ncerr, ncid)
171       ! Latitudes are in ascending order in the input file while
172       ! "rlatu" is in descending order so we need to invert order:
173       v1(1, :, 1:) = v1(1, jjm+1:1:-1, 1:)
174
175       ! Complete "v1" with the value at 0 pressure:
176       v1(1, :, 0) = top_value
177
178       ! Duplicate on all longitudes:
179       v1(2:, :, :) = spread(v1(1, :, :), dim=1, ncopies=iim-1)
180    end if
181    !$omp end master
182
183    call scatter2d(v1, v2)
184
185    ! Regrid in pressure at each horizontal position:
186    do i = 1, klon
187       v3(i, llm:1:-1) = regr1_lint(v2(i, :), (/0., plev/), pplay(i, llm:1:-1))
188       ! (invert order of indices because "pplay" is in descending order)
189    end do
190
191  end subroutine regr_pr_int
192
193end module regr_pr
Note: See TracBrowser for help on using the repository browser.