IODA
obsspace_mod.F90
Go to the documentation of this file.
1 !
2 ! (C) Copyright 2017 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 
7 !> Fortran interface to ObsSpace.
8 
10 
11 use, intrinsic :: iso_c_binding
12 use kinds
13 use string_f_c_mod
14 use datetime_mod
15 
16 implicit none
17 
18 private
19 public obsspace_construct
20 public obsspace_destruct
21 public obsspace_obsname
24 public obsspace_get_nlocs
26 public obsspace_get_nrecs
27 public obsspace_get_nvars
31 public obsspace_get_comm
33 public obsspace_get_index
34 public obsspace_get_db
35 public obsspace_put_db
36 public obsspace_has
39 
40 #include "ioda/obsspace_interface.f"
41 
42 !-------------------------------------------------------------------------------
43 
44 interface obsspace_get_db
45  module procedure obsspace_get_db_int32
46  module procedure obsspace_get_db_int64
47  module procedure obsspace_get_db_real32
48  module procedure obsspace_get_db_real64
49  module procedure obsspace_get_db_datetime
50 end interface
51 
52 interface obsspace_put_db
53  module procedure obsspace_put_db_int32
54  module procedure obsspace_put_db_int64
55  module procedure obsspace_put_db_real32
56  module procedure obsspace_put_db_real64
57 end interface
58 
59 !-------------------------------------------------------------------------------
60 contains
61 !-------------------------------------------------------------------------------
62 
63 type(c_ptr) function obsspace_construct(c_conf, tbegin, tend)
64  use fckit_configuration_module, only: fckit_configuration
65  use datetime_mod, only: datetime
66  implicit none
67  type(fckit_configuration), intent(in) :: c_conf
68  type(datetime), intent(in) :: tbegin, tend
69  type(c_ptr) :: c_tbegin, c_tend
70 
71  call f_c_datetime(tbegin, c_tbegin)
72  call f_c_datetime(tend, c_tend)
73  obsspace_construct = c_obsspace_construct(c_conf%c_ptr(), c_tbegin, c_tend)
74 end function obsspace_construct
75 
76 !-------------------------------------------------------------------------------
77 
78 subroutine obsspace_destruct(c_obss)
79  implicit none
80  type(c_ptr), intent(inout) :: c_obss
81 
82  call c_obsspace_destruct(c_obss)
83  c_obss = c_null_ptr
84 end subroutine obsspace_destruct
85 
86 !-------------------------------------------------------------------------------
87 
88 !> Get obsname from ObsSpace
89 
90 subroutine obsspace_obsname(obss, obsname)
91  use string_f_c_mod
92  implicit none
93 
94  type(c_ptr), value, intent(in) :: obss
95  character(*), intent(inout) :: obsname
96 
97  integer(c_size_t) :: lcname
98 
99  !< If changing the length of name and cname, need to change the ASSERT in
100  !obsspace_f.cc also
101  character(kind=c_char,len=1) :: cname(101)
102 
103  call c_obsspace_obsname(obss, lcname, cname)
104  call c_f_string(cname, obsname)
105  obsname = obsname(1:lcname)
106 
107 end subroutine obsspace_obsname
108 
109 !-------------------------------------------------------------------------------
110 
111 !> Get obsvariables from ObsSpace
112 
113 type(oops_variables) function obsspace_obsvariables(obss)
114  use oops_variables_mod
115  implicit none
116  type(c_ptr), value, intent(in) :: obss
117 
118  obsspace_obsvariables = oops_variables(c_obsspace_obsvariables(obss))
119 end function obsspace_obsvariables
120 
121 !-------------------------------------------------------------------------------
122 !> Return the number of observational locations in the input obs file
123 
124 integer function obsspace_get_gnlocs(c_obss)
125  implicit none
126  type(c_ptr), intent(in) :: c_obss
127 
128  ! Implicit conversion from c_size_t to integer which is safe in this case
130 end function obsspace_get_gnlocs
131 
132 !-------------------------------------------------------------------------------
133 !> Return the number of observational locations in the obs space
134 
135 integer function obsspace_get_nlocs(c_obss)
136  implicit none
137  type(c_ptr), intent(in) :: c_obss
138 
139  ! Implicit conversion from c_size_t to integer which is safe in this case
141 end function obsspace_get_nlocs
142 
143 !-------------------------------------------------------------------------------
144 !> Return the number of channels in obs space (zero if conventional obs type)
145 
146 integer function obsspace_get_nchans(c_obss)
147  implicit none
148  type(c_ptr), intent(in) :: c_obss
149 
150  ! Implicit conversion from c_size_t to integer which is safe in this case
152 end function obsspace_get_nchans
153 
154 !-------------------------------------------------------------------------------
155 !> Return the number of observational records (profiles)
156 
157 integer function obsspace_get_nrecs(c_obss)
158  implicit none
159  type(c_ptr), intent(in) :: c_obss
160 
161  ! Implicit conversion from c_size_t to integer which is safe in this case
163 end function obsspace_get_nrecs
164 
165 !-------------------------------------------------------------------------------
166 !> Return the number of observational variables
167 
168 integer function obsspace_get_nvars(c_obss)
169  implicit none
170  type(c_ptr), intent(in) :: c_obss
171 
172  ! Implicit conversion from c_size_t to integer which is safe in this case
174 end function obsspace_get_nvars
175 
176 !-------------------------------------------------------------------------------
177 !> Return the ObsSpace dimension name given the dimension id
178 
179 subroutine obsspace_get_dim_name(obss, dim_id, dim_name)
180  use string_f_c_mod
181  implicit none
182 
183  type(c_ptr), value, intent(in) :: obss
184  integer, intent(in) :: dim_id
185  character(*), intent(inout) :: dim_name
186 
187  integer(c_size_t) :: len_dim_name
188 
189  !< If changing the length of dim_name and c_dim_name, need to change the ASSERT in
190  !obsspace_f.cc also
191  character(kind=c_char,len=1) :: c_dim_name(101)
192 
193  call c_obsspace_get_dim_name(obss, dim_id, len_dim_name, c_dim_name)
194  call c_f_string(c_dim_name, dim_name)
195  dim_name = dim_name(1:len_dim_name)
196 
197 end subroutine obsspace_get_dim_name
198 
199 !-------------------------------------------------------------------------------
200 !> Return the size of the ObsSpace dimension given the dimension id
201 
202 integer function obsspace_get_dim_size(obss, dim_id)
203  implicit none
204  type(c_ptr), intent(in) :: obss
205  integer, intent(in) :: dim_id
206 
207  ! Implicit conversion from c_size_t to integer which is safe in this case
209 end function obsspace_get_dim_size
210 
211 !-------------------------------------------------------------------------------
212 !> Return the ObsSpace dimension id given the dimension name
213 
214 integer function obsspace_get_dim_id(obss, dim_name)
215  use string_f_c_mod
216  implicit none
217  type(c_ptr), intent(in) :: obss
218  character(len=*), intent(in) :: dim_name
219 
220  character(kind=c_char,len=1), allocatable :: c_dim_name(:)
221  call f_c_string(dim_name, c_dim_name)
222 
223  ! Implicit conversion from c_size_t to integer which is safe in this case
224  obsspace_get_dim_id = c_obsspace_get_dim_id(obss, c_dim_name)
225 end function obsspace_get_dim_id
226 
227 !-------------------------------------------------------------------------------
228 !> Return the name and name length of obsspace communicator
229 subroutine obsspace_get_comm(obss, f_comm)
230  use fckit_mpi_module, only: fckit_mpi_comm
231  use string_f_c_mod
232  implicit none
233 
234  type(c_ptr), intent(in) :: obss
235  type(fckit_mpi_comm),intent(out) :: f_comm
236 
237  integer :: lcname
238  !< If changing the length of name and cname, need to change the ASSERT in obsspace_f.cc also
239  character(kind=c_char,len=1) :: cname(101)
240  character(len=100) :: name
241  character(len=:), allocatable :: name_comm
242 
243  call c_obsspace_get_comm(obss, lcname, cname)
244  call c_f_string(cname, name)
245 
246  name_comm = name(1:lcname)
247  f_comm = fckit_mpi_comm(name_comm)
248 end subroutine obsspace_get_comm
249 
250 !-------------------------------------------------------------------------------
251 !> Get the record number vector
252 subroutine obsspace_get_recnum(obss, recnum)
253  implicit none
254  type(c_ptr), intent(in) :: obss
255  integer(c_size_t), intent(inout) :: recnum(:)
256 
257  integer(c_size_t) :: length
258 
259  length = size(recnum)
260  call c_obsspace_get_recnum(obss, length, recnum)
261 end subroutine obsspace_get_recnum
262 
263 !-------------------------------------------------------------------------------
264 !> Get the index vector
265 subroutine obsspace_get_index(obss, indx)
266  implicit none
267  type(c_ptr), intent(in) :: obss
268  integer(c_size_t), intent(inout) :: indx(:)
269 
270  integer(c_size_t) :: length
271 
272  length = size(indx)
273  call c_obsspace_get_index(obss, length, indx)
274 end subroutine obsspace_get_index
275 
276 !-------------------------------------------------------------------------------
277 
278 !> Return true if variable exists in database
279 
280 logical function obsspace_has(c_obss, group, vname)
281  implicit none
282  type(c_ptr), intent(in) :: c_obss
283  character(len=*), intent(in) :: group
284  character(len=*), intent(in) :: vname
285 
286  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
287 
288  call f_c_string(group, c_group)
289  call f_c_string(vname, c_vname)
290  obsspace_has = c_obsspace_has(c_obss, c_group, c_vname)
291 end function obsspace_has
292 
293 !-------------------------------------------------------------------------------
294 
295 !> Get a variable from the ObsSapce database
296 
297 subroutine obsspace_get_db_int32(obss, group, vname, vect, chan_select)
298  implicit none
299  type(c_ptr), value, intent(in) :: obss
300  character(len=*), intent(in) :: group
301  character(len=*), intent(in) :: vname
302  integer(c_int32_t), intent(inout) :: vect(:)
303  integer(c_int), intent(in), optional :: chan_select(:)
304 
305  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
306  integer(c_size_t) :: length
307  integer(c_size_t) :: len_cs
308  integer(c_int) :: dummy_chan_select(1) ! used as a fallback if chan_select is not present
309 
310  ! Translate query from Fortran string to C++ char[].
311  call f_c_string(group, c_group)
312  call f_c_string(vname, c_vname)
313  length = size(vect)
314  if (present(chan_select)) then
315  len_cs = size(chan_select)
316  call c_obsspace_get_int32(obss, c_group, c_vname, length, vect, len_cs, chan_select)
317  else
318  ! Note: we say that the number of channels is zero, which means that the contents of
319  ! dummy_chan_select won't actually be accessed.
320  call c_obsspace_get_int32(obss, c_group, c_vname, length, vect, 0_c_size_t, dummy_chan_select)
321  endif
322 
323  deallocate(c_group, c_vname)
324 end subroutine obsspace_get_db_int32
325 
326 
327 !-------------------------------------------------------------------------------
328 
329 !> Get a variable from the ObsSapce database
330 
331 subroutine obsspace_get_db_int64(obss, group, vname, vect, chan_select)
332  implicit none
333  type(c_ptr), value, intent(in) :: obss
334  character(len=*), intent(in) :: group
335  character(len=*), intent(in) :: vname
336  integer(c_int64_t), intent(inout) :: vect(:)
337  integer(c_int), intent(in), optional :: chan_select(:)
338 
339  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
340  integer(c_size_t) :: length
341  integer(c_size_t) :: len_cs
342  integer(c_int) :: dummy_chan_select(1) ! used as a fallback if chan_select is not present
343 
344  ! Translate query from Fortran string to C++ char[].
345  call f_c_string(group, c_group)
346  call f_c_string(vname, c_vname)
347  length = size(vect)
348  if (present(chan_select)) then
349  len_cs = size(chan_select)
350  call c_obsspace_get_int64(obss, c_group, c_vname, length, vect, len_cs, chan_select)
351  else
352  ! Note: we say that the number of channels is zero, which means that the contents of
353  ! dummy_chan_select won't actually be accessed.
354  call c_obsspace_get_int64(obss, c_group, c_vname, length, vect, 0_c_size_t, dummy_chan_select)
355  endif
356 
357  deallocate(c_group, c_vname)
358 end subroutine obsspace_get_db_int64
359 
360 !-------------------------------------------------------------------------------
361 
362 !> Get a variable from the ObsSapce database
363 
364 subroutine obsspace_get_db_real32(obss, group, vname, vect, chan_select)
365  implicit none
366  type(c_ptr), value, intent(in) :: obss
367  character(len=*), intent(in) :: group
368  character(len=*), intent(in) :: vname
369  real(c_float), intent(inout) :: vect(:)
370  integer(c_int), intent(in), optional :: chan_select(:)
371 
372  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
373  integer(c_size_t) :: length
374  integer(c_size_t) :: len_cs
375  integer(c_int) :: dummy_chan_select(1) ! used as a fallback if chan_select is not present
376 
377  ! Translate query from Fortran string to C++ char[].
378  call f_c_string(group, c_group)
379  call f_c_string(vname, c_vname)
380  length = size(vect)
381  if (present(chan_select)) then
382  len_cs = size(chan_select)
383  call c_obsspace_get_real32(obss, c_group, c_vname, length, vect, len_cs, chan_select)
384  else
385  ! Note: we say that the number of channels is zero, which means that the contents of
386  ! dummy_chan_select won't actually be accessed.
387  call c_obsspace_get_real32(obss, c_group, c_vname, length, vect, 0_c_size_t, dummy_chan_select)
388  endif
389 
390  deallocate(c_group, c_vname)
391 end subroutine obsspace_get_db_real32
392 
393 !-------------------------------------------------------------------------------
394 
395 !> Get a variable from the ObsSapce database
396 
397 subroutine obsspace_get_db_real64(obss, group, vname, vect, chan_select)
398  implicit none
399  type(c_ptr), value, intent(in) :: obss
400  character(len=*), intent(in) :: group
401  character(len=*), intent(in) :: vname
402  real(c_double), intent(inout) :: vect(:)
403  ! optional; if omitted, all channels will be retrieved
404  integer(c_int), intent(in), optional :: chan_select(:)
405 
406  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
407  integer(c_size_t) :: length
408  integer(c_size_t) :: len_cs
409  integer(c_int) :: dummy_chan_select(1) ! used as a fallback if chan_select is not present
410 
411  ! Translate query from Fortran string to C++ char[].
412  call f_c_string(group, c_group)
413  call f_c_string(vname, c_vname)
414  length = size(vect)
415 
416  if (present(chan_select)) then
417  len_cs = size(chan_select)
418  call c_obsspace_get_real64(obss, c_group, c_vname, length, vect, len_cs, chan_select)
419  else
420  ! Note: we say that the number of channels is zero, which means that the contents of
421  ! dummy_chan_select won't actually be accessed.
422  call c_obsspace_get_real64(obss, c_group, c_vname, length, vect, 0_c_size_t, dummy_chan_select)
423  endif
424 
425  deallocate(c_group, c_vname)
426 end subroutine obsspace_get_db_real64
427 
428 !-------------------------------------------------------------------------------
429 
430 !> Get datetime from the ObsSapce database
431 
432 subroutine obsspace_get_db_datetime(obss, group, vname, vect, chan_select)
433  implicit none
434  type(c_ptr), value, intent(in) :: obss
435  character(len=*), intent(in) :: group
436  character(len=*), intent(in) :: vname
437  type(datetime), intent(inout) :: vect(:)
438  integer(c_int), intent(in), optional :: chan_select(:)
439 
440  integer(c_size_t) :: length, i
441  integer(c_size_t) :: len_cs
442  integer(c_int) :: dummy_chan_select(1) ! used as a fallback if chan_select is not present
443  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
444  integer(c_int32_t), dimension(:), allocatable :: date
445  integer(c_int32_t), dimension(:), allocatable :: time
446  character(len=20) :: fstring
447 
448  call f_c_string(group, c_group)
449  call f_c_string(vname, c_vname)
450  length = size(vect)
451 
452  allocate(date(length), time(length))
453 
454  if (present(chan_select)) then
455  len_cs = size(chan_select)
456  call c_obsspace_get_datetime(obss, c_group, c_vname, length, date, time, len_cs, chan_select)
457  else
458  ! Note: we say that the number of channels is zero, which means that the contents of
459  ! dummy_chan_select won't actually be accessed.
460  call c_obsspace_get_datetime(obss, c_group, c_vname, length, date, time, &
461  0_c_size_t, dummy_chan_select)
462  endif
463 
464  ! Constrct datatime based on date and time
465  do i = 1, length
466  write(fstring, "(i4.4, a, i2.2, a, i2.2, a, i2.2, a, i2.2, a, i2.2, a)") &
467  date(i)/10000, '-', mod(date(i), 10000)/100, '-', mod(mod(date(i), 10000), 100), 'T', &
468  time(i)/10000, ':', mod(time(i), 10000)/100, ':', mod(mod(time(i), 10000), 100), 'Z'
469  call datetime_create(fstring, vect(i))
470  enddo
471 
472  ! Clean up
473  deallocate(date, time)
474 end subroutine obsspace_get_db_datetime
475 
476 !-------------------------------------------------------------------------------
477 
478 !> Store a vector in ObsSpace database
479 
480 subroutine obsspace_put_db_int32(obss, group, vname, vect, dim_ids)
481  implicit none
482  type(c_ptr), value, intent(in) :: obss
483  character(len=*), intent(in) :: group
484  character(len=*), intent(in) :: vname
485  integer(c_int32_t), intent(in) :: vect(:)
486  integer(c_int), intent(in), optional :: dim_ids(:) ! if not set, defaults to nlocs
487 
488  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
489  integer(c_size_t) :: length
490  integer(c_size_t) :: ndims
491  integer(c_int) :: fallback_dim_ids(1)
492 
493  ! Translate query from Fortran string to C++ char[].
494  call f_c_string(group, c_group)
495  call f_c_string(vname, c_vname)
496  length = size(vect)
497  if (present(dim_ids)) then
498  ndims = size(dim_ids)
499  call c_obsspace_put_int32(obss, c_group, c_vname, length, vect, ndims, dim_ids)
500  else
501  ndims = 1
502  fallback_dim_ids = (/ obsspace_get_nlocs_dim_id() /)
503  call c_obsspace_put_int32(obss, c_group, c_vname, length, vect, ndims, fallback_dim_ids)
504  endif
505 
506  deallocate(c_group, c_vname)
507 end subroutine obsspace_put_db_int32
508 
509 !-------------------------------------------------------------------------------
510 
511 !> Store a vector in ObsSpace database
512 
513 subroutine obsspace_put_db_int64(obss, group, vname, vect, dim_ids)
514  implicit none
515  type(c_ptr), value, intent(in) :: obss
516  character(len=*), intent(in) :: group
517  character(len=*), intent(in) :: vname
518  integer(c_int64_t), intent(in) :: vect(:)
519  integer(c_int), intent(in), optional :: dim_ids(:) ! if not set, defaults to nlocs
520 
521  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
522  integer(c_size_t) :: length
523  integer(c_size_t) :: ndims
524  integer(c_int) :: fallback_dim_ids(1)
525 
526  ! Translate query from Fortran string to C++ char[].
527  call f_c_string(group, c_group)
528  call f_c_string(vname, c_vname)
529  length = size(vect)
530 
531  if (present(dim_ids)) then
532  ndims = size(dim_ids)
533  call c_obsspace_put_int64(obss, c_group, c_vname, length, vect, ndims, dim_ids)
534  else
535  ndims = 1
536  fallback_dim_ids = (/ obsspace_get_nlocs_dim_id() /)
537  call c_obsspace_put_int64(obss, c_group, c_vname, length, vect, ndims, fallback_dim_ids)
538  endif
539 
540  deallocate(c_group, c_vname)
541 end subroutine obsspace_put_db_int64
542 
543 !-------------------------------------------------------------------------------
544 
545 !> Store a vector in ObsSpace database
546 
547 subroutine obsspace_put_db_real32(obss, group, vname, vect, dim_ids)
548  implicit none
549  type(c_ptr), value, intent(in) :: obss
550  character(len=*), intent(in) :: group
551  character(len=*), intent(in) :: vname
552  real(c_float), intent(in) :: vect(:)
553  integer(c_int), intent(in), optional :: dim_ids(:) ! if not set, defaults to nlocs
554 
555  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
556  integer(c_size_t) :: length
557  integer(c_size_t) :: ndims
558  integer(c_int) :: fallback_dim_ids(1)
559 
560  ! Translate query from Fortran string to C++ char[].
561  call f_c_string(group, c_group)
562  call f_c_string(vname, c_vname)
563  length = size(vect)
564  ndims = size(dim_ids)
565 
566  if (present(dim_ids)) then
567  ndims = size(dim_ids)
568  call c_obsspace_put_real32(obss, c_group, c_vname, length, vect, ndims, dim_ids)
569  else
570  ndims = 1
571  fallback_dim_ids = (/ obsspace_get_nlocs_dim_id() /)
572  call c_obsspace_put_real32(obss, c_group, c_vname, length, vect, ndims, fallback_dim_ids)
573  endif
574 
575  deallocate(c_group, c_vname)
576 end subroutine obsspace_put_db_real32
577 
578 !-------------------------------------------------------------------------------
579 
580 !> Store a vector in ObsSpace database
581 
582 subroutine obsspace_put_db_real64(obss, group, vname, vect, dim_ids)
583  implicit none
584  type(c_ptr), value, intent(in) :: obss
585  character(len=*), intent(in) :: group
586  character(len=*), intent(in) :: vname
587  real(c_double), intent(in) :: vect(:)
588  integer(c_int), intent(in), optional :: dim_ids(:) ! if not set, defaults to nlocs
589 
590  character(kind=c_char,len=1), allocatable :: c_group(:), c_vname(:)
591  integer(c_size_t) :: length
592  integer(c_size_t) :: ndims
593  integer(c_int) :: fallback_dim_ids(1)
594 
595  ! Translate query from Fortran string to C++ char[].
596  call f_c_string(group, c_group)
597  call f_c_string(vname, c_vname)
598  length = size(vect)
599 
600  if (present(dim_ids)) then
601  ndims = size(dim_ids)
602  call c_obsspace_put_real64(obss, c_group, c_vname, length, vect, ndims, dim_ids)
603  else
604  ndims = 1
605  fallback_dim_ids = (/ obsspace_get_nlocs_dim_id() /)
606  call c_obsspace_put_real64(obss, c_group, c_vname, length, vect, ndims, fallback_dim_ids)
607  endif
608 
609  deallocate(c_group, c_vname)
610 end subroutine obsspace_put_db_real64
611 
612 !-------------------------------------------------------------------------------
613 
614 !> Return the identifier of the nlocs dimension.
615 integer(c_int) function obsspace_get_nlocs_dim_id()
616  implicit none
618 end function obsspace_get_nlocs_dim_id
619 
620 !-------------------------------------------------------------------------------
621 
622 !> Return the identifier of the nchans dimension.
623 integer(c_int) function obsspace_get_nchans_dim_id()
624  implicit none
626 end function obsspace_get_nchans_dim_id
627 
628 end module obsspace_mod
Define interface for C++ ObsSpace code called from Fortran.
Fortran interface to ObsSpace.
Definition: obsspace_mod.F90:9
type(oops_variables) function, public obsspace_obsvariables(obss)
Get obsvariables from ObsSpace.
integer function, public obsspace_get_gnlocs(c_obss)
Return the number of observational locations in the input obs file.
subroutine obsspace_put_db_real64(obss, group, vname, vect, dim_ids)
Store a vector in ObsSpace database.
integer(c_int) function, public obsspace_get_nchans_dim_id()
Return the identifier of the nchans dimension.
subroutine, public obsspace_obsname(obss, obsname)
Get obsname from ObsSpace.
subroutine obsspace_get_db_int64(obss, group, vname, vect, chan_select)
Get a variable from the ObsSapce database.
subroutine obsspace_get_db_real64(obss, group, vname, vect, chan_select)
Get a variable from the ObsSapce database.
integer function, public obsspace_get_nlocs(c_obss)
Return the number of observational locations in the obs space.
integer(c_int) function, public obsspace_get_nlocs_dim_id()
Return the identifier of the nlocs dimension.
integer function, public obsspace_get_dim_id(obss, dim_name)
Return the ObsSpace dimension id given the dimension name.
subroutine, public obsspace_get_index(obss, indx)
Get the index vector.
subroutine, public obsspace_get_dim_name(obss, dim_id, dim_name)
Return the ObsSpace dimension name given the dimension id.
integer function, public obsspace_get_nchans(c_obss)
Return the number of channels in obs space (zero if conventional obs type)
subroutine obsspace_get_db_int32(obss, group, vname, vect, chan_select)
Get a variable from the ObsSapce database.
integer function, public obsspace_get_nrecs(c_obss)
Return the number of observational records (profiles)
subroutine obsspace_put_db_real32(obss, group, vname, vect, dim_ids)
Store a vector in ObsSpace database.
subroutine obsspace_get_db_datetime(obss, group, vname, vect, chan_select)
Get datetime from the ObsSapce database.
subroutine, public obsspace_get_recnum(obss, recnum)
Get the record number vector.
subroutine obsspace_put_db_int64(obss, group, vname, vect, dim_ids)
Store a vector in ObsSpace database.
subroutine obsspace_get_db_real32(obss, group, vname, vect, chan_select)
Get a variable from the ObsSapce database.
logical function, public obsspace_has(c_obss, group, vname)
Return true if variable exists in database.
subroutine obsspace_put_db_int32(obss, group, vname, vect, dim_ids)
Store a vector in ObsSpace database.
integer function, public obsspace_get_nvars(c_obss)
Return the number of observational variables.
integer function, public obsspace_get_dim_size(obss, dim_id)
Return the size of the ObsSpace dimension given the dimension id.
subroutine, public obsspace_destruct(c_obss)
subroutine, public obsspace_get_comm(obss, f_comm)
Return the name and name length of obsspace communicator.
type(c_ptr) function, public obsspace_construct(c_conf, tbegin, tend)