/[lmdze]/trunk/libf/phylmd/Mobidic/regr_pr_coefoz.f90
ViewVC logotype

Contents of /trunk/libf/phylmd/Mobidic/regr_pr_coefoz.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 52 - (show annotations)
Fri Sep 23 12:28:01 2011 UTC (12 years, 7 months ago) by guez
File size: 6092 byte(s)
Split "conflx.f" into single-procedure files in directory "Conflx".

Split "cv_routines.f" into single-procedure files in directory
"CV_routines". Made module "cvparam" from included file
"cvparam.h". No included file other than "netcdf.inc" left in LMDZE.

1 module regr_pr_coefoz
2
3 ! Both procedures of this module read a single Mobidic ozone
4 ! coefficient from "coefoz_LMDZ.nc", at the current day, regrid this
5 ! coefficient in pressure to the LMDZ vertical grid and pack it to the LMDZ
6 ! horizontal "physics" grid.
7 ! The input data is a 2D latitude -- pressure field.
8 ! The target horizontal LMDZ grid is the "scalar" grid: "rlonv", "rlatu".
9 ! We assume that the input data is already on the LMDZ "rlatu"
10 ! latitude grid.
11
12 implicit none
13
14 contains
15
16 subroutine regr_pr_av_coefoz(ncid, name, julien, v3)
17
18 ! "regr_pr_av_coefoz" stands for "regrid pressure averaging
19 ! coefficient ozone".
20 ! The target vertical LMDZ grid is the grid of layer boundaries.
21 ! The input data does not depend on longitude, but the pressure
22 ! at LMDZ layers does.
23 ! Therefore, the values on the LMDZ grid do depend on longitude.
24 ! Regridding in pressure is done by averaging a step function of pressure.
25
26 use dimens_m, only: iim, jjm, llm
27 use dimphy, only: klon
28 use grid_change, only: dyn_phy
29 use netcdf, only: nf90_get_var
30 use netcdf95, only: nf95_inq_varid, handle_err
31 use nr_util, only: assert
32 use numer_rec, only: regr1_step_av
33 use press_coefoz_m, only: press_in_edg
34 use pressure_var, only: p3d
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(out):: v3(:, :) ! (klon, llm)
41 ! (ozone coefficient from Mobidic on the "physics" grid)
42 ! ("v3(i, k)" is at longitude "xlon(i)", latitude
43 ! "xlat(i)", in layer "k".)
44
45 ! Variables local to the procedure:
46
47 integer varid, ncerr ! for NetCDF
48
49 real v1(jjm + 1, size(press_in_edg) - 1)
50 ! (ozone coefficient from "coefoz_LMDZ.nc" at day "julien")
51 ! ("v1(j, k)" is at latitude "rlatu(j)" and for
52 ! pressure interval "[press_in_edg(k), press_in_edg(k+1)]".)
53
54 real v2(iim + 1, jjm + 1, llm)
55 ! (ozone coefficient on the "dynamics" grid)
56 ! ("v2(i, j, k)" is at longitude "rlonv(i)", latitude
57 ! "rlatu(j)" and for pressure interval "[p3d(i, j, k+1), p3d(i, j, k)]".)
58
59 integer i, j, k
60
61 !--------------------------------------------
62
63 call assert(shape(v3) == (/klon, llm/), "regr_pr_av_coefoz")
64
65 call nf95_inq_varid(ncid, name, varid)
66
67 ! Get data at the right day from the input file:
68 ncerr = nf90_get_var(ncid, varid, v1, start=(/1, 1, julien/))
69 call handle_err("regr_pr_av_coefoz nf90_get_var " // name, ncerr, ncid)
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 = v1(jjm+1:1:-1, :)
73
74 ! Regrid in pressure at each horizontal position:
75 do j = 1, jjm + 1
76 do i = 1, iim
77 if (dyn_phy(i, j)) then
78 v2(i, j, llm:1:-1) &
79 = regr1_step_av(v1(j, :), press_in_edg, &
80 p3d(i, j, llm+1:1:-1))
81 ! (invert order of indices because "p3d" is in descending order)
82 end if
83 end do
84 end do
85
86 forall (k = 1:llm) v3(:, k) = pack(v2(:, :, k), dyn_phy)
87
88 end subroutine regr_pr_av_coefoz
89
90 !***************************************************************
91
92 subroutine regr_pr_int_coefoz(ncid, name, julien, top_value, v3)
93
94 ! "regr_pr_int_coefoz" stands for "regrid pressure interpolation
95 ! coefficient ozone".
96 ! The target vertical LMDZ grid is the grid of mid-layers.
97 ! The input data does not depend on longitude, but the pressure
98 ! at LMDZ mid-layers does.
99 ! Therefore, the values on the LMDZ grid do depend on longitude.
100 ! Regridding is by linear interpolation.
101
102 use dimens_m, only: iim, jjm, llm
103 use dimphy, only: klon
104 use grid_change, only: dyn_phy
105 use netcdf, only: nf90_get_var
106 use netcdf95, only: nf95_inq_varid, handle_err
107 use nr_util, only: assert
108 use numer_rec, only: regr1_lint
109 use press_coefoz_m, only: plev
110 use pressure_var, only: pls
111
112 integer, intent(in):: ncid ! NetCDF ID of the file
113 character(len=*), intent(in):: name ! of the NetCDF variable
114 integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
115
116 real, intent(in):: top_value
117 ! (extra value of ozone coefficient at 0 pressure)
118
119 real, intent(out):: v3(:, :) ! (klon, llm)
120 ! (ozone coefficient from Mobidic on the "physics" grid)
121 ! ("v3(i, k)" is at longitude "xlon(i)", latitude
122 ! "xlat(i)", middle of layer "k".)
123
124 ! Variables local to the procedure:
125
126 integer varid, ncerr ! for NetCDF
127
128 real v1(jjm + 1, 0:size(plev))
129 ! (ozone coefficient from "coefoz_LMDZ.nc" at day "julien")
130 ! ("v1(j, k >=1)" is at latitude "rlatu(j)" and pressure "plev(k)".)
131
132 real v2(iim + 1, jjm + 1, llm)
133 ! (ozone coefficient on the "dynamics" grid)
134 ! "v2(i, j, k)" is at longitude "rlonv(i)", latitude
135 ! "rlatu(j)" and pressure "pls(i, j, k)".)
136
137 integer i, j, k
138
139 !--------------------------------------------
140
141 call assert(shape(v3) == (/klon, llm/), "regr_pr_int_coefoz")
142
143 call nf95_inq_varid(ncid, name, varid)
144
145 ! Get data at the right day from the input file:
146 ncerr = nf90_get_var(ncid, varid, v1(:, 1:), start=(/1, 1, julien/))
147 call handle_err("regr_pr_int_coefoz nf90_get_var " // name, ncerr, ncid)
148 ! Latitudes are in ascending order in the input file while
149 ! "rlatu" is in descending order so we need to invert order:
150 v1(:, 1:) = v1(jjm+1:1:-1, 1:)
151
152 ! Complete "v1" with the value at 0 pressure:
153 v1(:, 0) = top_value
154
155 ! Regrid in pressure at each horizontal position:
156 do j = 1, jjm + 1
157 do i = 1, iim
158 if (dyn_phy(i, j)) then
159 v2(i, j, llm:1:-1) &
160 = regr1_lint(v1(j, :), (/0., plev/), pls(i, j, llm:1:-1))
161 ! (invert order of indices because "pls" is in descending order)
162 end if
163 end do
164 end do
165
166 forall (k = 1:llm) v3(:, k) = pack(v2(:, :, k), dyn_phy)
167
168 end subroutine regr_pr_int_coefoz
169
170 end module regr_pr_coefoz

  ViewVC Help
Powered by ViewVC 1.1.21