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

Annotation of /trunk/Sources/phylmd/Mobidic/regr_pr_av.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 90 - (hide annotations)
Wed Mar 12 21:16:36 2014 UTC (10 years, 2 months ago) by guez
Original Path: trunk/phylmd/Mobidic/regr_pr_av.f
File size: 2906 byte(s)
Removed procedures ini_histday, ini_histhf, write_histday and
write_histhf.

Divided file regr_pr_coefoz.f into regr_pr_av.f and
regr_pr_int.f. (Following LMDZ.) Divided module regr_pr_coefoz into
modules regr_pr_av_m and regr_pr_int_m. Renamed regr_pr_av_coefoz to
regr_pr_av and regr_pr_int_coefoz to regr_pr_int. The idea is that
those procedures are more general than Mobidic.

Removed argument dudyn of calfis and physiq. dudyn is not used either
in LMDZ. Removed computation in calfis of unused variable zpsrf (not
used either in LMDZ). Removed useless computation of dqfi in calfis
(part 62): the results were overwritten. (Same in LMDZ.)

1 guez 90 module regr_pr_av_m
2 guez 3
3     implicit none
4    
5     contains
6    
7 guez 90 subroutine regr_pr_av(ncid, name, julien, v3)
8 guez 3
9 guez 90 ! "regr_pr_av" stands for "regrid pressure averaging".
10    
11     ! This procedure reads a 2D latitude-pressure field from a NetCDF
12     ! file, at a given day, regrids this field in pressure to the LMDZ
13     ! vertical grid and packs it to the LMDZ horizontal "physics"
14     ! grid.
15    
16     ! We assume that, in the input file, the field has 3 dimensions:
17     ! latitude, pressure, julian day.
18    
19     ! We assume that the input field is already on the LMDZ "rlatu"
20     ! latitudes, except that latitudes are in ascending order in the
21     ! input file.
22    
23 guez 19 ! The target vertical LMDZ grid is the grid of layer boundaries.
24 guez 20 ! Regridding in pressure is done by averaging a step function of pressure.
25 guez 3
26     use dimens_m, only: iim, jjm, llm
27     use dimphy, only: klon
28 guez 52 use grid_change, only: dyn_phy
29 guez 90 use netcdf95, only: nf95_inq_varid, nf95_get_var
30 guez 36 use nr_util, only: assert
31 guez 61 use numer_rec_95, only: regr1_step_av
32 guez 17 use press_coefoz_m, only: press_in_edg
33 guez 19 use pressure_var, only: p3d
34 guez 3
35     integer, intent(in):: ncid ! NetCDF ID of the file
36     character(len=*), intent(in):: name ! of the NetCDF variable
37 guez 7 integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
38 guez 3
39 guez 7 real, intent(out):: v3(:, :) ! (klon, llm)
40 guez 90 ! regridded field on the partial "physics" grid
41     ! "v3(i, k)" is at longitude "xlon(i)", latitude "xlat(i)", in
42     ! layer "k".
43 guez 7
44     ! Variables local to the procedure:
45    
46 guez 20 integer varid, ncerr ! for NetCDF
47    
48 guez 7 real v1(jjm + 1, size(press_in_edg) - 1)
49 guez 90 ! input field at day "julien"
50     ! "v1(j, k)" is at latitude "rlatu(j)" and for
51     ! pressure interval "[press_in_edg(k), press_in_edg(k+1)]".
52 guez 7
53     real v2(iim + 1, jjm + 1, llm)
54 guez 90 ! regridded field on the "dynamics" grid
55 guez 20 ! ("v2(i, j, k)" is at longitude "rlonv(i)", latitude
56     ! "rlatu(j)" and for pressure interval "[p3d(i, j, k+1), p3d(i, j, k)]".)
57 guez 7
58 guez 20 integer i, j, k
59    
60 guez 7 !--------------------------------------------
61    
62 guez 90 call assert(shape(v3) == (/klon, llm/), "regr_pr_av")
63 guez 7
64     call nf95_inq_varid(ncid, name, varid)
65    
66     ! Get data at the right day from the input file:
67 guez 90 call nf95_get_var(ncid, varid, v1, start=(/1, 1, julien/))
68 guez 20 ! Latitudes are in ascending order in the input file while
69     ! "rlatu" is in descending order so we need to invert order:
70 guez 7 v1 = v1(jjm+1:1:-1, :)
71    
72     ! Regrid in pressure at each horizontal position:
73 guez 19 do j = 1, jjm + 1
74     do i = 1, iim
75     if (dyn_phy(i, j)) then
76     v2(i, j, llm:1:-1) &
77     = regr1_step_av(v1(j, :), press_in_edg, &
78     p3d(i, j, llm+1:1:-1))
79 guez 20 ! (invert order of indices because "p3d" is in descending order)
80 guez 19 end if
81     end do
82     end do
83 guez 7
84     forall (k = 1:llm) v3(:, k) = pack(v2(:, :, k), dyn_phy)
85    
86 guez 90 end subroutine regr_pr_av
87 guez 7
88 guez 90 end module regr_pr_av_m

  ViewVC Help
Powered by ViewVC 1.1.21