1 | ! $Id$ |
---|
2 | module regr_lat_time_climoz_m |
---|
3 | |
---|
4 | ! Author: Lionel GUEZ |
---|
5 | |
---|
6 | implicit none |
---|
7 | |
---|
8 | private |
---|
9 | public regr_lat_time_climoz |
---|
10 | |
---|
11 | contains |
---|
12 | |
---|
13 | subroutine regr_lat_time_climoz |
---|
14 | |
---|
15 | ! "regr_lat_time_climoz" stands for "regrid latitude time |
---|
16 | ! climatology ozone". |
---|
17 | |
---|
18 | ! This procedure reads a climatology of ozone from a NetCDF file, |
---|
19 | ! regrids it in latitude and time, and writes the regridded field |
---|
20 | ! to a new NetCDF file. |
---|
21 | |
---|
22 | ! The input field depends on time, pressure level and latitude. |
---|
23 | |
---|
24 | ! If the input field has missing values, they must be signaled by |
---|
25 | ! the "missing_value" attribute. At each latitude and each date, the |
---|
26 | ! missing values are replaced by the lowest valid value above missing |
---|
27 | ! values. |
---|
28 | |
---|
29 | ! We assume that the input field is a step function of latitude |
---|
30 | ! and that the input latitude coordinate gives the centers of steps. |
---|
31 | ! Regridding in latitude is made by averaging, with a cosine of |
---|
32 | ! latitude factor. |
---|
33 | ! The target LMDZ latitude grid is the "scalar" grid: "rlatu". |
---|
34 | ! The values of "rlatu" are taken to be the centers of intervals. |
---|
35 | |
---|
36 | ! Regridding in time is by linear interpolation. |
---|
37 | ! Monthly values are processed to get daily values, on the basis |
---|
38 | ! of a 360-day calendar. |
---|
39 | ! The input file may contain either values for 12 months or values |
---|
40 | ! for 14 months. |
---|
41 | ! If there are 14 months then we assume that we have (in that order): |
---|
42 | ! December, January, February, ..., November, December, January |
---|
43 | ! We use the first December value to interpolate values between January |
---|
44 | ! 1st and mid-January. |
---|
45 | ! We use the last January value to interpolate values between |
---|
46 | ! mid-December and end of December. |
---|
47 | ! If there are only 12 months in the input file than we assume |
---|
48 | ! periodicity for interpolation at the beginning and at the end of the |
---|
49 | ! year. |
---|
50 | |
---|
51 | ! We assume that in the input file: |
---|
52 | ! -- latitude is in degrees; |
---|
53 | ! -- pressure is in hPa (even though we do not use pressure values |
---|
54 | ! here, we write the unit of pressure in the NetCDF header, and we |
---|
55 | ! will use the assumption later, when we regrid in pressure). |
---|
56 | ! -- latitude and pressure are strictly monotonic (as all NetCDF |
---|
57 | ! coordinate variables should be); |
---|
58 | ! -- time increases (even though we do not use values of the input |
---|
59 | ! time coordinate); |
---|
60 | ! -- missing values are vertically contiguous, at the bottom of |
---|
61 | ! the vertical domain; |
---|
62 | ! -- there is no missing value at the top level of the atmosphere. |
---|
63 | |
---|
64 | use regr1_step_av_m, only: regr1_step_av |
---|
65 | use regr3_lint_m, only: regr3_lint |
---|
66 | use netcdf95, only: nf95_open, nf95_close, nf95_inq_varid, handle_err, & |
---|
67 | nf95_put_var, nf95_gw_var |
---|
68 | use netcdf, only: nf90_nowrite, nf90_get_att, nf90_noerr |
---|
69 | |
---|
70 | ! Variables local to the procedure: |
---|
71 | |
---|
72 | include "dimensions.h" |
---|
73 | ! (for "jjm") |
---|
74 | include "paramet.h" |
---|
75 | ! (for the other included files) |
---|
76 | include "comgeom2.h" |
---|
77 | ! (for "rlatv") |
---|
78 | include "comconst.h" |
---|
79 | ! (for "pi") |
---|
80 | |
---|
81 | integer ncid_in, ncid_out ! NetCDF IDs for input and output files |
---|
82 | integer n_plev ! number of pressure levels in the input data |
---|
83 | integer n_lat! number of latitudes in the input data |
---|
84 | |
---|
85 | real, pointer:: latitude(:) |
---|
86 | ! (of input data, converted to rad, sorted in strictly ascending order) |
---|
87 | |
---|
88 | real, allocatable:: lat_in_edg(:) |
---|
89 | ! (edges of latitude intervals for input data, in rad, in strictly |
---|
90 | ! ascending order) |
---|
91 | |
---|
92 | real, pointer:: plev(:) |
---|
93 | ! pressure levels of input data, sorted in strictly ascending order |
---|
94 | |
---|
95 | logical desc_lat ! latitude in descending order in the input file |
---|
96 | logical desc_plev ! pressure levels in descending order in the input file |
---|
97 | |
---|
98 | real, pointer:: o3_in(:, :, :) ! (n_lat, n_plev, 12 or 0:13) |
---|
99 | ! (ozone climatology from the input file) |
---|
100 | ! ("o3_in(j, k, l)" is at latitude "latitude(j)" and pressure |
---|
101 | ! level "plev(k)". "l" is between 1 and 12 or between 0 and 13) |
---|
102 | |
---|
103 | real missing_value |
---|
104 | |
---|
105 | real, allocatable:: o3_regr_lat(:, :, :) ! (jjm + 1, n_plev, 0:13) |
---|
106 | ! (mean of "o3_in" over a latitude interval of LMDZ) |
---|
107 | ! (First dimension is latitude interval. |
---|
108 | ! The latitude interval for "o3_regr_lat(j,:, :)" contains "rlatu(j)". |
---|
109 | ! If "j" is between 2 and "jjm" then the interval is: |
---|
110 | ! [rlatv(j), rlatv(j-1)] |
---|
111 | ! If "j" is 1 or "jjm + 1" then the interval is: |
---|
112 | ! [rlatv(1), pi / 2] |
---|
113 | ! or: |
---|
114 | ! [- pi / 2, rlatv(jjm)] |
---|
115 | ! respectively. |
---|
116 | ! "o3_regr_lat(:, k, :)" is for pressure level "plev(k)". |
---|
117 | ! Last dimension is month number.) |
---|
118 | |
---|
119 | real, allocatable:: o3_out(:, :, :) ! (jjm + 1, n_plev, 360) |
---|
120 | ! (regridded ozone climatology) |
---|
121 | ! ("o3_out(j, k, l)" is at latitude "rlatu(j)", pressure |
---|
122 | ! level "plev(k)" and date "January 1st 0h" + "tmidday(l)", in a |
---|
123 | ! 360-day calendar.) |
---|
124 | |
---|
125 | integer j, k, l |
---|
126 | |
---|
127 | ! For NetCDF: |
---|
128 | integer varid_in, varid_out, varid_plev, varid_time, varid, ncerr |
---|
129 | |
---|
130 | real, parameter:: tmidmonth(0:13) = (/(-15. + 30. * l, l = 0, 13)/) |
---|
131 | ! (time to middle of month, in days since January 1st 0h, in a |
---|
132 | ! 360-day calendar) |
---|
133 | ! (We add values -15 and 375 so that, for example, day 3 of the year is |
---|
134 | ! interpolated between the December and the January value.) |
---|
135 | |
---|
136 | real, parameter:: tmidday(360) = (/(l + 0.5, l = 0, 359)/) |
---|
137 | ! (time to middle of day, in days since January 1st 0h, in a |
---|
138 | ! 360-day calendar) |
---|
139 | |
---|
140 | !--------------------------------- |
---|
141 | |
---|
142 | print *, "Call sequence information: regr_lat_time_climoz" |
---|
143 | |
---|
144 | call nf95_open("climoz.nc", nf90_nowrite, ncid_in) |
---|
145 | |
---|
146 | ! Get coordinates from the input file: |
---|
147 | |
---|
148 | call nf95_inq_varid(ncid_in, "latitude", varid) |
---|
149 | call nf95_gw_var(ncid_in, varid, latitude) |
---|
150 | ! Convert from degrees to rad, because we will take the sine of latitude: |
---|
151 | latitude = latitude / 180. * pi |
---|
152 | n_lat = size(latitude) |
---|
153 | ! We need to supply the latitudes to "regr1_step_av" in |
---|
154 | ! ascending order, so invert order if necessary: |
---|
155 | desc_lat = latitude(1) > latitude(n_lat) |
---|
156 | if (desc_lat) latitude = latitude(n_lat:1:-1) |
---|
157 | |
---|
158 | ! Compute edges of latitude intervals: |
---|
159 | allocate(lat_in_edg(n_lat + 1)) |
---|
160 | lat_in_edg(1) = - pi / 2 |
---|
161 | forall (j = 2:n_lat) lat_in_edg(j) = (latitude(j - 1) + latitude(j)) / 2 |
---|
162 | lat_in_edg(n_lat + 1) = pi / 2 |
---|
163 | deallocate(latitude) ! pointer |
---|
164 | |
---|
165 | call nf95_inq_varid(ncid_in, "plev", varid) |
---|
166 | call nf95_gw_var(ncid_in, varid, plev) |
---|
167 | n_plev = size(plev) |
---|
168 | ! We only need the pressure coordinate to copy it to the output file. |
---|
169 | ! The program "gcm" will assume that pressure levels are in |
---|
170 | ! ascending order in the regridded climatology so invert order if |
---|
171 | ! necessary: |
---|
172 | desc_plev = plev(1) > plev(n_plev) |
---|
173 | if (desc_plev) plev = plev(n_plev:1:-1) |
---|
174 | |
---|
175 | ! Create the output file and get the variable IDs: |
---|
176 | call prepare_out(ncid_in, n_plev, ncid_out, varid_out, varid_plev, & |
---|
177 | varid_time) |
---|
178 | |
---|
179 | ! Write remaining coordinate variables: |
---|
180 | call nf95_put_var(ncid_out, varid_plev, plev) |
---|
181 | call nf95_put_var(ncid_out, varid_time, tmidday) |
---|
182 | |
---|
183 | deallocate(plev) ! pointer |
---|
184 | |
---|
185 | call nf95_inq_varid(ncid_in, "tro3", varid_in) |
---|
186 | call nf95_gw_var(ncid_in, varid_in, o3_in) |
---|
187 | if (desc_lat) o3_in = o3_in(n_lat:1:-1, :, :) |
---|
188 | if (desc_plev) o3_in = o3_in(:, n_plev:1:-1, :) |
---|
189 | |
---|
190 | ncerr = nf90_get_att(ncid_in, varid_in, "missing_value", missing_value) |
---|
191 | if (ncerr == nf90_noerr) then |
---|
192 | do l = 1, size(o3_in, 3) |
---|
193 | do j = 1, n_lat |
---|
194 | ! Find missing values, starting from top of atmosphere |
---|
195 | ! and going down. |
---|
196 | ! We assume that the highest level contains no missing value. |
---|
197 | k = 2 |
---|
198 | do |
---|
199 | if (o3_in(j, k, l) == missing_value .or. k == n_plev) exit |
---|
200 | k = k + 1 |
---|
201 | end do |
---|
202 | ! Replace missing values with the valid value at the |
---|
203 | ! lowest level above missing values: |
---|
204 | if (o3_in(j, k, l) == missing_value) & |
---|
205 | o3_in(j, k:n_plev, l) = o3_in(j, k-1, l) |
---|
206 | end do |
---|
207 | end do |
---|
208 | else |
---|
209 | print *, "regr_lat_time_climoz: no missing value attribute" |
---|
210 | end if |
---|
211 | |
---|
212 | call nf95_close(ncid_in) |
---|
213 | |
---|
214 | allocate(o3_regr_lat(jjm + 1, n_plev, 0:13)) |
---|
215 | allocate(o3_out(jjm + 1, n_plev, 360)) |
---|
216 | |
---|
217 | ! Regrid in latitude: |
---|
218 | ! We average with respect to sine of latitude, which is |
---|
219 | ! equivalent to weighting by cosine of latitude: |
---|
220 | if (size(o3_in, 3) == 12) then |
---|
221 | print *, "Found 12 months in ozone climatology, assuming periodicity..." |
---|
222 | o3_regr_lat(jjm+1:1:-1, :, 1:12) = regr1_step_av(o3_in, & |
---|
223 | xs=sin(lat_in_edg), xt=sin((/- pi / 2, rlatv(jjm:1:-1), pi / 2/))) |
---|
224 | ! (invert order of indices in "o3_regr_lat" because "rlatu" is |
---|
225 | ! in descending order) |
---|
226 | |
---|
227 | ! Duplicate January and December values, in preparation of time |
---|
228 | ! interpolation: |
---|
229 | o3_regr_lat(:, :, 0) = o3_regr_lat(:, :, 12) |
---|
230 | o3_regr_lat(:, :, 13) = o3_regr_lat(:, :, 1) |
---|
231 | else |
---|
232 | print *, "Using 14 months in ozone climatology..." |
---|
233 | o3_regr_lat(jjm+1:1:-1, :, :) = regr1_step_av(o3_in, & |
---|
234 | xs=sin(lat_in_edg), xt=sin((/- pi / 2, rlatv(jjm:1:-1), pi / 2/))) |
---|
235 | ! (invert order of indices in "o3_regr_lat" because "rlatu" is |
---|
236 | ! in descending order) |
---|
237 | end if |
---|
238 | |
---|
239 | deallocate(o3_in) ! pointer |
---|
240 | |
---|
241 | ! Regrid in time by linear interpolation: |
---|
242 | o3_out = regr3_lint(o3_regr_lat, tmidmonth, tmidday) |
---|
243 | |
---|
244 | ! Write to file: |
---|
245 | call nf95_put_var(ncid_out, varid_out, o3_out(jjm+1:1:-1, :, :)) |
---|
246 | ! (The order of "rlatu" is inverted in the output file) |
---|
247 | |
---|
248 | call nf95_close(ncid_out) |
---|
249 | |
---|
250 | end subroutine regr_lat_time_climoz |
---|
251 | |
---|
252 | !******************************************** |
---|
253 | |
---|
254 | subroutine prepare_out(ncid_in, n_plev, ncid_out, varid_out, varid_plev, & |
---|
255 | varid_time) |
---|
256 | |
---|
257 | ! This subroutine creates the NetCDF output file, defines |
---|
258 | ! dimensions and variables, and writes one of the coordinate variables. |
---|
259 | |
---|
260 | use netcdf95, only: nf95_create, nf95_def_dim, nf95_def_var, & |
---|
261 | nf95_put_att, nf95_enddef, nf95_copy_att, nf95_put_var |
---|
262 | use netcdf, only: nf90_clobber, nf90_float, nf90_global |
---|
263 | |
---|
264 | integer, intent(in):: ncid_in, n_plev |
---|
265 | integer, intent(out):: ncid_out, varid_out, varid_plev, varid_time |
---|
266 | |
---|
267 | ! Variables local to the procedure: |
---|
268 | |
---|
269 | include "dimensions.h" |
---|
270 | ! (for "jjm") |
---|
271 | include "paramet.h" |
---|
272 | ! (for the other included files) |
---|
273 | include "comgeom2.h" |
---|
274 | ! (for "rlatu") |
---|
275 | include "comconst.h" |
---|
276 | ! (for "pi") |
---|
277 | |
---|
278 | integer ncerr |
---|
279 | integer dimid_rlatu, dimid_plev, dimid_time |
---|
280 | integer varid_rlatu |
---|
281 | |
---|
282 | !--------------------------- |
---|
283 | |
---|
284 | print *, "Call sequence information: prepare_out" |
---|
285 | |
---|
286 | call nf95_create("climoz_LMDZ.nc", nf90_clobber, ncid_out) |
---|
287 | |
---|
288 | ! Dimensions: |
---|
289 | call nf95_def_dim(ncid_out, "time", 360, dimid_time) |
---|
290 | call nf95_def_dim(ncid_out, "plev", n_plev, dimid_plev) |
---|
291 | call nf95_def_dim(ncid_out, "rlatu", jjm + 1, dimid_rlatu) |
---|
292 | |
---|
293 | ! Define coordinate variables: |
---|
294 | |
---|
295 | call nf95_def_var(ncid_out, "time", nf90_float, dimid_time, varid_time) |
---|
296 | call nf95_put_att(ncid_out, varid_time, "units", "days since 2000-1-1") |
---|
297 | call nf95_put_att(ncid_out, varid_time, "calendar", "360_day") |
---|
298 | call nf95_put_att(ncid_out, varid_time, "standard_name", "time") |
---|
299 | |
---|
300 | call nf95_def_var(ncid_out, "plev", nf90_float, dimid_plev, varid_plev) |
---|
301 | call nf95_put_att(ncid_out, varid_plev, "units", "millibar") |
---|
302 | call nf95_put_att(ncid_out, varid_plev, "standard_name", "air_pressure") |
---|
303 | call nf95_put_att(ncid_out, varid_plev, "long_name", "air pressure") |
---|
304 | |
---|
305 | call nf95_def_var(ncid_out, "rlatu", nf90_float, dimid_rlatu, varid_rlatu) |
---|
306 | call nf95_put_att(ncid_out, varid_rlatu, "units", "degrees_north") |
---|
307 | call nf95_put_att(ncid_out, varid_rlatu, "standard_name", "latitude") |
---|
308 | |
---|
309 | ! Define the primary variable: |
---|
310 | |
---|
311 | call nf95_def_var(ncid_out, "tro3", nf90_float, & |
---|
312 | (/dimid_rlatu, dimid_plev, dimid_time/), varid_out) |
---|
313 | call nf95_put_att(ncid_out, varid_out, "long_name", "ozone mole fraction") |
---|
314 | call nf95_put_att(ncid_out, varid_out, "standard_name", & |
---|
315 | "mole_fraction_of_ozone_in_air") |
---|
316 | |
---|
317 | ! Global attributes: |
---|
318 | |
---|
319 | ! The following commands, copying attributes, may fail. |
---|
320 | ! That is OK. |
---|
321 | ! It should just mean that the attribute is not defined in the input file. |
---|
322 | |
---|
323 | call nf95_copy_att(ncid_in, nf90_global, "Conventions", ncid_out, & |
---|
324 | nf90_global, ncerr) |
---|
325 | call handle_err_copy_att("Conventions") |
---|
326 | |
---|
327 | call nf95_copy_att(ncid_in, nf90_global, "title", ncid_out, nf90_global, & |
---|
328 | ncerr) |
---|
329 | call handle_err_copy_att("title") |
---|
330 | |
---|
331 | call nf95_copy_att(ncid_in, nf90_global, "institution", ncid_out, & |
---|
332 | nf90_global, ncerr) |
---|
333 | call handle_err_copy_att("institution") |
---|
334 | |
---|
335 | call nf95_copy_att(ncid_in, nf90_global, "source", ncid_out, nf90_global, & |
---|
336 | ncerr) |
---|
337 | call handle_err_copy_att("source") |
---|
338 | |
---|
339 | call nf95_put_att(ncid_out, nf90_global, "comment", "Regridded for LMDZ") |
---|
340 | |
---|
341 | call nf95_enddef(ncid_out) |
---|
342 | |
---|
343 | ! Write one of the coordinate variables: |
---|
344 | call nf95_put_var(ncid_out, varid_rlatu, rlatu(jjm+1:1:-1) / pi * 180.) |
---|
345 | ! (convert from rad to degrees and sort in ascending order) |
---|
346 | |
---|
347 | contains |
---|
348 | |
---|
349 | subroutine handle_err_copy_att(att_name) |
---|
350 | |
---|
351 | use netcdf, only: nf90_noerr, nf90_strerror |
---|
352 | |
---|
353 | character(len=*), intent(in):: att_name |
---|
354 | |
---|
355 | !---------------------------------------- |
---|
356 | |
---|
357 | if (ncerr /= nf90_noerr) then |
---|
358 | print *, "regr_lat_time_climoz_m prepare_out nf95_copy_att " & |
---|
359 | // att_name // " -- " // trim(nf90_strerror(ncerr)) |
---|
360 | end if |
---|
361 | |
---|
362 | end subroutine handle_err_copy_att |
---|
363 | |
---|
364 | end subroutine prepare_out |
---|
365 | |
---|
366 | end module regr_lat_time_climoz_m |
---|