source: LMDZ6/trunk/libf/obsolete/regr_lat_time_climoz_m.F90 @ 4684

Last change on this file since 4684 was 2788, checked in by dcugnet, 8 years ago

Changes in ce0l about the way ozone forcing files are generated:

1) 3D raw input files "climoz.nc" are now handled.
2) Default behaviour is now to let the gcm interpolate in time online.

This helps to avoid huge forcing files (in particular for 3D fields).
In this case, the output files "climoz_LMDZ.nc" all have 14 records:

  • records 2-13 are obtained with records 1-12 of "climoz.nc".
  • records 1 and 14 are obtained respectively with:
    • record 12 of "climoz_m.nc" if available, of "climoz.nc" otherwise.
    • record 1 of "climoz_p.nc" if available, of "climoz.nc" otherwise.

3) If ok_daily_climoz key is TRUE, the time interpolation (one record

a day) is forced, using the 14 records described below.
This now depends on the calendar (it was on a 360 days basis only).

Changes in the gcm about the way zone forcing files are read/interpolated:

1) 3D horizontally interpolated "climoz_LMDZ.nc" files are now handled.
2) Daily files (already interpolated in time) are still handled, but their

number of records must match the expected number of days, that depends
on the calendar (records step is no longer 1/360 year).

3) 14 records monthly files are now handled (and prefered). This reduces

the I/O to a minimum and the aditional computational cost is low (simple
online linear time interpolation).

4) If adjust_tropopause key is TRUE, the input fields are stretched using

following method:

  • LMDZ dynamical tropopause is detected: Ptrop_lmdz = MAX ( P(Potential Vorticity==2PVU), P(theta==380K) )
  • file chemical tropopause is detected: Ptrop_file = P( tro3 == o3t ), where:

o3t = 91. + 28. * SIN(PI*(month-2)/6) (ppbV)

