source: LMDZ4/branches/LMDZ4-dev/libf/phylmd/regr_lat_time_climoz_m.F90 @ 1237

Last change on this file since 1237 was 1220, checked in by lguez, 15 years ago

-- Replaced "integer*4" declarations by "integer", "real*8" by

"real(kind=8)" and "real*4" by "real". Note that these are the only
modifications in the files "radiation_AR4.F" and "sw_aeroAR4.F90".

-- Corrected the kind of arguments to "max" and "min".

-- Replaced "nH" edit descriptors, which is a deleted feature in

Fortran 95, by character strings.

-- "regr_lat_time_climoz" now allows the pressure coordinate in the

input file to be in descending order.

-- Replaced call to not standard function "float" by call to intrinsic

function "real".

-- Included file "radepsi.h" in "physiq" was not used. Removed it.

The following set of modifications is related to the management of time.

-- In "gcm", "leapfrog" and "sortvarc0", "day_ini" was defined as 1

plus number of days between the reference date "(annee_ref,
day_ref)" and the first day of the current simulation. Changed
definition: "(annee_ref, day_ini)" is the first day of the current
simulation. There is an accompanying modification for "day_end".

-- Corrected bug in call to "ioconf_startdate" in "gcm".

-- Added call to "ioconf_calendar" in "create_etat0_limit".

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