OOPS
qg_obsdb_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
2 ! (C) Copyright 2017-2019 UCAR.
3 !
4 ! This software is licensed under the terms of the Apache Licence Version 2.0
5 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6 ! In applying this licence, ECMWF does not waive the privileges and immunities
7 ! granted to it by virtue of its status as an intergovernmental organisation nor
8 ! does it submit to any jurisdiction.
9 
11 
12 use atlas_module
13 use datetime_mod
14 use duration_mod
15 use fckit_configuration_module, only: fckit_configuration
16 use fckit_log_module, only: fckit_log
17 use iso_c_binding
18 use kinds
19 use netcdf
21 use qg_locs_mod
22 use qg_obsvec_mod
24 use qg_tools_mod
25 use random_mod
26 use string_f_c_mod
27 
28 implicit none
29 
30 private
31 public :: qg_obsdb
32 public :: qg_obsdb_registry
34 ! ------------------------------------------------------------------------------
35 integer,parameter :: rseed = 1 !< Random seed (for reproducibility)
36 
38  character(len=50) :: colname !< Column name
39  type(column_data),pointer :: next => null() !< Next column
40  integer :: nlev !< Number of levels
41  real(kind_real),allocatable :: values(:,:) !< Values
42 end type column_data
43 
45  character(len=50) :: grpname !< Group name
46  type(group_data),pointer :: next => null() !< Next group
47  integer :: nobs !< Number of observations
48  type(datetime),allocatable :: times(:) !< Time-slots
49  type(column_data),pointer :: colhead => null() !< Head column
50 end type group_data
51 
53  integer :: ngrp !< Number of groups
54  character(len=1024) :: filein !< Input filename
55  character(len=1024) :: fileout !< Output filename
56  type(group_data),pointer :: grphead => null() !< Head group
57 end type qg_obsdb
58 
59 #define LISTED_TYPE qg_obsdb
60 
61 !> Linked list interface - defines registry_t type
62 #include "oops/util/linkedList_i.f"
63 
64 !> Global registry
65 type(registry_t) :: qg_obsdb_registry
66 ! ------------------------------------------------------------------------------
67 contains
68 ! ------------------------------------------------------------------------------
69 ! Public
70 ! ------------------------------------------------------------------------------
71 !> Linked list implementation
72 #include "oops/util/linkedList_c.f"
73 ! ------------------------------------------------------------------------------
74 !> Setup observation data
75 subroutine qg_obsdb_setup(self,f_conf,winbgn,winend)
76 use string_utils
77 
78 implicit none
79 
80 ! Passed variables
81 type(qg_obsdb),intent(inout) :: self !< Observation data
82 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
83 type(datetime),intent(in) :: winbgn !< Start of window
84 type(datetime),intent(in) :: winend !< End of window
85 
86 ! Local variables
87 character(len=1024) :: fin,fout
88 character(len=:),allocatable :: str
89 
90 ! Input file
91 if (f_conf%has("obsdatain")) then
92  call f_conf%get_or_die("obsdatain.obsfile",str)
93  fin = str
94 else
95  fin = ''
96 endif
97 call fckit_log%info('qg_obsdb_setup: file in = '//trim(fin))
98 
99 ! Output file
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)
103 
104  fout = str
105  call fckit_log%info('qg_obsdb_setup: file out = '//trim(fout))
106 else
107  fout = ''
108 endif
109 
110 ! Set attributes
111 self%ngrp = 0
112 self%filein = fin
113 self%fileout = fout
114 
115 ! Read observation data
116 if (self%filein/='') call qg_obsdb_read(self,winbgn,winend)
117 
118 end subroutine qg_obsdb_setup
119 ! ------------------------------------------------------------------------------
120 !> Delete observation data
121 subroutine qg_obsdb_delete(self)
122 
123 implicit none
124 
125 ! Passed variables
126 type(qg_obsdb),intent(inout) :: self !< Observation data
127 
128 ! Local variables
129 type(group_data),pointer :: jgrp
130 type(column_data),pointer :: jcol
131 integer :: jobs
132 
133 ! Write observation data
134 if (self%fileout/='') call qg_obsdb_write(self)
135 
136 ! Release memory
137 do while (associated(self%grphead))
138  jgrp => self%grphead
139  self%grphead => jgrp%next
140  do jobs=1,jgrp%nobs
141  call datetime_delete(jgrp%times(jobs))
142  enddo
143  deallocate(jgrp%times)
144  do while (associated(jgrp%colhead))
145  jcol => jgrp%colhead
146  jgrp%colhead => jcol%next
147  deallocate(jcol%values)
148  deallocate(jcol)
149  enddo
150  deallocate(jgrp)
151 enddo
152 
153 end subroutine qg_obsdb_delete
154 ! ------------------------------------------------------------------------------
155 !> Get observation data
156 subroutine qg_obsdb_get(self,grp,col,ovec)
157 
158 implicit none
159 
160 ! Passed variables
161 type(qg_obsdb),intent(in) :: self !< Observation data
162 character(len=*),intent(in) :: grp !< Group
163 character(len=*),intent(in) :: col !< Column
164 type(qg_obsvec),intent(inout) :: ovec !< Observation vector
165 
166 ! Local variables
167 type(group_data),pointer :: jgrp
168 type(column_data),pointer :: jcol
169 integer :: jobs,jlev
170 
171 ! Print Group and column
172 call fckit_log%debug('qg_obsdb_get: grp = '//trim(grp))
173 call fckit_log%debug('qg_obsdb_get: col = '//trim(col))
174 
175 ! Find observation group
176 call qg_obsdb_find_group(self,grp,jgrp)
177 if (.not.associated(jgrp)) then
178  jgrp => self%grphead
179  do while (associated(jgrp))
180  call fckit_log%info('qg_obsdb_get: group '//trim(jgrp%grpname)//' exists')
181  jgrp => jgrp%next
182  enddo
183  call fckit_log%error('qg_obsdb_get: cannot find '//trim(grp))
184  call abor1_ftn('qg_obsdb_get: obs group not found')
185 endif
186 
187 ! Find observation column
188 call qg_obsdb_find_column(jgrp,col,jcol)
189 if (.not.associated(jcol)) call abor1_ftn('qg_obsdb_get: obs column not found')
190 
191 ! Get observation data
192 if (allocated(ovec%values)) deallocate(ovec%values)
193 ovec%nlev = jcol%nlev
194 ! get all the obs
195 ovec%nobs = jgrp%nobs
196 allocate(ovec%values(ovec%nlev,ovec%nobs))
197 do jobs=1,jgrp%nobs
198  do jlev=1,jcol%nlev
199  ovec%values(jlev,jobs) = jcol%values(jlev,jobs)
200  enddo
201 enddo
202 
203 end subroutine qg_obsdb_get
204 ! ------------------------------------------------------------------------------
205 !> Put observations data
206 subroutine qg_obsdb_put(self,grp,col,ovec)
207 
208 implicit none
209 
210 ! Passed variables
211 type(qg_obsdb),intent(inout) :: self !< Observation data
212 character(len=*),intent(in) :: grp !< Group
213 character(len=*),intent(in) :: col !< Column
214 type(qg_obsvec),intent(in) :: ovec !< Observation vector
215 
216 ! Local variables
217 type(group_data),pointer :: jgrp
218 type(column_data),pointer :: jcol
219 integer :: jobs,jlev
220 
221 ! Find observation group
222 call qg_obsdb_find_group(self,grp,jgrp)
223 if (.not.associated(jgrp)) then
224  jgrp => self%grphead
225  do while (associated(jgrp))
226  call fckit_log%info('qg_obsdb_put: group '//trim(jgrp%grpname)//' exists')
227  jgrp => jgrp%next
228  enddo
229  call fckit_log%error('qg_obsdb_put: cannot find '//trim(grp))
230  call abor1_ftn('qg_obsdb_put: obs group not found')
231 endif
232 
233 ! Find observation column (and add it if not there)
234 call qg_obsdb_find_column(jgrp,col,jcol)
235 if (.not.associated(jcol)) then
236  if (.not.associated(jgrp%colhead)) call abor1_ftn('qg_obsdb_put: no locations')
237  jcol => jgrp%colhead
238  do while (associated(jcol%next))
239  jcol => jcol%next
240  enddo
241  allocate(jcol%next)
242  jcol => jcol%next
243  jcol%colname = col
244  jcol%nlev = ovec%nlev
245  allocate(jcol%values(jcol%nlev,jgrp%nobs))
246 endif
247 
248 ! Put observation data
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')
251 do jobs=1,jgrp%nobs
252  do jlev=1,jcol%nlev
253  jcol%values(jlev,jobs) = ovec%values(jlev,jobs)
254  enddo
255 enddo
256 
257 end subroutine qg_obsdb_put
258 ! ------------------------------------------------------------------------------
259 !> Get locations from observation data
260 subroutine qg_obsdb_locations(self,grp,fields,c_times)
261 
262 implicit none
263 
264 ! Passed variables
265 type(qg_obsdb),intent(in) :: self !< Observation data
266 character(len=*),intent(in) :: grp !< Group
267 type(atlas_fieldset), intent(inout) :: fields !< Locations FieldSet
268 type(c_ptr), intent(in), value :: c_times !< pointer to times array in C++
269 
270 ! Local variables
271 integer :: nlocs, jo
272 character(len=8),parameter :: col = 'Location'
273 type(group_data),pointer :: jgrp
274 type(column_data),pointer :: jcol
275 type(atlas_field) :: field_z, field_lonlat
276 real(kind_real), pointer :: z(:), lonlat(:,:)
277 
278 ! Find observation group
279 call qg_obsdb_find_group(self,grp,jgrp)
280 if (.not.associated(jgrp)) call abor1_ftn('qg_obsdb_locations: obs group not found')
281 nlocs = jgrp%nobs
282 
283 ! Find observation column
284 call qg_obsdb_find_column(jgrp,col,jcol)
285 if (.not.associated(jcol)) call abor1_ftn('qg_obsdb_locations: obs column not found')
286 
287 ! Set number of observations
288 
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])
291 
292 call field_lonlat%data(lonlat)
293 call field_z%data(z)
294 
295 ! Copy coordinates
296 do jo = 1, nlocs
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))
301 enddo
302 
303 call fields%add(field_lonlat)
304 call fields%add(field_z)
305 
306 ! release pointers
307 call field_lonlat%final()
308 call field_z%final()
309 
310 end subroutine qg_obsdb_locations
311 ! ------------------------------------------------------------------------------
312 !> Generate observation data
313 subroutine qg_obsdb_generate(self,grp,f_conf,bgn,step,ktimes,kobs)
314 
315 implicit none
316 
317 ! Passed variables
318 type(qg_obsdb),intent(inout) :: self !< Observation data
319 character(len=*),intent(in) :: grp !< Group
320 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
321 type(datetime),intent(in) :: bgn !< Start time
322 type(duration),intent(in) :: step !< Time-step
323 integer,intent(in) :: ktimes !< Number of time-slots
324 integer,intent(inout) :: kobs !< Number of observations
325 
326 ! Local variables
327 integer :: nlev,nlocs
328 real(kind_real) :: err
329 type(datetime),allocatable :: times(:)
330 type(qg_obsvec) :: obsloc,obserr
331 
332 ! Get number of observations
333 call f_conf%get_or_die("obs_density",nlocs)
334 kobs = nlocs*ktimes
335 
336 ! Allocation
337 allocate(times(kobs))
338 
339 ! Generate locations
340 call qg_obsdb_generate_locations(nlocs,ktimes,bgn,step,times,obsloc)
341 
342 ! Create observations data
343 call qg_obsdb_create(self,trim(grp),times,obsloc)
344 
345 ! Create observation error
346 call f_conf%get_or_die("obs_error",err)
347 call f_conf%get_or_die("nval",nlev)
348 call qg_obsvec_setup(obserr,nlev,kobs)
349 obserr%values(:,:) = err
350 call qg_obsdb_put(self,trim(grp),'ObsError',obserr)
351 
352 ! Release memory
353 deallocate(times)
354 deallocate(obsloc%values)
355 deallocate(obserr%values)
356 
357 end subroutine qg_obsdb_generate
358 ! ------------------------------------------------------------------------------
359 !> Get observation data size
360 subroutine qg_obsdb_nobs(self,grp,kobs)
361 
362 implicit none
363 
364 ! Passed variables
365 type(qg_obsdb),intent(in) :: self !< Observation data
366 character(len=*),intent(in) :: grp !< Group
367 integer,intent(inout) :: kobs !< Number of observations
368 
369 ! Local variables
370 type(group_data),pointer :: jgrp
371 
372 ! Find group
373 call qg_obsdb_find_group(self,grp,jgrp)
374 
375 ! Get observation data size
376 if (associated(jgrp)) then
377  kobs = jgrp%nobs
378 else
379  kobs = 0
380 endif
381 
382 end subroutine qg_obsdb_nobs
383 ! ------------------------------------------------------------------------------
384 ! Private
385 ! ------------------------------------------------------------------------------
386 !> Read observation data
387 subroutine qg_obsdb_read(self,winbgn,winend)
388 
389 implicit none
390 
391 ! Passed variables
392 type(qg_obsdb),intent(inout) :: self !< Observation data
393 type(datetime),intent(in) :: winbgn !< Start of window
394 type(datetime),intent(in) :: winend !< End of window
395 
396 ! Local variables
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
401 type(group_data),pointer :: jgrp
402 type(column_data),pointer :: jcol
403 character(len=50) :: stime
404 logical, allocatable :: inwindow(:)
405 type(datetime) :: tobs
406 type(datetime), allocatable :: alltimes(:)
407 real(kind_real),allocatable :: readbuf(:,:)
408 
409 ! Open NetCDF file
410 call ncerr(nf90_open(trim(self%filein),nf90_nowrite,ncid))
411 
412 ! Get groups ids
413 call ncerr(nf90_inq_grps(ncid,self%ngrp,grp_ids))
414 
415 do igrp=1,self%ngrp
416  ! Allocation
417  if (igrp==1) then
418  allocate(self%grphead)
419  jgrp => self%grphead
420  else
421  allocate(jgrp%next)
422  jgrp => jgrp%next
423  endif
424 
425  ! Get group name
426  call ncerr(nf90_inq_grpname(grp_ids(igrp),jgrp%grpname))
427 
428  ! Get dimension id
429  call ncerr(nf90_inq_dimid(grp_ids(igrp),'nobs',nobs_id))
430 
431  ! Get dimension
432  call ncerr(nf90_inquire_dimension(grp_ids(igrp),nobs_id,len=nobs_in_grp))
433 
434  ! Get variable id
435  call ncerr(nf90_inq_varid(grp_ids(igrp),'times',times_id))
436 
437  ! Allocation
438  allocate(inwindow(nobs_in_grp))
439  allocate(alltimes(nobs_in_grp))
440 
441  ! Read in times
442  jgrp%nobs = 0
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
450  else
451  inwindow(iobs) = .false.
452  endif
453  enddo
454 
455  ! Allocation
456  allocate(jgrp%times(jgrp%nobs))
457 
458  ! Copy times
459  jobs=0
460  do iobs=1,nobs_in_grp
461  if (inwindow(iobs)) then
462  jobs = jobs+1
463  jgrp%times(jobs) = alltimes(iobs)
464  endif
465  end do
466 
467  ! Count columns
468  call ncerr(nf90_inq_grps(grp_ids(igrp),ncol,col_ids))
469 
470  ! Loop over columns
471  do icol=1,ncol
472  ! Allocation
473  if (icol==1) then
474  allocate(jgrp%colhead)
475  jcol => jgrp%colhead
476  else
477  allocate(jcol%next)
478  jcol => jcol%next
479  endif
480 
481  ! Get column name
482  call ncerr(nf90_inq_grpname(col_ids(icol),jcol%colname))
483 
484  ! Get dimension id
485  call ncerr(nf90_inq_dimid(col_ids(icol),'nlev',nlev_id))
486 
487  ! Get dimension
488  call ncerr(nf90_inquire_dimension(col_ids(icol),nlev_id,len=jcol%nlev))
489 
490  ! Get variable id
491  call ncerr(nf90_inq_varid(col_ids(icol),'values',values_id))
492 
493  ! Allocation
494  allocate(readbuf(jcol%nlev,nobs_in_grp))
495  allocate(jcol%values(jcol%nlev,jgrp%nobs))
496 
497  ! Get values
498  call ncerr(nf90_get_var(col_ids(icol),values_id,readbuf(1:jcol%nlev,:),(/1,1/),(/jcol%nlev,nobs_in_grp/)))
499 
500  ! Copy values
501  jobs = 0
502  do iobs=1,nobs_in_grp
503  if (inwindow(iobs)) then
504  jobs = jobs+1
505  jcol%values(:,jobs) = readbuf(:,iobs)
506  endif
507  enddo
508 
509  ! Release memory
510  deallocate(readbuf)
511  enddo
512 
513  ! Release memory
514  deallocate(alltimes)
515  deallocate(inwindow)
516 enddo
517 
518 ! Close NetCDF file
519 call ncerr(nf90_close(ncid))
520 
521 end subroutine qg_obsdb_read
522 ! ------------------------------------------------------------------------------
523 !> Write observation data
524 subroutine qg_obsdb_write(self)
525 
526 implicit none
527 
528 ! Passed variables
529 type(qg_obsdb),intent(in) :: self !< Observation data
530 
531 ! Local variables
532 integer :: iobs
533 integer :: ncid,nstrmax_id,grp_id,nobs_id,times_id,col_id,nlev_id,values_id
534 type(group_data),pointer :: jgrp
535 type(column_data),pointer :: jcol
536 character(len=50) :: stime
537 
538 ! Create NetCDF file
539 call ncerr(nf90_create(trim(self%fileout),or(nf90_clobber,nf90_netcdf4),ncid))
540 
541 ! Define dimensions
542 call ncerr(nf90_def_dim(ncid,'nstrmax',50,nstrmax_id))
543 
544 ! Loop over groups
545 jgrp => self%grphead
546 do while (associated(jgrp))
547  if (jgrp%nobs > 0) then
548  ! Create group
549  call ncerr(nf90_def_grp(ncid,jgrp%grpname,grp_id))
550 
551  ! Define dimension
552  call ncerr(nf90_def_dim(grp_id,'nobs',jgrp%nobs,nobs_id))
553 
554  ! Define variable
555  call ncerr(nf90_def_var(grp_id,'times',nf90_char,(/nstrmax_id,nobs_id/),times_id))
556 
557  ! Put variable
558  do iobs=1,jgrp%nobs
559  call datetime_to_string(jgrp%times(iobs),stime)
560  call ncerr(nf90_put_var(grp_id,times_id,stime,(/1,iobs/),(/50,1/)))
561  end do
562 
563  ! Loop over columns
564  jcol => jgrp%colhead
565  do while (associated(jcol))
566  ! Create subgroup
567  call ncerr(nf90_def_grp(grp_id,jcol%colname,col_id))
568 
569  ! Define dimension
570  call ncerr(nf90_def_dim(col_id,'nlev',jcol%nlev,nlev_id))
571 
572  ! Define variable
573  call ncerr(nf90_def_var(col_id,'values',nf90_double,(/nlev_id,nobs_id/),values_id))
574 
575  ! Put variable
576  call ncerr(nf90_put_var(col_id,values_id,jcol%values(1:jcol%nlev,:),(/1,1/),(/jcol%nlev,jgrp%nobs/)))
577 
578  ! Update
579  jcol => jcol%next
580  enddo
581  endif
582 
583  ! Update
584  jgrp=>jgrp%next
585 end do
586 
587 ! Close NetCDF file
588 call ncerr(nf90_close(ncid))
589 
590 end subroutine qg_obsdb_write
591 ! ------------------------------------------------------------------------------
592 !> Find observation data group
593 subroutine qg_obsdb_find_group(self,grp,find)
594 
595 implicit none
596 
597 ! Passed variables
598 type(qg_obsdb),intent(in) :: self !< Observation data
599 character(len=*),intent(in) :: grp !< Group
600 type(group_data),pointer,intent(inout) :: find !< Result
601 
602 ! Initialization
603 find => self%grphead
604 
605 ! Loop
606 do while (associated(find))
607  if (find%grpname==grp) exit
608  find => find%next
609 enddo
610 
611 end subroutine qg_obsdb_find_group
612 ! ------------------------------------------------------------------------------
613 !> Find observation data column
614 subroutine qg_obsdb_find_column(grp,col,find)
615 
616 implicit none
617 
618 ! Passed variables
619 type(group_data),intent(in) :: grp !< Observation data
620 character(len=*),intent(in) :: col !< Column
621 type(column_data),pointer,intent(inout) :: find !< Result
622 
623 ! Initialization
624 find=>grp%colhead
625 
626 ! Loop
627 do while (associated(find))
628  if (find%colname==col) exit
629  find => find%next
630 enddo
631 
632 end subroutine qg_obsdb_find_column
633 ! ------------------------------------------------------------------------------
634 !> Generate random locations
635 subroutine qg_obsdb_generate_locations(nlocs,ntimes,bgn,step,times,obsloc)
636 
637 implicit none
638 
639 ! Passed variables
640 integer,intent(in) :: nlocs !< Number of locations
641 integer,intent(in) :: ntimes !< Number of time-slots
642 type(datetime),intent(in) :: bgn !< Start time
643 type(duration),intent(in) :: step !< Time-step
644 type(datetime),intent(inout) :: times(nlocs*ntimes) !< Time-slots
645 type(qg_obsvec),intent(inout) :: obsloc !< Observation locations
646 
647 ! Local variables
648 integer :: jobs,iobs,jstep
649 real(kind_real) :: x(nlocs),y(nlocs),z(nlocs),lon(nlocs),lat(nlocs)
650 type(datetime) :: now
651 
652 ! Generate random locations
653 call uniform_distribution(x,0.0_kind_real,domain_zonal,rseed)
654 call uniform_distribution(y,0.0_kind_real,domain_meridional,rseed)
655 call uniform_distribution(z,0.0_kind_real,domain_depth,rseed)
656 
657 ! Convert to lon/lat
658 do jobs=1,nlocs
659  call xy_to_lonlat(x(jobs),y(jobs),lon(jobs),lat(jobs))
660 enddo
661 
662 ! Setup observation vector
663 call qg_obsvec_setup(obsloc,3,nlocs*ntimes)
664 
665 ! Set observation locations
666 now = bgn
667 iobs=0
668 do jstep=1,ntimes
669  do jobs=1,nlocs
670  iobs = iobs+1
671  times(iobs) = now
672  obsloc%values(:,iobs) = (/lon(jobs),lat(jobs),z(jobs)/)
673  enddo
674  call datetime_update(now,step)
675 enddo
676 
677 ! Release memory
678 call datetime_delete(now)
679 
680 end subroutine qg_obsdb_generate_locations
681 ! ------------------------------------------------------------------------------
682 !> Create observation data
683 subroutine qg_obsdb_create(self,grp,times,locs)
684 
685 implicit none
686 
687 ! Passed varaibles
688 type(qg_obsdb),intent(inout) :: self !< Observation data
689 character(len=*),intent(in) :: grp !< Group
690 type(datetime),intent(in) :: times(:) !< Time-slots
691 type(qg_obsvec),intent(in) :: locs !< Locations
692 
693 ! Local variables
694 type(group_data),pointer :: igrp
695 integer :: jobs,jlev
696 
697 ! Find observation group
698 call qg_obsdb_find_group(self,grp,igrp)
699 if (associated(igrp)) call abor1_ftn('qg_obsdb_create: obs group already exists')
700 if (associated(self%grphead)) then
701  igrp => self%grphead
702  do while (associated(igrp%next))
703  igrp => igrp%next
704  enddo
705  allocate(igrp%next)
706  igrp => igrp%next
707 else
708  allocate(self%grphead)
709  igrp => self%grphead
710 endif
711 
712 ! Create observation data
713 igrp%grpname = grp
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')
723 do jobs=1,igrp%nobs
724  do jlev=1,3
725  igrp%colhead%values(jlev,jobs) = locs%values(jlev,jobs)
726  enddo
727 enddo
728 self%ngrp = self%ngrp+1
729 
730 end subroutine qg_obsdb_create
731 ! ------------------------------------------------------------------------------
732 end module qg_obsdb_mod
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)
Definition: qg_locs_mod.F90:20
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.
subroutine, public ncerr(info)
Check NetCDF status.