UFO
GeoVaLs.interface.F90
Go to the documentation of this file.
1 !
2 ! (C) Copyright 2017-2018 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 !
8 
9 use fckit_configuration_module, only: fckit_configuration
10 use iso_c_binding
12 use kinds
13 
14 implicit none
15 
16 public :: ufo_geovals_registry
17 
18 private
19 integer, parameter :: max_string=800
20 
21 #define LISTED_TYPE ufo_geovals
22 
23 !> Linked list interface - defines registry_t type
24 #include "oops/util/linkedList_i.f"
25 
26 !> Global registry
27 type(registry_t) :: ufo_geovals_registry
28 
29 ! ------------------------------------------------------------------------------
30 contains
31 ! ------------------------------------------------------------------------------
32 !> Linked list implementation
33 #include "oops/util/linkedList_c.f"
34 ! ------------------------------------------------------------------------------
35 !> Setup GeoVaLs (don't store anything; don't do allocation yet)
36 subroutine ufo_geovals_default_constr_c(c_key_self) bind(c,name='ufo_geovals_default_constr_f90')
37 implicit none
38 integer(c_int), intent(inout) :: c_key_self
39 type(ufo_geovals), pointer :: self
40 
41 call ufo_geovals_registry%init()
42 call ufo_geovals_registry%add(c_key_self)
43 call ufo_geovals_registry%get(c_key_self, self)
45 
46 end subroutine ufo_geovals_default_constr_c
47 
48 subroutine ufo_geovals_setup_c(c_key_self, c_nlocs, c_vars, c_nvars, c_sizes) bind(c,name='ufo_geovals_setup_f90')
49 use oops_variables_mod
50 implicit none
51 integer(c_int), intent(inout) :: c_key_self
52 integer(c_int), intent(in) :: c_nlocs, c_nvars
53 type(c_ptr), value, intent(in) :: c_vars
54 integer(c_size_t), intent(in) :: c_sizes(c_nvars)
55 
56 type(ufo_geovals), pointer :: self
57 type(oops_variables) :: vars
58 
59 call ufo_geovals_registry%init()
60 call ufo_geovals_registry%add(c_key_self)
61 call ufo_geovals_registry%get(c_key_self, self)
62 
63 vars = oops_variables(c_vars)
64 call ufo_geovals_setup(self, vars, c_nlocs, c_nvars, c_sizes)
65 
66 end subroutine ufo_geovals_setup_c
67 
68 !> Setup GeoVaLs (store nlocs, variables; don't do allocation yet)
69 subroutine ufo_geovals_partial_setup_c(c_key_self, c_nlocs, c_vars) bind(c,name='ufo_geovals_partial_setup_f90')
70 use oops_variables_mod
71 implicit none
72 integer(c_int), intent(inout) :: c_key_self
73 integer(c_int), intent(in) :: c_nlocs
74 type(c_ptr), value, intent(in) :: c_vars
75 
76 type(ufo_geovals), pointer :: self
77 type(oops_variables) :: vars
78 
79 call ufo_geovals_registry%init()
80 call ufo_geovals_registry%add(c_key_self)
81 call ufo_geovals_registry%get(c_key_self, self)
82 
83 vars = oops_variables(c_vars)
84 call ufo_geovals_partial_setup(self, vars, c_nlocs)
85 
86 end subroutine ufo_geovals_partial_setup_c
87 
88 !> Allocate GeoVaLs
89 subroutine ufo_geovals_allocate_c(c_key_self, c_nlevels, c_vars) bind(c,name='ufo_geovals_allocate_f90')
90 use oops_variables_mod
91 implicit none
92 integer(c_int), intent(inout) :: c_key_self
93 integer(c_int), intent(in) :: c_nlevels
94 type(c_ptr), value, intent(in) :: c_vars
95 
96 type(ufo_geovals), pointer :: self
97 type(oops_variables) :: vars
98 
99 call ufo_geovals_registry%get(c_key_self, self)
100 
101 vars = oops_variables(c_vars)
102 call ufo_geovals_allocate(self, vars, c_nlevels)
103 
104 end subroutine ufo_geovals_allocate_c
105 
106 ! ------------------------------------------------------------------------------
107 !> Copy one GeoVaLs object into another
108 
109 subroutine ufo_geovals_copy_c(c_key_self, c_key_other) bind(c,name='ufo_geovals_copy_f90')
110 implicit none
111 integer(c_int), intent(in) :: c_key_self
112 integer(c_int), intent(inout) :: c_key_other
113 type(ufo_geovals), pointer :: self
114 type(ufo_geovals), pointer :: other
115 
116 call ufo_geovals_registry%init()
117 call ufo_geovals_registry%add(c_key_other)
118 call ufo_geovals_registry%get(c_key_self, self)
119 call ufo_geovals_registry%get(c_key_other, other)
120 
121 call ufo_geovals_copy(self, other)
122 
123 end subroutine ufo_geovals_copy_c
124 
125 ! ------------------------------------------------------------------------------
126 !> Copy one GeoVaLs location into another object
127 
128 subroutine ufo_geovals_copy_one_c(c_key_self, c_key_other, c_ind) bind(c,name='ufo_geovals_copy_one_f90')
129 implicit none
130 integer(c_int), intent(inout) :: c_key_self
131 integer(c_int), intent(in) :: c_key_other
132 integer(c_int), intent(in) :: c_ind
133 type(ufo_geovals), pointer :: self
134 type(ufo_geovals), pointer :: other
135 integer :: ind
136 
137 ! Convert location index from the C++ to the Fortran convention.
138 ind = c_ind + 1
139 
140 call ufo_geovals_registry%init()
141 call ufo_geovals_registry%add(c_key_self)
142 call ufo_geovals_registry%get(c_key_self, self)
143 call ufo_geovals_registry%get(c_key_other, other)
144 
145 call ufo_geovals_copy_one(self, other, ind)
146 
147 end subroutine ufo_geovals_copy_one_c
148 
149 ! ------------------------------------------------------------------------------
150 !> Analytic init
151 
152 subroutine ufo_geovals_analytic_init_c(c_key_self, c_locs, c_conf) bind(c,name='ufo_geovals_analytic_init_f90')
154 implicit none
155 integer(c_int), intent(in) :: c_key_self
156 type(c_ptr), value, intent(in) :: c_locs
157 type(c_ptr), value, intent(in) :: c_conf
158 
159 type(ufo_geovals), pointer :: self
160 type(ufo_locations) ::locs
161 character(len=30) :: ic
162 character(len=:), allocatable :: str
163 type(fckit_configuration) :: f_conf
164 
165 call ufo_geovals_registry%get(c_key_self, self)
166 
167 f_conf = fckit_configuration(c_conf)
168 call f_conf%get_or_die("analytic_init",str)
169 ic = str
170 locs = ufo_locations(c_locs)
171 call ufo_geovals_analytic_init(self,locs,ic)
172 
173 end subroutine ufo_geovals_analytic_init_c
174 
175 ! ------------------------------------------------------------------------------
176 
177 subroutine ufo_geovals_delete_c(c_key_self) bind(c,name='ufo_geovals_delete_f90')
178 
179 implicit none
180 integer(c_int), intent(inout) :: c_key_self
181 
182 type(ufo_geovals), pointer :: self
183 
184 call ufo_geovals_registry%get(c_key_self, self)
185 
186 call ufo_geovals_delete(self)
187 
188 call ufo_geovals_registry%remove(c_key_self)
189 
190 end subroutine ufo_geovals_delete_c
191 
192 ! ------------------------------------------------------------------------------
193 
194 subroutine ufo_geovals_zero_c(c_key_self) bind(c,name='ufo_geovals_zero_f90')
195 implicit none
196 integer(c_int), intent(in) :: c_key_self
197 type(ufo_geovals), pointer :: self
198 
199 call ufo_geovals_registry%get(c_key_self, self)
200 
201 call ufo_geovals_zero(self)
202 
203 end subroutine ufo_geovals_zero_c
204 
205 ! ------------------------------------------------------------------------------
206 
207 subroutine ufo_geovals_reorderzdir_c(c_key_self, lvar, c_var, lvar1, c_var1) bind(c,name='ufo_geovals_reorderzdir_f90')
208 use ufo_vars_mod, only: maxvarlen
209 use string_f_c_mod
210 implicit none
211 integer(c_int), intent(in) :: c_key_self
212 integer(c_int), intent(in) :: lvar
213 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
214 integer(c_int), intent(in) :: lvar1
215 character(kind=c_char, len=1), intent(in) :: c_var1(lvar1+1)
216 character(len=MAXVARLEN) :: varname
217 character(len=MAXVARLEN) :: vardir
218 type(ufo_geovals), pointer :: self
219 
220 call c_f_string(c_var, varname)
221 call c_f_string(c_var1, vardir)
222 call ufo_geovals_registry%get(c_key_self, self)
223 
224 call ufo_geovals_reorderzdir(self, varname, vardir)
225 
226 end subroutine ufo_geovals_reorderzdir_c
227 
228 ! ------------------------------------------------------------------------------
229 
230 subroutine ufo_geovals_abs_c(c_key_self) bind(c,name='ufo_geovals_abs_f90')
231 implicit none
232 integer(c_int), intent(in) :: c_key_self
233 type(ufo_geovals), pointer :: self
234 
235 call ufo_geovals_registry%get(c_key_self, self)
236 
237 call ufo_geovals_abs(self)
238 
239 end subroutine ufo_geovals_abs_c
240 
241 ! ------------------------------------------------------------------------------
242 
243 subroutine ufo_geovals_rms_c(c_key_self,vrms) bind(c,name='ufo_geovals_rms_f90')
244 implicit none
245 integer(c_int), intent(in) :: c_key_self
246 real(c_double), intent(inout) :: vrms
247 type(ufo_geovals), pointer :: self
248 
249 call ufo_geovals_registry%get(c_key_self, self)
250 
251 call ufo_geovals_rms(self,vrms)
252 
253 end subroutine ufo_geovals_rms_c
254 
255 ! ------------------------------------------------------------------------------
256 
257 subroutine ufo_geovals_random_c(c_key_self) bind(c,name='ufo_geovals_random_f90')
258 implicit none
259 integer(c_int), intent(in) :: c_key_self
260 type(ufo_geovals), pointer :: self
261 
262 call ufo_geovals_registry%get(c_key_self, self)
263 
264 call ufo_geovals_random(self)
265 
266 end subroutine ufo_geovals_random_c
267 
268 ! ------------------------------------------------------------------------------
269 
270 subroutine ufo_geovals_scalmult_c(c_key_self, zz) bind(c,name='ufo_geovals_scalmult_f90')
271 implicit none
272 integer(c_int), intent(in) :: c_key_self
273 real(c_double), intent(in) :: zz
274 type(ufo_geovals), pointer :: self
275 
276 call ufo_geovals_registry%get(c_key_self, self)
277 
278 call ufo_geovals_scalmult(self, zz)
279 
280 end subroutine ufo_geovals_scalmult_c
281 
282 ! ------------------------------------------------------------------------------
283 
284 subroutine ufo_geovals_profmult_c(c_key_self, nlocs, values) bind(c,name='ufo_geovals_profmult_f90')
285 implicit none
286 integer(c_int), intent(in) :: c_key_self
287 integer(c_int), intent(in) :: nlocs
288 real(c_float), intent(in) :: values(nlocs)
289 type(ufo_geovals), pointer :: self
290 
291 call ufo_geovals_registry%get(c_key_self, self)
292 
293 call ufo_geovals_profmult(self, nlocs, values)
294 
295 end subroutine ufo_geovals_profmult_c
296 
297 ! ------------------------------------------------------------------------------
298 
299 subroutine ufo_geovals_assign_c(c_key_self, c_key_rhs) bind(c,name='ufo_geovals_assign_f90')
300 implicit none
301 integer(c_int), intent(in) :: c_key_self
302 integer(c_int), intent(in) :: c_key_rhs
303 type(ufo_geovals), pointer :: self
304 type(ufo_geovals), pointer :: rhs
305 
306 call ufo_geovals_registry%get(c_key_self, self)
307 call ufo_geovals_registry%get(c_key_rhs, rhs)
308 
309 call ufo_geovals_assign(self, rhs)
310 
311 end subroutine ufo_geovals_assign_c
312 
313 ! ------------------------------------------------------------------------------
314 
315 subroutine ufo_geovals_add_c(c_key_self, c_key_other) bind(c,name='ufo_geovals_add_f90')
316 implicit none
317 integer(c_int), intent(in) :: c_key_self
318 integer(c_int), intent(in) :: c_key_other
319 type(ufo_geovals), pointer :: self
320 type(ufo_geovals), pointer :: other
321 
322 call ufo_geovals_registry%get(c_key_self, self)
323 call ufo_geovals_registry%get(c_key_other, other)
324 
325 call ufo_geovals_add(self, other)
326 
327 end subroutine ufo_geovals_add_c
328 
329 ! ------------------------------------------------------------------------------
330 
331 subroutine ufo_geovals_diff_c(c_key_self, c_key_other) bind(c,name='ufo_geovals_diff_f90')
332 implicit none
333 integer(c_int), intent(in) :: c_key_self
334 integer(c_int), intent(in) :: c_key_other
335 type(ufo_geovals), pointer :: self
336 type(ufo_geovals), pointer :: other
337 
338 call ufo_geovals_registry%get(c_key_self, self)
339 call ufo_geovals_registry%get(c_key_other, other)
340 
341 call ufo_geovals_diff(self, other)
342 
343 end subroutine ufo_geovals_diff_c
344 
345 ! ------------------------------------------------------------------------------
346 
347 subroutine ufo_geovals_schurmult_c(c_key_self, c_key_other) bind(c,name='ufo_geovals_schurmult_f90')
348 implicit none
349 integer(c_int), intent(in) :: c_key_self
350 integer(c_int), intent(in) :: c_key_other
351 type(ufo_geovals), pointer :: self
352 type(ufo_geovals), pointer :: other
353 
354 call ufo_geovals_registry%get(c_key_self, self)
355 call ufo_geovals_registry%get(c_key_other, other)
356 
357 call ufo_geovals_schurmult(self, other)
358 
359 end subroutine ufo_geovals_schurmult_c
360 
361 ! ------------------------------------------------------------------------------
362 
363 subroutine ufo_geovals_normalize_c(c_key_self, c_key_other) bind(c,name='ufo_geovals_normalize_f90')
364 implicit none
365 integer(c_int), intent(in) :: c_key_self
366 integer(c_int), intent(in) :: c_key_other
367 type(ufo_geovals), pointer :: self
368 type(ufo_geovals), pointer :: other
369 
370 call ufo_geovals_registry%get(c_key_self, self)
371 call ufo_geovals_registry%get(c_key_other, other)
372 
373 call ufo_geovals_normalize(self, other)
374 
375 end subroutine ufo_geovals_normalize_c
376 
377 ! ------------------------------------------------------------------------------
378 
379 subroutine ufo_geovals_split_c(c_key_self, c_key_other1, c_key_other2) bind(c,name='ufo_geovals_split_f90')
380 implicit none
381 integer(c_int), intent(in) :: c_key_self, c_key_other1, c_key_other2
382 type(ufo_geovals), pointer :: self, other1, other2
383 
384 call ufo_geovals_registry%get(c_key_self, self)
385 call ufo_geovals_registry%get(c_key_other1, other1)
386 call ufo_geovals_registry%get(c_key_other2, other2)
387 
388 call ufo_geovals_split(self, other1, other2)
389 
390 end subroutine ufo_geovals_split_c
391 
392 ! ------------------------------------------------------------------------------
393 
394 subroutine ufo_geovals_merge_c(c_key_self, c_key_other1, c_key_other2) bind(c,name='ufo_geovals_merge_f90')
395 implicit none
396 integer(c_int), intent(in) :: c_key_self, c_key_other1, c_key_other2
397 type(ufo_geovals), pointer :: self, other1, other2
398 
399 call ufo_geovals_registry%get(c_key_self, self)
400 call ufo_geovals_registry%get(c_key_other1, other1)
401 call ufo_geovals_registry%get(c_key_other2, other2)
402 
403 call ufo_geovals_merge(self, other1, other2)
404 
405 end subroutine ufo_geovals_merge_c
406 
407 ! ------------------------------------------------------------------------------
408 
409 subroutine ufo_geovals_minmaxavg_c(c_key_self, kobs, kvar, pmin, pmax, prms) bind(c,name='ufo_geovals_minmaxavg_f90')
410 implicit none
411 integer(c_int), intent(in) :: c_key_self
412 integer(c_int), intent(inout) :: kobs
413 integer(c_int), intent(in) :: kvar
414 real(c_double), intent(inout) :: pmin, pmax, prms
415 type(ufo_geovals), pointer :: self
416 
417 call ufo_geovals_registry%get(c_key_self, self)
418 
419 call ufo_geovals_minmaxavg(self, kobs, kvar, pmin, pmax, prms)
420 
421 end subroutine ufo_geovals_minmaxavg_c
422 
423 ! ------------------------------------------------------------------------------
424 
425 subroutine ufo_geovals_nlocs_c(c_key_self, kobs) bind(c, name='ufo_geovals_nlocs_f90')
426 implicit none
427 integer(c_int), intent(in) :: c_key_self
428 integer(c_size_t), intent(inout) :: kobs
429 type(ufo_geovals), pointer :: self
430 
431 call ufo_geovals_registry%get(c_key_self, self)
432 kobs = self%nlocs
433 
434 end subroutine ufo_geovals_nlocs_c
435 
436 ! ------------------------------------------------------------------------------
437 
438 subroutine ufo_geovals_nlevs_c(c_key_self, lvar, c_var, nlevs) bind(c, name='ufo_geovals_nlevs_f90')
439 use ufo_vars_mod, only: maxvarlen
440 use string_f_c_mod
441 implicit none
442 integer(c_int), intent(in) :: c_key_self
443 integer(c_int), intent(in) :: lvar
444 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
445 integer(c_int), intent(out) :: nlevs
446 
447 type(ufo_geoval), pointer :: geoval
448 character(len=MAXVARLEN) :: varname
449 type(ufo_geovals), pointer :: self
450 
451 call c_f_string(c_var, varname)
452 call ufo_geovals_registry%get(c_key_self, self)
453 
454 call ufo_geovals_get_var(self, varname, geoval)
455 
456 nlevs = geoval%nval
457 
458 end subroutine ufo_geovals_nlevs_c
459 
460 ! ------------------------------------------------------------------------------
461 
462 subroutine ufo_geovals_get2d_c(c_key_self, lvar, c_var, nlocs, values) bind(c, name='ufo_geovals_get2d_f90')
463 use ufo_vars_mod, only: maxvarlen
464 use string_f_c_mod
465 implicit none
466 integer(c_int), intent(in) :: c_key_self
467 integer(c_int), intent(in) :: lvar
468 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
469 integer(c_int), intent(in) :: nlocs
470 real(c_double), intent(inout) :: values(nlocs)
471 
472 character(max_string) :: err_msg
473 type(ufo_geoval), pointer :: geoval
474 character(len=MAXVARLEN) :: varname
475 type(ufo_geovals), pointer :: self
476 
477 call c_f_string(c_var, varname)
478 call ufo_geovals_registry%get(c_key_self, self)
479 
480 call ufo_geovals_get_var(self, varname, geoval)
481 
482 if (size(geoval%vals,1) /= 1) then
483  write(err_msg,*)'ufo_geovals_get2d_f90',trim(varname),'is not a 2D var:',size(geoval%vals,1), ' levels'
484  call abor1_ftn(err_msg)
485 endif
486 if (nlocs /= size(geoval%vals,2)) then
487  write(err_msg,*)'ufo_geovals_get2d_f90',trim(varname),'error locs number:',nlocs,size(geoval%vals,2)
488  call abor1_ftn(err_msg)
489 endif
490 
491 values(:) = geoval%vals(1,:)
492 
493 end subroutine ufo_geovals_get2d_c
494 
495 ! ------------------------------------------------------------------------------
496 
497 subroutine ufo_geovals_get_c(c_key_self, lvar, c_var, c_lev, nlocs, values) bind(c, name='ufo_geovals_get_f90')
498 use ufo_vars_mod, only: maxvarlen
499 use string_f_c_mod
500 implicit none
501 integer(c_int), intent(in) :: c_key_self
502 integer(c_int), intent(in) :: lvar
503 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
504 integer(c_int), intent(in) :: c_lev
505 integer(c_int), intent(in) :: nlocs
506 real(c_float), intent(inout) :: values(nlocs)
507 
508 character(max_string) :: err_msg
509 type(ufo_geoval), pointer :: geoval
510 character(len=MAXVARLEN) :: varname
511 type(ufo_geovals), pointer :: self
512 integer(c_int) :: lev
513 
514 ! Convert level index from the C++ to the Fortran convention.
515 lev = c_lev + 1
516 
517 call c_f_string(c_var, varname)
518 call ufo_geovals_registry%get(c_key_self, self)
519 
520 call ufo_geovals_get_var(self, varname, geoval)
521 
522 if (lev<1 .or. lev>size(geoval%vals,1)) then
523  write(err_msg,*)'ufo_geovals_get_f90 "',trim(varname),'" level out of range: 1~', &
524  size(geoval%vals,1), ', lev=', lev
525  call abor1_ftn(err_msg)
526 endif
527 if (nlocs /= size(geoval%vals,2)) then
528  write(err_msg,*)'ufo_geovals_get_f90 "',trim(varname),'" error locs number:',nlocs,&
529  ' /= ',size(geoval%vals,2)
530  call abor1_ftn(err_msg)
531 endif
532 
533 values(:) = geoval%vals(lev,:)
534 
535 end subroutine ufo_geovals_get_c
536 
537 ! ------------------------------------------------------------------------------
538 
539 subroutine ufo_geovals_get_loc_c(c_key_self, lvar, c_var, c_loc, nlevs, values) bind(c, name='ufo_geovals_get_loc_f90')
540 use ufo_vars_mod, only: maxvarlen
541 use string_f_c_mod
542 implicit none
543 integer(c_int), intent(in) :: c_key_self
544 integer(c_int), intent(in) :: lvar
545 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
546 integer(c_int), intent(in) :: c_loc
547 integer(c_int), intent(in) :: nlevs
548 real(c_double), intent(inout) :: values(nlevs)
549 
550 character(max_string) :: err_msg
551 type(ufo_geoval), pointer :: geoval
552 character(len=MAXVARLEN) :: varname
553 type(ufo_geovals), pointer :: self
554 integer(c_int) :: loc
555 
556 call c_f_string(c_var, varname)
557 call ufo_geovals_registry%get(c_key_self, self)
558 
559 call ufo_geovals_get_var(self, varname, geoval)
560 
561 ! Convert location index from the C++ to the Fortran convention.
562 loc = c_loc + 1
563 
564 if (loc<1 .or. loc>size(geoval%vals,2)) then
565  write(err_msg,*)'ufo_geovals_get_loc_f90',trim(varname),'location out of range:',loc,size(geoval%vals,2)
566  call abor1_ftn(err_msg)
567 endif
568 if (nlevs /= size(geoval%vals,1)) then
569  write(err_msg,*)'ufo_geovals_get_loc_f90',trim(varname),'incorrect number of levels:',nlevs,size(geoval%vals,1)
570  call abor1_ftn(err_msg)
571 endif
572 
573 values(:) = geoval%vals(:,loc)
574 
575 end subroutine ufo_geovals_get_loc_c
576 
577 ! ------------------------------------------------------------------------------
578 
579 subroutine ufo_geovals_getdouble_c(c_key_self, lvar, c_var, c_lev, nlocs, values)&
580  bind(c, name='ufo_geovals_getdouble_f90')
581 use ufo_vars_mod, only: maxvarlen
582 use string_f_c_mod
583 implicit none
584 integer(c_int), intent(in) :: c_key_self
585 integer(c_int), intent(in) :: lvar
586 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
587 integer(c_int), intent(in) :: c_lev
588 integer(c_int), intent(in) :: nlocs
589 real(c_double), intent(inout) :: values(nlocs)
590 
591 type(ufo_geoval), pointer :: geoval
592 character(len=MAXVARLEN) :: varname
593 type(ufo_geovals), pointer :: self
594 integer(c_int) :: lev
595 
596 ! Convert level index from the C++ to the Fortran convention.
597 lev = c_lev + 1
598 
599 call c_f_string(c_var, varname)
600 call ufo_geovals_registry%get(c_key_self, self)
601 call ufo_geovals_get_var(self, varname, geoval)
602 values(:) = geoval%vals(lev,:)
603 
604 end subroutine ufo_geovals_getdouble_c
605 
606 ! ------------------------------------------------------------------------------
607 
608 subroutine ufo_geovals_putdouble_c(c_key_self, lvar, c_var, c_lev, nlocs, values) bind(c, name='ufo_geovals_putdouble_f90')
609 use ufo_vars_mod, only: maxvarlen
610 use string_f_c_mod
611 integer(c_int), intent(in) :: c_key_self
612 integer(c_int), intent(in) :: lvar
613 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
614 integer(c_int), intent(in) :: c_lev
615 integer(c_int), intent(in) :: nlocs
616 real(c_double), intent(in) :: values(nlocs)
617 
618 type(ufo_geoval), pointer :: geoval
619 character(len=MAXVARLEN) :: varname
620 type(ufo_geovals), pointer :: self
621 integer(c_int) :: lev
622 
623 ! Convert level index from the C++ to the Fortran convention.
624 lev = c_lev + 1
625 
626 call c_f_string(c_var, varname)
627 call ufo_geovals_registry%get(c_key_self, self)
628 call ufo_geovals_get_var(self, varname, geoval)
629 geoval%vals(lev,:) = values(:)
630 
631 end subroutine ufo_geovals_putdouble_c
632 
633 ! ------------------------------------------------------------------------------
634 subroutine ufo_geovals_put_loc_c(c_key_self, lvar, c_var, c_loc, nlevs, values) bind(c, name='ufo_geovals_put_loc_f90')
635 use ufo_vars_mod, only: maxvarlen
636 use string_f_c_mod
637 
638 implicit none
639 integer(c_int), intent(in) :: c_key_self
640 integer(c_int), intent(in) :: lvar
641 character(kind=c_char, len=1), intent(in) :: c_var(lvar+1)
642 integer(c_int), intent(in) :: c_loc
643 integer(c_int), intent(in) :: nlevs
644 real(c_double), intent(in) :: values(nlevs)
645 
646 character(max_string) :: err_msg
647 type(ufo_geoval), pointer :: geoval
648 character(len=MAXVARLEN) :: varname
649 type(ufo_geovals), pointer :: self
650 integer(c_int) :: loc
651 
652 call c_f_string(c_var, varname)
653 call ufo_geovals_registry%get(c_key_self, self)
654 call ufo_geovals_get_var(self, varname, geoval)
655 
656 ! Convert location index from the C++ to the Fortran convention.
657 loc = c_loc + 1
658 
659 if (loc<1 .or. loc>size(geoval%vals,2)) then
660  write(err_msg,*)'ufo_geovals_put_loc_f90',trim(varname),'location out of range:',loc,size(geoval%vals,2)
661  call abor1_ftn(err_msg)
662 endif
663 if (nlevs /= size(geoval%vals,1)) then
664  write(err_msg,*)'ufo_geovals_put_loc_f90',trim(varname),'incorrect number of levels:',nlevs,size(geoval%vals,1)
665  call abor1_ftn(err_msg)
666 endif
667 
668 geoval%vals(:,loc) = values(:)
669 
670 end subroutine ufo_geovals_put_loc_c
671 
672 ! ------------------------------------------------------------------------------
673 
674 subroutine ufo_geovals_maxloc_c(c_key_self, mxval, iloc, ivar) bind(c,name='ufo_geovals_maxloc_f90')
675 implicit none
676 integer(c_int), intent(in) :: c_key_self
677 real(c_double), intent(inout) :: mxval
678 integer(c_int), intent(inout) :: iloc, ivar
679 type(ufo_geovals), pointer :: self
680 
681 call ufo_geovals_registry%get(c_key_self, self)
682 
683 call ufo_geovals_maxloc(self, mxval, iloc, ivar)
684 
685 end subroutine ufo_geovals_maxloc_c
686 
687 ! ------------------------------------------------------------------------------
688 
689 subroutine ufo_geovals_read_file_c(c_key_self, c_conf, c_obspace, c_vars) bind(c,name='ufo_geovals_read_file_f90')
690 use oops_variables_mod
691 use datetime_mod
692 
693 implicit none
694 integer(c_int), intent(inout) :: c_key_self
695 type(c_ptr), value, intent(in) :: c_conf
696 type(c_ptr), value, intent(in) :: c_obspace
697 type(c_ptr), value, intent(in) :: c_vars
698 
699 type(ufo_geovals), pointer :: self
700 character(max_string) :: filename
701 integer :: loc_multiplier
702 character(len=:), allocatable :: str
703 type(fckit_configuration) :: f_conf
704 type(oops_variables) :: vars
705 
706 call ufo_geovals_registry%init()
707 call ufo_geovals_registry%add(c_key_self)
708 call ufo_geovals_registry%get(c_key_self, self)
709 
710 ! read filename for config
711 f_conf = fckit_configuration(c_conf)
712 call f_conf%get_or_die("filename",str)
713 filename = str
714 
715 if (f_conf%has("loc_multiplier")) then
716  call f_conf%get_or_die("loc_multiplier", loc_multiplier)
717 else
718  loc_multiplier = 1
719 endif
720 
721 vars = oops_variables(c_vars)
722 ! read geovals
723 call ufo_geovals_read_netcdf(self, filename, loc_multiplier, c_obspace, vars)
724 
725 end subroutine ufo_geovals_read_file_c
726 
727 ! ------------------------------------------------------------------------------
728 
729 subroutine ufo_geovals_write_file_c(c_key_self, c_conf, c_rank) bind(c,name='ufo_geovals_write_file_f90')
730 implicit none
731 integer(c_int), intent(in) :: c_key_self
732 type(c_ptr), value, intent(in) :: c_conf
733 integer(c_size_t), intent(in) :: c_rank ! mpi rank (to be added to filename)
734 
735 type(ufo_geovals), pointer :: self
736 character(max_string) :: fout, filename
737 
738 character(len=10) :: cproc
739 integer :: ppos
740 character(len=:), allocatable :: str
741 type(fckit_configuration) :: f_conf
742 
743 ! read filename for config
744 f_conf = fckit_configuration(c_conf)
745 call f_conf%get_or_die("filename",str)
746 filename = str
747 
748 write(cproc,fmt='(i4.4)') c_rank
749 
750 ! Find the left-most dot in the file name, and use that to pick off the file name
751 ! and file extension.
752 ppos = scan(trim(filename), '.', back=.true.)
753 if (ppos > 0) then
754  ! found a file extension
755  fout = filename(1:ppos-1) // '_' // trim(adjustl(cproc)) // trim(filename(ppos:))
756 else
757  ! no file extension
758  fout = trim(filename) // '_' // trim(adjustl(cproc))
759 endif
760 
761 call ufo_geovals_registry%get(c_key_self, self)
762 call ufo_geovals_write_netcdf(self, fout)
763 
764 end subroutine ufo_geovals_write_file_c
765 
766 ! ------------------------------------------------------------------------------
767 
768 end module ufo_geovals_mod_c
subroutine ufo_geovals_put_loc_c(c_key_self, lvar, c_var, c_loc, nlevs, values)
subroutine ufo_geovals_nlocs_c(c_key_self, kobs)
subroutine ufo_geovals_normalize_c(c_key_self, c_key_other)
subroutine ufo_geovals_setup_c(c_key_self, c_nlocs, c_vars, c_nvars, c_sizes)
subroutine ufo_geovals_abs_c(c_key_self)
subroutine ufo_geovals_merge_c(c_key_self, c_key_other1, c_key_other2)
subroutine ufo_geovals_get_loc_c(c_key_self, lvar, c_var, c_loc, nlevs, values)
subroutine ufo_geovals_getdouble_c(c_key_self, lvar, c_var, c_lev, nlocs, values)
subroutine ufo_geovals_nlevs_c(c_key_self, lvar, c_var, nlevs)
subroutine ufo_geovals_split_c(c_key_self, c_key_other1, c_key_other2)
subroutine ufo_geovals_write_file_c(c_key_self, c_conf, c_rank)
subroutine ufo_geovals_schurmult_c(c_key_self, c_key_other)
subroutine ufo_geovals_random_c(c_key_self)
subroutine ufo_geovals_allocate_c(c_key_self, c_nlevels, c_vars)
Allocate GeoVaLs.
subroutine ufo_geovals_read_file_c(c_key_self, c_conf, c_obspace, c_vars)
integer, parameter max_string
subroutine ufo_geovals_putdouble_c(c_key_self, lvar, c_var, c_lev, nlocs, values)
subroutine ufo_geovals_get2d_c(c_key_self, lvar, c_var, nlocs, values)
subroutine ufo_geovals_default_constr_c(c_key_self)
Linked list implementation.
subroutine ufo_geovals_copy_c(c_key_self, c_key_other)
Copy one GeoVaLs object into another.
subroutine ufo_geovals_reorderzdir_c(c_key_self, lvar, c_var, lvar1, c_var1)
subroutine ufo_geovals_copy_one_c(c_key_self, c_key_other, c_ind)
Copy one GeoVaLs location into another object.
subroutine ufo_geovals_get_c(c_key_self, lvar, c_var, c_lev, nlocs, values)
subroutine ufo_geovals_partial_setup_c(c_key_self, c_nlocs, c_vars)
Setup GeoVaLs (store nlocs, variables; don't do allocation yet)
type(registry_t), public ufo_geovals_registry
Linked list interface - defines registry_t type.
subroutine ufo_geovals_add_c(c_key_self, c_key_other)
subroutine ufo_geovals_diff_c(c_key_self, c_key_other)
subroutine ufo_geovals_zero_c(c_key_self)
subroutine ufo_geovals_delete_c(c_key_self)
subroutine ufo_geovals_minmaxavg_c(c_key_self, kobs, kvar, pmin, pmax, prms)
subroutine ufo_geovals_analytic_init_c(c_key_self, c_locs, c_conf)
Analytic init.
subroutine ufo_geovals_rms_c(c_key_self, vrms)
subroutine ufo_geovals_maxloc_c(c_key_self, mxval, iloc, ivar)
subroutine ufo_geovals_profmult_c(c_key_self, nlocs, values)
subroutine ufo_geovals_assign_c(c_key_self, c_key_rhs)
subroutine ufo_geovals_scalmult_c(c_key_self, zz)
subroutine, public ufo_geovals_reorderzdir(self, varname, zdir)
subroutine, public ufo_geovals_maxloc(self, mxval, iobs, ivar)
Location where the summed geovals value is maximum.
subroutine, public ufo_geovals_default_constr(self)
subroutine, public ufo_geovals_allocate(self, vars, nlevels)
Deprecated. Rely on ufo_geovals_setup to allocate GeoVaLs instead. Allocates GeoVaLs for vars variabl...
subroutine, public ufo_geovals_random(self)
subroutine, public ufo_geovals_setup(self, vars, nlocs, nvars, nvals)
Initializes and allocates self GeoVaLs with nlocs number of locations for vars variables....
subroutine, public ufo_geovals_write_netcdf(self, filename)
subroutine, public ufo_geovals_normalize(self, other)
Normalization of one GeoVaLs object by another.
subroutine, public ufo_geovals_split(self, other1, other2)
subroutine, public ufo_geovals_read_netcdf(self, filename, loc_multiplier, c_obspace, vars)
subroutine, public ufo_geovals_rms(self, vrms)
subroutine, public ufo_geovals_copy(self, other)
Copy one GeoVaLs object into another.
subroutine, public ufo_geovals_delete(self)
subroutine, public ufo_geovals_analytic_init(self, locs, ic)
Initialize a GeoVaLs object based on an analytic state.
subroutine, public ufo_geovals_diff(self, other)
Difference between two GeoVaLs objects.
subroutine, public ufo_geovals_profmult(self, nlocs, values)
subroutine, public ufo_geovals_copy_one(self, other, loc_index)
Copy one location from GeoVaLs into a new object.
subroutine, public ufo_geovals_minmaxavg(self, kobs, kvar, pmin, pmax, prms)
subroutine, public ufo_geovals_scalmult(self, zz)
subroutine, public ufo_geovals_assign(self, rhs)
subroutine, public ufo_geovals_add(self, other)
Sum of two GeoVaLs objects.
subroutine, public ufo_geovals_get_var(self, varname, geoval)
subroutine, public ufo_geovals_partial_setup(self, vars, nlocs)
Deprecated, use ufo_geovals_setup instead. Partially initializes self GeoVaLs with nlocs number of lo...
subroutine, public ufo_geovals_zero(self)
subroutine, public ufo_geovals_abs(self)
subroutine, public ufo_geovals_schurmult(self, other)
Schur product of two GeoVaLs objects.
subroutine, public ufo_geovals_merge(self, other1, other2)
Fortran interface to ufo::Locations.
integer, parameter, public maxvarlen
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators