/[lmdze]/trunk/libf/phylmd/Mobidic/regr_lat_time_coefoz.f90
ViewVC logotype

Annotation of /trunk/libf/phylmd/Mobidic/regr_lat_time_coefoz.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 47 - (hide annotations)
Fri Jul 1 15:00:48 2011 UTC (12 years, 11 months ago) by guez
File size: 11939 byte(s)
Split "thermcell.f" and "cv3_routines.f".
Removed copies of files that are now in "L_util".
Moved "mva9" and "diagetpq" to their own files.
Unified variable names across procedures.

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

  ViewVC Help
Powered by ViewVC 1.1.21