orog_mask_tools  1.13.0
netcdf_io.F90
Go to the documentation of this file.
1 
4 
21  subroutine write_netcdf(im, jm, slm, land_frac, oro, orf, hprime, ntiles, tile, geolon, geolat, lon, lat)
22  implicit none
23  integer, intent(in):: im, jm, ntiles, tile
24  real, intent(in) :: lon(im), lat(jm)
25  real, intent(in), dimension(im,jm) :: slm, oro, orf, geolon, geolat, land_frac
26  real, intent(in), dimension(im,jm,14):: hprime
27  character(len=128) :: outfile
28  integer :: error, ncid
29  integer :: header_buffer_val = 16384
30  integer :: fsize=65536, inital = 0
31  integer :: dim1, dim2
32  integer :: dim_lon, dim_lat
33  integer :: id_geolon,id_geolat
34  integer :: id_slmsk,id_orog_raw,id_orog_filt,id_land_frac
35  integer :: id_stddev,id_convex
36  integer :: id_oa1,id_oa2,id_oa3,id_oa4
37  integer :: id_ol1,id_ol2,id_ol3,id_ol4
38  integer :: id_theta,id_gamma,id_sigma,id_elvmax
39  include "netcdf.inc"
40 
41  if(ntiles > 1) then
42  write(outfile, '(a,i4.4,a)') 'out.oro.tile', tile, '.nc'
43  else
44  outfile = "out.oro.nc"
45  endif
46 
47  dim1=size(lon,1)
48  dim2=size(lat,1)
49  write(6,*) ' netcdf dims are: ',dim1, dim2
50 
51  !--- open the file
52  error = nf__create(outfile, ior(nf_netcdf4,nf_classic_model), inital, fsize, ncid)
53  call netcdf_err(error, 'Creating file '//trim(outfile) )
54  !--- define dimension
55  error = nf_def_dim(ncid, 'lon', dim1, dim_lon)
56  call netcdf_err(error, 'define dimension lon for file='//trim(outfile) )
57  error = nf_def_dim(ncid, 'lat', dim2, dim_lat)
58  call netcdf_err(error, 'define dimension lat for file='//trim(outfile) )
59 
60  !--- define field
61 !---geolon
62  error = nf_def_var(ncid, 'geolon', nf_float, 2, (/dim_lon,dim_lat/), id_geolon)
63  call netcdf_err(error, 'define var geolon for file='//trim(outfile) )
64  error = nf_put_att_text(ncid, id_geolon, "long_name", 9, "Longitude")
65  call netcdf_err(error, 'define geolon name for file='//trim(outfile) )
66  error = nf_put_att_text(ncid, id_geolon, "units", 12, "degrees_east")
67  call netcdf_err(error, 'define geolon units for file='//trim(outfile) )
68 !---geolat
69  error = nf_def_var(ncid, 'geolat', nf_float, 2, (/dim_lon,dim_lat/), id_geolat)
70  call netcdf_err(error, 'define var geolat for file='//trim(outfile) )
71  error = nf_put_att_text(ncid, id_geolat, "long_name", 8, "Latitude")
72  call netcdf_err(error, 'define geolat name for file='//trim(outfile) )
73  error = nf_put_att_text(ncid, id_geolat, "units", 13, "degrees_north")
74  call netcdf_err(error, 'define geolat units for file='//trim(outfile) )
75 !---slmsk
76  error = nf_def_var(ncid, 'slmsk', nf_float, 2, (/dim_lon,dim_lat/), id_slmsk)
77  call netcdf_err(error, 'define var slmsk for file='//trim(outfile) )
78  error = nf_put_att_text(ncid, id_slmsk, "coordinates", 13, "geolon geolat")
79  call netcdf_err(error, 'define slmsk coordinates for file='//trim(outfile) )
80 !--- land_frac
81  error = nf_def_var(ncid, 'land_frac', nf_float, 2, (/dim_lon,dim_lat/), id_land_frac)
82  call netcdf_err(error, 'define var land_frac for file='//trim(outfile) )
83  error = nf_put_att_text(ncid, id_land_frac, "coordinates", 13, "geolon geolat")
84  call netcdf_err(error, 'define land_frac coordinates for file='//trim(outfile) )
85 !---orography - raw
86  error = nf_def_var(ncid, 'orog_raw', nf_float, 2, (/dim_lon,dim_lat/), id_orog_raw)
87  call netcdf_err(error, 'define var orog_raw for file='//trim(outfile) )
88  error = nf_put_att_text(ncid, id_orog_raw, "coordinates", 13, "geolon geolat")
89  call netcdf_err(error, 'define orog_raw coordinates for file='//trim(outfile) )
90 !---orography - filtered
91  error = nf_def_var(ncid, 'orog_filt', nf_float, 2, (/dim_lon,dim_lat/), id_orog_filt)
92  call netcdf_err(error, 'define var orog_filt for file='//trim(outfile) )
93  error = nf_put_att_text(ncid, id_orog_filt, "coordinates", 13, "geolon geolat")
94  call netcdf_err(error, 'define orog_filt coordinates for file='//trim(outfile) )
95 !---stddev
96  error = nf_def_var(ncid, 'stddev', nf_float, 2, (/dim_lon,dim_lat/), id_stddev)
97  call netcdf_err(error, 'define var stddev for file='//trim(outfile) )
98  error = nf_put_att_text(ncid, id_stddev, "coordinates", 13, "geolon geolat")
99  call netcdf_err(error, 'define stddev coordinates for file='//trim(outfile) )
100 !---convexity
101  error = nf_def_var(ncid, 'convexity', nf_float, 2, (/dim_lon,dim_lat/), id_convex)
102  call netcdf_err(error, 'define var convexity for file='//trim(outfile) )
103  error = nf_put_att_text(ncid, id_convex, "coordinates", 13, "geolon geolat")
104  call netcdf_err(error, 'define convexity coordinates for file='//trim(outfile) )
105 !---oa1 -> oa4
106  error = nf_def_var(ncid, 'oa1', nf_float, 2, (/dim_lon,dim_lat/), id_oa1)
107  call netcdf_err(error, 'define var oa1 for file='//trim(outfile) )
108  error = nf_put_att_text(ncid, id_oa1, "coordinates", 13, "geolon geolat")
109  call netcdf_err(error, 'define oa1 coordinates for file='//trim(outfile) )
110  error = nf_def_var(ncid, 'oa2', nf_float, 2, (/dim_lon,dim_lat/), id_oa2)
111  call netcdf_err(error, 'define var oa2 for file='//trim(outfile) )
112  error = nf_put_att_text(ncid, id_oa2, "coordinates", 13, "geolon geolat")
113  call netcdf_err(error, 'define oa2 coordinates for file='//trim(outfile) )
114  error = nf_def_var(ncid, 'oa3', nf_float, 2, (/dim_lon,dim_lat/), id_oa3)
115  call netcdf_err(error, 'define var oa3 for file='//trim(outfile) )
116  error = nf_put_att_text(ncid, id_oa3, "coordinates", 13, "geolon geolat")
117  call netcdf_err(error, 'define oa3 coordinates for file='//trim(outfile) )
118  error = nf_def_var(ncid, 'oa4', nf_float, 2, (/dim_lon,dim_lat/), id_oa4)
119  call netcdf_err(error, 'define var oa4 for file='//trim(outfile) )
120  error = nf_put_att_text(ncid, id_oa4, "coordinates", 13, "geolon geolat")
121  call netcdf_err(error, 'define oa4 coordinates for file='//trim(outfile) )
122 !---ol1 -> ol4
123  error = nf_def_var(ncid, 'ol1', nf_float, 2, (/dim_lon,dim_lat/), id_ol1)
124  call netcdf_err(error, 'define var ol1 for file='//trim(outfile) )
125  error = nf_put_att_text(ncid, id_ol1, "coordinates", 13, "geolon geolat")
126  call netcdf_err(error, 'define ol1 coordinates for file='//trim(outfile) )
127  error = nf_def_var(ncid, 'ol2', nf_float, 2, (/dim_lon,dim_lat/), id_ol2)
128  call netcdf_err(error, 'define var ol2 for file='//trim(outfile) )
129  error = nf_put_att_text(ncid, id_ol2, "coordinates", 13, "geolon geolat")
130  call netcdf_err(error, 'define ol2 coordinates for file='//trim(outfile) )
131  error = nf_def_var(ncid, 'ol3', nf_float, 2, (/dim_lon,dim_lat/), id_ol3)
132  call netcdf_err(error, 'define var ol3 for file='//trim(outfile) )
133  error = nf_put_att_text(ncid, id_ol3, "coordinates", 13, "geolon geolat")
134  call netcdf_err(error, 'define ol3 coordinates for file='//trim(outfile) )
135  error = nf_def_var(ncid, 'ol4', nf_float, 2, (/dim_lon,dim_lat/), id_ol4)
136  call netcdf_err(error, 'define var ol4 for file='//trim(outfile) )
137  error = nf_put_att_text(ncid, id_ol4, "coordinates", 13, "geolon geolat")
138  call netcdf_err(error, 'define ol4 coordinates for file='//trim(outfile) )
139 !---theta gamma sigma elvmax
140  error = nf_def_var(ncid, 'theta', nf_float, 2, (/dim_lon,dim_lat/), id_theta)
141  call netcdf_err(error, 'define var theta for file='//trim(outfile) )
142  error = nf_put_att_text(ncid, id_theta, "coordinates", 13, "geolon geolat")
143  call netcdf_err(error, 'define theta coordinates for file='//trim(outfile) )
144  error = nf_def_var(ncid, 'gamma', nf_float, 2, (/dim_lon,dim_lat/), id_gamma)
145  call netcdf_err(error, 'define var gamma for file='//trim(outfile) )
146  error = nf_put_att_text(ncid, id_gamma, "coordinates", 13, "geolon geolat")
147  call netcdf_err(error, 'define gamma coordinates for file='//trim(outfile) )
148  error = nf_def_var(ncid, 'sigma', nf_float, 2, (/dim_lon,dim_lat/), id_sigma)
149  call netcdf_err(error, 'define var sigma for file='//trim(outfile) )
150  error = nf_put_att_text(ncid, id_sigma, "coordinates", 13, "geolon geolat")
151  call netcdf_err(error, 'define sigma coordinates for file='//trim(outfile) )
152  error = nf_def_var(ncid, 'elvmax', nf_float, 2, (/dim_lon,dim_lat/), id_elvmax)
153  call netcdf_err(error, 'define var elvmax for file='//trim(outfile) )
154  error = nf_put_att_text(ncid, id_elvmax, "coordinates", 13, "geolon geolat")
155  call netcdf_err(error, 'define elvmax coordinates for file='//trim(outfile) )
156 
157  error = nf__enddef(ncid, header_buffer_val,4,0,4)
158  call netcdf_err(error, 'end meta define for file='//trim(outfile) )
159 
160  !--- write out data
161  error = nf_put_var_double( ncid, id_geolon, geolon(:dim1,:dim2))
162  call netcdf_err(error, 'write var geolon for file='//trim(outfile) )
163  error = nf_put_var_double( ncid, id_geolat, geolat(:dim1,:dim2))
164  call netcdf_err(error, 'write var geolat for file='//trim(outfile) )
165 
166  error = nf_put_var_double( ncid, id_slmsk, slm(:dim1,:dim2))
167  call netcdf_err(error, 'write var slmsk for file='//trim(outfile) )
168  error = nf_put_var_double( ncid, id_land_frac, land_frac(:dim1,:dim2))
169  call netcdf_err(error, 'write var land_frac for file='//trim(outfile) )
170 
171  error = nf_put_var_double( ncid, id_orog_raw, oro(:dim1,:dim2))
172  call netcdf_err(error, 'write var orog_raw for file='//trim(outfile) )
173  error = nf_put_var_double( ncid, id_orog_filt, orf(:dim1,:dim2))
174  call netcdf_err(error, 'write var orog_filt for file='//trim(outfile) )
175 
176  error = nf_put_var_double( ncid, id_stddev, hprime(:dim1,:dim2,1))
177  call netcdf_err(error, 'write var stddev for file='//trim(outfile) )
178  error = nf_put_var_double( ncid, id_convex, hprime(:dim1,:dim2,2))
179  call netcdf_err(error, 'write var convex for file='//trim(outfile) )
180 
181  error = nf_put_var_double( ncid, id_oa1, hprime(:dim1,:dim2,3))
182  call netcdf_err(error, 'write var oa1 for file='//trim(outfile) )
183  error = nf_put_var_double( ncid, id_oa2, hprime(:dim1,:dim2,4))
184  call netcdf_err(error, 'write var oa2 for file='//trim(outfile) )
185  error = nf_put_var_double( ncid, id_oa3, hprime(:dim1,:dim2,5))
186  call netcdf_err(error, 'write var oa3 for file='//trim(outfile) )
187  error = nf_put_var_double( ncid, id_oa4, hprime(:dim1,:dim2,6))
188  call netcdf_err(error, 'write var oa4 for file='//trim(outfile) )
189 
190  error = nf_put_var_double( ncid, id_ol1, hprime(:dim1,:dim2,7))
191  call netcdf_err(error, 'write var ol1 for file='//trim(outfile) )
192  error = nf_put_var_double( ncid, id_ol2, hprime(:dim1,:dim2,8))
193  call netcdf_err(error, 'write var ol2 for file='//trim(outfile) )
194  error = nf_put_var_double( ncid, id_ol3, hprime(:dim1,:dim2,9))
195  call netcdf_err(error, 'write var ol3 for file='//trim(outfile) )
196  error = nf_put_var_double( ncid, id_ol4, hprime(:dim1,:dim2,10))
197  call netcdf_err(error, 'write var ol4 for file='//trim(outfile) )
198 
199  error = nf_put_var_double( ncid, id_theta, hprime(:dim1,:dim2,11))
200  call netcdf_err(error, 'write var theta for file='//trim(outfile) )
201  error = nf_put_var_double( ncid, id_gamma, hprime(:dim1,:dim2,12))
202  call netcdf_err(error, 'write var gamma for file='//trim(outfile) )
203  error = nf_put_var_double( ncid, id_sigma, hprime(:dim1,:dim2,13))
204  call netcdf_err(error, 'write var sigma for file='//trim(outfile) )
205  error = nf_put_var_double( ncid, id_elvmax, hprime(:dim1,:dim2,14))
206  call netcdf_err(error, 'write var elvmax for file='//trim(outfile) )
207 
208  error = nf_close(ncid)
209  call netcdf_err(error, 'close file='//trim(outfile) )
210 
211  end subroutine
212 
218  subroutine netcdf_err( err, string )
219  integer, intent(in) :: err
220  character(len=*), intent(in) :: string
221  character(len=256) :: errmsg
222  include "netcdf.inc"
223 
224  if( err.EQ.nf_noerr )return
225  errmsg = nf_strerror(err)
226  print*, 'FATAL ERROR: ', trim(string), ': ', trim(errmsg)
227  call abort
228 
229  return
230  end subroutine netcdf_err
231 
243 
244  subroutine write_mask_netcdf(im, jm, slm, land_frac, ntiles, tile, geolon, geolat)
245  implicit none
246  integer, intent(in):: im, jm, ntiles, tile
247  real, intent(in), dimension(im,jm) :: slm, geolon, geolat, land_frac
248  character(len=128) :: outfile
249  integer :: error, ncid
250  integer :: header_buffer_val = 16384
251  integer :: fsize=65536, inital = 0
252  integer :: dim1, dim2
253  integer :: dim_lon, dim_lat
254  integer :: id_geolon,id_geolat
255  integer :: id_slmsk,id_land_frac
256  include "netcdf.inc"
257 
258  if(ntiles > 1) then
259  write(outfile, '(a,i4.4,a)') 'out.oro.tile', tile, '.nc'
260  else
261  outfile = "out.oro.nc"
262  endif
263 
264  dim1=im
265  dim2=jm
266  write(6,*) ' netcdf dims are: ',dim1, dim2
267 
268  !--- open the file
269  error = nf__create(outfile, ior(nf_netcdf4,nf_classic_model), inital, fsize, ncid)
270  call netcdf_err(error, 'Creating file '//trim(outfile) )
271  !--- define dimension
272  error = nf_def_dim(ncid, 'lon', dim1, dim_lon)
273  call netcdf_err(error, 'define dimension lon for file='//trim(outfile) )
274  error = nf_def_dim(ncid, 'lat', dim2, dim_lat)
275  call netcdf_err(error, 'define dimension lat for file='//trim(outfile) )
276 
277  !--- define field
278 !---geolon
279  error = nf_def_var(ncid, 'geolon', nf_float, 2, (/dim_lon,dim_lat/), id_geolon)
280  call netcdf_err(error, 'define var geolon for file='//trim(outfile) )
281  error = nf_put_att_text(ncid, id_geolon, "long_name", 9, "Longitude")
282  call netcdf_err(error, 'define geolon name for file='//trim(outfile) )
283  error = nf_put_att_text(ncid, id_geolon, "units", 12, "degrees_east")
284  call netcdf_err(error, 'define geolon units for file='//trim(outfile) )
285 !---geolat
286  error = nf_def_var(ncid, 'geolat', nf_float, 2, (/dim_lon,dim_lat/), id_geolat)
287  call netcdf_err(error, 'define var geolat for file='//trim(outfile) )
288  error = nf_put_att_text(ncid, id_geolat, "long_name", 8, "Latitude")
289  call netcdf_err(error, 'define geolat name for file='//trim(outfile) )
290  error = nf_put_att_text(ncid, id_geolat, "units", 13, "degrees_north")
291  call netcdf_err(error, 'define geolat units for file='//trim(outfile) )
292 !---slmsk
293  error = nf_def_var(ncid, 'slmsk', nf_float, 2, (/dim_lon,dim_lat/), id_slmsk)
294  call netcdf_err(error, 'define var slmsk for file='//trim(outfile) )
295  error = nf_put_att_text(ncid, id_slmsk, "coordinates", 13, "geolon geolat")
296  call netcdf_err(error, 'define slmsk coordinates for file='//trim(outfile) )
297 !--- land_frac
298  error = nf_def_var(ncid, 'land_frac', nf_float, 2, (/dim_lon,dim_lat/), id_land_frac)
299  call netcdf_err(error, 'define var land_frac for file='//trim(outfile) )
300  error = nf_put_att_text(ncid, id_land_frac, "coordinates", 13, "geolon geolat")
301  call netcdf_err(error, 'define land_frac coordinates for file='//trim(outfile) )
302 
303  error = nf__enddef(ncid, header_buffer_val,4,0,4)
304  call netcdf_err(error, 'end meta define for file='//trim(outfile) )
305 
306  !--- write out data
307  error = nf_put_var_double( ncid, id_geolon, geolon(:dim1,:dim2))
308  call netcdf_err(error, 'write var geolon for file='//trim(outfile) )
309 
310  error = nf_put_var_double( ncid, id_geolat, geolat(:dim1,:dim2))
311  call netcdf_err(error, 'write var geolat for file='//trim(outfile) )
312 
313  error = nf_put_var_double( ncid, id_slmsk, slm(:dim1,:dim2))
314  call netcdf_err(error, 'write var slmsk for file='//trim(outfile) )
315 
316  error = nf_put_var_double( ncid, id_land_frac, land_frac(:dim1,:dim2))
317  call netcdf_err(error, 'write var land_frac for file='//trim(outfile) )
318 
319  error = nf_close(ncid)
320  call netcdf_err(error, 'close file='//trim(outfile) )
321 
322  end subroutine
323 
324 
334 
335  subroutine read_mask(merge_file,slm,land_frac,lake_frac,im,jm)
337  implicit none
338  include "netcdf.inc"
339 
340  character(len=*), intent(in) :: merge_file
341 
342  integer, intent(in) :: im, jm
343 
344  real, intent(out) :: land_frac(im,jm)
345  real, intent(out) :: lake_frac(im,jm)
346  real, intent(out) :: slm(im,jm)
347 
348  integer ncid, error, fsize, id_var
349 
350  fsize = 66536
351 
352  print*, "merge_file=", trim(merge_file)
353  error=nf__open(merge_file,nf_nowrite,fsize,ncid)
354  call netcdf_err(error, 'Open file '//trim(merge_file) )
355 
356  error=nf_inq_varid(ncid, 'land_frac', id_var)
357  call netcdf_err(error, 'inquire varid of land_frac')
358  error=nf_get_var_double(ncid, id_var, land_frac)
359  call netcdf_err(error, 'inquire data of land_frac')
360 
361  print*,'land_frac ',maxval(land_frac),minval(land_frac)
362 
363  error=nf_inq_varid(ncid, 'slmsk', id_var)
364  call netcdf_err(error, 'inquire varid of slmsk')
365  error=nf_get_var_double(ncid, id_var, slm)
366  call netcdf_err(error, 'inquire data of slmsk')
367 
368  print*,'slmsk ',maxval(slm),minval(slm)
369 
370  error=nf_inq_varid(ncid, 'lake_frac', id_var)
371  call netcdf_err(error, 'inquire varid of lake_frac')
372  error=nf_get_var_double(ncid, id_var, lake_frac)
373  call netcdf_err(error, 'inquire data of lake_frac')
374 
375  print*,'lake_frac ',maxval(lake_frac),minval(lake_frac)
376 
377  error = nf_close(ncid)
378  print*,'bot of read_mask'
379 
380  end subroutine
subroutine write_mask_netcdf(im, jm, slm, land_frac, ntiles, tile, geolon, geolat)
Write the land mask file.
Definition: netcdf_io.F90:245
subroutine read_mask(merge_file, slm, land_frac, lake_frac, im, jm)
Read the land mask file.
Definition: netcdf_io.F90:336
subroutine netcdf_err(err, string)
Check NetCDF error code and output the error message.
Definition: netcdf_io.F90:219
subroutine write_netcdf(im, jm, slm, land_frac, oro, orf, hprime, ntiles, tile, geolon, geolat, lon, lat)
Write out orography file in netcdf format.
Definition: netcdf_io.F90:22