IODA Bundle
ncio_mod.f90
Go to the documentation of this file.
1 module ncio_mod
2 
3 use netcdf
4 use kinds, only: i_kind, r_single, r_kind
14 
15 implicit none
16 
17 private
18 public :: write_obs
19 
20 contains
21 
22 subroutine write_obs (filedate, write_opt, outdir)
23 
24  implicit none
25 
26  character(len=*), intent(in) :: filedate
27  integer(i_kind), intent(in) :: write_opt
28  character(len=*), intent(in) :: outdir
29 
30  character(len=512) :: ncfname ! netcdf file name
31  integer(i_kind), dimension(n_ncdim) :: ncid_ncdim
32  integer(i_kind), dimension(n_ncdim) :: val_ncdim
33  character(len=nstring) :: ncname
34  integer(i_kind) :: ncfileid
35  integer(i_kind) :: ntype, nvar
36  integer(i_kind) :: i, ityp, ivar, ii
37  integer(i_kind) :: idim, dim1, dim2
38  character(len=nstring), allocatable :: str_nstring(:)
39  character(len=ndatetime), allocatable :: str_ndatetime(:)
40  character(len=nstring), allocatable :: name_var_tb(:)
41  character(len=nstring) :: str_tmp
42  character(len=4) :: c4
43 
44  if ( write_opt == write_nc_conv ) then
45  ntype = nobtype
46  else if ( write_opt == write_nc_radiance ) then
47  ntype = ninst
48  else
49  write(*,*) ' Error: unknwon write_opt = ', write_opt
50  return
51  end if
52 
53  val_ncdim(4) = nstring
54  val_ncdim(5) = ndatetime
55 
56  obtype_loop: do ityp = 1, ntype
57 
58  if ( xdata(ityp)%nlocs == 0 ) cycle obtype_loop
59 
60  if ( write_opt == write_nc_conv ) then
61  ncfname = trim(outdir)//trim(obtype_list(ityp))//'_obs_'//trim(filedate)//'.nc4'
62  else if ( write_opt == write_nc_radiance ) then
63  ncfname = trim(outdir)//trim(inst_list(ityp))//'_obs_'//trim(filedate)//'.nc4'
64  allocate (name_var_tb(xdata(ityp)%nvars))
65  end if
66  write(*,*) '--- writing ', trim(ncfname)
67  call open_netcdf_for_write(trim(ncfname),ncfileid)
68 
69  val_ncdim(1) = xdata(ityp)%nvars
70  val_ncdim(2) = xdata(ityp)%nlocs
71  val_ncdim(3) = xdata(ityp)%nrecs
72 
73  ! define netcdf dimensions
74  do i = 1, n_ncdim
75  call def_netcdf_dims(ncfileid,trim(name_ncdim(i)),val_ncdim(i),ncid_ncdim(i))
76  end do
77 
78  ! define netcdf variables
79  do i = 1, xdata(ityp) % nvars
80  if ( write_opt == write_nc_conv ) then
81  ivar = xdata(ityp) % var_idx(i)
82  ncname = trim(name_var_met(ivar))//'@ObsValue'
83  call def_netcdf_var(ncfileid,ncname,(/ncid_ncdim(2)/),nf90_float,'units',unit_var_met(ivar))
84  ncname = trim(name_var_met(ivar))//'@ObsError'
85  call def_netcdf_var(ncfileid,ncname,(/ncid_ncdim(2)/),nf90_float)
86  ncname = trim(name_var_met(ivar))//'@PreQC'
87  call def_netcdf_var(ncfileid,ncname,(/ncid_ncdim(2)/),nf90_int)
88  ncname = trim(name_var_met(ivar))//'@ObsType'
89  call def_netcdf_var(ncfileid,ncname,(/ncid_ncdim(2)/),nf90_int)
90  else if ( write_opt == write_nc_radiance ) then
91  write(unit=c4, fmt='(i4)') i
92  name_var_tb(i) = trim(var_tb)//'_'//trim(adjustl(c4))
93  ncname = trim(name_var_tb(i))//'@ObsValue'
94  call def_netcdf_var(ncfileid,ncname,(/ncid_ncdim(2)/),nf90_float,'units','K')
95  ncname = trim(name_var_tb(i))//'@ObsError'
96  call def_netcdf_var(ncfileid,ncname,(/ncid_ncdim(2)/),nf90_float)
97  ncname = trim(name_var_tb(i))//'@PreQC'
98  call def_netcdf_var(ncfileid,ncname,(/ncid_ncdim(2)/),nf90_int)
99  end if
100  end do
101 
102  do i = 1, nvar_info
103  ncname = trim(name_var_info(i))//'@MetaData'
104  if ( trim(name_var_info(i)) == 'variable_names' ) then
105  ncname = trim(name_var_info(i))//'@VarMetaData'
106  end if
108  dim1 = ncid_ncdim(idim)
109  if ( ufo_vars_getindex(name_ncdim, dim_var_info(2,i)) > 0 ) then
111  dim2 = ncid_ncdim(idim)
112  call def_netcdf_var(ncfileid,ncname,(/dim1,dim2/),type_var_info(i))
113  else
114  call def_netcdf_var(ncfileid,ncname,(/dim1/),type_var_info(i))
115  end if
116  end do ! nvar_info
117 
118  if ( write_opt == write_nc_radiance ) then
119  do i = 1, nsen_info
120  ncname = trim(name_sen_info(i))//'@MetaData'
121  if ( trim(name_sen_info(i)) == 'sensor_channel' ) then
122  ncname = trim(name_sen_info(i))//'@VarMetaData'
123  end if
125  dim1 = ncid_ncdim(idim)
126  if ( ufo_vars_getindex(name_ncdim, dim_sen_info(2,i)) > 0 ) then
128  dim2 = ncid_ncdim(idim)
129  call def_netcdf_var(ncfileid,ncname,(/dim1,dim2/),type_sen_info(i))
130  else
131  call def_netcdf_var(ncfileid,ncname,(/dim1/),type_sen_info(i))
132  end if
133  end do ! nsen_info
134  end if ! write_nc_radiance
135 
136  call def_netcdf_end(ncfileid)
137 
138  ! writing netcdf variables
139  var_loop: do i = 1, xdata(ityp) % nvars
140  if ( write_opt == write_nc_conv ) then
141  ivar = xdata(ityp) % var_idx(i)
142  if ( vflag(ivar,ityp) == itrue ) then
143  ncname = trim(name_var_met(ivar))//'@ObsValue'
144  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xfield(:,i)%val)
145  ncname = trim(name_var_met(ivar))//'@ObsError'
146  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xfield(:,i)%err)
147  ncname = trim(name_var_met(ivar))//'@PreQC'
148  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xfield(:,i)%qm)
149  ncname = trim(name_var_met(ivar))//'@ObsType'
150  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xfield(:,i)%rptype)
151  end if
152  else if ( write_opt == write_nc_radiance ) then
153  ncname = trim(name_var_tb(i))//'@ObsValue'
154  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xfield(:,i)%val)
155  ncname = trim(name_var_tb(i))//'@ObsError'
156  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xfield(:,i)%err)
157  ncname = trim(name_var_tb(i))//'@PreQC'
158  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xfield(:,i)%qm)
159  end if
160  end do var_loop
161 
162  do i = 1, nvar_info
163  ncname = trim(name_var_info(i))//'@MetaData'
164  if ( trim(name_var_info(i)) == 'variable_names' ) then
165  ncname = trim(name_var_info(i))//'@VarMetaData'
166  end if
167  if ( type_var_info(i) == nf90_int ) then
168  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xinfo_int(:,i))
169  else if ( type_var_info(i) == nf90_float ) then
170  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xinfo_float(:,i))
171  else if ( type_var_info(i) == nf90_char ) then
172  if ( trim(name_var_info(i)) == 'variable_names' ) then
173  if ( write_opt == write_nc_conv ) then
174  call put_netcdf_var(ncfileid,ncname,name_var_met(xdata(ityp)%var_idx(:)))
175  else if ( write_opt == write_nc_radiance ) then
176  call put_netcdf_var(ncfileid,ncname,name_var_tb(:))
177  end if
178  else if ( trim(name_var_info(i)) == 'station_id' ) then
179  allocate(str_nstring(val_ncdim(2))) ! nlocs
180  str_nstring(:) = xdata(ityp)%xinfo_char(:,i)
181  call put_netcdf_var(ncfileid,ncname,str_nstring)
182  deallocate(str_nstring)
183  else if ( trim(name_var_info(i)) == 'datetime' ) then
184  allocate(str_ndatetime(val_ncdim(2)))
185  do ii = 1, val_ncdim(2)
186  str_tmp = xdata(ityp)%xinfo_char(ii,i)
187  str_ndatetime(ii) = str_tmp(1:ndatetime)
188  end do
189  call put_netcdf_var(ncfileid,ncname,str_ndatetime)
190  deallocate(str_ndatetime)
191  end if
192  end if
193  end do
194 
195  if ( write_opt == write_nc_radiance ) then
196  do i = 1, nsen_info
197  ncname = trim(name_sen_info(i))//'@MetaData'
198  if ( type_sen_info(i) == nf90_int ) then
199  if ( trim(name_sen_info(i)) == 'sensor_channel' ) then
200  ncname = trim(name_sen_info(i))//'@VarMetaData'
201  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xseninfo_int(:,i))
202  end if
203  !call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xseninfo_int(:,i))
204  else if ( type_sen_info(i) == nf90_float ) then
205  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xseninfo_float(:,i))
206  else if ( type_sen_info(i) == nf90_char ) then
207  call put_netcdf_var(ncfileid,ncname,xdata(ityp)%xseninfo_char(:,i))
208  end if
209  end do
210  deallocate (name_var_tb)
211  end if ! write_nc_radiance
212 
213  call close_netcdf(trim(ncfname),ncfileid)
214 
215  end do obtype_loop
216 
217  ! deallocate xdata
218  do i = 1, ntype
219  if ( allocated(xdata(i)%var_idx) ) deallocate(xdata(i)%var_idx)
220  if ( allocated(xdata(i)%xfield) ) deallocate(xdata(i)%xfield)
221  if ( allocated(xdata(i)%xinfo_int) ) deallocate(xdata(i)%xinfo_int)
222  if ( allocated(xdata(i)%xinfo_float) ) deallocate(xdata(i)%xinfo_float)
223  if ( allocated(xdata(i)%xinfo_char) ) deallocate(xdata(i)%xinfo_char)
224  if ( allocated(xdata(i)%xseninfo_int) ) deallocate(xdata(i)%xseninfo_int)
225  if ( allocated(xdata(i)%xseninfo_float) ) deallocate(xdata(i)%xseninfo_float)
226  if ( allocated(xdata(i)%xseninfo_char) ) deallocate(xdata(i)%xseninfo_char)
227  end do
228  deallocate(xdata)
229 
230 end subroutine write_obs
231 
232 end module ncio_mod
integer(i_kind), parameter nstring
Definition: define_mod.f90:15
integer(i_kind), parameter write_nc_radiance
Definition: define_mod.f90:25
character(len=nstring), dimension(nobtype) obtype_list
Definition: define_mod.f90:28
character(len=nstring), dimension(2, nsen_info) dim_sen_info
Definition: define_mod.f90:152
character(len=nstring), dimension(nvar_info) name_var_info
Definition: define_mod.f90:96
integer(i_kind), parameter nobtype
Definition: define_mod.f90:17
integer(i_kind), dimension(nsen_info) type_sen_info
Definition: define_mod.f90:142
integer(i_kind), parameter itrue
Definition: define_mod.f90:13
character(len=nstring), dimension(n_ncdim) name_ncdim
Definition: define_mod.f90:88
integer(i_kind), parameter ifalse
Definition: define_mod.f90:14
integer(i_kind), parameter write_nc_conv
Definition: define_mod.f90:24
integer(i_kind), parameter ninst
Definition: define_mod.f90:22
integer(i_kind), parameter nsen_info
Definition: define_mod.f90:21
integer(i_kind), parameter ndatetime
Definition: define_mod.f90:16
character(len=nstring), dimension(nsen_info) name_sen_info
Definition: define_mod.f90:132
integer(i_kind), parameter nvar_info
Definition: define_mod.f90:20
character(len=nstring), dimension(nvar_met) name_var_met
Definition: define_mod.f90:38
character(len=nstring), dimension(ninst) inst_list
Definition: define_mod.f90:70
type(xdata_type), dimension(:), allocatable xdata
Definition: define_mod.f90:185
integer(i_kind), dimension(nvar_met, nobtype) vflag
Definition: define_mod.f90:49
integer(i_kind), parameter n_ncdim
Definition: define_mod.f90:18
character(len=nstring), dimension(nvar_met) unit_var_met
Definition: define_mod.f90:59
integer(i_kind), dimension(nvar_info) type_var_info
Definition: define_mod.f90:108
character(len=nstring), dimension(2, nvar_info) dim_var_info
Definition: define_mod.f90:120
Definition: kinds.f90:1
integer, parameter, public r_single
Definition: kinds.f90:79
integer, parameter, public i_kind
Definition: kinds.f90:71
integer, parameter, public r_kind
Definition: kinds.f90:100
subroutine, public write_obs(filedate, write_opt, outdir)
Definition: ncio_mod.f90:23
subroutine, public close_netcdf(fname, ncfileid)
Definition: netcdf_mod.f90:48
subroutine, public get_netcdf_dims(fileid, variable, output)
Definition: netcdf_mod.f90:59
subroutine, public def_netcdf_dims(fileid, variable, input, output)
Definition: netcdf_mod.f90:162
subroutine, public open_netcdf_for_write(fname, ncfileid)
Definition: netcdf_mod.f90:146
subroutine, public def_netcdf_end(fileid)
Definition: netcdf_mod.f90:229
subroutine, public def_netcdf_var(fileid, variable, dimids, nctype, attrib_name, attrib)
Definition: netcdf_mod.f90:187
integer function, public ufo_vars_getindex(vars, varname)