OOPS
qg_fields_interface.F90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 ! In applying this licence, ECMWF does not waive the privileges and immunities
6 ! granted to it by virtue of its status as an intergovernmental organisation nor
7 ! does it submit to any jurisdiction.
8 
10 
11 use atlas_module, only: atlas_fieldset
12 use datetime_mod
13 use fckit_configuration_module, only: fckit_configuration
14 use iso_c_binding
15 use kinds
17 use qg_fields_mod
18 use qg_geom_mod
20 use qg_locs_mod
21 
22 implicit none
23 
24 private
25 ! ------------------------------------------------------------------------------
26 contains
27 ! ------------------------------------------------------------------------------
28 !> Create fields from geometry and variables
29 subroutine qg_fields_create_c(c_key_self,c_key_geom,c_vars,c_lbc) bind(c,name='qg_fields_create_f90')
30 
31 implicit none
32 
33 ! Passed variables
34 integer(c_int),intent(inout) :: c_key_self !< Fields
35 integer(c_int),intent(in) :: c_key_geom !< Geometry
36 type(c_ptr),value,intent(in) :: c_vars !< List of variables
37 logical(c_bool),intent(in) :: c_lbc !< Boundaries flag
38 
39 ! Local variables
40 type(qg_fields),pointer :: self
41 type(qg_geom),pointer :: geom
42 type(oops_variables) :: vars
43 logical :: lbc
44 
45 ! Interface
46 call qg_fields_registry%init()
47 call qg_fields_registry%add(c_key_self)
48 call qg_fields_registry%get(c_key_self,self)
49 call qg_geom_registry%get(c_key_geom,geom)
50 vars = oops_variables(c_vars)
51 lbc = c_lbc
52 
53 ! Call Fortran
54 call qg_fields_create(self,geom,vars,lbc)
55 
56 end subroutine qg_fields_create_c
57 ! ------------------------------------------------------------------------------
58 !> Create fields from another one
59 subroutine qg_fields_create_from_other_c(c_key_self,c_key_other,c_key_geom) bind(c,name='qg_fields_create_from_other_f90')
60 
61 implicit none
62 
63 ! Passed variables
64 integer(c_int),intent(inout) :: c_key_self !< Fields
65 integer(c_int),intent(in) :: c_key_other !< Other fields
66 integer(c_int),intent(in) :: c_key_geom !< Geometry
67 
68 ! Local variables
69 type(qg_fields),pointer :: self
70 type(qg_fields),pointer :: other
71 type(qg_geom),pointer :: geom
72 
73 ! Interface
74 call qg_fields_registry%get(c_key_other,other)
75 call qg_fields_registry%init()
76 call qg_fields_registry%add(c_key_self)
77 call qg_fields_registry%get(c_key_self,self)
78 call qg_geom_registry%get(c_key_geom,geom)
79 
80 ! Call Fortran
81 call qg_fields_create_from_other(self,other,geom)
82 
83 end subroutine qg_fields_create_from_other_c
84 ! ------------------------------------------------------------------------------
85 !> Delete fields
86 subroutine qg_fields_delete_c(c_key_self) bind(c,name='qg_fields_delete_f90')
87 
88 implicit none
89 
90 ! Passed variables
91 integer(c_int),intent(inout) :: c_key_self !< Fields
92 
93 ! Local variables
94 type(qg_fields),pointer :: self
95 
96 ! Interface
97 call qg_fields_registry%get(c_key_self,self)
98 
99 ! Call Fortran
100 call qg_fields_delete(self)
101 
102 ! Clear interface
103 call qg_fields_registry%remove(c_key_self)
104 
105 end subroutine qg_fields_delete_c
106 ! ------------------------------------------------------------------------------
107 !> Set fields to zero
108 subroutine qg_fields_zero_c(c_key_self) bind(c,name='qg_fields_zero_f90')
109 
110 implicit none
111 
112 ! Passed variables
113 integer(c_int),intent(in) :: c_key_self !< Fields
114 
115 ! Local variables
116 type(qg_fields),pointer :: self
117 
118 ! Interface
119 call qg_fields_registry%get(c_key_self,self)
120 
121 ! Call Fortran
122 call qg_fields_zero(self)
123 
124 end subroutine qg_fields_zero_c
125 ! ------------------------------------------------------------------------------
126 !> Set fields to ones
127 subroutine qg_fields_ones_c(c_key_self) bind(c,name='qg_fields_ones_f90')
128 
129 implicit none
130 
131 ! Passed variables
132 integer(c_int),intent(in) :: c_key_self !< Fields
133 
134 ! Local variables
135 type(qg_fields),pointer :: self
136 
137 ! Interface
138 call qg_fields_registry%get(c_key_self,self)
139 
140 ! Call Fortran
141 call qg_fields_ones(self)
142 
143 end subroutine qg_fields_ones_c
144 ! ------------------------------------------------------------------------------
145 !> Set fields to Diracs
146 subroutine qg_fields_dirac_c(c_key_self,c_conf) bind(c,name='qg_fields_dirac_f90')
147 
148 implicit none
149 
150 ! Passed variables
151 integer(c_int),intent(in) :: c_key_self !< Fields
152 type(c_ptr),value,intent(in) :: c_conf !< Configuration
153 
154 ! Local variables
155 type(fckit_configuration) :: f_conf
156 type(qg_fields),pointer :: self
157 
158 ! Interface
159 f_conf = fckit_configuration(c_conf)
160 call qg_fields_registry%get(c_key_self,self)
161 
162 ! Call Fortran
163 call qg_fields_dirac(self,f_conf)
164 
165 end subroutine qg_fields_dirac_c
166 ! ------------------------------------------------------------------------------
167 !> Generate random fields
168 subroutine qg_fields_random_c(c_key_self,c_vars) bind(c,name='qg_fields_random_f90')
169 
170 implicit none
171 
172 ! Passed variables
173 integer(c_int),intent(in) :: c_key_self !< Fields
174 type(c_ptr),value,intent(in) :: c_vars !< List of variables
175 
176 ! Local variables
177 type(qg_fields),pointer :: self
178 type(oops_variables) :: vars
179 
180 ! Interface
181 call qg_fields_registry%get(c_key_self,self)
182 vars = oops_variables(c_vars)
183 
184 ! Call Fortran
185 if (vars%has('x')) then
186  call qg_fields_random(self,'x')
187 elseif (vars%has('q')) then
188  call qg_fields_random(self,'q')
189 else
190  call abor1_ftn('qg_fields_random_c: x or q required in output field')
191 endif
192 
193 end subroutine qg_fields_random_c
194 ! ------------------------------------------------------------------------------
195 !> Copy fields
196 subroutine qg_fields_copy_c(c_key_self,c_key_other) bind(c,name='qg_fields_copy_f90')
197 
198 implicit none
199 
200 ! Passed variables
201 integer(c_int),intent(in) :: c_key_self !< Fields
202 integer(c_int),intent(in) :: c_key_other !< Other fields
203 
204 ! Local variables
205 type(qg_fields),pointer :: self
206 type(qg_fields),pointer :: other
207 
208 ! Interface
209 call qg_fields_registry%get(c_key_self,self)
210 call qg_fields_registry%get(c_key_other,other)
211 
212 ! Call Fortran
213 call qg_fields_copy(self,other)
214 
215 end subroutine qg_fields_copy_c
216 ! ------------------------------------------------------------------------------
217 !> Add fields
218 subroutine qg_fields_self_add_c(c_key_self,c_key_rhs) bind(c,name='qg_fields_self_add_f90')
219 
220 implicit none
221 
222 ! Passed variables
223 integer(c_int),intent(in) :: c_key_self !< Fields
224 integer(c_int),intent(in) :: c_key_rhs !< Right-hand side
225 
226 ! Local variables
227 type(qg_fields),pointer :: self
228 type(qg_fields),pointer :: rhs
229 
230 ! Interface
231 call qg_fields_registry%get(c_key_self,self)
232 call qg_fields_registry%get(c_key_rhs,rhs)
233 
234 ! Call Fortran
235 call qg_fields_self_add(self,rhs)
236 
237 end subroutine qg_fields_self_add_c
238 ! ------------------------------------------------------------------------------
239 !> Subtract fields
240 subroutine qg_fields_self_sub_c(c_key_self,c_key_rhs) bind(c,name='qg_fields_self_sub_f90')
241 
242 implicit none
243 
244 ! Passed variables
245 integer(c_int),intent(in) :: c_key_self !< Fields
246 integer(c_int),intent(in) :: c_key_rhs !< Right-hand side
247 
248 ! Local variables
249 type(qg_fields),pointer :: self
250 type(qg_fields),pointer :: rhs
251 
252 ! Interface
253 call qg_fields_registry%get(c_key_self,self)
254 call qg_fields_registry%get(c_key_rhs,rhs)
255 
256 ! Call Fortran
257 call qg_fields_self_sub(self,rhs)
258 
259 end subroutine qg_fields_self_sub_c
260 ! ------------------------------------------------------------------------------
261 !> Multiply fields by a scalar
262 subroutine qg_fields_self_mul_c(c_key_self,c_zz) bind(c,name='qg_fields_self_mul_f90')
263 
264 implicit none
265 
266 ! Passed variables
267 integer(c_int),intent(in) :: c_key_self !< Fields
268 real(c_double),intent(in) :: c_zz !< Multiplier
269 
270 ! Local variables
271 type(qg_fields),pointer :: self
272 
273 ! Interface
274 call qg_fields_registry%get(c_key_self,self)
275 
276 ! Call Fortran
277 call qg_fields_self_mul(self,c_zz)
278 
279 end subroutine qg_fields_self_mul_c
280 ! ------------------------------------------------------------------------------
281 !> Apply axpy operator to fields
282 subroutine qg_fields_axpy_c(c_key_self,c_zz,c_key_rhs) bind(c,name='qg_fields_axpy_f90')
283 
284 implicit none
285 
286 ! Passed variables
287 integer(c_int),intent(in) :: c_key_self !< Fields
288 real(c_double),intent(in) :: c_zz !< Multiplier
289 integer(c_int),intent(in) :: c_key_rhs !< Right-hand side
290 
291 ! Local variables
292 type(qg_fields),pointer :: self
293 type(qg_fields),pointer :: rhs
294 
295 ! Interface
296 call qg_fields_registry%get(c_key_self,self)
297 call qg_fields_registry%get(c_key_rhs,rhs)
298 
299 ! Call Fortran
300 call qg_fields_axpy(self,c_zz,rhs)
301 
302 end subroutine qg_fields_axpy_c
303 ! ------------------------------------------------------------------------------
304 !> Schur product of fields
305 subroutine qg_fields_self_schur_c(c_key_self,c_key_rhs) bind(c,name='qg_fields_self_schur_f90')
306 
307 implicit none
308 
309 ! Passed variables
310 integer(c_int),intent(in) :: c_key_self !< Fields
311 integer(c_int),intent(in) :: c_key_rhs !< Right-hand side
312 
313 ! Local variables
314 type(qg_fields),pointer :: self
315 type(qg_fields),pointer :: rhs
316 
317 ! Interface
318 call qg_fields_registry%get(c_key_self,self)
319 call qg_fields_registry%get(c_key_rhs,rhs)
320 
321 ! Call Fortran
322 call qg_fields_self_schur(self,rhs)
323 
324 end subroutine qg_fields_self_schur_c
325 ! ------------------------------------------------------------------------------
326 !> Compute dot product for fields
327 subroutine qg_fields_dot_prod_c(c_key_fld1,c_key_fld2,c_prod) bind(c,name='qg_fields_dot_prod_f90')
328 
329 implicit none
330 
331 ! Passed variables
332 integer(c_int),intent(in) :: c_key_fld1 !< First fields
333 integer(c_int),intent(in) :: c_key_fld2 !< Second fields
334 real(c_double),intent(inout) :: c_prod !< Dot product
335 
336 ! Local variables
337 type(qg_fields),pointer :: fld1,fld2
338 
339 ! Interface
340 call qg_fields_registry%get(c_key_fld1,fld1)
341 call qg_fields_registry%get(c_key_fld2,fld2)
342 
343 ! Call Fortran
344 call qg_fields_dot_prod(fld1,fld2,c_prod)
345 
346 end subroutine qg_fields_dot_prod_c
347 ! ------------------------------------------------------------------------------
348 !> Add increment to fields
349 subroutine qg_fields_add_incr_c(c_key_self,c_key_rhs) bind(c,name='qg_fields_add_incr_f90')
350 
351 implicit none
352 
353 ! Passed variables
354 integer(c_int),intent(in) :: c_key_self !< Fields
355 integer(c_int),intent(in) :: c_key_rhs !< Right-hand side
356 
357 ! Local variables
358 type(qg_fields),pointer :: self
359 type(qg_fields),pointer :: rhs
360 
361 ! Interface
362 call qg_fields_registry%get(c_key_self,self)
363 call qg_fields_registry%get(c_key_rhs,rhs)
364 
365 ! Call Fortran
366 call qg_fields_add_incr(self,rhs)
367 
368 end subroutine qg_fields_add_incr_c
369 ! ------------------------------------------------------------------------------
370 !> Compute increment from the difference of two fields
371 subroutine qg_fields_diff_incr_c(c_key_lhs,c_key_fld1,c_key_fld2) bind(c,name='qg_fields_diff_incr_f90')
372 
373 implicit none
374 
375 ! Passed variables
376 integer(c_int),intent(in) :: c_key_lhs !< Left-hand side
377 integer(c_int),intent(in) :: c_key_fld1 !< First fields
378 integer(c_int),intent(in) :: c_key_fld2 !< Second fields
379 
380 ! Local variables
381 type(qg_fields),pointer :: lhs
382 type(qg_fields),pointer :: fld1
383 type(qg_fields),pointer :: fld2
384 
385 ! Interface
386 call qg_fields_registry%get(c_key_lhs,lhs)
387 call qg_fields_registry%get(c_key_fld1,fld1)
388 call qg_fields_registry%get(c_key_fld2,fld2)
389 
390 ! Call Fortran
391 call qg_fields_diff_incr(lhs,fld1,fld2)
392 
393 end subroutine qg_fields_diff_incr_c
394 ! ------------------------------------------------------------------------------
395 !> Change fields resolution
396 subroutine qg_fields_change_resol_c(c_key_fld,c_key_rhs) bind(c,name='qg_fields_change_resol_f90')
397 
398 implicit none
399 
400 ! Passed variables
401 integer(c_int),intent(in) :: c_key_fld !< Fields
402 integer(c_int),intent(in) :: c_key_rhs !< Right-hand side
403 
404 ! Local variables
405 type(qg_fields),pointer :: fld,rhs
406 
407 ! Interface
408 call qg_fields_registry%get(c_key_fld,fld)
409 call qg_fields_registry%get(c_key_rhs,rhs)
410 
411 ! Call Fortran
412 call qg_fields_change_resol(fld,rhs)
413 
414 end subroutine qg_fields_change_resol_c
415 ! ------------------------------------------------------------------------------
416 !> Read fields from file
417 subroutine qg_fields_read_file_c(c_key_fld,c_conf,c_dt) bind(c,name='qg_fields_read_file_f90')
418 
419 implicit none
420 
421 ! Passed variables
422 integer(c_int),intent(in) :: c_key_fld !< Fields
423 type(c_ptr),value,intent(in) :: c_conf !< Configuration
424 type(c_ptr),value,intent(in) :: c_dt !< Date and time
425 
426 ! Local variables
427 type(fckit_configuration) :: f_conf
428 type(qg_fields),pointer :: fld
429 type(datetime) :: fdate
430 
431 ! Interface
432 f_conf = fckit_configuration(c_conf)
433 call qg_fields_registry%get(c_key_fld,fld)
434 call c_f_datetime(c_dt,fdate)
435 
436 ! Call Fortran
437 call qg_fields_read_file(fld,f_conf,fdate)
438 
439 end subroutine qg_fields_read_file_c
440 ! ------------------------------------------------------------------------------
441 !> Write fields to file
442 subroutine qg_fields_write_file_c(c_key_fld,c_conf,c_dt) bind(c,name='qg_fields_write_file_f90')
443 
444 implicit none
445 
446 ! Passed variables
447 integer(c_int),intent(in) :: c_key_fld !< Fields
448 type(c_ptr),value,intent(in) :: c_conf !< Configuration
449 type(c_ptr),value,intent(in) :: c_dt !< Date and time
450 
451 ! Local variables
452 type(fckit_configuration) :: f_conf
453 type(qg_fields),pointer :: fld
454 type(datetime) :: fdate
455 
456 ! Interface
457 f_conf = fckit_configuration(c_conf)
458 call qg_fields_registry%get(c_key_fld,fld)
459 call c_f_datetime(c_dt,fdate)
460 
461 ! Call Fortran
462 call qg_fields_write_file(fld,f_conf,fdate)
463 
464 end subroutine qg_fields_write_file_c
465 ! ------------------------------------------------------------------------------
466 !> Analytic initialization of fields
467 subroutine qg_fields_analytic_init_c(c_key_fld,c_conf,c_dt) bind(c,name='qg_fields_analytic_init_f90')
468 
469 implicit none
470 
471 ! Passed variables
472 integer(c_int),intent(in) :: c_key_fld !< Fields
473 type(c_ptr),value,intent(in) :: c_conf !< Configuration
474 type(c_ptr),value,intent(in) :: c_dt !< Date and time
475 
476 ! Local variables
477 type(fckit_configuration) :: f_conf
478 type(qg_fields),pointer :: fld
479 type(datetime) :: fdate
480 
481 ! Interface
482 f_conf = fckit_configuration(c_conf)
483 call qg_fields_registry%get(c_key_fld,fld)
484 call c_f_datetime(c_dt,fdate)
485 
486 ! Call Fortran
487  call qg_fields_analytic_init(fld,f_conf,fdate)
488 
489 end subroutine qg_fields_analytic_init_c
490 ! ------------------------------------------------------------------------------
491 !> Fields statistics
492 subroutine qg_fields_gpnorm_c(c_key_fld,vpresent,vmin,vmax,vrms) bind(c,name='qg_fields_gpnorm_f90')
493 
494 implicit none
495 
496 ! Passed variables
497 integer(c_int),intent(in) :: c_key_fld !< Fields
498 integer(c_int),intent(inout) :: vpresent(6) !< Variables presence flag
499 real(c_double),intent(inout) :: vmin(6) !< Variables minimum
500 real(c_double),intent(inout) :: vmax(6) !< Variables maximum
501 real(c_double),intent(inout) :: vrms(6) !< Variables RMS
502 
503 ! Local variables
504 type(qg_fields),pointer :: fld
505 
506 ! Interface
507 call qg_fields_registry%get(c_key_fld,fld)
508 
509 ! Call Fortran
510 call qg_fields_gpnorm(fld,vpresent,vmin,vmax,vrms)
511 
512 end subroutine qg_fields_gpnorm_c
513 ! ------------------------------------------------------------------------------
514 !> Fields RMS
515 subroutine qg_fields_rms_c(c_key_fld,prms) bind(c,name='qg_fields_rms_f90')
516 
517 implicit none
518 
519 ! Passed variables
520 integer(c_int),intent(in) :: c_key_fld !< Fields
521 real(c_double),intent(inout) :: prms
522 
523 ! Local variables
524 type(qg_fields),pointer :: fld
525 
526 ! Interface
527 call qg_fields_registry%get(c_key_fld,fld)
528 
529 ! Call Fortran
530 call qg_fields_rms(fld,prms)
531 
532 end subroutine qg_fields_rms_c
533 ! ------------------------------------------------------------------------------
534 !> Get fields geometry
535 subroutine qg_fields_sizes_c(c_key_fld,c_nx,c_ny,c_nz) bind(c,name='qg_fields_sizes_f90')
536 
537 implicit none
538 
539 ! Passed variables
540 integer(c_int),intent(in) :: c_key_fld !< Fields
541 integer(c_int),intent(inout) :: c_nx !< X size
542 integer(c_int),intent(inout) :: c_ny !< Y size
543 integer(c_int),intent(inout) :: c_nz !< Z size
544 
545 ! Local variables
546 type(qg_fields),pointer :: fld
547 
548 ! Interface
549 call qg_fields_registry%get(c_key_fld,fld)
550 
551 ! Call Fortran
552 call qg_fields_sizes(fld,c_nx,c_ny,c_nz)
553 
554 end subroutine qg_fields_sizes_c
555 ! ------------------------------------------------------------------------------
556 !> Get fields geometry
557 subroutine qg_fields_lbc_c(c_key_fld,c_lbc) bind(c,name='qg_fields_lbc_f90')
558 
559 implicit none
560 
561 ! Passed variables
562 integer(c_int),intent(in) :: c_key_fld !< Fields
563 integer(c_int),intent(inout) :: c_lbc !< LBC presence
564 
565 ! Local variables
566 type(qg_fields),pointer :: fld
567 
568 ! Interface
569 call qg_fields_registry%get(c_key_fld,fld)
570 
571 ! Call Fortran
572 call qg_fields_lbc(fld,c_lbc)
573 
574 end subroutine qg_fields_lbc_c
575 ! ------------------------------------------------------------------------------
576 !> Create ATLAS fields
577 subroutine qg_fields_set_atlas_c(c_key_fld,c_vars,c_afieldset) bind (c,name='qg_fields_set_atlas_f90')
578 
579 implicit none
580 
581 ! Passed variables
582 integer(c_int),intent(in) :: c_key_fld !< Fields
583 type(c_ptr),value,intent(in) :: c_vars !< List of variables
584 type(c_ptr),intent(in),value :: c_afieldset !< ATLAS fieldset pointer
585 
586 ! Local variables
587 type(qg_fields),pointer :: fld
588 type(oops_variables) :: vars
589 type(atlas_fieldset) :: afieldset
590 
591 ! Interface
592 call qg_fields_registry%get(c_key_fld,fld)
593 vars = oops_variables(c_vars)
594 afieldset = atlas_fieldset(c_afieldset)
595 
596 ! Call Fortran
597 call qg_fields_set_atlas(fld,vars,afieldset)
598 
599 end subroutine qg_fields_set_atlas_c
600 ! ------------------------------------------------------------------------------
601 !> Convert fields to ATLAS
602 subroutine qg_fields_to_atlas_c(c_key_fld,c_vars,c_afieldset) bind (c,name='qg_fields_to_atlas_f90')
603 
604 implicit none
605 
606 ! Passed variables
607 integer(c_int),intent(in) :: c_key_fld !< Fields
608 type(c_ptr),value,intent(in) :: c_vars !< List of variables
609 type(c_ptr),intent(in),value :: c_afieldset !< ATLAS fieldset pointer
610 
611 ! Local variables
612 type(qg_fields),pointer :: fld
613 type(oops_variables) :: vars
614 type(atlas_fieldset) :: afieldset
615 
616 ! Interface
617 call qg_fields_registry%get(c_key_fld,fld)
618 vars = oops_variables(c_vars)
619 afieldset = atlas_fieldset(c_afieldset)
620 
621 ! Call Fortran
622 call qg_fields_to_atlas(fld,vars,afieldset)
623 
624 end subroutine qg_fields_to_atlas_c
625 ! ------------------------------------------------------------------------------
626 !> Get fields from ATLAS
627 subroutine qg_fields_from_atlas_c(c_key_fld,c_vars,c_afieldset) bind (c,name='qg_fields_from_atlas_f90')
628 
629 implicit none
630 
631 ! Passed variables
632 integer(c_int),intent(in) :: c_key_fld !< Fields
633 type(c_ptr),value,intent(in) :: c_vars !< List of variables
634 type(c_ptr),intent(in),value :: c_afieldset !< ATLAS fieldset pointer
635 
636 ! Local variables
637 type(qg_fields),pointer :: fld
638 type(oops_variables) :: vars
639 type(atlas_fieldset) :: afieldset
640 
641 ! Interface
642 call qg_fields_registry%get(c_key_fld,fld)
643 vars = oops_variables(c_vars)
644 afieldset = atlas_fieldset(c_afieldset)
645 
646 ! Call Fortran
647 call qg_fields_from_atlas(fld,vars,afieldset)
648 
649 end subroutine qg_fields_from_atlas_c
650 ! ------------------------------------------------------------------------------
651 !> Get points from fields
652 subroutine qg_fields_getpoint_c(c_key_fld,c_key_iter,c_nval,c_vals) bind(c,name='qg_fields_getpoint_f90')
653 
654 implicit none
655 
656 ! Passed variables
657 integer(c_int),intent(in) :: c_key_fld !< Fields
658 integer(c_int),intent(in) :: c_key_iter !< Geometry iterator
659 integer(c_int),intent(in) :: c_nval !< Number of values
660 real(c_double),intent(inout) :: c_vals(c_nval) !< Values
661 
662 ! Local variables
663 type(qg_fields),pointer :: fld
664 type(qg_geom_iter),pointer :: iter
665 
666 ! Interface
667 call qg_fields_registry%get(c_key_fld,fld)
668 call qg_geom_iter_registry%get(c_key_iter,iter)
669 
670 ! Call Fortran
671 call qg_fields_getpoint(fld,iter,c_nval,c_vals)
672 
673 end subroutine qg_fields_getpoint_c
674 ! ------------------------------------------------------------------------------
675 !> Set points for the fields
676 subroutine qg_fields_setpoint_c(c_key_fld,c_key_iter,c_nval,c_vals) bind(c,name='qg_fields_setpoint_f90')
677 
678 implicit none
679 
680 ! Passed variables
681 integer(c_int),intent(in) :: c_key_fld !< Fields
682 integer(c_int),intent(in) :: c_key_iter !< Geometry iterator
683 integer(c_int),intent(in) :: c_nval !< Number of values
684 real(c_double),intent(in) :: c_vals(c_nval) !< Values
685 
686 ! Local variables
687 type(qg_fields),pointer :: fld
688 type(qg_geom_iter),pointer :: iter
689 
690 ! Interface
691 call qg_fields_registry%get(c_key_fld,fld)
692 call qg_geom_iter_registry%get(c_key_iter,iter)
693 
694 ! Call Fortran
695 call qg_fields_setpoint(fld,iter,c_nval,c_vals)
696 
697 end subroutine qg_fields_setpoint_c
698 ! ------------------------------------------------------------------------------
699 !> Serialize fields
700 subroutine qg_fields_serialize_c(c_key_fld,c_vsize,c_vect_fld) bind(c,name='qg_fields_serialize_f90')
701 
702 implicit none
703 
704 ! Passed variables
705 integer(c_int),intent(in) :: c_key_fld !< Fields
706 integer(c_int),intent(in) :: c_vsize !< Size
707 real(c_double),intent(out) :: c_vect_fld(c_vsize) !< Vector
708 
709 ! Local variables
710 type(qg_fields),pointer :: fld
711 
712 ! Interface
713 call qg_fields_registry%get(c_key_fld,fld)
714 
715 ! Call Fortran
716 call qg_fields_serialize(fld,c_vsize,c_vect_fld)
717 
718 end subroutine qg_fields_serialize_c
719 ! ------------------------------------------------------------------------------
720 !> Deserialize fields
721 subroutine qg_fields_deserialize_c(c_key_self,c_vsize,c_vect_fld,c_index) bind(c,name='qg_fields_deserialize_f90')
722 
723 implicit none
724 
725 ! Passed variables
726 integer(c_int),intent(in) :: c_key_self !< Fields
727 integer(c_int),intent(in) :: c_vsize !< Size
728 real(c_double),intent(in) :: c_vect_fld(c_vsize) !< Vector
729 integer(c_int), intent(inout):: c_index !< Index
730 
731 ! Local variables
732 type(qg_fields),pointer :: self
733 
734 ! Interface
735 call qg_fields_registry%get(c_key_self,self)
736 
737 ! Call Fortran
738 call qg_fields_deserialize(self,c_vsize,c_vect_fld,c_index)
739 
740 end subroutine qg_fields_deserialize_c
741 ! ------------------------------------------------------------------------------
742 end module qg_fields_interface
Fortran interface to Variables.
subroutine qg_fields_copy_c(c_key_self, c_key_other)
Copy fields.
subroutine qg_fields_from_atlas_c(c_key_fld, c_vars, c_afieldset)
Get fields from ATLAS.
subroutine qg_fields_create_c(c_key_self, c_key_geom, c_vars, c_lbc)
Create fields from geometry and variables.
subroutine qg_fields_dirac_c(c_key_self, c_conf)
Set fields to Diracs.
subroutine qg_fields_delete_c(c_key_self)
Delete fields.
subroutine qg_fields_analytic_init_c(c_key_fld, c_conf, c_dt)
Analytic initialization of fields.
subroutine qg_fields_self_add_c(c_key_self, c_key_rhs)
Add fields.
subroutine qg_fields_self_schur_c(c_key_self, c_key_rhs)
Schur product of fields.
subroutine qg_fields_write_file_c(c_key_fld, c_conf, c_dt)
Write fields to file.
subroutine qg_fields_create_from_other_c(c_key_self, c_key_other, c_key_geom)
Create fields from another one.
subroutine qg_fields_lbc_c(c_key_fld, c_lbc)
Get fields geometry.
subroutine qg_fields_read_file_c(c_key_fld, c_conf, c_dt)
Read fields from file.
subroutine qg_fields_axpy_c(c_key_self, c_zz, c_key_rhs)
Apply axpy operator to fields.
subroutine qg_fields_getpoint_c(c_key_fld, c_key_iter, c_nval, c_vals)
Get points from fields.
subroutine qg_fields_serialize_c(c_key_fld, c_vsize, c_vect_fld)
Serialize fields.
subroutine qg_fields_rms_c(c_key_fld, prms)
Fields RMS.
subroutine qg_fields_gpnorm_c(c_key_fld, vpresent, vmin, vmax, vrms)
Fields statistics.
subroutine qg_fields_to_atlas_c(c_key_fld, c_vars, c_afieldset)
Convert fields to ATLAS.
subroutine qg_fields_change_resol_c(c_key_fld, c_key_rhs)
Change fields resolution.
subroutine qg_fields_setpoint_c(c_key_fld, c_key_iter, c_nval, c_vals)
Set points for the fields.
subroutine qg_fields_diff_incr_c(c_key_lhs, c_key_fld1, c_key_fld2)
Compute increment from the difference of two fields.
subroutine qg_fields_self_mul_c(c_key_self, c_zz)
Multiply fields by a scalar.
subroutine qg_fields_add_incr_c(c_key_self, c_key_rhs)
Add increment to fields.
subroutine qg_fields_ones_c(c_key_self)
Set fields to ones.
subroutine qg_fields_dot_prod_c(c_key_fld1, c_key_fld2, c_prod)
Compute dot product for fields.
subroutine qg_fields_self_sub_c(c_key_self, c_key_rhs)
Subtract fields.
subroutine qg_fields_sizes_c(c_key_fld, c_nx, c_ny, c_nz)
Get fields geometry.
subroutine qg_fields_set_atlas_c(c_key_fld, c_vars, c_afieldset)
Create ATLAS fields.
subroutine qg_fields_zero_c(c_key_self)
Set fields to zero.
subroutine qg_fields_random_c(c_key_self, c_vars)
Generate random fields.
subroutine qg_fields_deserialize_c(c_key_self, c_vsize, c_vect_fld, c_index)
Deserialize fields.
subroutine, public qg_fields_set_atlas(self, vars, afieldset)
Set ATLAS field.
subroutine, public qg_fields_deserialize(self, vsize, vect_fld, index)
Deserialize fields.
subroutine, public qg_fields_zero(self)
Set fields to zero.
subroutine, public qg_fields_axpy(self, zz, rhs)
Apply axpy operator to fields.
subroutine, public qg_fields_to_atlas(self, vars, afieldset)
Convert fields to ATLAS.
subroutine, public qg_fields_delete(self)
Delete fields.
subroutine, public qg_fields_dirac(self, f_conf)
Set fields to Diracs.
subroutine, public qg_fields_getpoint(fld, iter, nval, vals)
Get points from fields.
subroutine, public qg_fields_gpnorm(fld, vpresent, vmin, vmax, vrms)
Fields statistics.
subroutine, public qg_fields_rms(fld, prms)
Fields RMS.
subroutine, public qg_fields_self_schur(self, rhs)
Schur product of fields.
subroutine, public qg_fields_self_mul(self, zz)
Multiply fields by a scalar.
type(registry_t), public qg_fields_registry
Linked list interface - defines registry_t type.
subroutine, public qg_fields_create(self, geom, vars, lbc)
Linked list implementation.
subroutine, public qg_fields_from_atlas(self, vars, afieldset)
Get fields from ATLAS.
subroutine, public qg_fields_analytic_init(fld, f_conf, vdate)
Analytic initialization of fields.
subroutine, public qg_fields_self_sub(self, rhs)
Subtract fields.
subroutine, public qg_fields_copy(self, other)
Copy fields.
subroutine, public qg_fields_dot_prod(fld1, fld2, zprod)
Compute dot product for fields.
subroutine, public qg_fields_write_file(fld, f_conf, vdate)
Write fields to file.
subroutine, public qg_fields_add_incr(self, rhs)
Add increment to fields.
subroutine, public qg_fields_self_add(self, rhs)
Add fields.
subroutine, public qg_fields_serialize(fld, vsize, vect_fld)
Serialize fields.
subroutine, public qg_fields_read_file(fld, f_conf, vdate)
Read fields from file.
subroutine, public qg_fields_sizes(fld, nx, ny, nz)
Get fields geometry.
subroutine, public qg_fields_setpoint(fld, iter, nval, vals)
Set fields values at a specified gridpoint.
subroutine, public qg_fields_lbc(fld, lbc)
Get LBC presence.
subroutine, public qg_fields_random(self, var)
Generate random fields.
subroutine, public qg_fields_ones(self)
Set fields to ones.
subroutine, public qg_fields_diff_incr(lhs, fld1, fld2)
Compute increment from the difference of two fields.
subroutine, public qg_fields_create_from_other(self, other, geom)
Create fields from another one.
subroutine, public qg_fields_change_resol(fld, rhs)
Change fields resolution.
type(registry_t), public qg_geom_iter_registry
Linked list interface - defines registry_t type.
type(registry_t), public qg_geom_registry
Linked list interface - defines registry_t type.
Definition: qg_geom_mod.F90:53