source: LMDZ5/trunk/libf/phylmd/regr_lat_time_coefoz_m.F90 @ 1907

Last change on this file since 1907 was 1907, checked in by lguez, 10 years ago

Added a copyright property to every file of the distribution, except
for the fcm files (which have their own copyright). Use svn propget on
a file to see the copyright. For instance:

$ svn propget copyright libf/phylmd/physiq.F90
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

Also added the files defining the CeCILL version 2 license, in French
and English, at the top of the LMDZ tree.

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