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

Last change on this file since 1905 was 1263, checked in by lguez, 15 years ago

1) Reactivated ability to read ozone (that was deactivated because of
dependency on version of IOIPSL). Added ability to read a pressure
coordinate in Pa in "regr_lat_time_climoz".

2) Added the ability to read a second ozone climatology, corresponding to
daylight ozone:

-- "read_climoz" is now an integer variable, instead of a logical
variable.

-- Added argument "read_climoz" to "phys_state_var_init",
"phys_output_open" and "regr_lat_time_climoz".

-- Created new variable "ozone_daylight" for "hist*.nc" output files.

-- Added a third dimension to variable "wo" in module
"phys_state_var_mod" and variable "POZON" in "radlwsw": index 1 for
average day-night ozone, index 2 for daylight ozone.

-- Added a fourth dimension to variables "o3_in", "o3_regr_lat" and
"o3_out" in "regr_lat_time_climoz": index 1 for average day-night
ozone, index 2 for daylight ozone.

-- In "physiq", moved call to "conf_phys" before call to
"phys_state_var_init". Thus, "conf_phys" is now inside the block "if
(first)" instead of "IF (debut)". There were definitions of "bl95_b0"
and "bl95_b1" that were useless because the variables were overwritten
by "conf_phys". Removed those definitions.

-- In "radlwsw", we pass the average day-night ozone to "LW_LMDAR4"
and the daylight ozone, if we have it, to "SW_LMDAR4" or
"SW_AEROAR4". If we do not have a specific field for daylight ozone
then "SW_LMDAR4" or "SW_AEROAR4" just get the average day-night ozone.

-- "regr_lat_time_climoz" now manages latitudes where the input ozone
field is missing at all levels (polar night).

-- Encapsulated "radlwsw" in a module.

3) Modifications to make sequential and parallel versions of
"create_etat0_limit" almost identical:

-- In "dyn3dpar/create_etat0_limit.F". No need to call
"phys_state_var_init", removed "use phys_state_var_mod" statement. No
need for "clesphys.h", removed "include" statement.

-- In "dyn3dpar/etat0_netcdf.F". Added argument "tau_ratqs" in call to
"conf_phys" (this bug was already corrected in "dyn3d"). Moved call to
"inifilr" after call to "infotrac_init" (as in "dyn3d").

4) Other peripheral modifications:

-- Added procedures "nf95_get_att" and "nf95_def_var_scalar" in
NetCDF95 interface. Overloaded "nf95_put_var" with three more
procedures: "nf95_put_var_FourByteReal", "nf95_put_var_FourByteInt",
"nf95_put_var_1D_FourByteInt".

-- Overloaded "regr1_step_av" with one more procedure:
"regr14_step_av". Overloaded "regr3_lint" with one more procedure:
"regr34_lint".

-- Corrected call to "Init_Phys_lmdz" in "dyn3d/create_etat0_limit.F":
the last argument should be an array, not a scalar.

-- Encapsulated "conf_phys" in a module.

-- Splitted module "regr_pr" into "regr_pr_av_m" and "regr_pr_int_m".

5) Tests:

This revision was compared to revision 1259, with optimization options
"debug" and "dev", parallelization options "none", "mpi", "omp" and
"mpi_omp", 1 and 2 MPI processes, 1 and 2 OpenMP threads, with the
compiler "FORTRAN90/SX Version 2.0 for SX-8". Both programs
"create_etat0_limit" and "gcm" were tested. In all cases,
parallelization does not change the results. With "read_climoz = 0" in
the ".def" files, the results of revision 1259 and of this revision
are the same.

