OOPS
qg_obsdb_interface.F90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
2 ! (C) Copyright 2017-2019 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 ! In applying this licence, ECMWF does not waive the privileges and immunities
7 ! granted to it by virtue of its status as an intergovernmental organisation nor
8 ! does it submit to any jurisdiction.
9 
11 
12 use atlas_module
13 use fckit_configuration_module, only: fckit_configuration
14 use datetime_mod
15 use duration_mod
16 use fckit_log_module, only: fckit_log
17 use iso_c_binding
18 use string_f_c_mod
19 use qg_locs_mod
20 use qg_obsdb_mod
21 use qg_obsvec_mod
22 use kinds
23 
24 private
25 ! ------------------------------------------------------------------------------
26 contains
27 ! ------------------------------------------------------------------------------
28 !> Setup observation data
29 subroutine qg_obsdb_setup_c(c_key_self,c_conf,c_winbgn,c_winend) bind(c,name='qg_obsdb_setup_f90')
30 
31 implicit none
32 
33 ! Passed variables
34 integer(c_int),intent(inout) :: c_key_self !< Observation data
35 type(c_ptr),value,intent(in) :: c_conf !< Configuration
36 type(c_ptr),value,intent(in) :: c_winbgn !< Start of window
37 type(c_ptr),value,intent(in) :: c_winend !< End of window
38 
39 ! Local variables
40 type(fckit_configuration) :: f_conf
41 type(qg_obsdb),pointer :: self
42 type(datetime) :: winbgn
43 type(datetime) :: winend
44 
45 ! Interface
46 f_conf = fckit_configuration(c_conf)
47 call qg_obsdb_registry%init()
48 call qg_obsdb_registry%add(c_key_self)
49 call qg_obsdb_registry%get(c_key_self,self)
50 call c_f_datetime(c_winbgn,winbgn)
51 call c_f_datetime(c_winend,winend)
52 
53 ! Call Fortran
54 call qg_obsdb_setup(self,f_conf,winbgn,winend)
55 
56 end subroutine qg_obsdb_setup_c
57 ! ------------------------------------------------------------------------------
58 !> Delete observation data
59 subroutine qg_obsdb_delete_c(c_key_self) bind(c,name='qg_obsdb_delete_f90')
60 
61 implicit none
62 
63 ! Passed variables
64 integer(c_int),intent(inout) :: c_key_self !< Observation data
65 
66 ! Local variables
67 type(qg_obsdb),pointer :: self
68 
69 ! Interface
70 call qg_obsdb_registry%get(c_key_self,self)
71 
72 ! Call Fortran
73 call qg_obsdb_delete(self)
74 
75 ! Clear interface
76 call qg_obsdb_registry%remove(c_key_self)
77 
78 end subroutine qg_obsdb_delete_c
79 ! ------------------------------------------------------------------------------
80 !> Get observation data
81 subroutine qg_obsdb_get_c(c_key_self,lgrp,c_grp,lcol,c_col,c_key_ovec) bind(c,name='qg_obsdb_get_f90')
82 
83 implicit none
84 
85 ! Passed variables
86 integer(c_int),intent(in) :: c_key_self !< Observation data
87 integer(c_int),intent(in) :: lgrp !< Group size
88 character(kind=c_char,len=1),intent(in) :: c_grp(lgrp+1) !< Group name
89 integer(c_int),intent(in) :: lcol !< Column size
90 character(kind=c_char,len=1),intent(in) :: c_col(lcol+1) !< Column name
91 integer(c_int),intent(in) :: c_key_ovec !< Observation vector
92 
93 ! Local variables
94 type(qg_obsdb),pointer :: self
95 type(qg_obsvec),pointer :: ovec
96 character(len=lgrp) :: grp
97 character(len=lcol) :: col
98 
99 ! Interface
100 call qg_obsdb_registry%get(c_key_self,self)
101 call c_f_string(c_grp,grp)
102 call c_f_string(c_col,col)
103 call qg_obsvec_registry%get(c_key_ovec,ovec)
104 
105 ! Call Fortran
106 call qg_obsdb_get(self,trim(grp),trim(col),ovec)
107 
108 end subroutine qg_obsdb_get_c
109 ! ------------------------------------------------------------------------------
110 !> Put observation data
111 subroutine qg_obsdb_put_c(c_key_self,lgrp,c_grp,lcol,c_col,c_key_ovec) bind(c,name='qg_obsdb_put_f90')
112 
113 implicit none
114 
115 ! Passed variables
116 integer(c_int),intent(in) :: c_key_self !< Observation data
117 integer(c_int),intent(in) :: lgrp !< Group size
118 character(kind=c_char,len=1),intent(in) :: c_grp(lgrp+1) !< Group name
119 integer(c_int),intent(in) :: lcol !< Column size
120 character(kind=c_char,len=1),intent(in) :: c_col(lcol+1) !< Column name
121 integer(c_int),intent(in) :: c_key_ovec !< Observation vector
122 
123 ! Local variables
124 type(qg_obsdb),pointer :: self
125 type(qg_obsvec),pointer :: ovec
126 character(len=lgrp) :: grp
127 character(len=lcol) :: col
128 
129 ! Interface
130 call qg_obsdb_registry%get(c_key_self,self)
131 call c_f_string(c_grp,grp)
132 call c_f_string(c_col,col)
133 call qg_obsvec_registry%get(c_key_ovec,ovec)
134 
135 ! Call Fortran
136 call qg_obsdb_put(self,trim(grp),trim(col),ovec)
137 
138 end subroutine qg_obsdb_put_c
139 ! ------------------------------------------------------------------------------
140 !> Get locations from observation data
141 subroutine qg_obsdb_locations_c(c_key_self,lgrp,c_grp,c_fields,c_times) bind(c,name='qg_obsdb_locations_f90')
142 
143 implicit none
144 
145 ! Passed variables
146 integer(c_int),intent(in) :: c_key_self !< Observation data
147 integer(c_int),intent(in) :: lgrp !< Group size
148 character(kind=c_char,len=1),intent(in) :: c_grp(lgrp+1) !< Group name
149 type(c_ptr), intent(in), value :: c_fields !< Locations fieldset
150 type(c_ptr), intent(in), value :: c_times !< times
151 
152 ! Local variables
153 type(qg_obsdb),pointer :: self
154 character(len=lgrp) :: grp
155 type(atlas_fieldset) :: fields
156 
157 ! Interface
158 call qg_obsdb_registry%get(c_key_self,self)
159 call c_f_string(c_grp,grp)
160 fields = atlas_fieldset(c_fields)
161 
162 ! Call Fortran
163 call qg_obsdb_locations(self,grp,fields,c_times)
164 
165 call fields%final()
166 
167 end subroutine qg_obsdb_locations_c
168 ! ------------------------------------------------------------------------------
169 !> Generate observation data
170 subroutine qg_obsdb_generate_c(c_key_self,lgrp,c_grp,c_conf,c_bgn,c_step,ktimes,kobs) bind(c,name='qg_obsdb_generate_f90')
171 
172 implicit none
173 
174 ! Passed variables
175 integer(c_int),intent(in) :: c_key_self !< Observation data
176 integer(c_int),intent(in) :: lgrp !< Group size
177 character(kind=c_char,len=1),intent(in) :: c_grp(lgrp+1) !< Group name
178 type(c_ptr),value,intent(in) :: c_conf !< Configuration
179 type(c_ptr),value,intent(in) :: c_bgn !< Start time
180 type(c_ptr),value,intent(in) :: c_step !< Time-step
181 integer(c_int),intent(in) :: ktimes !< Number of time-slots
182 integer(c_int),intent(inout) :: kobs !< Number of observations
183 
184 ! Mocal variables
185 type(fckit_configuration) :: f_conf
186 type(qg_obsdb),pointer :: self
187 character(len=lgrp) :: grp
188 type(datetime) :: bgn
189 type(duration) :: step
190 
191 ! Interface
192 f_conf = fckit_configuration(c_conf)
193 call qg_obsdb_registry%get(c_key_self,self)
194 call c_f_string(c_grp,grp)
195 call c_f_datetime(c_bgn,bgn)
196 call c_f_duration(c_step,step)
197 
198 ! Call Fortran
199 call qg_obsdb_generate(self,grp,f_conf,bgn,step,ktimes,kobs)
200 
201 end subroutine qg_obsdb_generate_c
202 ! ------------------------------------------------------------------------------
203 !> Get observation data size
204 subroutine qg_obsdb_nobs_c(c_key_self,lgrp,c_grp,kobs) bind(c,name='qg_obsdb_nobs_f90')
205 
206 implicit none
207 
208 ! Passed variables
209 integer(c_int),intent(in) :: c_key_self !< Observation data
210 integer(c_int),intent(in) :: lgrp !< Group size
211 character(kind=c_char,len=1),intent(in) :: c_grp(lgrp+1) !< Group name
212 integer(c_int),intent(inout) :: kobs !< Number of observations
213 
214 ! Local variables
215 type(qg_obsdb),pointer :: self
216 character(len=lgrp) :: grp
217 
218 ! Interface
219 call qg_obsdb_registry%get(c_key_self,self)
220 call c_f_string(c_grp,grp)
221 
222 ! Call Fortran
223 call qg_obsdb_nobs(self,grp,kobs)
224 
225 end subroutine qg_obsdb_nobs_c
226 ! ------------------------------------------------------------------------------
227 end module qg_obsdb_interface
subroutine qg_obsdb_nobs_c(c_key_self, lgrp, c_grp, kobs)
Get observation data size.
subroutine qg_obsdb_setup_c(c_key_self, c_conf, c_winbgn, c_winend)
Setup observation data.
subroutine qg_obsdb_delete_c(c_key_self)
Delete observation data.
subroutine qg_obsdb_locations_c(c_key_self, lgrp, c_grp, c_fields, c_times)
Get locations from observation data.
subroutine qg_obsdb_put_c(c_key_self, lgrp, c_grp, lcol, c_col, c_key_ovec)
Put observation data.
subroutine qg_obsdb_generate_c(c_key_self, lgrp, c_grp, c_conf, c_bgn, c_step, ktimes, kobs)
Generate observation data.
subroutine qg_obsdb_get_c(c_key_self, lgrp, c_grp, lcol, c_col, c_key_ovec)
Get observation data.
subroutine, public qg_obsdb_generate(self, grp, f_conf, bgn, step, ktimes, kobs)
Generate observation data.
subroutine, public qg_obsdb_put(self, grp, col, ovec)
Put observations data.
subroutine, public qg_obsdb_locations(self, grp, fields, c_times)
Get locations from observation data.
subroutine, public qg_obsdb_delete(self)
Delete observation data.
subroutine, public qg_obsdb_get(self, grp, col, ovec)
Get observation data.
type(registry_t), public qg_obsdb_registry
Linked list interface - defines registry_t type.
subroutine, public qg_obsdb_setup(self, f_conf, winbgn, winend)
Linked list implementation.
subroutine, public qg_obsdb_nobs(self, grp, kobs)
Get observation data size.
type(registry_t), public qg_obsvec_registry
Linked list interface - defines registry_t type.