OOPS
qg_obsvec_interface.F90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
2 ! (C) Copyright 2017-2021 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 iso_c_binding
13 use qg_obsvec_mod
14 
15 implicit none
16 
17 private
18 ! ------------------------------------------------------------------------------
19 contains
20 ! ------------------------------------------------------------------------------
21 !> Setup observation vector
22 subroutine qg_obsvec_setup_c(c_key_self,nlev,nobs) bind(c,name='qg_obsvec_setup_f90')
23 
24 implicit none
25 
26 ! Passed variables
27 integer(c_int),intent(inout) :: c_key_self !< Observation vector
28 integer(c_int),intent(in) :: nlev !< Number of levels
29 integer(c_int),intent(in) :: nobs !< Number of observations
30 
31 ! Local variables
32 type(qg_obsvec),pointer :: self
33 
34 ! Interface
35 call qg_obsvec_registry%init()
36 call qg_obsvec_registry%add(c_key_self)
37 call qg_obsvec_registry%get(c_key_self,self)
38 
39 ! Call Fortran
40 call qg_obsvec_setup(self,nlev,nobs)
41 
42 end subroutine qg_obsvec_setup_c
43 ! ------------------------------------------------------------------------------
44 !> Clone observation vector
45 subroutine qg_obsvec_clone_c(c_key_self,c_key_other) bind(c,name='qg_obsvec_clone_f90')
46 
47 implicit none
48 
49 ! Passed variables
50 integer(c_int),intent(inout) :: c_key_self !< Observation vector
51 integer(c_int),intent(in) :: c_key_other !< Other observation vector
52 
53 ! Local variables
54 type(qg_obsvec),pointer :: self,other
55 
56 ! Interface
57 call qg_obsvec_registry%get(c_key_other,other)
58 call qg_obsvec_registry%init()
59 call qg_obsvec_registry%add(c_key_self)
60 call qg_obsvec_registry%get(c_key_self,self)
61 
62 ! Call Fortran
63 call qg_obsvec_clone(self,other)
64 
65 end subroutine qg_obsvec_clone_c
66 ! ------------------------------------------------------------------------------
67 !> Delete observation vector
68 subroutine qg_obsvec_delete_c(c_key_self) bind(c,name='qg_obsvec_delete_f90')
69 
70 implicit none
71 
72 ! Passed variables
73 integer(c_int),intent(inout) :: c_key_self !< Observation vector
74 
75 ! Local variables
76 type(qg_obsvec),pointer :: self
77 
78 ! Interface
79 call qg_obsvec_registry%get(c_key_self,self)
80 
81 ! Call Fortran
82 call qg_obsvec_delete(self)
83 
84 ! Clear interface
85 call qg_obsvec_registry%remove(c_key_self)
86 
87 end subroutine qg_obsvec_delete_c
88 ! ------------------------------------------------------------------------------
89 !> Copy observation vector
90 subroutine qg_obsvec_copy_c(c_key_self,c_key_other) bind(c,name='qg_obsvec_copy_f90')
91 
92 implicit none
93 
94 ! Passed variables
95 integer(c_int),intent(in) :: c_key_self !< Observation vector
96 integer(c_int),intent(in) :: c_key_other !< Other observation vector
97 
98 ! Local variables
99 type(qg_obsvec),pointer :: self,other
100 
101 ! Interface
102 call qg_obsvec_registry%get(c_key_self,self)
103 call qg_obsvec_registry%get(c_key_other,other)
104 
105 ! Call Fortran
106 call qg_obsvec_copy(self,other)
107 
108 end subroutine qg_obsvec_copy_c
109 ! ------------------------------------------------------------------------------
110 !> Set observation vector to zero
111 subroutine qg_obsvec_zero_c(c_key_self) bind(c,name='qg_obsvec_zero_f90')
112 
113 implicit none
114 
115 ! Passed variables
116 integer(c_int),intent(in) :: c_key_self !< Observation vector
117 
118 ! Local variables
119 type(qg_obsvec),pointer :: self
120 
121 ! Interface
122 call qg_obsvec_registry%get(c_key_self,self)
123 
124 ! Call Fortran
125 call qg_obsvec_zero(self)
126 
127 end subroutine qg_obsvec_zero_c
128 ! ------------------------------------------------------------------------------
129 !> Set i-th value of the observation vector to missing value
130 subroutine qg_obsvec_settomissing_ith_c(c_key_self, i) bind(c,name='qg_obsvec_settomissing_ith_f90')
131 
132 implicit none
133 
134 ! Passed variables
135 integer(c_int),intent(in) :: c_key_self !< Observation vector
136 integer(c_int),intent(in) :: i !< index of value to be set to missing value
137 ! Local variables
138 type(qg_obsvec),pointer :: self
139 
140 ! Interface
141 call qg_obsvec_registry%get(c_key_self,self)
142 
143 ! Call Fortran
144 ! increase index by 1 (C indices start with 0; Fortran indices start with 1)
145 call qg_obsvec_settomissing_ith(self, i+1)
146 
147 end subroutine qg_obsvec_settomissing_ith_c
148 ! ------------------------------------------------------------------------------
149 !> Set observation vector to ones
150 subroutine qg_obsvec_ones_c(c_key_self) bind(c,name='qg_obsvec_ones_f90')
151 
152 implicit none
153 
154 ! Passed variables
155 integer(c_int),intent(in) :: c_key_self !< Observation vector
156 
157 ! Local variables
158 type(qg_obsvec),pointer :: self
159 
160 ! Interface
161 call qg_obsvec_registry%get(c_key_self,self)
162 
163 ! Call Fortran
164 call qg_obsvec_ones(self)
165 
166 end subroutine qg_obsvec_ones_c
167 ! ------------------------------------------------------------------------------
168 !> Mask self observation vector (set values to missing where mask is set)
169 subroutine qg_obsvec_mask_c(c_key_self,c_key_mask) bind(c,name='qg_obsvec_mask_f90')
170 
171 implicit none
172 
173 ! Passed variables
174 integer(c_int),intent(in) :: c_key_self !< Observation vector
175 integer(c_int),intent(in) :: c_key_mask !< Mask
176 
177 ! Local variables
178 type(qg_obsvec),pointer :: self,mask
179 
180 ! Interface
181 call qg_obsvec_registry%get(c_key_self,self)
182 call qg_obsvec_registry%get(c_key_mask,mask)
183 
184 ! Call Fortran
185 call qg_obsvec_mask(self,mask)
186 
187 end subroutine qg_obsvec_mask_c
188 ! ------------------------------------------------------------------------------
189 !> Mask self observation vector (set values to missing where mask is a missing value)
190 subroutine qg_obsvec_mask_with_missing_c(c_key_self,c_key_mask) bind(c,name='qg_obsvec_mask_with_missing_f90')
191 
192 implicit none
193 
194 ! Passed variables
195 integer(c_int),intent(in) :: c_key_self !< Observation vector
196 integer(c_int),intent(in) :: c_key_mask !< Mask
197 
198 ! Local variables
199 type(qg_obsvec),pointer :: self,mask
200 
201 ! Interface
202 call qg_obsvec_registry%get(c_key_self,self)
203 call qg_obsvec_registry%get(c_key_mask,mask)
204 
205 ! Call Fortran
206 call qg_obsvec_mask_with_missing(self,mask)
207 
208 end subroutine qg_obsvec_mask_with_missing_c
209 ! ------------------------------------------------------------------------------
210 !> Multiply observation vector with a scalar
211 subroutine qg_obsvec_mul_scal_c(c_key_self,zz) bind(c,name='qg_obsvec_mul_scal_f90')
212 
213 implicit none
214 
215 ! Passed variables
216 integer(c_int),intent(in) :: c_key_self !< Observation vector
217 real(c_double),intent(in) :: zz !< Multiplier
218 
219 ! Local variables
220 type(qg_obsvec),pointer :: self
221 
222 ! Interface
223 call qg_obsvec_registry%get(c_key_self,self)
224 
225 ! Call Fortran
226 call qg_obsvec_mul_scal(self,zz)
227 
228 end subroutine qg_obsvec_mul_scal_c
229 ! ------------------------------------------------------------------------------
230 !> Add observation vector
231 subroutine qg_obsvec_add_c(c_key_self,c_key_other) bind(c,name='qg_obsvec_add_f90')
232 
233 implicit none
234 
235 ! Passed variables
236 integer(c_int),intent(in) :: c_key_self !< Observation vector
237 integer(c_int),intent(in) :: c_key_other !< Other observation vector
238 
239 ! Local variables
240 type(qg_obsvec),pointer :: self,other
241 
242 ! Interface
243 call qg_obsvec_registry%get(c_key_self,self)
244 call qg_obsvec_registry%get(c_key_other,other)
245 
246 ! Call Fortran
247 call qg_obsvec_add(self,other)
248 
249 end subroutine qg_obsvec_add_c
250 ! ------------------------------------------------------------------------------
251 !> Subtract observation vector
252 subroutine qg_obsvec_sub_c(c_key_self,c_key_other) bind(c,name='qg_obsvec_sub_f90')
253 
254 implicit none
255 
256 ! Passed variables
257 integer(c_int),intent(in) :: c_key_self !< Observation vector
258 integer(c_int),intent(in) :: c_key_other !< Other observation vector
259 
260 ! Local variables
261 type(qg_obsvec),pointer :: self,other
262 
263 ! Interface
264 call qg_obsvec_registry%get(c_key_self,self)
265 call qg_obsvec_registry%get(c_key_other,other)
266 
267 ! Call Fortran
268 call qg_obsvec_sub(self,other)
269 
270 end subroutine qg_obsvec_sub_c
271 ! ------------------------------------------------------------------------------
272 !> Multiply observation vector
273 subroutine qg_obsvec_mul_c(c_key_self,c_key_other) bind(c,name='qg_obsvec_mul_f90')
274 
275 implicit none
276 
277 ! Passed variables
278 integer(c_int),intent(in) :: c_key_self !< Observation vector
279 integer(c_int),intent(in) :: c_key_other !< Other observation vector
280 
281 ! Local variables
282 type(qg_obsvec),pointer :: self,other
283 
284 ! Interface
285 call qg_obsvec_registry%get(c_key_self,self)
286 call qg_obsvec_registry%get(c_key_other,other)
287 
288 ! Call Fortran
289 call qg_obsvec_mul(self,other)
290 
291 end subroutine qg_obsvec_mul_c
292 ! ------------------------------------------------------------------------------
293 !> Divide observation vector
294 subroutine qg_obsvec_div_c(c_key_self,c_key_other) bind(c,name='qg_obsvec_div_f90')
295 
296 implicit none
297 
298 ! Passed variables
299 integer(c_int),intent(in) :: c_key_self !< Observation vector
300 integer(c_int),intent(in) :: c_key_other !< Other observation vector
301 
302 ! Local variables
303 type(qg_obsvec),pointer :: self,other
304 
305 ! Interface
306 call qg_obsvec_registry%get(c_key_self,self)
307 call qg_obsvec_registry%get(c_key_other,other)
308 
309 ! Call Fortran
310 call qg_obsvec_div(self,other)
311 
312 end subroutine qg_obsvec_div_c
313 ! ------------------------------------------------------------------------------
314 !> Apply axpy on observation vector
315 subroutine qg_obsvec_axpy_c(c_key_self,zz,c_key_other) bind(c,name='qg_obsvec_axpy_f90')
316 
317 implicit none
318 
319 ! Passed variables
320 integer(c_int),intent(in) :: c_key_self !< Observation vector
321 real(c_double),intent(in) :: zz !< Multiplier
322 integer(c_int),intent(in) :: c_key_other !< Other observation vector
323 
324 ! Local variables
325 type(qg_obsvec),pointer :: self,other
326 
327 ! Interface
328 call qg_obsvec_registry%get(c_key_self,self)
329 call qg_obsvec_registry%get(c_key_other,other)
330 
331 ! Call Fortran
332 call qg_obsvec_axpy(self,zz,other)
333 
334 end subroutine qg_obsvec_axpy_c
335 ! ------------------------------------------------------------------------------
336 !> Invert observation vector
337 subroutine qg_obsvec_invert_c(c_key_self) bind(c,name='qg_obsvec_invert_f90')
338 
339 implicit none
340 
341 ! Passed variables
342 integer(c_int),intent(in) :: c_key_self !< Observation vector
343 
344 ! Local variables
345 type(qg_obsvec),pointer :: self
346 
347 ! Interface
348 call qg_obsvec_registry%get(c_key_self,self)
349 
350 ! Call Fortran
351 call qg_obsvec_invert(self)
352 
353 end subroutine qg_obsvec_invert_c
354 ! ------------------------------------------------------------------------------
355 !> Generate random observation vector
356 subroutine qg_obsvec_random_c(c_odb,c_self) bind(c,name='qg_obsvec_random_f90')
357 
358 implicit none
359 
360 ! Passed variables
361 type(c_ptr),intent(in) :: c_odb !< Observation data base
362 integer(c_int),intent(in) :: c_self !< Observation vector
363 
364 ! Local variables
365 type(qg_obsvec),pointer :: self
366 
367 ! Interface
368 call qg_obsvec_registry%get(c_self,self)
369 
370 ! Call Fortran
371 call qg_obsvec_random(c_odb,self)
372 
373 end subroutine qg_obsvec_random_c
374 ! ------------------------------------------------------------------------------
375 !> Compute dot product between observation vectors
376 subroutine qg_obsvec_dotprod_c(c_key_obsvec1,c_key_obsvec2,zz) bind(c,name='qg_obsvec_dotprod_f90')
377 
378 implicit none
379 
380 ! Passed variables
381 integer(c_int),intent(in) :: c_key_obsvec1 !< Observation vector 1
382 integer(c_int),intent(in) :: c_key_obsvec2 !< Observation vector 2
383 real(c_double),intent(inout) :: zz !< Dot product
384 
385 ! Local variables
386 type(qg_obsvec),pointer :: obsvec1,obsvec2
387 
388 ! Interface
389 call qg_obsvec_registry%get(c_key_obsvec1,obsvec1)
390 call qg_obsvec_registry%get(c_key_obsvec2,obsvec2)
391 
392 ! Call Fortran
393 call qg_obsvec_dotprod(obsvec1,obsvec2,zz)
394 
395 end subroutine qg_obsvec_dotprod_c
396 ! ------------------------------------------------------------------------------
397 !> Compute observation vector statistics
398 subroutine qg_obsvec_stats_c(c_key_self,zmin,zmax,zavg) bind(c,name='qg_obsvec_stats_f90')
399 
400 implicit none
401 
402 ! Passed variables
403 integer(c_int),intent(in) :: c_key_self !< Observation vector
404 real(c_double),intent(inout) :: zmin !< Minimum
405 real(c_double),intent(inout) :: zmax !< Maximum
406 real(c_double),intent(inout) :: zavg !< Average
407 
408 ! Local variables
409 type(qg_obsvec),pointer :: self
410 
411 ! Interface
412 call qg_obsvec_registry%get(c_key_self,self)
413 
414 ! Call Fortran
415 call qg_obsvec_stats(self,zmin,zmax,zavg)
416 
417 end subroutine qg_obsvec_stats_c
418 ! ------------------------------------------------------------------------------
419 !> Get number of observations (not missing)
420 subroutine qg_obsvec_nobs_c(c_key_self,kobs) bind(c,name='qg_obsvec_nobs_f90')
421 
422 implicit none
423 
424 ! Passed variables
425 integer(c_int),intent(in) :: c_key_self !< Observation vector
426 integer(c_int),intent(inout) :: kobs !< Observation vector size
427 
428 ! Local vector
429 type(qg_obsvec),pointer :: self
430 
431 ! Interface
432 call qg_obsvec_registry%get(c_key_self,self)
433 
434 ! Call Fortran
435 call qg_obsvec_nobs(self,kobs)
436 
437 end subroutine qg_obsvec_nobs_c
438 ! ------------------------------------------------------------------------------
439 !> Get observation vector size
440 subroutine qg_obsvec_size_c(c_key_self,kobs) bind(c,name='qg_obsvec_size_f90')
441 
442 implicit none
443 
444 ! Passed variables
445 integer(c_int),intent(in) :: c_key_self !< Observation vector
446 integer(c_int),intent(inout) :: kobs !< Observation vector size
447 
448 ! Local vector
449 type(qg_obsvec),pointer :: self
450 
451 ! Interface
452 call qg_obsvec_registry%get(c_key_self,self)
453 
454 ! Call Fortran
455 call qg_obsvec_size(self,kobs)
456 
457 end subroutine qg_obsvec_size_c
458 ! ------------------------------------------------------------------------------
459 !> Get observation vector size (only non-masked observations)
460 subroutine qg_obsvec_nobs_withmask_c(c_key_self,c_key_mask,kobs) bind(c,name='qg_obsvec_nobs_withmask_f90')
461 
462 implicit none
463 
464 ! Passed variables
465 integer(c_int),intent(in) :: c_key_self !< Observation vector
466 integer(c_int),intent(in) :: c_key_mask !< Mask
467 integer(c_int),intent(inout) :: kobs !< Observation vector size
468 
469 ! Local vector
470 type(qg_obsvec),pointer :: self, mask
471 
472 ! Interface
473 call qg_obsvec_registry%get(c_key_self,self)
474 call qg_obsvec_registry%get(c_key_mask,mask)
475 
476 ! Call Fortran
477 call qg_obsvec_nobs_withmask(self,mask,kobs)
478 
479 end subroutine qg_obsvec_nobs_withmask_c
480 
481 ! ------------------------------------------------------------------------------
482 !> Get all non-masked out observation values
483 subroutine qg_obsvec_get_withmask_c(c_key_self,c_key_mask,vals,nvals) bind(c,name='qg_obsvec_get_withmask_f90')
484 
485 implicit none
486 
487 ! Passed variables
488 integer(c_int),intent(in) :: c_key_self !< Observation vector
489 integer(c_int),intent(in) :: c_key_mask !< Mask
490 integer(c_int),intent(in) :: nvals !< number of obs
491 real(c_double),intent(out),dimension(nvals) :: vals !< ob. values
492 
493 ! Local vector
494 type(qg_obsvec),pointer :: self, mask
495 
496 ! Interface
497 call qg_obsvec_registry%get(c_key_self,self)
498 call qg_obsvec_registry%get(c_key_mask,mask)
499 
500 ! Call Fortran
501 call qg_obsvec_get_withmask(self,mask,vals,nvals)
502 
503 end subroutine qg_obsvec_get_withmask_c
504 
505 ! ------------------------------------------------------------------------------
506 end module qg_obsvec_interface
subroutine qg_obsvec_invert_c(c_key_self)
Invert observation vector.
subroutine qg_obsvec_ones_c(c_key_self)
Set observation vector to ones.
subroutine qg_obsvec_clone_c(c_key_self, c_key_other)
Clone observation vector.
subroutine qg_obsvec_mask_with_missing_c(c_key_self, c_key_mask)
Mask self observation vector (set values to missing where mask is a missing value)
subroutine qg_obsvec_get_withmask_c(c_key_self, c_key_mask, vals, nvals)
Get all non-masked out observation values.
subroutine qg_obsvec_stats_c(c_key_self, zmin, zmax, zavg)
Compute observation vector statistics.
subroutine qg_obsvec_sub_c(c_key_self, c_key_other)
Subtract observation vector.
subroutine qg_obsvec_random_c(c_odb, c_self)
Generate random observation vector.
subroutine qg_obsvec_copy_c(c_key_self, c_key_other)
Copy observation vector.
subroutine qg_obsvec_dotprod_c(c_key_obsvec1, c_key_obsvec2, zz)
Compute dot product between observation vectors.
subroutine qg_obsvec_nobs_c(c_key_self, kobs)
Get number of observations (not missing)
subroutine qg_obsvec_mul_scal_c(c_key_self, zz)
Multiply observation vector with a scalar.
subroutine qg_obsvec_nobs_withmask_c(c_key_self, c_key_mask, kobs)
Get observation vector size (only non-masked observations)
subroutine qg_obsvec_zero_c(c_key_self)
Set observation vector to zero.
subroutine qg_obsvec_mul_c(c_key_self, c_key_other)
Multiply observation vector.
subroutine qg_obsvec_size_c(c_key_self, kobs)
Get observation vector size.
subroutine qg_obsvec_settomissing_ith_c(c_key_self, i)
Set i-th value of the observation vector to missing value.
subroutine qg_obsvec_mask_c(c_key_self, c_key_mask)
Mask self observation vector (set values to missing where mask is set)
subroutine qg_obsvec_add_c(c_key_self, c_key_other)
Add observation vector.
subroutine qg_obsvec_setup_c(c_key_self, nlev, nobs)
Setup observation vector.
subroutine qg_obsvec_div_c(c_key_self, c_key_other)
Divide observation vector.
subroutine qg_obsvec_delete_c(c_key_self)
Delete observation vector.
subroutine qg_obsvec_axpy_c(c_key_self, zz, c_key_other)
Apply axpy on observation vector.
subroutine, public qg_obsvec_dotprod(obsvec1, obsvec2, zz)
Compute dot product between observation vectors.
subroutine, public qg_obsvec_random(c_odb, self)
Generate random observation vector.
subroutine, public qg_obsvec_ones(self)
Set observation vector to ones.
subroutine, public qg_obsvec_mask(self, mask)
Mask observation vector (set values to missing values where mask == 1)
subroutine, public qg_obsvec_mul(self, other)
Multiply observation vector.
subroutine, public qg_obsvec_zero(self)
Set observation vector to zero.
subroutine, public qg_obsvec_settomissing_ith(self, i)
Set i-th value of observation vector to missing value.
subroutine, public qg_obsvec_delete(self)
Delete observation vector.
subroutine, public qg_obsvec_mul_scal(self, zz)
Multiply observation vector with a scalar.
type(registry_t), public qg_obsvec_registry
Linked list interface - defines registry_t type.
subroutine, public qg_obsvec_size(self, kobs)
Get observation vector size.
subroutine, public qg_obsvec_mask_with_missing(self, mask)
Mask observation vector (set values to missing values where mask == missing value)
subroutine, public qg_obsvec_clone(self, other)
Clone observation vector.
subroutine, public qg_obsvec_get_withmask(self, obsmask, vals, nvals)
Get non-missing values from observation vector into vals array.
subroutine, public qg_obsvec_div(self, other)
Divide observation vector.
subroutine, public qg_obsvec_copy(self, other)
Copy observation vector.
subroutine, public qg_obsvec_setup(self, nlev, nobs)
Linked list implementation.
subroutine, public qg_obsvec_invert(self)
Invert observation vector.
subroutine, public qg_obsvec_nobs_withmask(self, obsmask, kobs)
Get observation vector size (only non-masked observations)
subroutine, public qg_obsvec_add(self, other)
Add observation vector.
subroutine, public qg_obsvec_stats(self, zmin, zmax, zavg)
Compute observation vector statistics.
subroutine, public qg_obsvec_axpy(self, zz, other)
Apply axpy on observation vector.
subroutine, public qg_obsvec_sub(self, other)
Subtract observation vector.
subroutine, public qg_obsvec_nobs(self, kobs)
Get observation vector size.