Go to the documentation of this file.
11 use atlas_module,
only: atlas_fieldset
13 use fckit_configuration_module,
only: fckit_configuration
30 subroutine qg_fields_create_c(c_key_self,c_key_geom,c_vars,c_lbc) bind(c,name='qg_fields_create_f90')
35 integer(c_int),
intent(inout) :: c_key_self
36 integer(c_int),
intent(in) :: c_key_geom
37 type(c_ptr),
value,
intent(in) :: c_vars
38 logical(c_bool),
intent(in) :: c_lbc
65 integer(c_int),
intent(inout) :: c_key_self
66 integer(c_int),
intent(in) :: c_key_other
89 integer(c_int),
intent(inout) :: c_key_self
111 integer(c_int),
intent(in) :: c_key_self
130 integer(c_int),
intent(in) :: c_key_self
149 integer(c_int),
intent(in) :: c_key_self
150 type(c_ptr),
value,
intent(in) :: c_conf
153 type(fckit_configuration) :: f_conf
157 f_conf = fckit_configuration(c_conf)
171 integer(c_int),
intent(in) :: c_key_self
190 integer(c_int),
intent(in) :: c_key_self
191 integer(c_int),
intent(in) :: c_key_other
212 integer(c_int),
intent(in) :: c_key_self
213 integer(c_int),
intent(in) :: c_key_rhs
234 integer(c_int),
intent(in) :: c_key_self
235 integer(c_int),
intent(in) :: c_key_rhs
256 integer(c_int),
intent(in) :: c_key_self
257 real(c_double),
intent(in) :: c_zz
276 integer(c_int),
intent(in) :: c_key_self
277 real(c_double),
intent(in) :: c_zz
278 integer(c_int),
intent(in) :: c_key_rhs
299 integer(c_int),
intent(in) :: c_key_self
300 integer(c_int),
intent(in) :: c_key_rhs
321 integer(c_int),
intent(in) :: c_key_fld1
322 integer(c_int),
intent(in) :: c_key_fld2
323 real(c_double),
intent(inout) :: c_prod
343 integer(c_int),
intent(in) :: c_key_self
344 integer(c_int),
intent(in) :: c_key_rhs
365 integer(c_int),
intent(in) :: c_key_lhs
366 integer(c_int),
intent(in) :: c_key_fld1
367 integer(c_int),
intent(in) :: c_key_fld2
390 integer(c_int),
intent(in) :: c_key_fld
391 integer(c_int),
intent(in) :: c_key_rhs
411 integer(c_int),
intent(in) :: c_key_fld
412 type(c_ptr),
value,
intent(in) :: c_conf
413 type(c_ptr),
value,
intent(in) :: c_dt
416 type(fckit_configuration) :: f_conf
418 type(datetime) :: fdate
421 f_conf = fckit_configuration(c_conf)
423 call c_f_datetime(c_dt,fdate)
436 integer(c_int),
intent(in) :: c_key_fld
437 type(c_ptr),
value,
intent(in) :: c_conf
438 type(c_ptr),
value,
intent(in) :: c_dt
441 type(fckit_configuration) :: f_conf
443 type(datetime) :: fdate
446 f_conf = fckit_configuration(c_conf)
448 call c_f_datetime(c_dt,fdate)
461 integer(c_int),
intent(in) :: c_key_fld
462 type(c_ptr),
value,
intent(in) :: c_conf
463 type(c_ptr),
value,
intent(in) :: c_dt
466 type(fckit_configuration) :: f_conf
468 type(datetime) :: fdate
471 f_conf = fckit_configuration(c_conf)
473 call c_f_datetime(c_dt,fdate)
486 integer(c_int),
intent(in) :: c_key_fld
487 integer(c_int),
intent(in) :: nb
488 real(c_double),
intent(inout) :: pstat(4*(1+nb))
507 integer(c_int),
intent(in) :: c_key_fld
508 real(c_double),
intent(inout) :: prms
527 integer(c_int),
intent(in) :: c_key_fld
528 integer(c_int),
intent(inout) :: c_nx
529 integer(c_int),
intent(inout) :: c_ny
530 integer(c_int),
intent(inout) :: c_nz
531 integer(c_int),
intent(inout) :: c_nb
550 integer(c_int),
intent(in) :: c_key_fld
551 integer(c_int),
intent(inout) :: c_lq
552 integer(c_int),
intent(inout) :: c_lbc
571 integer(c_int),
intent(in) :: c_key_fld
572 type(c_ptr),
value,
intent(in) :: c_vars
573 type(c_ptr),
intent(in),
value :: c_afieldset
578 type(atlas_fieldset) :: afieldset
583 afieldset = atlas_fieldset(c_afieldset)
596 integer(c_int),
intent(in) :: c_key_fld
597 type(c_ptr),
value,
intent(in) :: c_vars
598 type(c_ptr),
intent(in),
value :: c_afieldset
603 type(atlas_fieldset) :: afieldset
608 afieldset = atlas_fieldset(c_afieldset)
621 integer(c_int),
intent(in) :: c_key_fld
622 type(c_ptr),
value,
intent(in) :: c_vars
623 type(c_ptr),
intent(in),
value :: c_afieldset
628 type(atlas_fieldset) :: afieldset
633 afieldset = atlas_fieldset(c_afieldset)
646 integer(c_int),
intent(in) :: c_key_fld
647 integer(c_int),
intent(in) :: c_key_iter
648 integer(c_int),
intent(in) :: c_nval
649 real(c_double),
intent(inout) :: c_vals(c_nval)
670 integer(c_int),
intent(in) :: c_key_fld
671 integer(c_int),
intent(in) :: c_key_iter
672 integer(c_int),
intent(in) :: c_nval
673 real(c_double),
intent(in) :: c_vals(c_nval)
694 integer(c_int),
intent(in) :: c_key_fld
695 integer(c_int),
intent(in) :: c_vsize
696 real(c_double),
intent(out) :: c_vect_fld(c_vsize)
716 integer(c_int),
intent(in) :: c_key_self
717 integer(c_int),
intent(in) :: c_vsize
718 real(c_double),
intent(in) :: c_vect_fld(c_vsize)
719 integer(c_int),
intent(inout):: c_index
subroutine, public qg_fields_rms(fld, prms)
Fields RMS.
subroutine qg_fields_self_sub_c(c_key_self, c_key_rhs)
Subtract fields.
subroutine, public qg_fields_read_file(fld, f_conf, vdate)
Read fields from file.
subroutine qg_fields_self_mul_c(c_key_self, c_zz)
Multiply fields by a scalar.
subroutine qg_fields_diff_incr_c(c_key_lhs, c_key_fld1, c_key_fld2)
Compute increment from the difference of two fields.
subroutine, public qg_fields_axpy(self, zz, rhs)
Apply axpy operator to fields.
subroutine, public qg_fields_sizes(fld, nx, ny, nz, nb)
Get fields geometry.
subroutine, public qg_fields_ones(self)
Set fields to ones.
subroutine qg_fields_gpnorm_c(c_key_fld, nb, pstat)
Fields statistics.
subroutine qg_fields_delete_c(c_key_self)
Delete fields.
subroutine qg_fields_change_resol_c(c_key_fld, c_key_rhs)
Change fields resolution.
subroutine, public qg_fields_set_atlas(self, vars, afieldset)
Set ATLAS field.
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_zero(self)
Set fields to zero.
subroutine qg_fields_random_c(c_key_self)
Generate random fields.
subroutine, public qg_fields_self_schur(self, rhs)
Schur product of fields.
subroutine, public qg_fields_vars(fld, lq, lbc)
Get fields variables.
subroutine qg_fields_ones_c(c_key_self)
Set fields to ones.
subroutine qg_fields_to_atlas_c(c_key_fld, c_vars, c_afieldset)
Convert fields to ATLAS.
subroutine, public qg_fields_create(self, geom, vars, lbc)
Linked list implementation.
Fortran interface to Variables.
subroutine qg_fields_zero_c(c_key_self)
Set fields to zero.
subroutine, public qg_fields_from_atlas(self, vars, afieldset)
Get fields from ATLAS.
type(registry_t), public qg_geom_registry
Linked list interface - defines registry_t type.
subroutine, public qg_fields_analytic_init(fld, f_conf, vdate)
Analytic initialization of fields.
subroutine qg_fields_dot_prod_c(c_key_fld1, c_key_fld2, c_prod)
Compute dot product for fields.
subroutine, public qg_fields_self_add(self, rhs)
Add fields.
subroutine, public qg_fields_diff_incr(lhs, fld1, fld2)
Compute increment from the difference of two fields.
subroutine qg_fields_rms_c(c_key_fld, prms)
Fields RMS.
subroutine, public qg_fields_gpnorm(fld, nb, pstat)
Fields statistics.
subroutine qg_fields_dirac_c(c_key_self, c_conf)
Set fields to Diracs.
subroutine qg_fields_copy_c(c_key_self, c_key_other)
Copy fields.
subroutine, public qg_fields_self_mul(self, zz)
Multiply fields by a scalar.
subroutine qg_fields_from_atlas_c(c_key_fld, c_vars, c_afieldset)
Get fields from ATLAS.
subroutine qg_fields_getpoint_c(c_key_fld, c_key_iter, c_nval, c_vals)
Get points from fields.
type(registry_t), public qg_geom_iter_registry
Linked list interface - defines registry_t type.
subroutine qg_fields_set_atlas_c(c_key_fld, c_vars, c_afieldset)
Create ATLAS fields.
subroutine, public qg_fields_serialize(fld, vsize, vect_fld)
Serialize fields.
subroutine, public qg_fields_setpoint(fld, iter, nval, vals)
Set fields values at a specified gridpoint.
subroutine qg_fields_vars_c(c_key_fld, c_lq, c_lbc)
Get fields variables.
subroutine, public qg_fields_getpoint(fld, iter, nval, vals)
Get points from fields.
subroutine qg_fields_create_c(c_key_self, c_key_geom, c_vars, c_lbc)
Create fields from geometry and variables.
subroutine, public qg_fields_self_sub(self, rhs)
Subtract fields.
subroutine qg_fields_serialize_c(c_key_fld, c_vsize, c_vect_fld)
Serialize fields.
subroutine, public qg_fields_create_from_other(self, other)
Create fields from another one.
subroutine qg_fields_deserialize_c(c_key_self, c_vsize, c_vect_fld, c_index)
Deserialize fields.
subroutine qg_fields_add_incr_c(c_key_self, c_key_rhs)
Add increment to fields.
subroutine qg_fields_analytic_init_c(c_key_fld, c_conf, c_dt)
Analytic initialization of fields.
subroutine, public qg_fields_random(self)
Generate random fields.
subroutine qg_fields_self_schur_c(c_key_self, c_key_rhs)
Schur product of fields.
subroutine qg_fields_read_file_c(c_key_fld, c_conf, c_dt)
Read fields from file.
subroutine, public qg_fields_add_incr(self, rhs)
Add increment to fields.
subroutine, public qg_fields_deserialize(self, vsize, vect_fld, index)
Deserialize fields.
type(registry_t), public qg_fields_registry
Linked list interface - defines registry_t type.
subroutine qg_fields_self_add_c(c_key_self, c_key_rhs)
Add fields.
subroutine qg_fields_setpoint_c(c_key_fld, c_key_iter, c_nval, c_vals)
Set points for the fields.
subroutine qg_fields_sizes_c(c_key_fld, c_nx, c_ny, c_nz, c_nb)
Get fields geometry.
subroutine, public qg_fields_copy(self, other, bconly)
Copy fields.
subroutine qg_fields_write_file_c(c_key_fld, c_conf, c_dt)
Write fields to file.
subroutine, public qg_fields_dirac(self, f_conf)
Set fields to Diracs.
subroutine, public qg_fields_to_atlas(self, vars, afieldset)
Convert fields to ATLAS.
subroutine qg_fields_axpy_c(c_key_self, c_zz, c_key_rhs)
Apply axpy operator to fields.
subroutine qg_fields_create_from_other_c(c_key_self, c_key_other)
Create fields from another one.
subroutine, public qg_fields_change_resol(fld, rhs)
Change fields resolution.
subroutine, public qg_fields_delete(self)
Delete fields.