This formula comes from Thouret & al., ACP 6, 1033-1051, 2006.
The second term of the expression is multiplied by TANH(lat_deg/20.)
to account for latitude dependency.

  • File profile is streched in a +/- 5kms zone around the mean tropopause to ensure resulting tropopause matches the one of LMDZ. See procedure regr_pr_time_av for more details.
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 16.7 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 mod_grid_phy_lmdz, ONLY : nbp_lat
69    use regr1_conserv_m, only: regr1_conserv
70    use regr3_lint_m, only: regr3_lint
71    use netcdf95, only: handle_err, nf95_close, nf95_get_att, nf95_gw_var, &
72         nf95_inq_dimid, nf95_inq_varid, nf95_inquire_dimension, nf95_open, &
73         nf95_put_var
74    use netcdf, only: nf90_get_att, nf90_get_var, nf90_noerr, nf90_nowrite
75    use assert_m, only: assert
76    use regular_lonlat_mod, only : boundslat_reg, south
77    use nrtype, only: pi
78    use slopes_m, only: slopes
79
80    integer, intent(in):: read_climoz ! read ozone climatology
81    ! Allowed values are 1 and 2
82    ! 1: read a single ozone climatology that will be used day and night
83    ! 2: read two ozone climatologies, the average day and night
84    ! climatology and the daylight climatology
85
86    ! Variables local to the procedure:
87
88    integer n_plev ! number of pressure levels in the input data
89    integer n_lat ! number of latitudes in the input data
90    integer n_month ! number of months in the input data
91
92    real, pointer:: latitude(:)
93    ! (of input data, converted to rad, sorted in strictly ascending order)
94
95    real, allocatable:: sin_lat_in_edg(:)
96    ! (sine of edges of latitude intervals for input data, in rad, in strictly
97    ! ascending order)
98
99    real, pointer:: plev(:)
100    ! pressure levels of input data, sorted in strictly ascending
101    ! order, converted to hPa
102
103    logical desc_lat ! latitude in descending order in the input file
104    logical desc_plev ! pressure levels in descending order in the input file
105
106    real, allocatable:: o3_in(:, :, :, :)
107    ! (n_lat, n_plev, n_month, read_climoz)
108    ! ozone climatologies from the input file
109    ! "o3_in(j, k, :, :)" is at latitude "latitude(j)" and pressure
110    ! level "plev(k)".
111    ! Third dimension is month index, first value may be December or January.
112    ! "o3_in(:, :, :, 1)" is for the day- night average, "o3_in(:, :, :, 2)"
113    ! is for daylight.
114
115    real missing_value
116
117    real, allocatable:: o3_regr_lat(:, :, :, :)
118    ! (nbp_lat, n_plev, 0:13, read_climoz)
119    ! mean of "o3_in" over a latitude interval of LMDZ
120    ! First dimension is latitude interval.
121    ! The latitude interval for "o3_regr_lat(j,:, :, :)" contains "rlatu(j)".
122    ! If "j" is between 2 and "nbp_lat - 1" then the interval is:
123    ! [rlatv(j), rlatv(j-1)]
124    ! If "j" is 1 or "nbp_lat" then the interval is:
125    ! [rlatv(1), pi / 2]
126    ! or:
127    ! [- pi / 2, rlatv(nbp_lat - 1)]
128    ! respectively.
129    ! "o3_regr_lat(:, k, :, :)" is for pressure level "plev(k)".
130    ! Third dimension is month number, 1 for January.
131    ! "o3_regr_lat(:, :, :, 1)" is average day and night,
132    ! "o3_regr_lat(:, :, :, 2)" is for daylight.
133
134    real, allocatable:: o3_out(:, :, :, :)
135    ! (nbp_lat, n_plev, 360, read_climoz)
136    ! regridded ozone climatology
137    ! "o3_out(j, k, l, :)" is at latitude "rlatu(j)", pressure
138    ! level "plev(k)" and date "January 1st 0h" + "tmidday(l)", in a
139    ! 360-day calendar.
140    ! "o3_out(:, :, :, 1)" is average day and night,
141    ! "o3_out(:, :, :, 2)" is for daylight.
142
143    integer j, k, l,m
144
145    ! For NetCDF:
146    integer ncid_in, ncid_out ! IDs for input and output files
147    integer varid_plev, varid_time, varid, ncerr, dimid
148    character(len=80) press_unit ! pressure unit
149
150    integer varid_in(read_climoz), varid_out(read_climoz)
151    ! index 1 is for average ozone day and night, index 2 is for
152    ! daylight ozone.
153
154    real, parameter:: tmidmonth(0:13) = (/(-15. + 30. * l, l = 0, 13)/)
155    ! (time to middle of month, in days since January 1st 0h, in a
156    ! 360-day calendar)
157    ! (We add values -15 and 375 so that, for example, day 3 of the year is
158    ! interpolated between the December and the January value.)
159
160    real, parameter:: tmidday(360) = (/(l + 0.5, l = 0, 359)/)
161    ! (time to middle of day, in days since January 1st 0h, in a
162    ! 360-day calendar)
163
164    !---------------------------------
165
166    print *, "Call sequence information: regr_lat_time_climoz"
167    call assert(read_climoz == 1 .or. read_climoz == 2, "regr_lat_time_climoz")
168
169    call nf95_open("climoz.nc", nf90_nowrite, ncid_in)
170
171    ! Get coordinates from the input file:
172
173    call nf95_inq_varid(ncid_in, "latitude", varid)
174    call nf95_gw_var(ncid_in, varid, latitude)
175    ! Convert from degrees to rad, because we will take the sine of latitude:
176    latitude = latitude / 180. * pi
177    n_lat = size(latitude)
178    ! We need to supply the latitudes to "regr1_conserv" in
179    ! ascending order, so invert order if necessary:
180    desc_lat = latitude(1) > latitude(n_lat)
181    if (desc_lat) latitude = latitude(n_lat:1:-1)
182
183    ! Compute edges of latitude intervals:
184    allocate(sin_lat_in_edg(n_lat + 1))
185    sin_lat_in_edg(1) = - 1.
186    forall (j = 2:n_lat) sin_lat_in_edg(j) = sin((latitude(j - 1) &
187         + latitude(j)) / 2.)
188    sin_lat_in_edg(n_lat + 1) = 1.
189    deallocate(latitude) ! pointer
190
191    call nf95_inq_varid(ncid_in, "plev", varid)
192    call nf95_gw_var(ncid_in, varid, plev)
193    n_plev = size(plev)
194    ! We only need the pressure coordinate to copy it to the output file.
195    ! The program "gcm" will assume that pressure levels are in
196    ! ascending order in the regridded climatology so invert order if
197    ! necessary:
198    desc_plev = plev(1) > plev(n_plev)
199    if (desc_plev) plev = plev(n_plev:1:-1)
200    call nf95_get_att(ncid_in, varid, "units", press_unit)
201    if (press_unit == "Pa") then
202       ! Convert to hPa:
203       plev = plev / 100.
204    elseif (press_unit /= "hPa") then
205       print *, "regr_lat_time_climoz: the only recognized units are Pa " &
206            // "and hPa."
207       stop 1
208    end if
209
210    ! Create the output file and get the variable IDs:
211    call prepare_out(ncid_in, n_plev, ncid_out, varid_out, varid_plev, &
212         varid_time)
213
214    ! Write remaining coordinate variables:
215    call nf95_put_var(ncid_out, varid_plev, plev)
216    call nf95_put_var(ncid_out, varid_time, tmidday)
217
218    deallocate(plev) ! pointer
219
220    ! Get the  number of months:
221    call nf95_inq_dimid(ncid_in, "time", dimid)
222    call nf95_inquire_dimension(ncid_in, dimid, nclen=n_month)
223
224    allocate(o3_in(n_lat, n_plev, n_month, read_climoz))
225
226    call nf95_inq_varid(ncid_in, "tro3", varid_in(1))
227    ncerr = nf90_get_var(ncid_in, varid_in(1), o3_in(:, :, :, 1))
228    call handle_err("regr_lat_time_climoz nf90_get_var tro3", ncerr, ncid_in)
229
230    if (read_climoz == 2) then
231       call nf95_inq_varid(ncid_in, "tro3_daylight", varid_in(2))
232       ncerr = nf90_get_var(ncid_in, varid_in(2), o3_in(:, :, :, 2))
233       call handle_err("regr_lat_time_climoz nf90_get_var tro3_daylight", &
234            ncerr, ncid_in, varid_in(2))
235    end if
236
237    if (desc_lat) o3_in = o3_in(n_lat:1:-1, :, :, :)
238    if (desc_plev) o3_in = o3_in(:, n_plev:1:-1, :, :)
239
240    do m = 1, read_climoz
241       ncerr = nf90_get_att(ncid_in, varid_in(m), "missing_value", &
242            missing_value)
243       if (ncerr == nf90_noerr) then
244          do l = 1, n_month
245             ! Take care of latitudes where values are all missing:
246
247             ! Next to the south pole:
248             j = 1
249             do while (o3_in(j, 1, l, m) == missing_value)
250                j = j + 1
251             end do
252             if (j > 1) o3_in(:j-1, :, l, m) = &
253                  spread(o3_in(j, :, l, m), dim=1, ncopies=j-1)
254             
255             ! Next to the north pole:
256             j = n_lat
257             do while (o3_in(j, 1, l, m) == missing_value)
258                j = j - 1
259             end do
260             if (j < n_lat) o3_in(j+1:, :, l, m) = &
261                  spread(o3_in(j, :, l, m), dim=1, ncopies=n_lat-j)
262
263             ! Take care of missing values at high pressure:
264             do j = 1, n_lat
265                ! Find missing values, starting from top of atmosphere
266                ! and going down.
267                ! We have already taken care of latitudes full of
268                ! missing values so the highest level has a valid value.
269                k = 2
270                do while  (o3_in(j, k, l, m) /= missing_value .and. k < n_plev)
271                   k = k + 1
272                end do
273                ! Replace missing values with the valid value at the
274                ! lowest level above missing values:
275                if (o3_in(j, k, l, m) == missing_value) &
276                     o3_in(j, k:n_plev, l, m) = o3_in(j, k-1, l, m)
277             end do
278          end do
279       else
280          print *, "regr_lat_time_climoz: field ", m, &
281               ", no missing value attribute"
282       end if
283    end do
284
285    call nf95_close(ncid_in)
286
287    allocate(o3_regr_lat(nbp_lat, n_plev, 0:13, read_climoz))
288    allocate(o3_out(nbp_lat, n_plev, 360, read_climoz))
289
290    ! Regrid in latitude:
291    ! We average with respect to sine of latitude, which is
292    ! equivalent to weighting by cosine of latitude:
293    if (n_month == 12) then
294       print *, &
295            "Found 12 months in ozone climatologies, assuming periodicity..."
296       call regr1_conserv(o3_in, xs = sin_lat_in_edg, &
297            xt = (/- 1., sin(boundslat_reg(nbp_lat - 1:1:- 1, south)), 1./), &
298            vt = o3_regr_lat(nbp_lat:1:- 1, :, 1:12, :), &
299            slope = slopes(o3_in, sin_lat_in_edg))
300       ! (invert order of indices in "o3_regr_lat" because "rlatu" is
301       ! in descending order)
302
303       ! Duplicate January and December values, in preparation of time
304       ! interpolation:
305       o3_regr_lat(:, :, 0, :) = o3_regr_lat(:, :, 12, :)
306       o3_regr_lat(:, :, 13, :) = o3_regr_lat(:, :, 1, :)
307    else
308       print *, "Using 14 months in ozone climatologies..."
309       call regr1_conserv(o3_in, xs = sin_lat_in_edg, &
310            xt = (/- 1., sin(boundslat_reg(nbp_lat - 1:1:- 1, south)), 1./), &
311            vt = o3_regr_lat(nbp_lat:1:- 1, :, :, :), &
312            slope = slopes(o3_in, sin_lat_in_edg))
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(nbp_lat: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 mod_grid_phy_lmdz, ONLY : nbp_lat
339    use netcdf95, only: nf95_create, nf95_def_dim, nf95_def_var, &
340         nf95_put_att, nf95_enddef, nf95_copy_att, nf95_put_var
341    use netcdf, only: nf90_clobber, nf90_float, nf90_global
342    use nrtype, only: pi
343    use regular_lonlat_mod, only : lat_reg
344
345    integer, intent(in):: ncid_in, n_plev
346    integer, intent(out):: ncid_out, varid_plev, varid_time
347
348    integer, intent(out):: varid_out(:) ! dim(1 or 2)
349    ! "varid_out(1)" is for average ozone day and night,
350    ! "varid_out(2)" is for daylight ozone.
351
352    ! Variables local to the procedure:
353
354    integer ncerr
355    integer dimid_rlatu, dimid_plev, dimid_time
356    integer varid_rlatu
357
358    !---------------------------
359
360    print *, "Call sequence information: prepare_out"
361
362    call nf95_create("climoz_LMDZ.nc", nf90_clobber, ncid_out)
363
364    ! Dimensions:
365    call nf95_def_dim(ncid_out, "time", 360, dimid_time)
366    call nf95_def_dim(ncid_out, "plev", n_plev, dimid_plev)
367    call nf95_def_dim(ncid_out, "rlatu", nbp_lat, dimid_rlatu)
368
369    ! Define coordinate variables:
370
371    call nf95_def_var(ncid_out, "time", nf90_float, dimid_time, varid_time)
372    call nf95_put_att(ncid_out, varid_time, "units", "days since 2000-1-1")
373    call nf95_put_att(ncid_out, varid_time, "calendar", "360_day")
374    call nf95_put_att(ncid_out, varid_time, "standard_name", "time")
375
376    call nf95_def_var(ncid_out, "plev", nf90_float, dimid_plev, varid_plev)
377    call nf95_put_att(ncid_out, varid_plev, "units", "millibar")
378    call nf95_put_att(ncid_out, varid_plev, "standard_name", "air_pressure")
379    call nf95_put_att(ncid_out, varid_plev, "long_name", "air pressure")
380
381    call nf95_def_var(ncid_out, "rlatu", nf90_float, dimid_rlatu, varid_rlatu)
382    call nf95_put_att(ncid_out, varid_rlatu, "units", "degrees_north")
383    call nf95_put_att(ncid_out, varid_rlatu, "standard_name", "latitude")
384
385    ! Define the primary variables:
386
387    call nf95_def_var(ncid_out, "tro3", nf90_float, &
388         (/dimid_rlatu, dimid_plev, dimid_time/), varid_out(1))
389    call nf95_put_att(ncid_out, varid_out(1), "long_name", &
390         "ozone mole fraction")
391    call nf95_put_att(ncid_out, varid_out(1), "standard_name", &
392         "mole_fraction_of_ozone_in_air")
393
394    if (size(varid_out) == 2) then
395       call nf95_def_var(ncid_out, "tro3_daylight", nf90_float, &
396            (/dimid_rlatu, dimid_plev, dimid_time/), varid_out(2))
397       call nf95_put_att(ncid_out, varid_out(2), "long_name", &
398            "ozone mole fraction in daylight")
399    end if
400
401    ! Global attributes:
402
403    ! The following commands, copying attributes, may fail.
404    ! That is OK.
405    ! It should just mean that the attribute is not defined in the input file.
406
407    call nf95_copy_att(ncid_in, nf90_global, "Conventions", ncid_out, &
408         nf90_global, ncerr)
409    call handle_err_copy_att("Conventions")
410
411    call nf95_copy_att(ncid_in, nf90_global, "title", ncid_out, nf90_global, &
412         ncerr)
413    call handle_err_copy_att("title")
414
415    call nf95_copy_att(ncid_in, nf90_global, "institution", ncid_out, &
416         nf90_global, ncerr)
417    call handle_err_copy_att("institution")
418
419    call nf95_copy_att(ncid_in, nf90_global, "source", ncid_out, nf90_global, &
420         ncerr)
421    call handle_err_copy_att("source")
422
423    call nf95_put_att(ncid_out, nf90_global, "comment", "Regridded for LMDZ")
424
425    call nf95_enddef(ncid_out)
426
427    ! Write one of the coordinate variables:
428    call nf95_put_var(ncid_out, varid_rlatu, lat_reg(nbp_lat:1:-1) / pi * 180.)
429    ! (convert from rad to degrees and sort in ascending order)
430
431  contains
432
433    subroutine handle_err_copy_att(att_name)
434
435      use netcdf, only: nf90_noerr, nf90_strerror
436
437      character(len=*), intent(in):: att_name
438
439      !----------------------------------------
440
441      if (ncerr /= nf90_noerr) then
442         print *, "regr_lat_time_climoz_m prepare_out nf95_copy_att " &
443              // att_name // " -- " // trim(nf90_strerror(ncerr))
444      end if
445
446    end subroutine handle_err_copy_att
447
448  end subroutine prepare_out
449
450end module regr_lat_time_climoz_m
Note: See TracBrowser for help on using the repository browser.