15 use fckit_configuration_module,
only: fckit_configuration
16 use fckit_log_module,
only: fckit_log
38 character(len=50) :: colname
41 real(kind_real),
allocatable :: values(:,:)
45 character(len=50) :: grpname
48 type(datetime),
allocatable :: times(:)
54 character(len=1024) :: filein
55 character(len=1024) :: fileout
59 #define LISTED_TYPE qg_obsdb
62 #include "oops/util/linkedList_i.f"
72 #include "oops/util/linkedList_c.f"
82 type(fckit_configuration),
intent(in) :: f_conf
83 type(datetime),
intent(in) :: winbgn
84 type(datetime),
intent(in) :: winend
87 character(len=1024) :: fin,fout
88 character(len=:),
allocatable :: str
91 if (f_conf%has(
"obsdatain"))
then
92 call f_conf%get_or_die(
"obsdatain.obsfile",str)
97 call fckit_log%info(
'qg_obsdb_setup: file in = '//trim(fin))
100 if (f_conf%has(
"obsdataout"))
then
101 call f_conf%get_or_die(
"obsdataout.obsfile",str)
102 call swap_name_member(f_conf, str)
105 call fckit_log%info(
'qg_obsdb_setup: file out = '//trim(fout))
126 type(
qg_obsdb),
intent(inout) :: self
137 do while (
associated(self%grphead))
139 self%grphead => jgrp%next
141 call datetime_delete(jgrp%times(jobs))
143 deallocate(jgrp%times)
144 do while (
associated(jgrp%colhead))
146 jgrp%colhead => jcol%next
147 deallocate(jcol%values)
162 character(len=*),
intent(in) :: grp
163 character(len=*),
intent(in) :: col
172 call fckit_log%debug(
'qg_obsdb_get: grp = '//trim(grp))
173 call fckit_log%debug(
'qg_obsdb_get: col = '//trim(col))
177 if (.not.
associated(jgrp))
then
179 do while (
associated(jgrp))
180 call fckit_log%info(
'qg_obsdb_get: group '//trim(jgrp%grpname)//
' exists')
183 call fckit_log%error(
'qg_obsdb_get: cannot find '//trim(grp))
184 call abor1_ftn(
'qg_obsdb_get: obs group not found')
189 if (.not.
associated(jcol))
call abor1_ftn(
'qg_obsdb_get: obs column not found')
192 if (
allocated(ovec%values))
deallocate(ovec%values)
193 ovec%nlev = jcol%nlev
195 ovec%nobs = jgrp%nobs
196 allocate(ovec%values(ovec%nlev,ovec%nobs))
199 ovec%values(jlev,jobs) = jcol%values(jlev,jobs)
211 type(
qg_obsdb),
intent(inout) :: self
212 character(len=*),
intent(in) :: grp
213 character(len=*),
intent(in) :: col
223 if (.not.
associated(jgrp))
then
225 do while (
associated(jgrp))
226 call fckit_log%info(
'qg_obsdb_put: group '//trim(jgrp%grpname)//
' exists')
229 call fckit_log%error(
'qg_obsdb_put: cannot find '//trim(grp))
230 call abor1_ftn(
'qg_obsdb_put: obs group not found')
235 if (.not.
associated(jcol))
then
236 if (.not.
associated(jgrp%colhead))
call abor1_ftn(
'qg_obsdb_put: no locations')
238 do while (
associated(jcol%next))
244 jcol%nlev = ovec%nlev
245 allocate(jcol%values(jcol%nlev,jgrp%nobs))
249 if (ovec%nobs/=jgrp%nobs)
call abor1_ftn(
'qg_obsdb_put: error obs number')
250 if (ovec%nlev/=jcol%nlev)
call abor1_ftn(
'qg_obsdb_put: error col number')
253 jcol%values(jlev,jobs) = ovec%values(jlev,jobs)
266 character(len=*),
intent(in) :: grp
267 type(atlas_fieldset),
intent(inout) :: fields
268 type(c_ptr),
intent(in),
value :: c_times
272 character(len=8),
parameter :: col =
'Location'
275 type(atlas_field) :: field_z, field_lonlat
276 real(kind_real),
pointer :: z(:), lonlat(:,:)
280 if (.not.
associated(jgrp))
call abor1_ftn(
'qg_obsdb_locations: obs group not found')
285 if (.not.
associated(jcol))
call abor1_ftn(
'qg_obsdb_locations: obs column not found')
289 field_lonlat = atlas_field(name=
"lonlat", kind=atlas_real(kind_real), shape=[2,nlocs])
290 field_z = atlas_field(name=
"altitude", kind=atlas_real(kind_real), shape=[nlocs])
292 call field_lonlat%data(lonlat)
297 lonlat(1,jo) = jcol%values(1,jo)
298 lonlat(2,jo) = jcol%values(2,jo)
299 z(jo) = jcol%values(3,jo)
300 call f_c_push_to_datetime_vector(c_times, jgrp%times(jo))
303 call fields%add(field_lonlat)
304 call fields%add(field_z)
307 call field_lonlat%final()
318 type(
qg_obsdb),
intent(inout) :: self
319 character(len=*),
intent(in) :: grp
320 type(fckit_configuration),
intent(in) :: f_conf
321 type(datetime),
intent(in) :: bgn
322 type(duration),
intent(in) :: step
323 integer,
intent(in) :: ktimes
324 integer,
intent(inout) :: kobs
327 integer :: nlev,nlocs
328 real(kind_real) :: err
329 type(datetime),
allocatable :: times(:)
333 call f_conf%get_or_die(
"obs_density",nlocs)
337 allocate(times(kobs))
346 call f_conf%get_or_die(
"obs_error",err)
347 call f_conf%get_or_die(
"nval",nlev)
349 obserr%values(:,:) = err
354 deallocate(obsloc%values)
355 deallocate(obserr%values)
366 character(len=*),
intent(in) :: grp
367 integer,
intent(inout) :: kobs
376 if (
associated(jgrp))
then
392 type(
qg_obsdb),
intent(inout) :: self
393 type(datetime),
intent(in) :: winbgn
394 type(datetime),
intent(in) :: winend
397 integer,
parameter :: ngrpmax = 3
398 integer,
parameter :: ncolmax = 7
399 integer :: grp_ids(ngrpmax),igrp,nobs_in_grp,iobs,jobs,ncol,col_ids(ncolmax),icol,nlev_id,values_id
400 integer :: ncid,nobs_id,times_id
403 character(len=50) :: stime
404 logical,
allocatable :: inwindow(:)
405 type(datetime) :: tobs
406 type(datetime),
allocatable :: alltimes(:)
407 real(kind_real),
allocatable :: readbuf(:,:)
410 call ncerr(nf90_open(trim(self%filein),nf90_nowrite,ncid))
413 call ncerr(nf90_inq_grps(ncid,self%ngrp,grp_ids))
418 allocate(self%grphead)
426 call ncerr(nf90_inq_grpname(grp_ids(igrp),jgrp%grpname))
429 call ncerr(nf90_inq_dimid(grp_ids(igrp),
'nobs',nobs_id))
432 call ncerr(nf90_inquire_dimension(grp_ids(igrp),nobs_id,len=nobs_in_grp))
435 call ncerr(nf90_inq_varid(grp_ids(igrp),
'times',times_id))
438 allocate(inwindow(nobs_in_grp))
439 allocate(alltimes(nobs_in_grp))
443 do iobs=1,nobs_in_grp
444 call ncerr(nf90_get_var(grp_ids(igrp),times_id,stime,(/1,iobs/),(/50,1/)))
445 call datetime_create(stime,tobs)
446 if ((tobs > winbgn).and.(tobs <= winend))
then
447 inwindow(iobs) = .true.
448 alltimes(iobs) = tobs
449 jgrp%nobs = jgrp%nobs+1
451 inwindow(iobs) = .false.
456 allocate(jgrp%times(jgrp%nobs))
460 do iobs=1,nobs_in_grp
461 if (inwindow(iobs))
then
463 jgrp%times(jobs) = alltimes(iobs)
468 call ncerr(nf90_inq_grps(grp_ids(igrp),ncol,col_ids))
474 allocate(jgrp%colhead)
482 call ncerr(nf90_inq_grpname(col_ids(icol),jcol%colname))
485 call ncerr(nf90_inq_dimid(col_ids(icol),
'nlev',nlev_id))
488 call ncerr(nf90_inquire_dimension(col_ids(icol),nlev_id,len=jcol%nlev))
491 call ncerr(nf90_inq_varid(col_ids(icol),
'values',values_id))
494 allocate(readbuf(jcol%nlev,nobs_in_grp))
495 allocate(jcol%values(jcol%nlev,jgrp%nobs))
498 call ncerr(nf90_get_var(col_ids(icol),values_id,readbuf(1:jcol%nlev,:),(/1,1/),(/jcol%nlev,nobs_in_grp/)))
502 do iobs=1,nobs_in_grp
503 if (inwindow(iobs))
then
505 jcol%values(:,jobs) = readbuf(:,iobs)
519 call ncerr(nf90_close(ncid))
533 integer :: ncid,nstrmax_id,grp_id,nobs_id,times_id,col_id,nlev_id,values_id
536 character(len=50) :: stime
539 call ncerr(nf90_create(trim(self%fileout),or(nf90_clobber,nf90_netcdf4),ncid))
542 call ncerr(nf90_def_dim(ncid,
'nstrmax',50,nstrmax_id))
546 do while (
associated(jgrp))
547 if (jgrp%nobs > 0)
then
549 call ncerr(nf90_def_grp(ncid,jgrp%grpname,grp_id))
552 call ncerr(nf90_def_dim(grp_id,
'nobs',jgrp%nobs,nobs_id))
555 call ncerr(nf90_def_var(grp_id,
'times',nf90_char,(/nstrmax_id,nobs_id/),times_id))
559 call datetime_to_string(jgrp%times(iobs),stime)
560 call ncerr(nf90_put_var(grp_id,times_id,stime,(/1,iobs/),(/50,1/)))
565 do while (
associated(jcol))
567 call ncerr(nf90_def_grp(grp_id,jcol%colname,col_id))
570 call ncerr(nf90_def_dim(col_id,
'nlev',jcol%nlev,nlev_id))
573 call ncerr(nf90_def_var(col_id,
'values',nf90_double,(/nlev_id,nobs_id/),values_id))
576 call ncerr(nf90_put_var(col_id,values_id,jcol%values(1:jcol%nlev,:),(/1,1/),(/jcol%nlev,jgrp%nobs/)))
588 call ncerr(nf90_close(ncid))
599 character(len=*),
intent(in) :: grp
600 type(
group_data),
pointer,
intent(inout) :: find
606 do while (
associated(find))
607 if (find%grpname==grp)
exit
620 character(len=*),
intent(in) :: col
627 do while (
associated(find))
628 if (find%colname==col)
exit
640 integer,
intent(in) :: nlocs
641 integer,
intent(in) :: ntimes
642 type(datetime),
intent(in) :: bgn
643 type(duration),
intent(in) :: step
644 type(datetime),
intent(inout) :: times(nlocs*ntimes)
648 integer :: jobs,iobs,jstep
649 real(kind_real) :: x(nlocs),y(nlocs),z(nlocs),lon(nlocs),lat(nlocs)
650 type(datetime) :: now
672 obsloc%values(:,iobs) = (/lon(jobs),lat(jobs),z(jobs)/)
674 call datetime_update(now,step)
678 call datetime_delete(now)
688 type(
qg_obsdb),
intent(inout) :: self
689 character(len=*),
intent(in) :: grp
690 type(datetime),
intent(in) :: times(:)
699 if (
associated(igrp))
call abor1_ftn(
'qg_obsdb_create: obs group already exists')
700 if (
associated(self%grphead))
then
702 do while (
associated(igrp%next))
708 allocate(self%grphead)
714 igrp%nobs =
size(times)
715 allocate(igrp%times(igrp%nobs))
716 igrp%times(:) = times(:)
717 allocate(igrp%colhead)
718 igrp%colhead%colname =
'Location'
719 igrp%colhead%nlev = 3
720 allocate(igrp%colhead%values(3,igrp%nobs))
721 if (locs%nlev/=3)
call abor1_ftn(
'qg_obsdb_create: error locations not 3D')
722 if (locs%nobs/=igrp%nobs)
call abor1_ftn(
'qg_obsdb_create: error locations number')
725 igrp%colhead%values(jlev,jobs) = locs%values(jlev,jobs)
728 self%ngrp = self%ngrp+1
real(kind_real), parameter, public domain_depth
Model depth (m)
real(kind_real), parameter, public domain_meridional
Model domain in meridional direction (m)
real(kind_real), parameter, public domain_zonal
Model domain in zonal direction (m)
integer, parameter rseed
Random seed (for reproducibility)
subroutine, public qg_obsdb_generate(self, grp, f_conf, bgn, step, ktimes, kobs)
Generate observation data.
subroutine qg_obsdb_create(self, grp, times, locs)
Create observation data.
subroutine, public qg_obsdb_put(self, grp, col, ovec)
Put observations data.
subroutine, public qg_obsdb_locations(self, grp, fields, c_times)
Get locations from observation data.
subroutine, public qg_obsdb_delete(self)
Delete observation data.
subroutine qg_obsdb_generate_locations(nlocs, ntimes, bgn, step, times, obsloc)
Generate random locations.
subroutine qg_obsdb_find_group(self, grp, find)
Find observation data group.
subroutine qg_obsdb_write(self)
Write observation data.
subroutine, public qg_obsdb_get(self, grp, col, ovec)
Get observation data.
subroutine qg_obsdb_find_column(grp, col, find)
Find observation data column.
type(registry_t), public qg_obsdb_registry
Linked list interface - defines registry_t type.
subroutine, public qg_obsdb_setup(self, f_conf, winbgn, winend)
Linked list implementation.
subroutine, public qg_obsdb_nobs(self, grp, kobs)
Get observation data size.
subroutine qg_obsdb_read(self, winbgn, winend)
Read observation data.
subroutine, public qg_obsvec_setup(self, nlev, nobs)
Linked list implementation.
subroutine, public xy_to_lonlat(x, y, lon, lat, mapfac)
Convert x/y to lon/lat.