FV3 Bundle
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 
9 !> Interfaces to be called from C++ for Fortran handling of QG model fields
10 
11 ! ------------------------------------------------------------------------------
12 
13 subroutine qg_field_create_c(c_key_self, c_key_geom, c_vars) bind(c,name='qg_field_create_f90')
14 use iso_c_binding
15 use qg_fields
16 use qg_geom_mod
17 use qg_vars_mod
18 implicit none
19 integer(c_int), intent(inout) :: c_key_self
20 integer(c_int), intent(in) :: c_key_geom !< Geometry
21 integer(c_int), dimension(*), intent(in) :: c_vars !< List of variables
22 
23 type(qg_field), pointer :: self
24 type(qg_geom), pointer :: geom
25 type(qg_vars) :: vars
26 
27 call qg_geom_registry%get(c_key_geom, geom)
28 call qg_field_registry%init()
29 call qg_field_registry%add(c_key_self)
30 call qg_field_registry%get(c_key_self,self)
31 
32 call qg_vars_create(vars, c_vars)
33 
34 call create(self, geom, vars)
35 
36 end subroutine qg_field_create_c
37 
38 ! ------------------------------------------------------------------------------
39 
40 subroutine qg_field_delete_c(c_key_self) bind(c,name='qg_field_delete_f90')
41 use iso_c_binding
42 use qg_fields
43 implicit none
44 integer(c_int), intent(inout) :: c_key_self
45 type(qg_field), pointer :: self
46 
47 call qg_field_registry%get(c_key_self,self)
48 
49 call delete(self)
50 
51 call qg_field_registry%remove(c_key_self)
52 
53 end subroutine qg_field_delete_c
54 
55 ! ------------------------------------------------------------------------------
56 
57 subroutine qg_field_zero_c(c_key_self) bind(c,name='qg_field_zero_f90')
58 use iso_c_binding
59 use qg_fields
60 implicit none
61 integer(c_int), intent(in) :: c_key_self
62 type(qg_field), pointer :: self
63 
64 call qg_field_registry%get(c_key_self,self)
65 call zeros(self)
66 
67 end subroutine qg_field_zero_c
68 
69 ! ------------------------------------------------------------------------------
70 
71 subroutine qg_field_dirac_c(c_key_self,c_conf) bind(c,name='qg_field_dirac_f90')
72 use iso_c_binding
73 use qg_fields
74 implicit none
75 integer(c_int), intent(in) :: c_key_self
76 type(c_ptr), intent(in) :: c_conf !< Configuration
77 type(qg_field), pointer :: self
78 
79 call qg_field_registry%get(c_key_self,self)
80 call dirac(self,c_conf)
81 
82 end subroutine qg_field_dirac_c
83 
84 ! ------------------------------------------------------------------------------
85 
86 subroutine qg_field_random_c(c_key_self) bind(c,name='qg_field_random_f90')
87 use iso_c_binding
88 use qg_fields
89 implicit none
90 integer(c_int), intent(in) :: c_key_self
91 type(qg_field), pointer :: self
92 
93 call qg_field_registry%get(c_key_self,self)
94 call random(self)
95 
96 end subroutine qg_field_random_c
97 
98 ! ------------------------------------------------------------------------------
99 
100 subroutine qg_field_copy_c(c_key_self,c_key_rhs) bind(c,name='qg_field_copy_f90')
101 use iso_c_binding
102 use qg_fields
103 implicit none
104 integer(c_int), intent(in) :: c_key_self
105 integer(c_int), intent(in) :: c_key_rhs
106 
107 type(qg_field), pointer :: self
108 type(qg_field), pointer :: rhs
109 call qg_field_registry%get(c_key_self,self)
110 call qg_field_registry%get(c_key_rhs,rhs)
111 
112 call copy(self, rhs)
113 
114 end subroutine qg_field_copy_c
115 
116 ! ------------------------------------------------------------------------------
117 
118 subroutine qg_field_self_add_c(c_key_self,c_key_rhs) bind(c,name='qg_field_self_add_f90')
119 use iso_c_binding
120 use qg_fields
121 implicit none
122 integer(c_int), intent(in) :: c_key_self
123 integer(c_int), intent(in) :: c_key_rhs
124 
125 type(qg_field), pointer :: self
126 type(qg_field), pointer :: rhs
127 call qg_field_registry%get(c_key_self,self)
128 call qg_field_registry%get(c_key_rhs,rhs)
129 
130 call self_add(self,rhs)
131 
132 end subroutine qg_field_self_add_c
133 
134 ! ------------------------------------------------------------------------------
135 
136 subroutine qg_field_self_schur_c(c_key_self,c_key_rhs) bind(c,name='qg_field_self_schur_f90')
137 use iso_c_binding
138 use qg_fields
139 implicit none
140 integer(c_int), intent(in) :: c_key_self
141 integer(c_int), intent(in) :: c_key_rhs
142 
143 type(qg_field), pointer :: self
144 type(qg_field), pointer :: rhs
145 call qg_field_registry%get(c_key_self,self)
146 call qg_field_registry%get(c_key_rhs,rhs)
147 
148 call self_schur(self,rhs)
149 
150 end subroutine qg_field_self_schur_c
151 
152 ! ------------------------------------------------------------------------------
153 
154 subroutine qg_field_self_sub_c(c_key_self,c_key_rhs) bind(c,name='qg_field_self_sub_f90')
155 use iso_c_binding
156 use qg_fields
157 implicit none
158 integer(c_int), intent(in) :: c_key_self
159 integer(c_int), intent(in) :: c_key_rhs
160 
161 type(qg_field), pointer :: self
162 type(qg_field), pointer :: rhs
163 call qg_field_registry%get(c_key_self,self)
164 call qg_field_registry%get(c_key_rhs,rhs)
165 
166 call self_sub(self,rhs)
167 
168 end subroutine qg_field_self_sub_c
169 
170 ! ------------------------------------------------------------------------------
171 
172 subroutine qg_field_self_mul_c(c_key_self,c_zz) bind(c,name='qg_field_self_mul_f90')
173 use iso_c_binding
174 use qg_fields
175 use kinds
176 implicit none
177 integer(c_int), intent(in) :: c_key_self
178 real(c_double), intent(in) :: c_zz
179 type(qg_field), pointer :: self
180 real(kind=kind_real) :: zz
181 
182 call qg_field_registry%get(c_key_self,self)
183 zz = c_zz
184 
185 call self_mul(self,zz)
186 
187 end subroutine qg_field_self_mul_c
188 
189 ! ------------------------------------------------------------------------------
190 
191 subroutine qg_field_axpy_c(c_key_self,c_zz,c_key_rhs) bind(c,name='qg_field_axpy_f90')
192 use iso_c_binding
193 use qg_fields
194 use kinds
195 implicit none
196 integer(c_int), intent(in) :: c_key_self
197 real(c_double), intent(in) :: c_zz
198 integer(c_int), intent(in) :: c_key_rhs
199 
200 type(qg_field), pointer :: self
201 type(qg_field), pointer :: rhs
202 real(kind=kind_real) :: zz
203 
204 call qg_field_registry%get(c_key_self,self)
205 call qg_field_registry%get(c_key_rhs,rhs)
206 zz = c_zz
207 
208 call axpy(self,zz,rhs)
209 
210 end subroutine qg_field_axpy_c
211 
212 ! ------------------------------------------------------------------------------
213 
214 subroutine qg_field_dot_prod_c(c_key_fld1,c_key_fld2,c_prod) bind(c,name='qg_field_dot_prod_f90')
215 use iso_c_binding
216 use qg_fields
217 use kinds
218 implicit none
219 integer(c_int), intent(in) :: c_key_fld1, c_key_fld2
220 real(c_double), intent(inout) :: c_prod
221 real(kind=kind_real) :: zz
222 type(qg_field), pointer :: fld1, fld2
223 
224 call qg_field_registry%get(c_key_fld1,fld1)
225 call qg_field_registry%get(c_key_fld2,fld2)
226 
227 call dot_prod(fld1,fld2,zz)
228 
229 c_prod = zz
230 
231 end subroutine qg_field_dot_prod_c
232 
233 ! ------------------------------------------------------------------------------
234 
235 subroutine qg_field_add_incr_c(c_key_self,c_key_rhs) bind(c,name='qg_field_add_incr_f90')
236 use iso_c_binding
237 use qg_fields
238 implicit none
239 integer(c_int), intent(in) :: c_key_self
240 integer(c_int), intent(in) :: c_key_rhs
241 type(qg_field), pointer :: self
242 type(qg_field), pointer :: rhs
243 
244 call qg_field_registry%get(c_key_self,self)
245 call qg_field_registry%get(c_key_rhs,rhs)
246 
247 call add_incr(self,rhs)
248 
249 end subroutine qg_field_add_incr_c
250 
251 ! ------------------------------------------------------------------------------
252 
253 subroutine qg_field_diff_incr_c(c_key_lhs,c_key_x1,c_key_x2) bind(c,name='qg_field_diff_incr_f90')
254 use iso_c_binding
255 use qg_fields
256 implicit none
257 integer(c_int), intent(in) :: c_key_lhs
258 integer(c_int), intent(in) :: c_key_x1
259 integer(c_int), intent(in) :: c_key_x2
260 type(qg_field), pointer :: lhs
261 type(qg_field), pointer :: x1
262 type(qg_field), pointer :: x2
263 
264 call qg_field_registry%get(c_key_lhs,lhs)
265 call qg_field_registry%get(c_key_x1,x1)
266 call qg_field_registry%get(c_key_x2,x2)
267 
268 call diff_incr(lhs,x1,x2)
269 
270 end subroutine qg_field_diff_incr_c
271 
272 ! ------------------------------------------------------------------------------
273 
274 subroutine qg_field_change_resol_c(c_key_fld,c_key_rhs) bind(c,name='qg_field_change_resol_f90')
275 use iso_c_binding
276 use qg_fields
277 implicit none
278 integer(c_int), intent(in) :: c_key_fld
279 integer(c_int), intent(in) :: c_key_rhs
280 type(qg_field), pointer :: fld, rhs
281 
282 call qg_field_registry%get(c_key_fld,fld)
283 call qg_field_registry%get(c_key_rhs,rhs)
284 
285 call change_resol(fld,rhs)
286 
287 end subroutine qg_field_change_resol_c
288 
289 ! ------------------------------------------------------------------------------
290 
291 subroutine qg_field_ug_coord_c(c_key_fld, c_key_ug, c_colocated) bind (c,name='qg_field_ug_coord_f90')
292 use iso_c_binding
293 use qg_fields
295 implicit none
296 integer(c_int), intent(in) :: c_key_fld
297 integer(c_int), intent(in) :: c_key_ug
298 integer(c_int), intent(in) :: c_colocated
299 type(qg_field), pointer :: fld
300 type(unstructured_grid), pointer :: ug
301 integer :: colocated
302 
303 call qg_field_registry%get(c_key_fld,fld)
304 call unstructured_grid_registry%get(c_key_ug,ug)
305 colocated = c_colocated
306 
307 call ug_coord(fld, ug, colocated)
308 
309 end subroutine qg_field_ug_coord_c
310 
311 ! ------------------------------------------------------------------------------
312 
313 subroutine qg_field_field_to_ug_c(c_key_fld, c_key_ug, c_colocated) bind (c,name='qg_field_field_to_ug_f90')
314 use iso_c_binding
315 use qg_fields
317 implicit none
318 integer(c_int), intent(in) :: c_key_fld
319 integer(c_int), intent(in) :: c_key_ug
320 integer(c_int), intent(in) :: c_colocated
321 type(qg_field), pointer :: fld
322 type(unstructured_grid), pointer :: ug
323 integer :: colocated
324 
325 call qg_field_registry%get(c_key_fld,fld)
326 call unstructured_grid_registry%get(c_key_ug,ug)
327 colocated = c_colocated
328 
329 call field_to_ug(fld, ug, colocated)
330 
331 end subroutine qg_field_field_to_ug_c
332 
333 ! ------------------------------------------------------------------------------
334 
335 subroutine qg_field_field_from_ug_c(c_key_fld, c_key_ug) bind (c,name='qg_field_field_from_ug_f90')
336 use iso_c_binding
337 use qg_fields
339 implicit none
340 integer(c_int), intent(in) :: c_key_fld
341 integer(c_int), intent(in) :: c_key_ug
342 type(qg_field), pointer :: fld
343 type(unstructured_grid), pointer :: ug
344 
345 call qg_field_registry%get(c_key_fld,fld)
346 call unstructured_grid_registry%get(c_key_ug,ug)
347 
348 call field_from_ug(fld, ug)
349 
350 end subroutine qg_field_field_from_ug_c
351 
352 ! ------------------------------------------------------------------------------
353 
354 subroutine qg_field_read_file_c(c_key_fld, c_conf, c_dt) bind(c,name='qg_field_read_file_f90')
355 use iso_c_binding
356 use qg_fields
357 use datetime_mod
358 
359 implicit none
360 integer(c_int), intent(in) :: c_key_fld !< Fields
361 type(c_ptr), intent(in) :: c_conf !< Configuration
362 type(c_ptr), intent(inout) :: c_dt !< DateTime
363 
364 type(qg_field), pointer :: fld
365 type(datetime) :: fdate
366 
367 call qg_field_registry%get(c_key_fld,fld)
368 call c_f_datetime(c_dt, fdate)
369 call read_file(fld, c_conf, fdate)
370 
371 end subroutine qg_field_read_file_c
372 
373 ! ------------------------------------------------------------------------------
374 subroutine qg_field_analytic_init_c(c_key_fld, c_key_geom, c_conf, c_dt) bind(c,name='qg_field_analytic_init_f90')
375 
376  use iso_c_binding
377  use qg_fields
378  use qg_geom_mod
379  use datetime_mod
380 
381  implicit none
382  integer(c_int), intent(in) :: c_key_fld !< Fields
383  integer(c_int), intent(in) :: c_key_geom !< Grid information
384  type(c_ptr), intent(in) :: c_conf !< Configuration
385  type(c_ptr), intent(inout) :: c_dt !< DateTime
386 
387  type(qg_field), pointer :: fld
388  type(qg_geom), pointer :: geom
389  type(datetime) :: fdate
390 
391  call qg_field_registry%get(c_key_fld,fld)
392  call qg_geom_registry%get(c_key_geom,geom)
393  call c_f_datetime(c_dt, fdate)
394 
395  call analytic_init(fld,geom,c_conf,fdate)
396 
397 end subroutine qg_field_analytic_init_c
398 
399 ! ------------------------------------------------------------------------------
400 
401 subroutine qg_field_write_file_c(c_key_fld, c_conf, c_dt) bind(c,name='qg_field_write_file_f90')
402 use iso_c_binding
403 use qg_fields
404 use datetime_mod
405 
406 implicit none
407 integer(c_int), intent(in) :: c_key_fld !< Fields
408 type(c_ptr), intent(in) :: c_conf !< Configuration
409 type(c_ptr), intent(in) :: c_dt !< DateTime
410 
411 type(qg_field), pointer :: fld
412 type(datetime) :: fdate
413 
414 call qg_field_registry%get(c_key_fld,fld)
415 call c_f_datetime(c_dt, fdate)
416 call write_file(fld, c_conf, fdate)
417 
418 end subroutine qg_field_write_file_c
419 
420 ! ------------------------------------------------------------------------------
421 
422 subroutine qg_field_gpnorm_c(c_key_fld, kf, pstat) bind(c,name='qg_field_gpnorm_f90')
423 use iso_c_binding
424 use qg_fields
425 use kinds
426 implicit none
427 integer(c_int), intent(in) :: c_key_fld
428 integer(c_int), intent(in) :: kf
429 real(c_double), intent(inout) :: pstat(3*kf)
430 
431 type(qg_field), pointer :: fld
432 real(kind=kind_real) :: zstat(3, kf)
433 integer :: jj, js, jf
434 
435 call qg_field_registry%get(c_key_fld,fld)
436 
437 call gpnorm(fld, kf, zstat)
438 jj=0
439 do jf = 1, kf
440  do js = 1, 3
441  jj=jj+1
442  pstat(jj) = zstat(js,jf)
443  enddo
444 enddo
445 
446 end subroutine qg_field_gpnorm_c
447 
448 ! ------------------------------------------------------------------------------
449 
450 subroutine qg_field_getpoint_c(c_key_fld, c_key_iter, c_nval, c_vals) bind(c,name='qg_field_getpoint_f90')
451 use iso_c_binding
452 use qg_fields
454 use kinds
455 implicit none
456 integer(c_int), intent(in) :: c_key_fld
457 integer(c_int), intent(in) :: c_key_iter
458 integer(c_int), intent(in) :: c_nval
459 real(c_double), intent(inout) :: c_vals(c_nval)
460 
461 type(qg_field), pointer :: fld
462 type(qg_geom_iter), pointer :: iter
463 
464 call qg_field_registry%get(c_key_fld, fld)
465 call qg_geom_iter_registry%get(c_key_iter, iter)
466 
467 !AS NOTE: this check fails; for some reason the pointers aren't pointing to
468 ! the same thing (???)
469 !if (.not. associated(fld%geom, iter%geom)) then
470 ! print *, 'qg_field_getpoint ERROR: geometries are different!'
471 ! stop
472 !endif
473 
474 if (fld%nf*fld%nl /= c_nval) then
475  print *, 'qg_field_getpoint ERROR: array sizes are different: ', fld%nf*fld%nl, c_nval
476  stop
477 endif
478 
479 c_vals = fld%gfld3d(iter%ilon, iter%ilat, :)
480 
481 end subroutine qg_field_getpoint_c
482 
483 ! ------------------------------------------------------------------------------
484 
485 subroutine qg_field_rms_c(c_key_fld, prms) bind(c,name='qg_field_rms_f90')
486 use iso_c_binding
487 use qg_fields
488 use kinds
489 implicit none
490 integer(c_int), intent(in) :: c_key_fld
491 real(c_double), intent(inout) :: prms
492 
493 type(qg_field), pointer :: fld
494 real(kind=kind_real) :: zz
495 
496 call qg_field_registry%get(c_key_fld,fld)
497 
498 call fldrms(fld, zz)
499 
500 prms = zz
501 
502 end subroutine qg_field_rms_c
503 
504 ! ------------------------------------------------------------------------------
505 
506 subroutine qg_field_interp_tl_c(c_key_fld,c_key_loc,c_vars,c_key_gom) bind(c,name='qg_field_interp_tl_f90')
507 use iso_c_binding
508 use qg_fields
509 use qg_locs_mod
510 use qg_vars_mod
511 use qg_goms_mod
512 implicit none
513 integer(c_int), intent(in) :: c_key_fld !< Fields to be interpolated
514 integer(c_int), intent(in) :: c_key_loc !< List of requested locations
515 integer(c_int), dimension(*), intent(in) :: c_vars !< List of variables
516 integer(c_int), intent(in) :: c_key_gom !< Interpolated values
517 type(qg_field), pointer :: fld
518 type(qg_locs), pointer :: locs
519 type(qg_goms), pointer :: gom
520 type(qg_vars) :: vars
521 
522 call qg_field_registry%get(c_key_fld, fld)
523 call qg_locs_registry%get(c_key_loc, locs)
524 call qg_goms_registry%get(c_key_gom, gom)
525 call qg_vars_create(vars, c_vars)
526 
527 call interp_tl(fld, locs, vars, gom)
528 
529 end subroutine qg_field_interp_tl_c
530 
531 ! ------------------------------------------------------------------------------
532 
533 subroutine qg_field_interp_ad_c(c_key_fld,c_key_loc,c_vars,c_key_gom) bind(c,name='qg_field_interp_ad_f90')
534 use iso_c_binding
535 use qg_fields
536 use qg_locs_mod
537 use qg_vars_mod
538 use qg_goms_mod
539 implicit none
540 integer(c_int), intent(in) :: c_key_fld !< Fields to be interpolated
541 integer(c_int), intent(in) :: c_key_loc !< List of requested locations
542 integer(c_int), dimension(*), intent(in) :: c_vars !< List of variables
543 integer(c_int), intent(in) :: c_key_gom !< Interpolated values
544 type(qg_field), pointer :: fld
545 type(qg_locs), pointer :: locs
546 type(qg_goms), pointer :: gom
547 type(qg_vars) :: vars
548 
549 call qg_field_registry%get(c_key_fld, fld)
550 call qg_locs_registry%get(c_key_loc, locs)
551 call qg_goms_registry%get(c_key_gom, gom)
552 call qg_vars_create(vars, c_vars)
553 
554 call interp_ad(fld, locs, vars, gom)
555 
556 end subroutine qg_field_interp_ad_c
557 
558 ! ------------------------------------------------------------------------------
559 
560 subroutine qg_fieldnum_c(c_key_fld, c_nx, c_ny, c_nf, c_nb, c_nl) bind(c,name='qg_field_sizes_f90')
561 use iso_c_binding
562 use qg_fields
563 implicit none
564 integer(c_int), intent(in) :: c_key_fld
565 integer(kind=c_int), intent(inout) :: c_nx, c_ny, c_nf, c_nb, c_nl
566 type(qg_field), pointer :: fld
567 
568 call qg_field_registry%get(c_key_fld,fld)
569 
570 c_nx = fld%geom%nx
571 c_ny = fld%geom%ny
572 c_nf = fld%nf
573 c_nl = fld%nl
574 c_nb =0
575 if (fld%lbc) c_nb = 2
576 
577 end subroutine qg_fieldnum_c
578 
579 ! ------------------------------------------------------------------------------
subroutine, public self_mul(self, zz)
Definition: qg_fields.F90:315
subroutine qg_field_field_from_ug_c(c_key_fld, c_key_ug)
subroutine, public axpy(self, zz, rhs)
Definition: qg_fields.F90:333
Fortran module handling geometry iterator for the QG model.
type(registry_t), public qg_field_registry
Linked list interface - defines registry_t type.
Definition: qg_fields.F90:65
subroutine, public read_file(fld, c_conf, vdate)
Definition: qg_fields.F90:524
subroutine qg_field_diff_incr_c(c_key_lhs, c_key_x1, c_key_x2)
subroutine qg_field_ug_coord_c(c_key_fld, c_key_ug, c_colocated)
subroutine qg_field_read_file_c(c_key_fld, c_conf, c_dt)
subroutine qg_field_getpoint_c(c_key_fld, c_key_iter, c_nval, c_vals)
subroutine qg_field_gpnorm_c(c_key_fld, kf, pstat)
type(registry_t), public qg_locs_registry
Linked list interface - defines registry_t type.
Definition: qg_locs_mod.F90:37
subroutine, public copy(self, rhs)
Definition: qg_fields.F90:225
subroutine, public create(self, geom, vars)
Linked list implementation.
Definition: qg_fields.F90:76
subroutine qg_field_change_resol_c(c_key_fld, c_key_rhs)
subroutine qg_field_random_c(c_key_self)
subroutine qg_field_interp_ad_c(c_key_fld, c_key_loc, c_vars, c_key_gom)
subroutine qg_field_add_incr_c(c_key_self, c_key_rhs)
subroutine, public random(self)
Definition: qg_fields.F90:208
type(registry_t), public qg_geom_registry
Linked list interface - defines registry_t type.
Definition: qg_geom_mod.F90:40
subroutine, public self_schur(self, rhs)
Definition: qg_fields.F90:273
subroutine qg_field_self_add_c(c_key_self, c_key_rhs)
subroutine, public self_add(self, rhs)
Definition: qg_fields.F90:252
subroutine qg_field_rms_c(c_key_fld, prms)
subroutine qg_field_self_schur_c(c_key_self, c_key_rhs)
subroutine qg_field_axpy_c(c_key_self, c_zz, c_key_rhs)
subroutine, public add_incr(self, rhs)
Definition: qg_fields.F90:392
subroutine qg_field_create_c(c_key_self, c_key_geom, c_vars)
Interfaces to be called from C++ for Fortran handling of QG model fields.
subroutine, public interp_tl(fld, locs, vars, gom)
Definition: qg_fields.F90:1166
subroutine qg_field_interp_tl_c(c_key_fld, c_key_loc, c_vars, c_key_gom)
Fortran module handling geometry for the QG model.
Definition: qg_geom_mod.F90:11
subroutine qg_fieldnum_c(c_key_fld, c_nx, c_ny, c_nf, c_nb, c_nl)
subroutine, public diff_incr(lhs, x1, x2)
Definition: qg_fields.F90:424
subroutine, public dirac(self, c_conf)
Definition: qg_fields.F90:167
subroutine, public ug_coord(self, ug, colocated)
Definition: qg_fields.F90:1422
subroutine, public field_to_ug(self, ug, colocated)
Definition: qg_fields.F90:1479
subroutine, public qg_vars_create(self, kvars)
Linked list implementation.
Definition: qg_vars_mod.F90:46
Fortran module handling observation locations.
Definition: qg_locs_mod.F90:11
subroutine, public write_file(fld, c_conf, vdate)
Definition: qg_fields.F90:994
type(registry_t), public qg_geom_iter_registry
Linked list interface - defines registry_t type.
subroutine qg_field_zero_c(c_key_self)
subroutine qg_field_dot_prod_c(c_key_fld1, c_key_fld2, c_prod)
type(registry_t), public qg_goms_registry
Linked list interface - defines registry_t type.
Definition: qg_goms_mod.F90:42
subroutine qg_field_self_mul_c(c_key_self, c_zz)
subroutine, public field_from_ug(self, ug)
Definition: qg_fields.F90:1534
subroutine qg_field_self_sub_c(c_key_self, c_key_rhs)
subroutine, public interp_ad(fld, locs, vars, gom)
Definition: qg_fields.F90:1253
type(registry_t), public unstructured_grid_registry
Linked list interface - defines registry_t type.
subroutine, public dot_prod(fld1, fld2, zprod)
Definition: qg_fields.F90:355
Fortran module to handle variables for the QG model.
Definition: qg_vars_mod.F90:11
Handle fields for the QG model.
Definition: qg_fields.F90:11
subroutine, public analytic_init(fld, geom, config, vdate)
Analytic Initialization for the QG model.
Definition: qg_fields.F90:707
subroutine, public zeros(self)
Definition: qg_fields.F90:151
subroutine qg_field_write_file_c(c_key_fld, c_conf, c_dt)
subroutine, public fldrms(fld, prms)
Definition: qg_fields.F90:1126
subroutine, public self_sub(self, rhs)
Definition: qg_fields.F90:294
subroutine qg_field_dirac_c(c_key_self, c_conf)
subroutine, public change_resol(fld, rhs)
Definition: qg_fields.F90:463
Fortran module handling interpolated (to obs locations) model variables.
Definition: qg_goms_mod.F90:11
subroutine qg_field_copy_c(c_key_self, c_key_rhs)
Fortran module for handling generic unstructured grid.
subroutine, public gpnorm(fld, nf, pstat)
Definition: qg_fields.F90:1089
subroutine qg_field_delete_c(c_key_self)
subroutine qg_field_field_to_ug_c(c_key_fld, c_key_ug, c_colocated)
subroutine qg_field_analytic_init_c(c_key_fld, c_key_geom, c_conf, c_dt)
subroutine, public delete(self)
Definition: qg_fields.F90:136