OOPS
qg_obsvec_mod.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 iso_c_binding
13 use kinds
14 use missing_values_mod
15 use random_mod
16 
17 implicit none
18 
19 private
20 public :: qg_obsvec
21 public :: qg_obsvec_registry
27 ! ------------------------------------------------------------------------------
28 interface
29  subroutine qg_obsvec_random_i(odb,nn,zz) bind(c,name='qg_obsvec_random_f')
30  use iso_c_binding
31  implicit none
32  type(c_ptr),intent(in) :: odb
33  integer(c_int),intent(in) :: nn
34  real(c_double),intent(inout) :: zz
35  end subroutine qg_obsvec_random_i
36 end interface
37 ! ------------------------------------------------------------------------------
39  integer :: nobs = 0 !< Number of observations
40  integer :: nlev = 0 !< Number of levels
41  real(kind_real),allocatable :: values(:,:) !< Values
42  real(kind_real) :: missing !< Missing value
43 end type qg_obsvec
44 
45 #define LISTED_TYPE qg_obsvec
46 
47 !> Linked list interface - defines registry_t type
48 #include "oops/util/linkedList_i.f"
49 
50 !> Global registry
51 type(registry_t) :: qg_obsvec_registry
52 ! ------------------------------------------------------------------------------
53 contains
54 ! ------------------------------------------------------------------------------
55 ! Public
56 ! ------------------------------------------------------------------------------
57 !> Linked list implementation
58 #include "oops/util/linkedList_c.f"
59 ! ------------------------------------------------------------------------------
60 !> Setup observation vector
61 subroutine qg_obsvec_setup(self,nlev,nobs)
62 
63 implicit none
64 
65 ! Passed variables
66 type(qg_obsvec),intent(inout) :: self !< Observation vector
67 integer,intent(in) :: nlev !< Number of levels
68 integer,intent(in) :: nobs !< Number of observations
69 
70 ! Set sizes
71 self%nlev = nlev
72 self%nobs = nobs
73 
74 ! Release memory
75 if (allocated(self%values)) deallocate(self%values)
76 
77 ! Allocation
78 allocate(self%values(self%nlev,self%nobs))
79 
80 ! Initialization
81 self%values = 0.0_kind_real
82 self%missing = missing_value(self%missing)
83 
84 end subroutine qg_obsvec_setup
85 ! ------------------------------------------------------------------------------
86 !> Clone observation vector
87 subroutine qg_obsvec_clone(self,other)
88 
89 implicit none
90 
91 ! Passed variables
92 type(qg_obsvec),intent(inout) :: self !< Observation vector
93 type(qg_obsvec),intent(in) :: other !< Other observation vector
94 
95 ! Set sizes
96 self%nlev = other%nlev
97 self%nobs = other%nobs
98 self%missing = other%missing
99 
100 ! Allocation
101 allocate(self%values(self%nlev,self%nobs))
102 
103 end subroutine qg_obsvec_clone
104 ! ------------------------------------------------------------------------------
105 !> Delete observation vector
106 subroutine qg_obsvec_delete(self)
107 
108 implicit none
109 
110 ! Passed variables
111 type(qg_obsvec),intent(inout) :: self !< Observation vector
112 
113 ! Release memory
114 deallocate(self%values)
115 
116 end subroutine qg_obsvec_delete
117 ! ------------------------------------------------------------------------------
118 !> Copy observation vector
119 subroutine qg_obsvec_copy(self,other)
120 
121 implicit none
122 
123 ! Passed variables
124 type(qg_obsvec),intent(inout) :: self !< Observation vector
125 type(qg_obsvec),intent(in) :: other !< Other observation vector
126 
127 if ((other%nlev/=self%nlev).or.(other%nobs/=self%nobs)) then
128  ! Release memory
129  deallocate(self%values)
130 
131  ! Set sizes
132  self%nlev = other%nlev
133  self%nobs = other%nobs
134  self%missing = other%missing
135 
136  ! Allocation
137  allocate(self%values(self%nlev,self%nobs))
138 endif
139 
140 ! Copy data
141 self%values = other%values
142 
143 end subroutine qg_obsvec_copy
144 ! ------------------------------------------------------------------------------
145 !> Set observation vector to zero
146 subroutine qg_obsvec_zero(self)
147 
148 implicit none
149 
150 ! Passed variables
151 type(qg_obsvec),intent(inout) :: self !< Observation vector
152 
153 ! Set observation vector to zero
154 self%values = 0.0
155 
156 end subroutine qg_obsvec_zero
157 ! ------------------------------------------------------------------------------
158 !> Set i-th value of observation vector to missing value
159 subroutine qg_obsvec_settomissing_ith(self, i)
160 
161 implicit none
162 
163 ! Passed variables
164 type(qg_obsvec),intent(inout) :: self !< Observation vector
165 integer, intent(in) :: i
166 
167 ! Set observation vector to zero
168 self%values(:,i) = self%missing
169 
170 end subroutine qg_obsvec_settomissing_ith
171 ! ------------------------------------------------------------------------------
172 !> Set observation vector to ones
173 subroutine qg_obsvec_ones(self)
174 
175 implicit none
176 
177 ! Passed variables
178 type(qg_obsvec),intent(inout) :: self !< Observation vector
179 
180 ! Set observation vector to ones
181 self%values = 1.0
182 
183 end subroutine qg_obsvec_ones
184 ! ------------------------------------------------------------------------------
185 !> Mask observation vector (set values to missing values where mask == 1)
186 subroutine qg_obsvec_mask(self,mask)
187 implicit none
188 
189 ! Passed variables
190 type(qg_obsvec),intent(inout) :: self !< Observation vector
191 type(qg_obsvec),intent(in) :: mask !< mask
192 
193 if ((self%nobs/=mask%nobs).or.(self%nlev/=mask%nlev)) call abor1_ftn('qg_obsvec_mask: inconsistent sizes')
194 
195 where(mask%values == 1) self%values = self%missing
196 
197 end subroutine qg_obsvec_mask
198 ! ------------------------------------------------------------------------------
199 !> Mask observation vector (set values to missing values where mask == missing value)
200 subroutine qg_obsvec_mask_with_missing(self,mask)
201 implicit none
202 
203 ! Passed variables
204 type(qg_obsvec),intent(inout) :: self !< Observation vector
205 type(qg_obsvec),intent(in) :: mask !< mask
206 
207 if ((self%nobs/=mask%nobs).or.(self%nlev/=mask%nlev)) call abor1_ftn('qg_obsvec_mask: inconsistent sizes')
208 
209 where(mask%values == mask%missing) self%values = self%missing
210 
211 end subroutine qg_obsvec_mask_with_missing
212 ! ------------------------------------------------------------------------------
213 !> Multiply observation vector with a scalar
214 subroutine qg_obsvec_mul_scal(self,zz)
215 
216 implicit none
217 
218 ! Passed variables
219 type(qg_obsvec),intent(inout) :: self !< Observation vector
220 real(kind_real),intent(in) :: zz !< Multiplier
221 
222 ! Multiply observation vector with a scalar
223 where(self%values /= self%missing) self%values = zz*self%values
224 
225 end subroutine qg_obsvec_mul_scal
226 ! ------------------------------------------------------------------------------
227 !> Add observation vector
228 subroutine qg_obsvec_add(self,other)
229 
230 implicit none
231 
232 ! Passed variables
233 type(qg_obsvec),intent(inout) :: self !< Observation vector
234 type(qg_obsvec),intent(in) :: other !< Other observation vector
235 
236 if ((self%nobs/=other%nobs).or.(self%nlev/=other%nlev)) call abor1_ftn('qg_obsvec_add: inconsistent sizes')
237 
238 ! Add observation vector
239 where(self%values /= self%missing .and. other%values /= other%missing)
240  self%values = self%values+other%values
241 elsewhere
242  self%values = self%missing
243 endwhere
244 
245 end subroutine qg_obsvec_add
246 ! ------------------------------------------------------------------------------
247 !> Subtract observation vector
248 subroutine qg_obsvec_sub(self,other)
249 
250 implicit none
251 
252 ! Passed variables
253 type(qg_obsvec),intent(inout) :: self !< Observation vector
254 type(qg_obsvec),intent(in) :: other !< Other observation vector
255 
256 if ((self%nobs/=other%nobs).or.(self%nlev/=other%nlev)) call abor1_ftn('qg_obsvec_sub: inconsistent sizes')
257 
258 ! Subtract observation vector
259 where(self%values /= self%missing .and. other%values /= other%missing)
260  self%values = self%values-other%values
261 elsewhere
262  self%values = self%missing
263 endwhere
264 
265 end subroutine qg_obsvec_sub
266 ! ------------------------------------------------------------------------------
267 !> Multiply observation vector
268 subroutine qg_obsvec_mul(self,other)
269 
270 implicit none
271 
272 ! Passed variables
273 type(qg_obsvec),intent(inout) :: self !< Observation vector
274 type(qg_obsvec),intent(in) :: other !< Other observation vector
275 
276 if ((self%nobs/=other%nobs).or.(self%nlev/=other%nlev)) call abor1_ftn('qg_obsvec_mul: inconsistent sizes')
277 
278 ! Multiply observation vector
279 where(self%values /= self%missing .and. other%values /= other%missing)
280  self%values = self%values*other%values
281 elsewhere
282  self%values = self%missing
283 endwhere
284 
285 end subroutine qg_obsvec_mul
286 ! ------------------------------------------------------------------------------
287 !> Divide observation vector
288 subroutine qg_obsvec_div(self,other)
289 
290 implicit none
291 
292 ! Passed variables
293 type(qg_obsvec),intent(inout) :: self !< Observation vector
294 type(qg_obsvec),intent(in) :: other !< Other observation vector
295 
296 if ((self%nobs/=other%nobs).or.(self%nlev/=other%nlev)) call abor1_ftn('qg_obsvec_div: inconsistent sizes')
297 
298 ! Divide observation vector
299 where(self%values /= self%missing .and. other%values /= other%missing)
300  self%values = self%values/other%values
301 elsewhere
302  self%values = self%missing
303 endwhere
304 
305 end subroutine qg_obsvec_div
306 ! ------------------------------------------------------------------------------
307 !> Apply axpy on observation vector
308 subroutine qg_obsvec_axpy(self,zz,other)
309 
310 implicit none
311 
312 ! Passed variables
313 type(qg_obsvec),intent(inout) :: self !< Observation vector
314 real(kind_real),intent(in) :: zz !< Multiplier
315 type(qg_obsvec),intent(in) :: other !< Other observation vector
316 
317 if ((self%nobs/=other%nobs).or.(self%nlev/=other%nlev)) call abor1_ftn('qg_obsvec_axpy: inconsistent sizes')
318 
319 ! Apply axpy on observation vector
320 where(self%values /= self%missing .and. other%values /= other%missing)
321  self%values = self%values+zz*other%values
322 elsewhere
323  self%values = self%missing
324 endwhere
325 
326 end subroutine qg_obsvec_axpy
327 ! ------------------------------------------------------------------------------
328 !> Invert observation vector
329 subroutine qg_obsvec_invert(self)
330 
331 implicit none
332 
333 ! Passed variables
334 type(qg_obsvec),intent(inout) :: self !< Observation vector
335 
336 ! Invert observation vector
337 where(self%values /= self%missing) self%values = 1.0/self%values
338 
339 end subroutine qg_obsvec_invert
340 ! ------------------------------------------------------------------------------
341 !> Generate random observation vector
342 subroutine qg_obsvec_random(c_odb,self)
343 
344 implicit none
345 
346 ! Passed variables
347 type(c_ptr),intent(in) :: c_odb !< Observation data base
348 type(qg_obsvec),intent(inout) :: self !< Observation vector
349 
350 ! Local variables
351 integer :: nval
352 
353 ! Compute total size
354 nval = self%nobs*self%nlev
355 
356 ! Get random values
357 call qg_obsvec_random_i(c_odb,nval,self%values(1,1))
358 
359 end subroutine qg_obsvec_random
360 ! ------------------------------------------------------------------------------
361 !> Compute dot product between observation vectors
362 subroutine qg_obsvec_dotprod(obsvec1,obsvec2,zz)
363 
364 implicit none
365 
366 ! Passed variables
367 type(qg_obsvec),intent(in) :: obsvec1 !< Observation vector 1
368 type(qg_obsvec),intent(in) :: obsvec2 !< Observation vector 2
369 real(kind_real),intent(inout) :: zz !< Dot product
370 
371 ! Local variables
372 integer :: jlev,jobs
373 
374 ! Check sizes
375 if ((obsvec1%nobs/=obsvec2%nobs).or.(obsvec1%nlev/=obsvec2%nlev)) call abor1_ftn('qg_obsvec_dotprod: inconsistent sizes')
376 
377 ! Initalization
378 zz = 0.0
379 
380 ! Loop over values
381 do jobs=1,obsvec1%nobs
382  do jlev=1,obsvec1%nlev
383  if (obsvec1%values(jlev, jobs) /= obsvec1%missing .and. &
384  obsvec2%values(jlev, jobs) /= obsvec2%missing) &
385  zz = zz+obsvec1%values(jlev,jobs)*obsvec2%values(jlev,jobs)
386  enddo
387 enddo
388 
389 end subroutine qg_obsvec_dotprod
390 ! ------------------------------------------------------------------------------
391 !> Compute observation vector statistics
392 subroutine qg_obsvec_stats(self,zmin,zmax,zavg)
393 
394 implicit none
395 
396 ! Passed variables
397 type(qg_obsvec),intent(in):: self !< Observation vector
398 real(kind_real),intent(inout) :: zmin !< Minimum
399 real(kind_real),intent(inout) :: zmax !< Maximum
400 real(kind_real),intent(inout) :: zavg !< Average
401 
402 if (self%nobs*self%nlev>0) then
403  ! Compute statistics
404  if (.not.allocated(self%values)) call abor1_ftn('qg_obsvec_stats: obs vector not allocated')
405  zmin = minval(self%values, mask = (self%values /= self%missing))
406  zmax = maxval(self%values, mask = (self%values /= self%missing))
407  zavg = sum(self%values, mask = (self%values /= self%missing)) / &
408  count(mask = (self%values /= self%missing))
409 else
410  ! Empty observation vector
411  zmin = 0.0
412  zmax = 0.0
413  zavg = 0.0
414 endif
415 
416 end subroutine qg_obsvec_stats
417 
418 ! ------------------------------------------------------------------------------
419 !> Get observation vector size
420 subroutine qg_obsvec_size(self,kobs)
421 implicit none
422 type(qg_obsvec),intent(in) :: self !< Observation vector
423 integer,intent(inout) :: kobs !< Observation vector size
424 kobs = size(self%values) + 2
425 end subroutine qg_obsvec_size
426 
427 ! ------------------------------------------------------------------------------
428 !> Get observation vector size
429 subroutine qg_obsvec_nobs(self,kobs)
430 
431 implicit none
432 
433 ! Passed variables
434 type(qg_obsvec),intent(in) :: self !< Observation vector
435 integer,intent(inout) :: kobs !< Observation vector size
436 
437 ! Get observation vector size
438 kobs = count(mask = (self%values /= self%missing))
439 
440 end subroutine qg_obsvec_nobs
441 
442 ! ------------------------------------------------------------------------------
443 !> Get observation vector size (only non-masked observations)
444 subroutine qg_obsvec_nobs_withmask(self,obsmask,kobs)
445 
446 implicit none
447 
448 ! Passed variables
449 type(qg_obsvec),intent(in) :: self !< Observation vector
450 type(qg_obsvec),intent(in) :: obsmask !< mask
451 integer,intent(inout) :: kobs !< Observation vector size
452 
453 ! Get observation vector size
454 kobs = count(mask = (self%values /= self%missing) .and. &
455  (obsmask%values /= obsmask%missing))
456 
457 end subroutine qg_obsvec_nobs_withmask
458 
459 ! ------------------------------------------------------------------------------
460 !> Get non-missing values from observation vector into vals array
461 subroutine qg_obsvec_get_withmask(self,obsmask,vals,nvals)
462 
463 implicit none
464 
465 ! Passed variables
466 type(qg_obsvec),intent(in) :: self !< Observation vector
467 type(qg_obsvec),intent(in) :: obsmask !< mask
468 integer,intent(in) :: nvals !< Number of non-missing values
469 real(kind_real), dimension(nvals), intent(out) :: vals!< returned value
470 
471 integer :: jobs, jlev, jval
472 
473 jval = 1
474 ! Loop over values
475 do jobs=1,self%nobs
476  do jlev=1,self%nlev
477  if ((self%values(jlev, jobs) /= self%missing) .and. &
478  (obsmask%values(jlev, jobs) /= obsmask%missing)) then
479  if (jval > nvals) call abor1_ftn('qg_obsvec_get: inconsistent vector size')
480  vals(jval) = self%values(jlev, jobs)
481  jval = jval + 1
482  endif
483  enddo
484 enddo
485 
486 end subroutine qg_obsvec_get_withmask
487 
488 ! ------------------------------------------------------------------------------
489 end module qg_obsvec_mod
subroutine, public qg_obsvec_dotprod(obsvec1, obsvec2, zz)
Compute dot product between observation vectors.
subroutine, public qg_obsvec_random(c_odb, self)
Generate random observation vector.
subroutine, public qg_obsvec_ones(self)
Set observation vector to ones.
subroutine, public qg_obsvec_mask(self, mask)
Mask observation vector (set values to missing values where mask == 1)
subroutine, public qg_obsvec_mul(self, other)
Multiply observation vector.
subroutine, public qg_obsvec_zero(self)
Set observation vector to zero.
subroutine, public qg_obsvec_settomissing_ith(self, i)
Set i-th value of observation vector to missing value.
subroutine, public qg_obsvec_delete(self)
Delete observation vector.
subroutine, public qg_obsvec_mul_scal(self, zz)
Multiply observation vector with a scalar.
type(registry_t), public qg_obsvec_registry
Linked list interface - defines registry_t type.
subroutine, public qg_obsvec_size(self, kobs)
Get observation vector size.
subroutine, public qg_obsvec_mask_with_missing(self, mask)
Mask observation vector (set values to missing values where mask == missing value)
subroutine, public qg_obsvec_clone(self, other)
Clone observation vector.
subroutine, public qg_obsvec_get_withmask(self, obsmask, vals, nvals)
Get non-missing values from observation vector into vals array.
subroutine, public qg_obsvec_div(self, other)
Divide observation vector.
subroutine, public qg_obsvec_copy(self, other)
Copy observation vector.
subroutine, public qg_obsvec_setup(self, nlev, nobs)
Linked list implementation.
subroutine, public qg_obsvec_invert(self)
Invert observation vector.
subroutine, public qg_obsvec_nobs_withmask(self, obsmask, kobs)
Get observation vector size (only non-masked observations)
subroutine, public qg_obsvec_add(self, other)
Add observation vector.
subroutine, public qg_obsvec_stats(self, zmin, zmax, zavg)
Compute observation vector statistics.
subroutine, public qg_obsvec_axpy(self, zz, other)
Apply axpy on observation vector.
subroutine, public qg_obsvec_sub(self, other)
Subtract observation vector.
subroutine, public qg_obsvec_nobs(self, kobs)
Get observation vector size.