Merge branch 'master' into devel
[wrffire.git] / standalone / util_netcdf.F
blob0cda1e8779fc22ecf54ea5ef884d8e1ebdec36cf
1 !*** A simple output module for wrf-fire standalone that produces output 
2 !    files and reads input files that are somewhat compatible with wrf files
5 module util_netcdf
6 use netcdf
7 implicit none
9 end type
11 contains
13 subroutine create_dim_attr(filename,grid)
15 ! Create an empty NetCDF file with proper dimension names 
16 ! (a.k.a. attributes) defined.
18 implicit none
20 !*** arguments
21 character(len=*),intent(in)::filename  ! filename to create
22 type(grid_info),intent(in)::grid       ! grid information structure
24 !*** local
25 integer::ncid,dimid,timeid,strlenid,varid
27 !*** executable
29 ! create an empty file
30 call check(nf90_create(filename,nf90_clobber,ncid))
32 ! define all dimensions
34 !!The function NF90_DEF_DIM adds a new dimension to an open netCDF dataset in 
35 !!define mode. It returns (as an argument) a dimension ID, given the netCDF ID,
36 !!the dimension name, and the dimension length. (From NetCDF docs)
38 call check(nf90_def_dim(ncid,dim_time,nf90_unlimited,dimid))
39 if(compat_fire_grid)then
40   call check(nf90_def_dim(ncid,dim_fire_x,(grid%natmx+1)*grid%sr_x,dimid))
41   call check(nf90_def_dim(ncid,dim_fire_y,(grid%natmy+1)*grid%sr_y,dimid))
42 else
43   call check(nf90_def_dim(ncid,dim_fire_x,grid%nfirex,dimid))
44   call check(nf90_def_dim(ncid,dim_fire_y,grid%nfirey,dimid))
45 endif
46 call check(nf90_def_dim(ncid,dim_atm_x,grid%natmx,dimid))
47 call check(nf90_def_dim(ncid,dim_atm_y,grid%natmy,dimid))
48 call check(nf90_def_dim(ncid,dim_atm_z,grid%natmz,dimid))
49 call check(nf90_def_dim(ncid,dim_atm_x_s,grid%natmx+1,dimid))
50 call check(nf90_def_dim(ncid,dim_atm_y_s,grid%natmy+1,dimid))
51 call check(nf90_def_dim(ncid,dim_atm_z_s,grid%natmz+1,dimid))
52 call check(nf90_def_dim(ncid,'DateStrLen',12,dimid))
54 ! define global attributes
55 call check(nf90_put_att(ncid,nf90_global,'DX',grid%dx))
56 call check(nf90_put_att(ncid,nf90_global,'DY',grid%dy))
57 call check(nf90_put_att(ncid,nf90_global,'DT',grid%dt))
58 call check(nf90_put_att(ncid,nf90_global,'STANDALONE_DRIVER',1))
60 call check(nf90_enddef(ncid))
61 call check(nf90_close(ncid))
63 end subroutine create_dim_attr
66 !***
69 subroutine create_output_file(filename,grid)
71 !*** Create a file containing meta-data suitable for output of this simulation.
72 !    Only creates dimensions, variables, and attributes... does not populate
73 !    data.  Squash file if it already exists.
75 implicit none
77 !*** arguments
78 character(len=*),intent(in)::filename  ! filename to create
79 type(grid_info),intent(in)::grid       ! grid information structure
81 !*** local
82 integer::ncid,dimid,timeid,strlenid,varid
84 !*** executable
86 call create_dim_attr(filename,grid)
88 call check(nf90_open(filename,nf90_write,ncid))
89 call check(nf90_inq_dimid(ncid,'DateStrLen',strlenid))
90 call check(nf90_inq_dimid(ncid,dim_time,timeid))
92 ! define a timekeeping variable
93 call check(nf90_redef(ncid))
94 call check(nf90_def_var(ncid,'Times',nf90_char,(/strlenid,timeid/),varid))
95 call check(nf90_enddef(ncid))
96 call check(nf90_close(ncid))
98 ! create all of the output variables
99 call define_var(filename,var_lfn,dim_fire,unit_lfn,desc_lfn)
100 call define_var(filename,var_tign,dim_fire,unit_tign,desc_tign)
101 call define_var(filename,var_vx,dim_fire,unit_vx,desc_vx)
102 call define_var(filename,var_vy,dim_fire,unit_vy,desc_vy)
103 call define_var(filename,var_grnhfx,dim_fire,unit_grnhfx,desc_grnhfx)
105 end subroutine create_output_file
108 !***
111 subroutine define_var(filename,varname,dims,units,description)
113 !*** define a variable in a netcdf data set, the file is assumed to exist and
114 !    have valid meta-data (as created by create_output_file)
116 implicit none
118 !*** arguments
119 character(len=*),intent(in)::filename,varname  ! create variable varname in filename
120 character(len=*),dimension(:),intent(in)::dims ! the dimension names of the variable
121 character(len=*),intent(in)::units,description ! attributes created by wrf (not used at the moment)
123 !*** local
124 integer::ncid,i,varid
125 integer,dimension(4)::dimids
126 character(len=*),parameter::memorder='XYZ'
127 character(len=3)::stag
129 !*** executable
130 ! open the file
131 call check(nf90_open(filename,nf90_write,ncid))
133 ! get dimension id's
134 do i=1,size(dims)
135   call check(nf90_inq_dimid(ncid,dims(i),dimids(i)))
136 enddo
138 ! enter define mode and define the variable
139 call check(nf90_redef(ncid))
140 call check(nf90_def_var(ncid,varname,vartype,dimids(1:size(dims)),varid))
142 ! add attributes
143 call check(nf90_put_att(ncid,varid,'FieldType',field_type))
144 call check(nf90_put_att(ncid,varid,'MemoryOrder',memorder(1:size(dims))))
145 call check(nf90_put_att(ncid,varid,'description',description))
146 call check(nf90_put_att(ncid,varid,'units',units))
147 if(size(dims).eq.3)then
148   stag='Z'
149 elseif(trim(dims(1)).eq.dim_atm_x_s)then
150   stag='X'
151 elseif(trim(dims(2)).eq.dim_atm_y_s)then
152   stag='Y'
153 elseif(trim(dims(3)).eq.dim_atm_z_s)then
154   stag='Z'
155 else
156   stag=''
157 endif
158 call check(nf90_put_att(ncid,varid,'stagger',stag))
159 call check(nf90_put_att(ncid,varid,'coordinates','XLONG XLAT'))
161 ! close file
162 call check(nf90_enddef(ncid))
163 call check(nf90_close(ncid))
164 end subroutine define_var
167 !***
170 subroutine write_inputs(filename,input,grid)
171 implicit none
173 ! write simulation input file in a way compatible 
174 ! with the wrf input file
176 !*** arguments
177 character(len=*),intent(in)::filename
178 type(input_vars),intent(in)::input
179 type(grid_info),intent(in)::grid
181 !*** local
182 integer::ncid,varid
183 integer,dimension(4)::s,c
185 !*** executable
186 call create_dim_attr(filename,grid)
187 call define_var(filename,var_nfuel_cat,dim_fire,'','')
188 call define_var(filename,var_ux,dim_atm_u,'','')
189 call define_var(filename,var_uy,dim_atm_v,'','')
190 call define_var(filename,var_dzdxf,dim_fire,'','')
191 call define_var(filename,var_dzdyf,dim_fire,'','')
192 call define_var(filename,var_zsf,dim_fire,'','')
194 call check(nf90_open(filename,nf90_write,ncid))
196 s(:)=1
197 c(:)=1
198 c(1)=grid%nfirex
199 c(2)=grid%nfirey
201 call check(nf90_inq_varid(ncid,var_nfuel_cat,varid))
202 call check(nf90_put_var(ncid,varid,input%nfuel_cat,start=s(1:3),count=c(1:3)))
203 call check(nf90_inq_varid(ncid,var_dzdxf,varid))
204 call check(nf90_put_var(ncid,varid,input%dzdxf,start=s(1:3),count=c(1:3)))
205 call check(nf90_inq_varid(ncid,var_dzdyf,varid))
206 call check(nf90_put_var(ncid,varid,input%dzdyf,start=s(1:3),count=c(1:3)))
207 call check(nf90_inq_varid(ncid,var_zsf,varid))
208 call check(nf90_put_var(ncid,varid,input%zsf,start=s(1:3),count=c(1:3)))
210 c(1)=grid%natmx+1
211 c(2)=grid%natmy
213 call check(nf90_inq_varid(ncid,var_ux,varid))
214 call check(nf90_put_var(ncid,varid,input%vx,start=s,count=c))
216 c(1)=grid%natmx
217 c(2)=grid%natmy+1
219 call check(nf90_inq_varid(ncid,var_uy,varid))
220 call check(nf90_put_var(ncid,varid,input%vy,start=s,count=c))
222 call check(nf90_close(ncid))
224 end subroutine write_inputs
227 !***
230 subroutine write_vars(filename,output,time)
231 !*** append variables to an output file (extending by the unlimited time
232 !    dimension)
233 implicit none
235 !*** arguments
236 character(len=*),intent(in)::filename
237 type(output_vars),intent(in)::output
238 real,intent(in)::time
240 !*** local
241 integer::ncid,tstep,dimid,varid
242 integer,dimension(3)::s,c
243 character(len=12)::timestr
245 !*** executable
247 ! open the file
248 call check(nf90_open(filename,nf90_write,ncid))
250 ! get the current number of timeslices already written
251 call check(nf90_inq_dimid(ncid,dim_time,dimid))
252 call check(nf90_inquire_dimension(ncid,dimid,len=tstep))
254 ! write out the current simulation time
255 call check(nf90_inq_varid(ncid,'Times',varid))
256 write(timestr,'(E12.3)')time
257 call check(nf90_put_var(ncid,varid,timestr,start=(/1,tstep+1/),count=(/12,1/)))
259 ! set up start and count variables for nf90_put_var so that it writes to the
260 ! correct time slice and write data for each variable
261 s(:)=1
262 c(:)=0
263 s(3)=tstep+1
264 c(3)=1
265 c(1:2)=(/size(output%lfn,1),size(output%lfn,2)/)
266 call check(nf90_inq_varid(ncid,var_lfn,varid))
267 call check(nf90_put_var(ncid,varid,output%lfn,start=s,count=c))
268 c(1:2)=(/size(output%tign,1),size(output%tign,2)/)
269 call check(nf90_inq_varid(ncid,var_tign,varid))
270 call check(nf90_put_var(ncid,varid,output%tign,start=s,count=c))
271 c(1:2)=(/size(output%vx,1),size(output%vx,2)/)
272 call check(nf90_inq_varid(ncid,var_vx,varid))
273 call check(nf90_put_var(ncid,varid,output%vx,start=s,count=c))
274 c(1:2)=(/size(output%vy,1),size(output%vy,2)/)
275 call check(nf90_inq_varid(ncid,var_vy,varid))
276 call check(nf90_put_var(ncid,varid,output%vy,start=s,count=c))
277 c(1:2)=(/size(output%grnhfx,1),size(output%grnhfx,2)/)
278 call check(nf90_inq_varid(ncid,var_grnhfx,varid))
279 call check(nf90_put_var(ncid,varid,output%grnhfx,start=s,count=c))
281 ! close file
282 call check(nf90_close(ncid))
283 end subroutine write_vars
286 !***
289 subroutine set_grid_from_file(filename,grid)
290 ! get grid sizes from input file
291 implicit none
293 !*** arguments
294 character(len=*),intent(in)::filename
295 type(grid_info),intent(out)::grid
297 !*** local
298 integer::ncid,it
299 integer,dimension(4)::fdimid,adimid
301 !*** executable
302 call check(nf90_open(filename,nf90_nowrite,ncid))
303 call check(nf90_inq_dimid(ncid,dim_time,fdimid(3)))
304 adimid(4)=fdimid(3)
305 call check(nf90_inq_dimid(ncid,dim_fire_x,fdimid(1)))
306 call check(nf90_inq_dimid(ncid,dim_fire_y,fdimid(2)))
307 call check(nf90_inq_dimid(ncid,dim_atm_x, adimid(1)))
308 call check(nf90_inq_dimid(ncid,dim_atm_y, adimid(2)))
309 call check(nf90_inq_dimid(ncid,dim_atm_z, adimid(3)))
311 call check(nf90_inquire_dimension(ncid,fdimid(3),len=it))
312 if(it.ne.1)then
313   call crash('invalid number of time steps in input file, must be 1')
314 endif
315 call check(nf90_inquire_dimension(ncid,fdimid(2),len=grid%nfirey))
316 call check(nf90_inquire_dimension(ncid,fdimid(1),len=grid%nfirex))
317 call check(nf90_inquire_dimension(ncid,adimid(2),len=grid%natmy))
318 call check(nf90_inquire_dimension(ncid,adimid(1),len=grid%natmx))
319 call check(nf90_inquire_dimension(ncid,adimid(3),len=grid%natmz))
321 call check(nf90_get_att(ncid,nf90_global,'DX',grid%dx))
322 call check(nf90_get_att(ncid,nf90_global,'DY',grid%dy))
323 call check(nf90_get_att(ncid,nf90_global,'DT',grid%dt))
325 if(compat_fire_grid)then
326   grid%sr_x=grid%nfirex/(grid%natmx+1)
327   grid%sr_y=grid%nfirey/(grid%natmy+1)
328   if( (grid%natmx+1)*grid%sr_x .ne. grid%nfirex .or. &
329       (grid%natmy+1)*grid%sr_y .ne. grid%nfirey)then
330       call crash('invalid dimensions in input file or compat_fire_grid set incorrectly')
331   endif
332   grid%nfirex=grid%natmx*grid%sr_x
333   grid%nfirey=grid%natmy*grid%sr_y
334 else
335   grid%sr_x=grid%nfirex/grid%natmx
336   grid%sr_y=grid%nfirey/grid%natmy
337   if( grid%natmx*grid%sr_x .ne. grid%nfirex .or. &
338       grid%natmy*grid%sr_y .ne. grid%nfirey)then
339       call crash('invalid dimensions in input file or compat_fire_grid set incorrectly')
340   endif
341 endif
342 call check(nf90_close(ncid))
343 end subroutine set_grid_from_file
346 !***
349 subroutine read_vars(filename,input,grid)
350 ! read all variables from input file
351 implicit none
353 !*** arguments
354 character(len=*),intent(in)::filename
355 type(input_vars),intent(inout)::input
356 type(grid_info),intent(in)::grid
358 !*** local
359 integer::ncid,varid,ierr
360 integer,dimension(4)::s,c
361 real,dimension(:,:),allocatable::u,v
363 !*** executable
364 call check(nf90_open(filename,nf90_nowrite,ncid))
365 s(:)=1
366 c(:)=1
367 c(1)=grid%nfirex
368 c(2)=grid%nfirey
369 call check(nf90_inq_varid(ncid,var_nfuel_cat,varid))
370 call check(nf90_get_var(ncid,varid,input%nfuel_cat,start=s(1:3),count=c(1:3)))
371 call check(nf90_inq_varid(ncid,var_dzdxf,varid))
372 call check(nf90_get_var(ncid,varid,input%dzdxf,start=s(1:3),count=c(1:3)))
373 call check(nf90_inq_varid(ncid,var_dzdyf,varid))
374 call check(nf90_get_var(ncid,varid,input%dzdyf,start=s(1:3),count=c(1:3)))
376 allocate(u(0:grid%natmx+2,0:grid%natmy+2),v(0:grid%natmx+2,0:grid%natmy+2),stat=ierr)
377 if(ierr.ne.0)call crash('memory allocation error')
379 c(1)=grid%natmx+1
380 c(2)=grid%natmy
381 c(3)=1
382 c(4)=1
384 call check(nf90_inq_varid(ncid,var_ux,varid))
385 call check(nf90_get_var(ncid,varid,u(1:grid%natmx+1,1:grid%natmy),start=s,count=c))
387 c(1)=grid%natmx
388 c(2)=grid%natmy+1
389 c(3)=1
390 c(4)=1
392 call check(nf90_inq_varid(ncid,var_uy,varid))
393 call check(nf90_get_var(ncid,varid,v(1:grid%natmx,1:grid%natmy+1),start=s,count=c))
395   call continue_at_boundary(1,0,0.,       &
396           0,grid%natmx+2,0,grid%natmy+2,  &
397           1,grid%natmx,1,grid%natmy+1,    &
398           1,grid%natmx,1,grid%natmy+1,    &
399           1,grid%natmx,1,grid%natmy+1,    &
400           v)
402   call continue_at_boundary(0,1,0.,       &
403           0,grid%natmx+2,0,grid%natmy+2,  &
404           1,grid%natmx+1,1,grid%natmy,    &
405           1,grid%natmx+1,1,grid%natmy,    &
406           1,grid%natmx+1,1,grid%natmy,    &
407           u)
410   call interpolate_2d(0,grid%natmx+2,0,grid%natmy+2,         &
411                       1,grid%natmx+1,0,grid%natmy+1,         &
412                       1,grid%nfirex,1,grid%nfirey,           &
413                       1,grid%nfirex,1,grid%nfirey,           &
414                       grid%sr_x,grid%sr_y,                   &
415                       1.,1.,.5,1-.5+(grid%sr_y)*.5,          &
416                       u,input%vx)
418   call interpolate_2d(0,grid%natmx+2,0,grid%natmy+2,         &
419                       0,grid%natmx+1,1,grid%natmy+1,         &
420                       1,grid%nfirex,1,grid%nfirey,           &
421                       1,grid%nfirex,1,grid%nfirey,           &
422                       grid%sr_x,grid%sr_y,                   &
423                       1.,1.,1-.5+(grid%sr_x)*.5,.5,          &
424                       v,input%vy)
426 deallocate(u,v,stat=ierr)
427 if(ierr.ne.0)call crash('deallocation error')
428 call check(nf90_close(ncid))
429 end subroutine read_vars
431 subroutine check(ncerr)
432 implicit none
433 integer,intent(in)::ncerr
434 if(ncerr.ne.nf90_noerr)then
435   print*,"Error calling NetCDF subroutine"
436   print*,trim(nf90_strerror(ncerr))
437   call crash("NETCDF ERROR")
438 endif
439 end subroutine 
441 end module wrf_netcdf