OOPS
qg_geom_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, atlas_functionspace_pointcloud
12 use fckit_configuration_module, only: fckit_configuration
13 use fckit_log_module,only: fckit_log
14 use kinds
15 use iso_c_binding
17 use qg_geom_mod
18 
19 implicit none
20 
21 private
22 ! ------------------------------------------------------------------------------
23 contains
24 ! ------------------------------------------------------------------------------
25 !> Setup geometry
26 subroutine qg_geom_setup_c(c_key_self,c_conf) bind(c,name='qg_geom_setup_f90')
27 
28 ! Passed variables
29 integer(c_int),intent(inout) :: c_key_self !< Geometry
30 type(c_ptr),value,intent(in) :: c_conf !< Configuration
31 
32 ! Local variables
33 type(fckit_configuration) :: f_conf
34 type(qg_geom),pointer :: self
35 
36 ! Interface
37 f_conf = fckit_configuration(c_conf)
38 call qg_geom_registry%init()
39 call qg_geom_registry%add(c_key_self)
40 call qg_geom_registry%get(c_key_self,self)
41 
42 ! Call Fortran
43 call qg_geom_setup(self,f_conf)
44 
45 end subroutine qg_geom_setup_c
46 ! ------------------------------------------------------------------------------
47 !> Set ATLAS lon/lat field
48 subroutine qg_geom_set_atlas_lonlat_c(c_key_self,c_afieldset) bind(c,name='qg_geom_set_atlas_lonlat_f90')
49 
50 ! Passed variables
51 integer(c_int),intent(in) :: c_key_self !< Geometry
52 type(c_ptr),intent(in),value :: c_afieldset !< ATLAS fieldset pointer
53 
54 ! Local variables
55 type(qg_geom),pointer :: self
56 type(atlas_fieldset) :: afieldset
57 
58 ! Interface
59 call qg_geom_registry%get(c_key_self,self)
60 afieldset = atlas_fieldset(c_afieldset)
61 
62 ! Call Fortran
63 call qg_geom_set_atlas_lonlat(self,afieldset)
64 
65 end subroutine qg_geom_set_atlas_lonlat_c
66 ! ------------------------------------------------------------------------------
67 !> Set ATLAS function space pointer
68 subroutine qg_geom_set_atlas_functionspace_pointer_c(c_key_self,c_afunctionspace) &
69  & bind(c,name='qg_geom_set_atlas_functionspace_pointer_f90')
70 
71 ! Passed variables
72 integer(c_int),intent(in) :: c_key_self !< Geometry
73 type(c_ptr),intent(in),value :: c_afunctionspace !< ATLAS function space pointer
74 
75 ! Local variables
76 type(qg_geom),pointer :: self
77 
78 ! Interface
79 call qg_geom_registry%get(c_key_self,self)
80 self%afunctionspace = atlas_functionspace_pointcloud(c_afunctionspace)
81 
83 ! ------------------------------------------------------------------------------
84 !> Fill ATLAS fieldset
85 subroutine qg_geom_fill_atlas_fieldset_c(c_key_self,c_afieldset) &
86  & bind(c,name='qg_geom_fill_atlas_fieldset_f90')
87 
88 ! Passed variables
89 integer(c_int),intent(in) :: c_key_self !< Geometry
90 type(c_ptr),intent(in),value :: c_afieldset !< ATLAS fieldset pointer
91 
92 ! Local variables
93 type(qg_geom),pointer :: self
94 type(atlas_fieldset) :: afieldset
95 
96 ! Interface
97 call qg_geom_registry%get(c_key_self,self)
98 afieldset = atlas_fieldset(c_afieldset)
99 
100 ! Call Fortran
101 call qg_geom_fill_atlas_fieldset(self,afieldset)
102 
103 end subroutine qg_geom_fill_atlas_fieldset_c
104 ! ------------------------------------------------------------------------------
105 !> Clone geometry
106 subroutine qg_geom_clone_c(c_key_self,c_key_other) bind(c,name='qg_geom_clone_f90')
107 
108 ! Passed variables
109 integer(c_int),intent(inout) :: c_key_self !< Geometry
110 integer(c_int),intent(in) :: c_key_other !< Other geometry
111 
112 ! Local variables
113 type(qg_geom),pointer :: self,other
114 
115 ! Interface
116 call qg_geom_registry%get(c_key_other,other)
117 call qg_geom_registry%init()
118 call qg_geom_registry%add(c_key_self)
119 call qg_geom_registry%get(c_key_self ,self )
120 
121 ! Call Fortran
122 call qg_geom_clone(self,other)
123 
124 end subroutine qg_geom_clone_c
125 ! ------------------------------------------------------------------------------
126 !> Delete geometry
127 subroutine qg_geom_delete_c(c_key_self) bind(c,name='qg_geom_delete_f90')
128 
129 ! Passed variables
130 integer(c_int),intent(inout) :: c_key_self !< Geometry
131 
132 ! Local variables
133 type(qg_geom),pointer :: self
134 
135 ! Interface
136 call qg_geom_registry%get(c_key_self,self)
137 
138 ! Call Fortran
139 call qg_geom_delete(self)
140 
141 ! Clear interface
142 call qg_geom_registry%remove(c_key_self)
143 
144 end subroutine qg_geom_delete_c
145 ! ------------------------------------------------------------------------------
146 !> Get geometry info
147 subroutine qg_geom_info_c(c_key_self,c_nx,c_ny,c_nz,c_deltax,c_deltay) bind(c,name='qg_geom_info_f90')
148 
149 ! Passed variables
150 integer(c_int),intent(in) :: c_key_self !< Geometry
151 integer(c_int),intent(inout) :: c_nx !< Number of points in the zonal direction
152 integer(c_int),intent(inout) :: c_ny !< Number of points in the meridional direction
153 integer(c_int),intent(inout) :: c_nz !< Number of vertical levels
154 real(c_double),intent(inout) :: c_deltax !< Zonal cell size
155 real(c_double),intent(inout) :: c_deltay !< Meridional cell size
156 
157 ! Local variables
158 type(qg_geom),pointer :: self
159 
160 ! Interface
161 call qg_geom_registry%get(c_key_self,self)
162 
163 ! Call Fortran
164 call qg_geom_info(self,c_nx,c_ny,c_nz,c_deltax,c_deltay)
165 
166 end subroutine qg_geom_info_c
167 ! ------------------------------------------------------------------------------
168 !> Get dimensions of computational domain
169 subroutine qg_geom_dimensions_c(lonmin, lonmax, latmin, latmax, zmax) bind(c,name='qg_geom_dimensions_f90')
170 
171  use qg_constants_mod
172  implicit none
173 
174  ! Passed variables
175  real(c_double), intent(inout) :: lonmin !< longitude min
176  real(c_double), intent(inout) :: lonmax !< longitude max
177  real(c_double), intent(inout) :: latmin !< latitude min
178  real(c_double), intent(inout) :: latmax !< latitude max
179  real(c_double), intent(inout) :: zmax !< vertical max
180 
181  call xy_to_lonlat(0.0_kind_real, 0.0_kind_real, lonmin, latmin)
182  call xy_to_lonlat(domain_zonal, domain_meridional, lonmax, latmax)
183 
184  zmax = domain_depth
185 
186 end subroutine qg_geom_dimensions_c
187 
188 ! ------------------------------------------------------------------------------
189 end module qg_geom_interface
real(kind_real), parameter, public domain_depth
Model depth (m)
real(kind_real), parameter, public domain_meridional
Model domain in meridional direction (m)
real(kind_real), parameter, public domain_zonal
Model domain in zonal direction (m)
subroutine qg_geom_delete_c(c_key_self)
Delete geometry.
subroutine qg_geom_info_c(c_key_self, c_nx, c_ny, c_nz, c_deltax, c_deltay)
Get geometry info.
subroutine qg_geom_dimensions_c(lonmin, lonmax, latmin, latmax, zmax)
Get dimensions of computational domain.
subroutine qg_geom_set_atlas_functionspace_pointer_c(c_key_self, c_afunctionspace)
Set ATLAS function space pointer.
subroutine qg_geom_setup_c(c_key_self, c_conf)
Setup geometry.
subroutine qg_geom_clone_c(c_key_self, c_key_other)
Clone geometry.
subroutine qg_geom_fill_atlas_fieldset_c(c_key_self, c_afieldset)
Fill ATLAS fieldset.
subroutine qg_geom_set_atlas_lonlat_c(c_key_self, c_afieldset)
Set ATLAS lon/lat field.
type(registry_t), public qg_geom_registry
Linked list interface - defines registry_t type.
Definition: qg_geom_mod.F90:53
subroutine, public qg_geom_clone(self, other)
Clone geometry.
subroutine, public qg_geom_delete(self)
Delete geometry.
subroutine, public qg_geom_setup(self, f_conf)
Linked list implementation.
Definition: qg_geom_mod.F90:64
subroutine, public qg_geom_fill_atlas_fieldset(self, afieldset)
Fill ATLAS fieldset.
subroutine, public qg_geom_set_atlas_lonlat(self, afieldset)
Set ATLAS lon/lat field.
subroutine, public qg_geom_info(self, nx, ny, nz, deltax, deltay)
Get geometry info.
subroutine, public xy_to_lonlat(x, y, lon, lat, mapfac)
Convert x/y to lon/lat.