source: LMDZ4/branches/LMDZ4-dev/libf/bibio/simple.F90 @ 1264

Last change on this file since 1264 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: 8.2 KB
Line 
1! $Id$
2module simple
3
4  implicit none
5
6contains
7
8  subroutine nf95_open(path, mode, ncid, chunksize, ncerr)
9
10    use netcdf, only: nf90_open
11    use handle_err_m, only: handle_err
12
13    character(len=*), intent(in):: path
14    integer, intent(in):: mode
15    integer, intent(out):: ncid
16    integer, intent(inout), optional:: chunksize
17    integer, intent(out), optional:: ncerr
18
19    ! Variable local to the procedure:
20    integer ncerr_not_opt
21
22    !-------------------
23
24    ncerr_not_opt = nf90_open(path, mode, ncid, chunksize)
25    if (present(ncerr)) then
26       ncerr = ncerr_not_opt
27    else
28       call handle_err("nf95_open " // path, ncerr_not_opt)
29    end if
30
31  end subroutine nf95_open
32
33  !************************
34
35  subroutine nf95_inq_dimid(ncid, name, dimid, ncerr)
36
37    use netcdf, only: nf90_inq_dimid
38    use handle_err_m, only: handle_err
39
40    integer,             intent( in) :: ncid
41    character (len = *), intent( in) :: name
42    integer,             intent(out) :: dimid
43    integer, intent(out), optional:: ncerr
44
45    ! Variable local to the procedure:
46    integer ncerr_not_opt
47
48    !-------------------
49
50    ncerr_not_opt = nf90_inq_dimid(ncid, name, dimid)
51    if (present(ncerr)) then
52       ncerr = ncerr_not_opt
53    else
54       call handle_err("nf95_inq_dimid", ncerr_not_opt, ncid)
55    end if
56
57  end subroutine nf95_inq_dimid
58
59  !************************
60
61  subroutine nf95_inquire_dimension(ncid, dimid, name, len, ncerr)
62
63    use netcdf, only: nf90_inquire_dimension
64    use handle_err_m, only: handle_err
65
66    integer,                       intent( in) :: ncid, dimid
67    character (len = *), optional, intent(out) :: name
68    integer,             optional, intent(out) :: len
69    integer, intent(out), optional:: ncerr
70
71    ! Variable local to the procedure:
72    integer ncerr_not_opt
73
74    !-------------------
75
76    ncerr_not_opt = nf90_inquire_dimension(ncid, dimid, name, len)
77    if (present(ncerr)) then
78       ncerr = ncerr_not_opt
79    else
80       call handle_err("nf95_inquire_dimension", ncerr_not_opt, ncid)
81    end if
82
83  end subroutine nf95_inquire_dimension
84
85  !************************
86
87  subroutine nf95_inq_varid(ncid, name, varid, ncerr)
88
89    use netcdf, only: nf90_inq_varid
90    use handle_err_m, only: handle_err
91
92    integer,             intent(in) :: ncid
93    character (len = *), intent(in) :: name
94    integer,             intent(out) :: varid
95    integer, intent(out), optional:: ncerr
96
97    ! Variable local to the procedure:
98    integer ncerr_not_opt
99
100    !-------------------
101
102    ncerr_not_opt = nf90_inq_varid(ncid, name, varid)
103    if (present(ncerr)) then
104       ncerr = ncerr_not_opt
105    else
106       call handle_err("nf95_inq_varid, name = " // name, ncerr_not_opt, ncid)
107    end if
108
109  end subroutine nf95_inq_varid
110
111  !************************
112
113  subroutine nf95_inquire_variable(ncid, varid, name, xtype, ndims, dimids, &
114       nAtts, ncerr)
115
116    ! In "nf90_inquire_variable", "dimids" is an assumed-size array.
117    ! This is the classical case of an array the size of which is
118    ! unknown in the calling procedure, before the call.
119    ! Here we use a better solution: a pointer argument array.
120    ! This procedure associates and defines "dimids" if it is present.
121
122    use netcdf, only: nf90_inquire_variable, nf90_max_var_dims
123    use handle_err_m, only: handle_err
124
125    integer, intent(in):: ncid, varid
126    character(len = *), optional, intent(out):: name
127    integer, optional, intent(out) :: xtype, ndims
128    integer, dimension(:), optional, pointer :: dimids
129    integer, optional, intent(out) :: nAtts
130    integer, intent(out), optional :: ncerr
131
132    ! Variable local to the procedure:
133    integer ncerr_not_opt
134    integer dimids_local(nf90_max_var_dims)
135    integer ndims_not_opt
136
137    !-------------------
138
139    if (present(dimids)) then
140       ncerr_not_opt = nf90_inquire_variable(ncid, varid, name, xtype, &
141            ndims_not_opt, dimids_local, nAtts)
142       allocate(dimids(ndims_not_opt)) ! also works if ndims_not_opt == 0
143       dimids = dimids_local(:ndims_not_opt)
144       if (present(ndims)) ndims = ndims_not_opt
145    else
146       ncerr_not_opt = nf90_inquire_variable(ncid, varid, name, xtype, ndims, &
147            nAtts=nAtts)
148    end if
149
150    if (present(ncerr)) then
151       ncerr = ncerr_not_opt
152    else
153       call handle_err("nf95_inquire_variable", ncerr_not_opt, ncid)
154    end if
155
156  end subroutine nf95_inquire_variable
157
158  !************************
159
160  subroutine nf95_create(path, cmode, ncid, initialsize, chunksize, ncerr)
161   
162    use netcdf, only: nf90_create
163    use handle_err_m, only: handle_err
164
165    character (len = *), intent(in   ) :: path
166    integer,             intent(in   ) :: cmode
167    integer,             intent(  out) :: ncid
168    integer, optional,   intent(in   ) :: initialsize
169    integer, optional,   intent(inout) :: chunksize
170    integer, intent(out), optional :: ncerr
171
172    ! Variable local to the procedure:
173    integer ncerr_not_opt
174
175    !-------------------
176
177    ncerr_not_opt = nf90_create(path, cmode, ncid, initialsize, chunksize)
178    if (present(ncerr)) then
179       ncerr = ncerr_not_opt
180    else
181       call handle_err("nf95_create " // path, ncerr_not_opt)
182    end if
183
184  end subroutine nf95_create
185
186  !************************
187
188  subroutine nf95_def_dim(ncid, name, len, dimid, ncerr)
189
190    use netcdf, only: nf90_def_dim
191    use handle_err_m, only: handle_err
192
193    integer,             intent( in) :: ncid
194    character (len = *), intent( in) :: name
195    integer,             intent( in) :: len
196    integer,             intent(out) :: dimid
197    integer, intent(out), optional :: ncerr
198
199    ! Variable local to the procedure:
200    integer ncerr_not_opt
201
202    !-------------------
203
204    ncerr_not_opt = nf90_def_dim(ncid, name, len, dimid)
205    if (present(ncerr)) then
206       ncerr = ncerr_not_opt
207    else
208       call handle_err("nf95_def_dim", ncerr_not_opt, ncid)
209    end if
210
211  end subroutine nf95_def_dim
212
213  !***********************
214
215  subroutine nf95_redef(ncid, ncerr)
216
217    use netcdf, only: nf90_redef
218    use handle_err_m, only: handle_err
219
220    integer, intent( in) :: ncid
221    integer, intent(out), optional :: ncerr
222
223    ! Variable local to the procedure:
224    integer ncerr_not_opt
225
226    !-------------------
227
228    ncerr_not_opt = nf90_redef(ncid)
229    if (present(ncerr)) then
230       ncerr = ncerr_not_opt
231    else
232       call handle_err("nf95_redef", ncerr_not_opt, ncid)
233    end if
234
235  end subroutine nf95_redef
236 
237  !***********************
238
239  subroutine nf95_enddef(ncid, h_minfree, v_align, v_minfree, r_align, ncerr)
240
241    use netcdf, only: nf90_enddef
242    use handle_err_m, only: handle_err
243
244    integer,           intent( in) :: ncid
245    integer, optional, intent( in) :: h_minfree, v_align, v_minfree, r_align
246    integer, intent(out), optional :: ncerr
247
248    ! Variable local to the procedure:
249    integer ncerr_not_opt
250
251    !-------------------
252
253    ncerr_not_opt = nf90_enddef(ncid, h_minfree, v_align, v_minfree, r_align)
254    if (present(ncerr)) then
255       ncerr = ncerr_not_opt
256    else
257       call handle_err("nf95_enddef", ncerr_not_opt, ncid)
258    end if
259
260  end subroutine nf95_enddef
261
262  !***********************
263
264  subroutine nf95_close(ncid, ncerr)
265
266    use netcdf, only: nf90_close
267    use handle_err_m, only: handle_err
268
269    integer, intent( in) :: ncid
270    integer, intent(out), optional :: ncerr
271
272    ! Variable local to the procedure:
273    integer ncerr_not_opt
274
275    !-------------------
276
277    ncerr_not_opt = nf90_close(ncid)
278    if (present(ncerr)) then
279       ncerr = ncerr_not_opt
280    else
281       call handle_err("nf95_close", ncerr_not_opt)
282    end if
283
284  end subroutine nf95_close
285
286  !***********************
287
288  subroutine nf95_copy_att(ncid_in, varid_in, name, ncid_out, varid_out, ncerr)
289
290    use netcdf, only: nf90_copy_att
291    use handle_err_m, only: handle_err
292
293    integer, intent( in):: ncid_in,  varid_in
294    character(len=*), intent( in):: name
295    integer, intent( in):: ncid_out, varid_out
296    integer, intent(out), optional:: ncerr
297
298    ! Variable local to the procedure:
299    integer ncerr_not_opt
300
301    !-------------------
302
303    ncerr_not_opt = nf90_copy_att(ncid_in, varid_in, name, ncid_out, varid_out)
304    if (present(ncerr)) then
305       ncerr = ncerr_not_opt
306    else
307       call handle_err("nf95_copy_att", ncerr_not_opt, ncid_out)
308    end if
309
310  end subroutine nf95_copy_att
311
312end module simple
Note: See TracBrowser for help on using the repository browser.