source: LMDZ4/branches/LMDZ4V5.0-LF/libf/phylmd/regr_pr_av_m.F90 @ 3603

Last change on this file since 3603 was 1279, checked in by Laurent Fairhead, 15 years ago

Merged LMDZ4-dev branch changes r1241:1278 into the trunk
Running trunk and LMDZ4-dev in LMDZOR configuration on local
machine (sequential) and SX8 (4-proc) yields identical results
(restart and restartphy are identical binarily)
Log history from r1241 to r1278 is available by switching to
source:LMDZ4/branches/LMDZ4-dev-20091210

File size: 4.4 KB
Line 
1! $Id$
2module regr_pr_av_m
3
4  ! Author: Lionel GUEZ
5
6  implicit none
7
8contains
9
10  subroutine regr_pr_av(ncid, name, julien, press_in_edg, paprs, v3)
11
12    ! "regr_pr_av" stands for "regrid pressure averaging".
13    ! In this procedure:
14    ! -- the root process reads 2D latitude-pressure fields from a
15    !    NetCDF file, at a given day.
16    ! -- the fields are packed to the LMDZ horizontal "physics"
17    !    grid and scattered to all threads of all processes;
18    ! -- in all the threads of all the processes, the fields are regridded in
19    !    pressure to the LMDZ vertical grid.
20    ! We assume that, in the input file, the fields have 3 dimensions:
21    ! latitude, pressure, julian day.
22    ! We assume that the input fields are already on the "rlatu"
23    ! latitudes, excepth that latitudes are in ascending order in the input
24    ! file.
25    ! We assume that the inputs fields have the same pressure coordinate.
26
27    ! The target vertical LMDZ grid is the grid of layer boundaries.
28    ! Regridding in pressure is done by averaging a step function of pressure.
29
30    ! All the fields are regridded as a single multi-dimensional array
31    ! so it saves CPU time to call this procedure once for several NetCDF
32    ! variables rather than several times, each time for a single
33    ! NetCDF variable.
34
35    use dimphy, only: klon
36    use netcdf95, only: nf95_inq_varid, handle_err
37    use netcdf, only: nf90_get_var
38    use assert_m, only: assert
39    use assert_eq_m, only: assert_eq
40    use regr1_step_av_m, only: regr1_step_av
41    use mod_phys_lmdz_mpi_data, only: is_mpi_root
42
43    use mod_phys_lmdz_transfert_para, only: scatter2d
44    ! (pack to the LMDZ horizontal "physics" grid and scatter)
45
46    integer, intent(in):: ncid ! NetCDF ID of the file
47    character(len=*), intent(in):: name(:) ! of the NetCDF variables
48    integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
49
50    real, intent(in):: press_in_edg(:)
51    ! edges of pressure intervals for input data, in Pa, in strictly
52    ! ascending order
53
54    real, intent(in):: paprs(:, :) ! (klon, llm + 1)
55    ! (pression pour chaque inter-couche, en Pa)
56
57    real, intent(out):: v3(:, :, :) ! (klon, llm, size(name))
58    ! regridded fields on the partial "physics" grid
59    ! "v3(i, k, l)" is at longitude "xlon(i)", latitude
60    ! "xlat(i)", in pressure interval "[paprs(i, k+1), paprs(i, k)]",
61    ! for NetCDF variable "name(l)".
62
63    ! Variables local to the procedure:
64
65    include "dimensions.h"
66    integer varid, ncerr ! for NetCDF
67
68    real  v1(iim, jjm + 1, size(press_in_edg) - 1, size(name))
69    ! input fields at day "julien", on the global "dynamics" horizontal grid
70    ! First dimension is for longitude.
71    ! The values are the same for all longitudes.
72    ! "v1(:, j, k, l)" is at latitude "rlatu(j)", for
73    ! pressure interval "[press_in_edg(k), press_in_edg(k+1)]" and
74    ! NetCDF variable "name(l)".
75
76    real v2(klon, size(press_in_edg) - 1, size(name))
77    ! fields scattered to the partial "physics" horizontal grid
78    ! "v2(i, k, l)" is at longitude "xlon(i)", latitude "xlat(i)",
79    ! for pressure interval "[press_in_edg(k), press_in_edg(k+1)]" and
80    ! NetCDF variable "name(l)".
81
82    integer i, n_var
83
84    !--------------------------------------------
85
86    call assert(size(v3, 1) == klon, size(v3, 2) == llm, "regr_pr_av v3 klon")
87    n_var = assert_eq(size(name), size(v3, 3), "regr_pr_av v3 n_var")
88    call assert(shape(paprs) == (/klon, llm+1/), "regr_pr_av paprs")
89
90    !$omp master
91    if (is_mpi_root) then
92       do i = 1, n_var
93          call nf95_inq_varid(ncid, name(i), varid)
94         
95          ! Get data at the right day from the input file:
96          ncerr = nf90_get_var(ncid, varid, v1(1, :, :, i), &
97               start=(/1, 1, julien/))
98          call handle_err("regr_pr_av nf90_get_var " // name(i), ncerr, ncid)
99       end do
100       
101       ! Latitudes are in ascending order in the input file while
102       ! "rlatu" is in descending order so we need to invert order:
103       v1(1, :, :, :) = v1(1, jjm+1:1:-1, :, :)
104
105       ! Duplicate on all longitudes:
106       v1(2:, :, :, :) = spread(v1(1, :, :, :), dim=1, ncopies=iim-1)
107    end if
108    !$omp end master
109
110    call scatter2d(v1, v2)
111
112    ! Regrid in pressure at each horizontal position:
113    do i = 1, klon
114       v3(i, llm:1:-1, :) = regr1_step_av(v2(i, :, :), press_in_edg, &
115            paprs(i, llm+1:1:-1))
116       ! (invert order of indices because "paprs" is in descending order)
117    end do
118
119  end subroutine regr_pr_av
120
121end module regr_pr_av_m
Note: See TracBrowser for help on using the repository browser.