File size: 16.6 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(read_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.
26
27    ! We assume that the input field is a step function of latitude
28    ! and that the input latitude coordinate gives the centers of steps.
29    ! Regridding in latitude is made by averaging, with a cosine of
30    ! latitude factor.
31    ! The target LMDZ latitude grid is the "scalar" grid: "rlatu".
32    ! The values of "rlatu" are taken to be the centers of intervals.
33
34    ! We assume that in the input file:
35
36    ! -- Latitude is in degrees.
37
38    ! -- Latitude and pressure are strictly monotonic (as all NetCDF
39    ! coordinate variables should be).
40
41    ! -- The time coordinate is in ascending order (even though we do
42    ! not use its values).
43    ! The input file may contain either values for 12 months or values
44    ! for 14 months.
45    ! If there are 14 months then we assume that we have (in that order):
46    ! December, January, February, ..., November, December, January
47
48    ! -- Missing values are contiguous, at the bottom of
49    ! the vertical domain and at the latitudinal boundaries.
50
51    ! If values are all missing at a given latitude and date, then we
52    ! replace those missing values by values at the closest latitude,
53    ! equatorward, with valid values.
54    ! Then, at each latitude and each date, the missing values are replaced
55    ! by the lowest valid value above missing values.
56
57    ! Regridding in time is by linear interpolation.
58    ! Monthly values are processed to get daily values, on the basis
59    ! of a 360-day calendar.
60    ! If there are 14 months, we use the first December value to
61    ! interpolate values between January 1st and mid-January.
62    ! We use the last January value to interpolate values between
63    ! mid-December and end of December.
64    ! If there are only 12 months in the input file then we assume
65    ! periodicity for interpolation at the beginning and at the end of the
66    ! year.
67
68    use regr1_step_av_m, only: regr1_step_av
69    use regr3_lint_m, only: regr3_lint
70    use netcdf95, only: handle_err, nf95_close, nf95_get_att, nf95_gw_var, &
71         nf95_inq_dimid, nf95_inq_varid, nf95_inquire_dimension, nf95_open, &
72         nf95_put_var
73    use netcdf, only: nf90_get_att, nf90_get_var, nf90_noerr, nf90_nowrite
74    use assert_m, only: assert
75
76    integer, intent(in):: read_climoz ! read ozone climatology
77    ! Allowed values are 1 and 2
78    ! 1: read a single ozone climatology that will be used day and night
79    ! 2: read two ozone climatologies, the average day and night
80    ! climatology and the daylight climatology
81
82    ! Variables local to the procedure:
83
84    include "dimensions.h"
85    ! (for "jjm")
86    include "paramet.h"
87    ! (for the other included files)
88    include "comgeom2.h"
89    ! (for "rlatv")
90    include "comconst.h"
91    ! (for "pi")
92
93    integer n_plev ! number of pressure levels in the input data
94    integer n_lat ! number of latitudes in the input data
95    integer n_month ! number of months in the input data
96
97    real, pointer:: latitude(:)
98    ! (of input data, converted to rad, sorted in strictly ascending order)
99
100    real, allocatable:: lat_in_edg(:)
101    ! (edges of latitude intervals for input data, in rad, in strictly
102    ! ascending order)
103
104    real, pointer:: plev(:)
105    ! pressure levels of input data, sorted in strictly ascending
106    ! order, converted to hPa
107
108    logical desc_lat ! latitude in descending order in the input file
109    logical desc_plev ! pressure levels in descending order in the input file
110
111    real, allocatable:: o3_in(:, :, :, :)
112    ! (n_lat, n_plev, n_month, read_climoz)
113    ! ozone climatologies from the input file
114    ! "o3_in(j, k, :, :)" is at latitude "latitude(j)" and pressure
115    ! level "plev(k)".
116    ! Third dimension is month index, first value may be December or January.
117    ! "o3_in(:, :, :, 1)" is for the day- night average, "o3_in(:, :, :, 2)"
118    ! is for daylight.
119
120    real missing_value
121
122    real, allocatable:: o3_regr_lat(:, :, :, :)
123    ! (jjm + 1, n_plev, 0:13, read_climoz)
124    ! mean of "o3_in" over a latitude interval of LMDZ
125    ! First dimension is latitude interval.
126    ! The latitude interval for "o3_regr_lat(j,:, :, :)" contains "rlatu(j)".
127    ! If "j" is between 2 and "jjm" then the interval is:
128    ! [rlatv(j), rlatv(j-1)]
129    ! If "j" is 1 or "jjm + 1" then the interval is:
130    ! [rlatv(1), pi / 2]
131    ! or:
132    ! [- pi / 2, rlatv(jjm)]
133    ! respectively.
134    ! "o3_regr_lat(:, k, :, :)" is for pressure level "plev(k)".
135    ! Third dimension is month number, 1 for January.
136    ! "o3_regr_lat(:, :, :, 1)" is average day and night,
137    ! "o3_regr_lat(:, :, :, 2)" is for daylight.
138
139    real, allocatable:: o3_out(:, :, :, :)
140    ! (jjm + 1, n_plev, 360, read_climoz)
141    ! regridded ozone climatology
142    ! "o3_out(j, k, l, :)" is at latitude "rlatu(j)", pressure
143    ! level "plev(k)" and date "January 1st 0h" + "tmidday(l)", in a
144    ! 360-day calendar.
145    ! "o3_out(:, :, :, 1)" is average day and night,
146    ! "o3_out(:, :, :, 2)" is for daylight.
147
148    integer j, k, l,m
149
150    ! For NetCDF:
151    integer ncid_in, ncid_out ! IDs for input and output files
152    integer varid_plev, varid_time, varid, ncerr, dimid
153    character(len=80) press_unit ! pressure unit
154
155    integer varid_in(read_climoz), varid_out(read_climoz)
156    ! index 1 is for average ozone day and night, index 2 is for
157    ! daylight ozone.
158
159    real, parameter:: tmidmonth(0:13) = (/(-15. + 30. * l, l = 0, 13)/)
160    ! (time to middle of month, in days since January 1st 0h, in a
161    ! 360-day calendar)
162    ! (We add values -15 and 375 so that, for example, day 3 of the year is
163    ! interpolated between the December and the January value.)
164
165    real, parameter:: tmidday(360) = (/(l + 0.5, l = 0, 359)/)
166    ! (time to middle of day, in days since January 1st 0h, in a
167    ! 360-day calendar)
168
169    !---------------------------------
170
171    print *, "Call sequence information: regr_lat_time_climoz"
172    call assert(read_climoz == 1 .or. read_climoz == 2, "regr_lat_time_climoz")
173
174    call nf95_open("climoz.nc", nf90_nowrite, ncid_in)
175
176    ! Get coordinates from the input file:
177
178    call nf95_inq_varid(ncid_in, "latitude", varid)
179    call nf95_gw_var(ncid_in, varid, latitude)
180    ! Convert from degrees to rad, because we will take the sine of latitude:
181    latitude = latitude / 180. * pi
182    n_lat = size(latitude)
183    ! We need to supply the latitudes to "regr1_step_av" in
184    ! ascending order, so invert order if necessary:
185    desc_lat = latitude(1) > latitude(n_lat)
186    if (desc_lat) latitude = latitude(n_lat:1:-1)
187
188    ! Compute edges of latitude intervals:
189    allocate(lat_in_edg(n_lat + 1))
190    lat_in_edg(1) = - pi / 2
191    forall (j = 2:n_lat) lat_in_edg(j) = (latitude(j - 1) + latitude(j)) / 2
192    lat_in_edg(n_lat + 1) = pi / 2
193    deallocate(latitude) ! pointer
194
195    call nf95_inq_varid(ncid_in, "plev", varid)
196    call nf95_gw_var(ncid_in, varid, plev)
197    n_plev = size(plev)
198    ! We only need the pressure coordinate to copy it to the output file.
199    ! The program "gcm" will assume that pressure levels are in
200    ! ascending order in the regridded climatology so invert order if
201    ! necessary:
202    desc_plev = plev(1) > plev(n_plev)
203    if (desc_plev) plev = plev(n_plev:1:-1)
204    call nf95_get_att(ncid_in, varid, "units", press_unit)
205    if (press_unit == "Pa") then
206       ! Convert to hPa:
207       plev = plev / 100.
208    elseif (press_unit /= "hPa") then
209       print *, "regr_lat_time_climoz: the only recognized units are Pa " &
210            // "and hPa."
211       stop 1
212    end if
213
214    ! Create the output file and get the variable IDs:
215    call prepare_out(ncid_in, n_plev, ncid_out, varid_out, varid_plev, &
216         varid_time)
217
218    ! Write remaining coordinate variables:
219    call nf95_put_var(ncid_out, varid_plev, plev)
220    call nf95_put_var(ncid_out, varid_time, tmidday)
221
222    deallocate(plev) ! pointer
223
224    ! Get the  number of months:
225    call nf95_inq_dimid(ncid_in, "time", dimid)
226    call nf95_inquire_dimension(ncid_in, dimid, len=n_month)
227
228    allocate(o3_in(n_lat, n_plev, n_month, read_climoz))
229
230    call nf95_inq_varid(ncid_in, "tro3", varid_in(1))
231    ncerr = nf90_get_var(ncid_in, varid_in(1), o3_in(:, :, :, 1))
232    call handle_err("regr_lat_time_climoz nf90_get_var tro3", ncerr, ncid_in)
233
234    if (read_climoz == 2) then
235       call nf95_inq_varid(ncid_in, "tro3_daylight", varid_in(2))
236       ncerr = nf90_get_var(ncid_in, varid_in(2), o3_in(:, :, :, 2))
237       call handle_err("regr_lat_time_climoz nf90_get_var tro3_daylight", &
238            ncerr, ncid_in, varid_in(2))
239    end if
240
241    if (desc_lat) o3_in = o3_in(n_lat:1:-1, :, :, :)
242    if (desc_plev) o3_in = o3_in(:, n_plev:1:-1, :, :)
243
244    do m = 1, read_climoz
245       ncerr = nf90_get_att(ncid_in, varid_in(m), "missing_value", &
246            missing_value)
247       if (ncerr == nf90_noerr) then
248          do l = 1, n_month
249             ! Take care of latitudes where values are all missing:
250
251             ! Next to the south pole:
252             j = 1
253             do while (o3_in(j, 1, l, m) == missing_value)
254                j = j + 1
255             end do
256             if (j > 1) o3_in(:j-1, :, l, m) = &
257                  spread(o3_in(j, :, l, m), dim=1, ncopies=j-1)
258             
259             ! Next to the north pole:
260             j = n_lat
261             do while (o3_in(j, 1, l, m) == missing_value)
262                j = j - 1
263             end do
264             if (j < n_lat) o3_in(j+1:, :, l, m) = &
265                  spread(o3_in(j, :, l, m), dim=1, ncopies=n_lat-j)
266
267             ! Take care of missing values at high pressure:
268             do j = 1, n_lat
269                ! Find missing values, starting from top of atmosphere
270                ! and going down.
271                ! We have already taken care of latitudes full of
272                ! missing values so the highest level has a valid value.
273                k = 2
274                do while  (o3_in(j, k, l, m) /= missing_value .and. k < n_plev)
275                   k = k + 1
276                end do
277                ! Replace missing values with the valid value at the
278                ! lowest level above missing values:
279                if (o3_in(j, k, l, m) == missing_value) &
280                     o3_in(j, k:n_plev, l, m) = o3_in(j, k-1, l, m)
281             end do
282          end do
283       else
284          print *, "regr_lat_time_climoz: field ", m, &
285               ", no missing value attribute"
286       end if
287    end do
288
289    call nf95_close(ncid_in)
290
291    allocate(o3_regr_lat(jjm + 1, n_plev, 0:13, read_climoz))
292    allocate(o3_out(jjm + 1, n_plev, 360, read_climoz))
293
294    ! Regrid in latitude:
295    ! We average with respect to sine of latitude, which is
296    ! equivalent to weighting by cosine of latitude:
297    if (n_month == 12) then
298       print *, &
299            "Found 12 months in ozone climatologies, assuming periodicity..."
300       o3_regr_lat(jjm+1:1:-1, :, 1:12, :) = regr1_step_av(o3_in, &
301            xs=sin(lat_in_edg), xt=sin((/- pi / 2, rlatv(jjm:1:-1), pi / 2/)))
302       ! (invert order of indices in "o3_regr_lat" because "rlatu" is
303       ! in descending order)
304
305       ! Duplicate January and December values, in preparation of time
306       ! interpolation:
307       o3_regr_lat(:, :, 0, :) = o3_regr_lat(:, :, 12, :)
308       o3_regr_lat(:, :, 13, :) = o3_regr_lat(:, :, 1, :)
309    else
310       print *, "Using 14 months in ozone climatologies..."
311       o3_regr_lat(jjm+1:1:-1, :, :, :) = regr1_step_av(o3_in, &
312            xs=sin(lat_in_edg), xt=sin((/- pi / 2, rlatv(jjm:1:-1), pi / 2/)))
313       ! (invert order of indices in "o3_regr_lat" because "rlatu" is
314       ! in descending order)
315    end if
316
317    ! Regrid in time by linear interpolation:
318    o3_out = regr3_lint(o3_regr_lat, tmidmonth, tmidday)
319
320    ! Write to file:
321    do m = 1, read_climoz
322       call nf95_put_var(ncid_out, varid_out(m), o3_out(jjm+1:1:-1, :, :, m))
323       ! (The order of "rlatu" is inverted in the output file)
324    end do
325
326    call nf95_close(ncid_out)
327
328  end subroutine regr_lat_time_climoz
329
330  !********************************************
331
332  subroutine prepare_out(ncid_in, n_plev, ncid_out, varid_out, varid_plev, &
333       varid_time)
334
335    ! This subroutine creates the NetCDF output file, defines
336    ! dimensions and variables, and writes one of the coordinate variables.
337
338    use netcdf95, only: nf95_create, nf95_def_dim, nf95_def_var, &
339         nf95_put_att, nf95_enddef, nf95_copy_att, nf95_put_var
340    use netcdf, only: nf90_clobber, nf90_float, nf90_global
341
342    integer, intent(in):: ncid_in, n_plev
343    integer, intent(out):: ncid_out, varid_plev, varid_time
344
345    integer, intent(out):: varid_out(:) ! dim(1 or 2)
346    ! "varid_out(1)" is for average ozone day and night,
347    ! "varid_out(2)" is for daylight ozone.
348
349    ! Variables local to the procedure:
350
351    include "dimensions.h"
352    ! (for "jjm")
353    include "paramet.h"
354    ! (for the other included files)
355    include "comgeom2.h"
356    ! (for "rlatu")
357    include "comconst.h"
358    ! (for "pi")
359
360    integer ncerr
361    integer dimid_rlatu, dimid_plev, dimid_time
362    integer varid_rlatu
363
364    !---------------------------
365
366    print *, "Call sequence information: prepare_out"
367
368    call nf95_create("climoz_LMDZ.nc", nf90_clobber, ncid_out)
369
370    ! Dimensions:
371    call nf95_def_dim(ncid_out, "time", 360, dimid_time)
372    call nf95_def_dim(ncid_out, "plev", n_plev, dimid_plev)
373    call nf95_def_dim(ncid_out, "rlatu", jjm + 1, dimid_rlatu)
374
375    ! Define coordinate variables:
376
377    call nf95_def_var(ncid_out, "time", nf90_float, dimid_time, varid_time)
378    call nf95_put_att(ncid_out, varid_time, "units", "days since 2000-1-1")
379    call nf95_put_att(ncid_out, varid_time, "calendar", "360_day")
380    call nf95_put_att(ncid_out, varid_time, "standard_name", "time")
381
382    call nf95_def_var(ncid_out, "plev", nf90_float, dimid_plev, varid_plev)
383    call nf95_put_att(ncid_out, varid_plev, "units", "millibar")
384    call nf95_put_att(ncid_out, varid_plev, "standard_name", "air_pressure")
385    call nf95_put_att(ncid_out, varid_plev, "long_name", "air pressure")
386
387    call nf95_def_var(ncid_out, "rlatu", nf90_float, dimid_rlatu, varid_rlatu)
388    call nf95_put_att(ncid_out, varid_rlatu, "units", "degrees_north")
389    call nf95_put_att(ncid_out, varid_rlatu, "standard_name", "latitude")
390
391    ! Define the primary variables:
392
393    call nf95_def_var(ncid_out, "tro3", nf90_float, &
394         (/dimid_rlatu, dimid_plev, dimid_time/), varid_out(1))
395    call nf95_put_att(ncid_out, varid_out(1), "long_name", &
396         "ozone mole fraction")
397    call nf95_put_att(ncid_out, varid_out(1), "standard_name", &
398         "mole_fraction_of_ozone_in_air")
399
400    if (size(varid_out) == 2) then
401       call nf95_def_var(ncid_out, "tro3_daylight", nf90_float, &
402            (/dimid_rlatu, dimid_plev, dimid_time/), varid_out(2))
403       call nf95_put_att(ncid_out, varid_out(2), "long_name", &
404            "ozone mole fraction in daylight")
405    end if
406
407    ! Global attributes:
408
409    ! The following commands, copying attributes, may fail.
410    ! That is OK.
411    ! It should just mean that the attribute is not defined in the input file.
412
413    call nf95_copy_att(ncid_in, nf90_global, "Conventions", ncid_out, &
414         nf90_global, ncerr)
415    call handle_err_copy_att("Conventions")
416
417    call nf95_copy_att(ncid_in, nf90_global, "title", ncid_out, nf90_global, &
418         ncerr)
419    call handle_err_copy_att("title")
420
421    call nf95_copy_att(ncid_in, nf90_global, "institution", ncid_out, &
422         nf90_global, ncerr)
423    call handle_err_copy_att("institution")
424
425    call nf95_copy_att(ncid_in, nf90_global, "source", ncid_out, nf90_global, &
426         ncerr)
427    call handle_err_copy_att("source")
428
429    call nf95_put_att(ncid_out, nf90_global, "comment", "Regridded for LMDZ")
430
431    call nf95_enddef(ncid_out)
432
433    ! Write one of the coordinate variables:
434    call nf95_put_var(ncid_out, varid_rlatu, rlatu(jjm+1:1:-1) / pi * 180.)
435    ! (convert from rad to degrees and sort in ascending order)
436
437  contains
438
439    subroutine handle_err_copy_att(att_name)
440
441      use netcdf, only: nf90_noerr, nf90_strerror
442
443      character(len=*), intent(in):: att_name
444
445      !----------------------------------------
446
447      if (ncerr /= nf90_noerr) then
448         print *, "regr_lat_time_climoz_m prepare_out nf95_copy_att " &
449              // att_name // " -- " // trim(nf90_strerror(ncerr))
450      end if
451
452    end subroutine handle_err_copy_att
453
454  end subroutine prepare_out
455
456end module regr_lat_time_climoz_m
Note: See TracBrowser for help on using the repository browser.