FV3 Bundle
fv3jedi_state_interface_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2018 UCAR
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 
6 ! ------------------------------------------------------------------------------
7 
9 
11 use config_mod
12 use datetime_mod
13 use duration_mod
14 use iso_c_binding
15 
22 
23 !GetValues
24 use ioda_locs_mod
26 use ufo_vars_mod
30 
31 private
32 public :: fv3jedi_state_registry
33 
34 ! ------------------------------------------------------------------------------
35 
36 contains
37 
38 ! ------------------------------------------------------------------------------
39 
40 subroutine fv3jedi_state_create_c(c_key_self, c_key_geom, c_vars) bind(c,name='fv3jedi_state_create_f90')
41 
42 implicit none
43 integer(c_int), intent(inout) :: c_key_self
44 integer(c_int), intent(in) :: c_key_geom !< Geometry
45 type(c_ptr), intent(in) :: c_vars !< List of variables
46 
47 type(fv3jedi_state), pointer :: self
48 type(fv3jedi_geom), pointer :: geom
49 type(fv3jedi_vars) :: vars
50 
51 call fv3jedi_geom_registry%get(c_key_geom, geom)
52 call fv3jedi_state_registry%init()
53 call fv3jedi_state_registry%add(c_key_self)
54 call fv3jedi_state_registry%get(c_key_self,self)
55 
56 call fv3jedi_vars_create(c_vars,vars)
57 call create(self, geom, vars)
58 
59 end subroutine fv3jedi_state_create_c
60 
61 ! ------------------------------------------------------------------------------
62 
63 subroutine fv3jedi_state_delete_c(c_key_self) bind(c,name='fv3jedi_state_delete_f90')
64 
65 implicit none
66 integer(c_int), intent(inout) :: c_key_self
67 type(fv3jedi_state), pointer :: self
68 
69 call fv3jedi_state_registry%get(c_key_self,self)
70 
71 call delete(self)
72 
73 call fv3jedi_state_registry%remove(c_key_self)
74 
75 end subroutine fv3jedi_state_delete_c
76 
77 ! ------------------------------------------------------------------------------
78 
79 subroutine fv3jedi_state_zero_c(c_key_self) bind(c,name='fv3jedi_state_zero_f90')
80 
81 implicit none
82 integer(c_int), intent(in) :: c_key_self
83 type(fv3jedi_state), pointer :: self
84 
85 call fv3jedi_state_registry%get(c_key_self,self)
86 call zeros(self)
87 
88 end subroutine fv3jedi_state_zero_c
89 
90 ! ------------------------------------------------------------------------------
91 
92 subroutine fv3jedi_state_copy_c(c_key_self,c_key_rhs) bind(c,name='fv3jedi_state_copy_f90')
93 
94 implicit none
95 integer(c_int), intent(in) :: c_key_self
96 integer(c_int), intent(in) :: c_key_rhs
97 
98 type(fv3jedi_state), pointer :: self
99 type(fv3jedi_state), pointer :: rhs
100 call fv3jedi_state_registry%get(c_key_self,self)
101 call fv3jedi_state_registry%get(c_key_rhs,rhs)
102 
103 call copy(self, rhs)
104 
105 end subroutine fv3jedi_state_copy_c
106 
107 ! ------------------------------------------------------------------------------
108 
109 subroutine fv3jedi_state_axpy_c(c_key_self,c_zz,c_key_rhs) bind(c,name='fv3jedi_state_axpy_f90')
111 implicit none
112 integer(c_int), intent(in) :: c_key_self
113 real(c_double), intent(in) :: c_zz
114 integer(c_int), intent(in) :: c_key_rhs
115 
116 type(fv3jedi_state), pointer :: self
117 type(fv3jedi_state), pointer :: rhs
118 real(kind=kind_real) :: zz
119 
120 call fv3jedi_state_registry%get(c_key_self,self)
121 call fv3jedi_state_registry%get(c_key_rhs,rhs)
122 zz = c_zz
123 
124 call axpy(self,zz,rhs)
125 
126 end subroutine fv3jedi_state_axpy_c
127 
128 ! ------------------------------------------------------------------------------
129 
130 subroutine fv3jedi_state_add_incr_c(c_key_geom,c_key_self,c_key_rhs) bind(c,name='fv3jedi_state_add_incr_f90')
132 implicit none
133 integer(c_int), intent(in) :: c_key_geom
134 integer(c_int), intent(in) :: c_key_self
135 integer(c_int), intent(in) :: c_key_rhs
136 type(fv3jedi_geom), pointer :: geom
137 type(fv3jedi_state), pointer :: self
138 type(fv3jedi_increment), pointer :: rhs
139 
140 call fv3jedi_state_registry%get(c_key_self,self)
141 call fv3jedi_increment_registry%get(c_key_rhs,rhs)
142 call fv3jedi_geom_registry%get(c_key_geom, geom)
143 
144 call add_incr(geom,self,rhs)
145 
146 end subroutine fv3jedi_state_add_incr_c
147 
148 ! ------------------------------------------------------------------------------
149 
150 subroutine fv3jedi_state_change_resol_c(c_key_state,c_key_rhs) bind(c,name='fv3jedi_state_change_resol_f90')
152 implicit none
153 integer(c_int), intent(in) :: c_key_state
154 integer(c_int), intent(in) :: c_key_rhs
155 type(fv3jedi_state), pointer :: state, rhs
156 
157 call fv3jedi_state_registry%get(c_key_state,state)
158 call fv3jedi_state_registry%get(c_key_rhs,rhs)
159 
160 call change_resol(state,rhs)
161 
162 end subroutine fv3jedi_state_change_resol_c
163 
164 ! ------------------------------------------------------------------------------
165 
166 subroutine fv3jedi_state_read_file_c(c_key_geom, c_key_state, c_conf, c_dt) bind(c,name='fv3jedi_state_read_file_f90')
168 implicit none
169 integer(c_int), intent(in) :: c_key_state !< State
170 type(c_ptr), intent(in) :: c_conf !< Configuration
171 type(c_ptr), intent(inout) :: c_dt !< DateTime
172 integer(c_int), intent(in) :: c_key_geom !< Geometry
173 
174 type(fv3jedi_state), pointer :: state
175 type(datetime) :: fdate
176 type(fv3jedi_geom), pointer :: geom
177 
178 call fv3jedi_geom_registry%get(c_key_geom, geom)
179 call fv3jedi_state_registry%get(c_key_state,state)
180 call c_f_datetime(c_dt, fdate)
181 call read_file(geom, state, c_conf, fdate)
182 
183 end subroutine fv3jedi_state_read_file_c
184 
185 ! ------------------------------------------------------------------------------
186 
187 subroutine fv3jedi_state_analytic_init_c(c_key_state, c_key_geom, c_conf, c_dt) bind(c,name='fv3jedi_state_analytic_init_f90')
189 implicit none
190 integer(c_int), intent(in) :: c_key_state !< State
191 integer(c_int), intent(in) :: c_key_geom !< Geometry
192 type(c_ptr), intent(in) :: c_conf !< Configuration
193 type(c_ptr), intent(inout) :: c_dt !< DateTime
194 
195 type(fv3jedi_state), pointer :: state
196 type(fv3jedi_geom), pointer :: geom
197 type(datetime) :: fdate
198 
199 call fv3jedi_geom_registry%get(c_key_geom, geom)
200 call fv3jedi_state_registry%get(c_key_state,state)
201 call c_f_datetime(c_dt, fdate)
202 call analytic_ic(state, geom, c_conf, fdate)
203 
204 end subroutine fv3jedi_state_analytic_init_c
205 
206 ! ------------------------------------------------------------------------------
207 
208 subroutine fv3jedi_state_write_file_c(c_key_geom, c_key_state, c_conf, c_dt) bind(c,name='fv3jedi_state_write_file_f90')
210 implicit none
211 integer(c_int), intent(in) :: c_key_state !< State
212 type(c_ptr), intent(in) :: c_conf !< Configuration
213 type(c_ptr), intent(in) :: c_dt !< DateTime
214 integer(c_int), intent(in) :: c_key_geom !< Geometry
215 
216 type(fv3jedi_state), pointer :: state
217 type(datetime) :: fdate
218 type(fv3jedi_geom), pointer :: geom
219 
220 call fv3jedi_geom_registry%get(c_key_geom, geom)
221 call fv3jedi_state_registry%get(c_key_state,state)
222 call c_f_datetime(c_dt, fdate)
223 call write_file(geom, state, c_conf, fdate)
224 
225 end subroutine fv3jedi_state_write_file_c
226 
227 ! ------------------------------------------------------------------------------
228 
229 subroutine fv3jedi_state_gpnorm_c(c_key_state, kf, pstat) bind(c,name='fv3jedi_state_gpnorm_f90')
231 implicit none
232 integer(c_int), intent(in) :: c_key_state
233 integer(c_int), intent(in) :: kf
234 real(c_double), intent(inout) :: pstat(3*kf)
235 
236 type(fv3jedi_state), pointer :: state
237 real(kind=kind_real) :: zstat(3, kf)
238 integer :: jj, js, jf
239 
240 call fv3jedi_state_registry%get(c_key_state,state)
241 
242 call gpnorm(state, kf, zstat)
243 jj=0
244 do jf = 1, kf
245  do js = 1, 3
246  jj=jj+1
247  pstat(jj) = zstat(js,jf)
248  enddo
249 enddo
250 
251 end subroutine fv3jedi_state_gpnorm_c
252 
253 ! ------------------------------------------------------------------------------
254 
255 subroutine fv3jedi_state_rms_c(c_key_state, prms) bind(c,name='fv3jedi_state_rms_f90')
257 implicit none
258 integer(c_int), intent(in) :: c_key_state
259 real(c_double), intent(inout) :: prms
260 
261 type(fv3jedi_state), pointer :: state
262 real(kind=kind_real) :: zz
263 
264 call fv3jedi_state_registry%get(c_key_state,state)
265 
266 call staterms(state, zz)
267 
268 prms = zz
269 
270 end subroutine fv3jedi_state_rms_c
271 
272 ! ------------------------------------------------------------------------------
273 
274 subroutine fv3jedi_state_getvalues_notraj_c(c_key_geom, c_key_state,c_key_loc,c_vars,c_key_gom) bind(c,name='fv3jedi_state_getvalues_notraj_f90')
276 implicit none
277 integer(c_int), intent(in) :: c_key_state !< State to be interpolated
278 integer(c_int), intent(in) :: c_key_loc !< List of requested locations
279 type(c_ptr), intent(in) :: c_vars !< List of requested variables
280 integer(c_int), intent(in) :: c_key_gom !< Interpolated values
281 integer(c_int), intent(in) :: c_key_geom !< Geometry
282 type(fv3jedi_state), pointer :: state
283 type(ioda_locs), pointer :: locs
284 type(ufo_geovals), pointer :: gom
285 type(ufo_vars) :: vars
286 type(fv3jedi_geom), pointer :: geom
287 
288 call ufo_vars_setup(vars, c_vars)
289 
290 call fv3jedi_geom_registry%get(c_key_geom, geom)
291 call fv3jedi_state_registry%get(c_key_state, state)
292 call ioda_locs_registry%get(c_key_loc, locs)
293 call ufo_geovals_registry%get(c_key_gom, gom)
294 
295 call getvalues(geom, state, locs, vars, gom)
296 
298 
299 ! ------------------------------------------------------------------------------
300 
301 subroutine fv3jedi_state_getvalues_c(c_key_geom, c_key_state,c_key_loc,c_vars,c_key_gom,c_key_traj) bind(c,name='fv3jedi_state_getvalues_f90')
303 implicit none
304 integer(c_int), intent(in) :: c_key_state !< State to be interpolated
305 integer(c_int), intent(in) :: c_key_loc !< List of requested locations
306 type(c_ptr), intent(in) :: c_vars !< List of requested variables
307 integer(c_int), intent(in) :: c_key_gom !< Interpolated values
308 integer(c_int), intent(in) :: c_key_traj !< Trajectory for interpolation/transforms
309 integer(c_int), intent(in) :: c_key_geom !< Geometry
310 
311 type(fv3jedi_state), pointer :: state
312 type(ioda_locs), pointer :: locs
313 type(ufo_geovals), pointer :: gom
314 type(ufo_vars) :: vars
315 type(fv3jedi_getvalues_traj), pointer :: traj
316 type(fv3jedi_geom), pointer :: geom
317 
318 call ufo_vars_setup(vars, c_vars)
319 
320 call fv3jedi_state_registry%get(c_key_state, state)
321 call fv3jedi_geom_registry%get(c_key_geom, geom)
322 call ioda_locs_registry%get(c_key_loc, locs)
323 call ufo_geovals_registry%get(c_key_gom, gom)
324 call fv3jedi_getvalues_traj_registry%get(c_key_traj, traj)
325 
326 call getvalues(geom, state, locs, vars, gom, traj)
327 
328 end subroutine fv3jedi_state_getvalues_c
329 
330 ! ------------------------------------------------------------------------------
331 
332 subroutine fv3jedi_state_sizes_c(c_key_self,nx,ny,nv) bind(c,name='fv3jedi_state_sizes_f90')
334 implicit none
335 integer(c_int), intent(in) :: c_key_self
336 integer(c_int), intent(inout) :: nx,ny,nv
337 type(fv3jedi_state), pointer :: self
338 
339 call fv3jedi_state_registry%get(c_key_self,self)
340 
341 nv = self%vars%nv
342 nx = self%npx
343 ny = self%npy
344 
345 end subroutine fv3jedi_state_sizes_c
346 
347 ! ------------------------------------------------------------------------------
348 
subroutine fv3jedi_state_getvalues_notraj_c(c_key_geom, c_key_state, c_key_loc, c_vars, c_key_gom)
subroutine, public gpnorm(inc, nf, pstat)
subroutine fv3jedi_state_copy_c(c_key_self, c_key_rhs)
subroutine, public staterms(state, prms)
type(registry_t), public fv3jedi_geom_registry
Linked list interface - defines registry_t type.
Fortran derived type to hold FV3JEDI increment.
subroutine fv3jedi_state_analytic_init_c(c_key_state, c_key_geom, c_conf, c_dt)
subroutine fv3jedi_state_add_incr_c(c_key_geom, c_key_self, c_key_rhs)
subroutine, public change_resol(inc, rhs)
subroutine, public write_file(geom, inc, c_conf, vdate)
subroutine, public copy(self, rhs)
Fortran module handling interpolation trajectory for the FV3 model.
type(registry_t), public fv3jedi_increment_registry
Linked list interface - defines registry_t type.
type(registry_t), public fv3jedi_getvalues_traj_registry
Linked list interface - defines registry_t type.
subroutine fv3jedi_state_create_c(c_key_self, c_key_geom, c_vars)
subroutine fv3jedi_state_axpy_c(c_key_self, c_zz, c_key_rhs)
Fortran derived type to hold geometry data for the FV3JEDI model.
subroutine fv3jedi_state_delete_c(c_key_self)
subroutine fv3jedi_state_change_resol_c(c_key_state, c_key_rhs)
Fortran module to handle variables for the FV3JEDI model.
subroutine, public delete(self)
type(registry_t), public ioda_locs_registry
Linked list interface - defines registry_t type.
Handle state for the FV3JEDI odel.
subroutine fv3jedi_state_gpnorm_c(c_key_state, kf, pstat)
subroutine fv3jedi_state_sizes_c(c_key_self, nx, ny, nv)
subroutine, public fv3jedi_vars_create(c_vars, self)
subroutine, public add_incr(self, rhs)
type(registry_t), public ufo_geovals_registry
Linked list interface - defines registry_t type.
Utilities for increment for the FV3JEDI model.
Utilities for state for the FV3JEDI model.
subroutine fv3jedi_state_read_file_c(c_key_geom, c_key_state, c_conf, c_dt)
subroutine, public read_file(geom, inc, c_conf, vdate)
subroutine, public analytic_ic(state, geom, c_conf, vdate)
Analytic Initialization for the FV3 Model.
Fortran module handling geometry for the FV3 model.
subroutine fv3jedi_state_getvalues_c(c_key_geom, c_key_state, c_key_loc, c_vars, c_key_gom, c_key_traj)
type(registry_t), public fv3jedi_state_registry
Linked list interface - defines registry_t type.
subroutine, public zeros(self)
subroutine fv3jedi_state_rms_c(c_key_state, prms)
Fortran module handling observation locations.
subroutine, public axpy(self, zz, rhs)
subroutine, public create(self, geom, vars)
Fortran module handling geometry for the FV3 model.
subroutine fv3jedi_state_write_file_c(c_key_geom, c_key_state, c_conf, c_dt)
subroutine, public ufo_vars_setup(self, c_vars)