UFO
ufo_geovals_mod.F90
Go to the documentation of this file.
1 !
2 ! (C) Copyright 2017-2018 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 !
8 
9 use iso_c_binding
10 use ufo_vars_mod
11 use kinds
12 use obsspace_mod
13 use missing_values_mod
14 
15 use fckit_mpi_module, only: fckit_mpi_comm, fckit_mpi_sum
16 
17 implicit none
18 private
19 integer, parameter :: max_string=800
20 
21 public :: ufo_geovals, ufo_geoval
22 public :: ufo_geovals_get_var
26 public :: ufo_geovals_profmult
34 
36 
37 ! ------------------------------------------------------------------------------
38 
39 !> type to hold interpolated field for one variable, one observation
40 type :: ufo_geoval
41  real(kind_real), allocatable :: vals(:,:) !< values (nval, nlocs)
42  integer :: nval = 0 !< number of values in profile
43  integer :: nlocs = 0 !< number of observations
44 end type ufo_geoval
45 
46 !> type to hold interpolated fields required by the obs operators
47 type :: ufo_geovals
48  integer :: nlocs = 0 !< number of observations
49  integer :: nvar = 0 !< number of variables (supposed to be the
50  ! same for same obs operator
51 
52  type(ufo_geoval), allocatable :: geovals(:) !< array of interpolated
53  ! vertical profiles for all obs (nvar)
54 
55  character(len=MAXVARLEN), allocatable :: variables(:) !< variable list
56 
57  real(c_double) :: missing_value !< obsspace missing value mark
58 
59  logical :: linit = .false. !< .true. if all the ufo_geoval arrays inside geovals
60  ! were allocated and have data
61 end type ufo_geovals
62 
63 ! ------------------------------------------------------------------------------
64 contains
65 ! ------------------------------------------------------------------------------
66 
68 implicit none
69 type(ufo_geovals), intent(inout) :: self
70 
71 self%nlocs = 0
72 self%missing_value = missing_value(0.0)
73 self%nvar = 0
74 self%linit = .false.
75 
76 end subroutine ufo_geovals_default_constr
77 
78 ! ------------------------------------------------------------------------------
79 !> Initializes and allocates \p self GeoVaLs with \p nlocs number of locations for
80 !> \p vars variables. \p nvals array contains number of values to allocate for
81 !> each of the variables
82 subroutine ufo_geovals_setup(self, vars, nlocs, nvars, nvals)
83 use oops_variables_mod
84 implicit none
85 type(ufo_geovals), intent(inout) :: self
86 type(oops_variables), intent(in) :: vars
87 integer, intent(in) :: nlocs, nvars
88 integer(c_size_t), intent(in) :: nvals(nvars)
89 
90 integer :: ivar
91 
92 call ufo_geovals_delete(self)
93 self%nlocs = nlocs
94 self%missing_value = missing_value(self%missing_value)
95 
96 self%nvar = vars%nvars()
97 allocate(self%geovals(self%nvar))
98 allocate(self%variables(self%nvar))
99 do ivar = 1, self%nvar
100  self%variables(ivar) = vars%variable(ivar)
101  self%geovals(ivar)%nlocs = nlocs
102  self%geovals(ivar)%nval = nvals(ivar)
103  allocate(self%geovals(ivar)%vals(nvals(ivar), nlocs))
104  self%geovals(ivar)%vals(:,:) = 0.0
105 enddo
106 self%linit = .true.
107 
108 end subroutine ufo_geovals_setup
109 
110 ! ------------------------------------------------------------------------------
111 !> Deprecated, use ufo_geovals_setup instead.
112 !> Partially initializes \p self GeoVaLs with \p nlocs number of locations
113 !> \p vars variables. Does not allocate geovals(i)%vals.
114 subroutine ufo_geovals_partial_setup(self, vars, nlocs)
115 use oops_variables_mod
116 implicit none
117 type(ufo_geovals), intent(inout) :: self
118 type(oops_variables), intent(in) :: vars
119 integer, intent(in) :: nlocs
120 
121 integer :: ivar
122 
123 call ufo_geovals_delete(self)
124 self%nlocs = nlocs
125 self%missing_value = missing_value(self%missing_value)
126 
127 self%nvar = vars%nvars()
128 allocate(self%geovals(self%nvar))
129 allocate(self%variables(self%nvar))
130 do ivar = 1, self%nvar
131  self%variables(ivar) = vars%variable(ivar)
132  self%geovals(ivar)%nlocs = nlocs
133  self%geovals(ivar)%nval = 0
134 enddo
135 
136 end subroutine ufo_geovals_partial_setup
137 
138 ! ------------------------------------------------------------------------------
139 !> Deprecated. Rely on ufo_geovals_setup to allocate GeoVaLs instead.
140 !> Allocates GeoVaLs for \p vars variables with \p nlevels number of levels.
141 !> If the GeoVaLs for this variable were allocated before with different size,
142 !> aborts.
143 subroutine ufo_geovals_allocate(self, vars, nlevels)
144 use oops_variables_mod
145 implicit none
146 type(ufo_geovals), intent(inout) :: self
147 type(oops_variables), intent(in) :: vars
148 integer, intent(in) :: nlevels
149 
150 integer :: ivar, ivar_gvals
151 character(max_string) :: err_msg
152 
153 do ivar = 1, vars%nvars()
154  ! find index of variable to be allocated
155  ivar_gvals = ufo_vars_getindex(self%variables, vars%variable(ivar))
156  ! abort if we are trying to allocate geovals for nonexistent variable
157  if (ivar_gvals < 0) then
158  write(err_msg,*) "ufo_geovals_allocate: ", trim(vars%variable(ivar)), " doesn't exist in geovals"
159  call abor1_ftn(err_msg)
160  endif
161  ! abort if we are trying to allocate geovals again, and with a different size
162  if (allocated(self%geovals(ivar_gvals)%vals) .and. (self%geovals(ivar_gvals)%nval /= nlevels)) then
163  write(err_msg,*) "ufo_geovals_allocate: attempting to allocate already allocated geovals for ", &
164  trim(vars%variable(ivar)), ". Previously allocated as ", self%geovals(ivar_gvals)%nval, &
165  " levels; now trying to allocate as ", nlevels, " levels."
166  call abor1_ftn(err_msg)
167  ! only allocate if not already allocated
168  elseif (.not. allocated(self%geovals(ivar_gvals)%vals)) then
169  self%geovals(ivar_gvals)%nval = nlevels
170  allocate(self%geovals(ivar_gvals)%vals(nlevels, self%nlocs))
171  endif
172 enddo
173 
174 ! check if all variables are now allocated, and set self%linit accordingly
175 self%linit = .true.
176 do ivar = 1, self%nvar
177  if (.not. allocated(self%geovals(ivar)%vals)) self%linit = .false.
178 enddo
179 
180 end subroutine ufo_geovals_allocate
181 
182 ! ------------------------------------------------------------------------------
183 
184 subroutine ufo_geovals_delete(self)
185 implicit none
186 type(ufo_geovals), intent(inout) :: self
187 
188 integer :: ivar
189 
190 if (allocated(self%geovals)) then
191  do ivar = 1, self%nvar
192  if (allocated(self%geovals(ivar)%vals)) deallocate(self%geovals(ivar)%vals)
193  enddo
194  deallocate(self%geovals)
195 endif
196 if (allocated(self%variables)) deallocate(self%variables)
197 self%nvar = 0
198 self%nlocs = 0
199 self%linit = .false.
200 
201 end subroutine ufo_geovals_delete
202 
203 ! ------------------------------------------------------------------------------
204 
205 subroutine ufo_geovals_get_var(self, varname, geoval)
206 implicit none
207 type(ufo_geovals), target, intent(in) :: self
208 character(len=*), intent(in) :: varname
209 type(ufo_geoval), pointer, intent(inout) :: geoval
210 
211 character(len=*), parameter :: myname_="ufo_geovals_get_var"
212 
213 character(max_string) :: err_msg
214 integer :: ivar, jv
215 
216 geoval => null()
217 if (.not. self%linit) then
218  !return
219 endif
220 
221 ivar = ufo_vars_getindex(self%variables, varname)
222 
223 if (ivar < 0) then
224  write(0,*)'ufo_geovals_get_var looking for ',trim(varname),' in:'
225  do jv=1,self%nvar
226  write(0,*)'ufo_geovals_get_var ',jv,trim(self%variables(jv))
227  enddo
228  write(err_msg,*) myname_, " ", trim(varname), ' doesnt exist'
229  call abor1_ftn(err_msg)
230 else
231  geoval => self%geovals(ivar)
232 endif
233 
234 end subroutine ufo_geovals_get_var
235 
236 ! ------------------------------------------------------------------------------
237 
238 subroutine ufo_geovals_zero(self)
239 implicit none
240 type(ufo_geovals), intent(inout) :: self
241 integer :: ivar
242 
243 if (.not. self%linit) then
244  call abor1_ftn("ufo_geovals_zero: geovals not initialized")
245 endif
246 do ivar = 1, self%nvar
247  self%geovals(ivar)%vals(:,:) = 0.0
248 enddo
249 
250 end subroutine ufo_geovals_zero
251 
252 ! ------------------------------------------------------------------------------
253 
254 subroutine ufo_geovals_abs(self)
255 implicit none
256 type(ufo_geovals), intent(inout) :: self
257 integer :: ivar
258 
259 if (.not. self%linit) then
260  call abor1_ftn("ufo_geovals_abs: geovals not initialized")
261 endif
262 do ivar = 1, self%nvar
263  self%geovals(ivar)%vals = abs(self%geovals(ivar)%vals)
264 enddo
265 
266 end subroutine ufo_geovals_abs
267 
268 ! ------------------------------------------------------------------------------
269 
270 subroutine ufo_geovals_rms(self,vrms)
271 implicit none
272 type(ufo_geovals), intent(in) :: self
273 real(kind_real), intent(inout) :: vrms
274 integer :: jv, jo
275 real(kind_real) :: n
276 
277 if (.not. self%linit) then
278  call abor1_ftn("ufo_geovals_rms: geovals not initialized")
279 endif
280 vrms=0.0_kind_real
281 n=0.0_kind_real
282 do jv = 1, self%nvar
283  do jo = 1, self%nlocs
284  vrms = vrms + sum(self%geovals(jv)%vals(:,jo)**2)
285  n=n+self%geovals(jv)%nval
286  enddo
287 enddo
288 
289 if ( n > 0) vrms = sqrt(vrms/n)
290 
291 end subroutine ufo_geovals_rms
292 
293 ! ------------------------------------------------------------------------------
294 
295 subroutine ufo_geovals_random(self)
296 use random_mod
297 implicit none
298 type(ufo_geovals), intent(inout) :: self
299 integer :: ivar
300 integer :: rseed = 7
301 
302 if (.not. self%linit) then
303  call abor1_ftn("ufo_geovals_random: geovals not initialized")
304 endif
305 do ivar = 1, self%nvar
306  call normal_distribution(self%geovals(ivar)%vals, 0.0_kind_real, 1.0_kind_real, rseed)
307 enddo
308 
309 end subroutine ufo_geovals_random
310 
311 ! ------------------------------------------------------------------------------
312 
313 subroutine ufo_geovals_scalmult(self, zz)
314 implicit none
315 type(ufo_geovals), intent(inout) :: self
316 real(kind_real), intent(in) :: zz
317 integer :: jv, jo, jz
318 
319 if (.not. self%linit) then
320  call abor1_ftn("ufo_geovals_scalmult: geovals not allocated")
321 endif
322 
323 do jv=1,self%nvar
324  do jo=1,self%nlocs
325  do jz = 1, self%geovals(jv)%nval
326  self%geovals(jv)%vals(jz,jo) = zz * self%geovals(jv)%vals(jz,jo)
327  enddo
328  enddo
329 enddo
330 
331 end subroutine ufo_geovals_scalmult
332 
333 ! ------------------------------------------------------------------------------
334 
335 subroutine ufo_geovals_profmult(self, nlocs, values)
336 implicit none
337 type(ufo_geovals), intent(inout) :: self
338 integer(c_int), intent(in) :: nlocs
339 real(c_float), intent(in) :: values(nlocs)
340 integer :: jv, jo
341 
342 if (.not. self%linit) then
343  call abor1_ftn("ufo_geovals_profmult: geovals not allocated")
344 endif
345 
346 do jv=1,self%nvar
347  do jo=1,self%nlocs
348  self%geovals(jv)%vals(:,jo) = values(jo) * self%geovals(jv)%vals(:,jo)
349  enddo
350 enddo
351 
352 end subroutine ufo_geovals_profmult
353 
354 ! ------------------------------------------------------------------------------
355 
356 
357 subroutine ufo_geovals_assign(self, rhs)
358 implicit none
359 type(ufo_geovals), intent(inout) :: self
360 type(ufo_geovals), intent(in) :: rhs
361 integer :: jv, jo, jz
362 integer :: iv
363 character(max_string) :: err_msg
364 
365 if (.not. self%linit) then
366  call abor1_ftn("ufo_geovals_scalmult: geovals not allocated")
367 endif
368 if (.not. rhs%linit) then
369  call abor1_ftn("ufo_geovals_scalmult: geovals not allocated")
370 endif
371 
372 if (self%nlocs /= rhs%nlocs) then
373  call abor1_ftn("ufo_geovals_assign: nlocs different between lhs and rhs")
374 endif
375 
376 do jv=1,self%nvar
377  iv = ufo_vars_getindex(rhs%variables, self%variables(jv))
378  if (iv < 0) then
379  write(err_msg,*) 'ufo_geovals_assign: var ', trim(self%variables(jv)), ' doesnt exist in rhs'
380  call abor1_ftn(trim(err_msg))
381  endif
382  if (self%geovals(jv)%nval /= rhs%geovals(iv)%nval) then
383  write(err_msg,*) 'ufo_geovals_assign: nvals for var ', trim(self%variables(jv)), ' are different in lhs and rhs'
384  call abor1_ftn(trim(err_msg))
385  endif
386  do jo=1,self%nlocs
387  do jz = 1, self%geovals(jv)%nval
388  self%geovals(jv)%vals(jz,jo) = rhs%geovals(iv)%vals(jz,jo)
389  enddo
390  enddo
391 enddo
392 
393 end subroutine ufo_geovals_assign
394 
395 ! ------------------------------------------------------------------------------
396 subroutine ufo_geovals_reorderzdir(self, varname, zdir)
397 implicit none
398 type(ufo_geovals),intent(inout) :: self
399 character(len=*), intent(in) :: varname
400 character(len=*), intent(in) :: zdir
401 
402 type(ufo_geovals) :: selfclone
403 type(ufo_geoval), pointer :: geoval
404 character(max_string) :: err_msg
405 integer:: iobs, ivar, ival, kval
406 logical :: do_flip !< .true. if all the ufo_geoval arrays inside geovals
407 
408 do_flip = .false.
409 
410 if (.not. self%linit) then
411  call abor1_ftn("ufo_geovals_reorderzdir: geovals not allocated")
412 endif
413 
414 ! Get vertical coordinate variable
415 call ufo_geovals_get_var(self, varname, geoval)
416 if (.not. associated(geoval)) then
417  write(err_msg, *) 'ufo_geovals_reorderzdir: geoval vertical coordinate variable ', trim(varname), ' doesnt exist'
418 endif
419 
420 ! Check if reorder variables is necessary based on the direction defined by zdir
421 if ((zdir == "bottom2top" .and. geoval%vals(1,1) < geoval%vals(geoval%nval,1)) .or. &
422  (zdir == "top2bottom" .and. geoval%vals(1,1) > geoval%vals(geoval%nval,1))) then
423  do_flip = .true.
424 else if (zdir /= "bottom2top" .or. zdir /= "top2bottom") then
425  write(err_msg, *) 'ufo_geovals_reorderzdir: z-coordinate direction ', trim(zdir), ' not defined'
426 else
427  return
428 endif
429 
430 call ufo_geovals_copy(self, selfclone)
431 
432 if (do_flip) then
433  do ivar = 1, self%nvar
434  do ival = 1, self%geovals(ivar)%nval
435  kval = self%geovals(ivar)%nval - ival + 1
436  self%geovals(ivar)%vals(ival,:) = selfclone%geovals(ivar)%vals(kval,:)
437  enddo
438  enddo
439 endif
440 
441 call ufo_geovals_delete(selfclone)
442 
443 end subroutine ufo_geovals_reorderzdir
444 
445 ! ------------------------------------------------------------------------------
446 !> Sum of two GeoVaLs objects
447 
448 subroutine ufo_geovals_add(self, other)
449 implicit none
450 type(ufo_geovals), intent(inout) :: self
451 type(ufo_geovals), intent(in) :: other
452 integer :: jv, jo, jz
453 integer :: iv
454 character(max_string) :: err_msg
455 
456 if (.not. self%linit) then
457  call abor1_ftn("ufo_geovals_add: geovals not allocated")
458 endif
459 if (.not. other%linit) then
460  call abor1_ftn("ufo_geovals_add: geovals not allocated")
461 endif
462 
463 if (self%nlocs /= other%nlocs) then
464  call abor1_ftn("ufo_geovals_add: nlocs different between lhs and rhs")
465 endif
466 
467 do jv=1,self%nvar
468  iv = ufo_vars_getindex(other%variables, self%variables(jv))
469  if (iv .ne. -1) then !Only add if exists in RHS
470  if (self%geovals(jv)%nval /= other%geovals(iv)%nval) then
471  write(err_msg,*) 'ufo_geovals_add: nvals for var ', trim(self%variables(jv)), ' are different in lhs and rhs'
472  call abor1_ftn(trim(err_msg))
473  endif
474  do jo=1,self%nlocs
475  do jz = 1, self%geovals(jv)%nval
476  self%geovals(jv)%vals(jz,jo) = self%geovals(jv)%vals(jz,jo) + other%geovals(iv)%vals(jz,jo)
477  enddo
478  enddo
479  endif
480 enddo
481 
482 end subroutine ufo_geovals_add
483 
484 ! ------------------------------------------------------------------------------
485 !> Difference between two GeoVaLs objects
486 
487 subroutine ufo_geovals_diff(self, other)
488 implicit none
489 type(ufo_geovals), intent(inout) :: self
490 type(ufo_geovals), intent(in) :: other
491 integer :: jv, jo, jz
492 integer :: iv
493 character(max_string) :: err_msg
494 
495 if (.not. self%linit) then
496  call abor1_ftn("ufo_geovals_diff: geovals not allocated")
497 endif
498 if (.not. other%linit) then
499  call abor1_ftn("ufo_geovals_diff: geovals not allocated")
500 endif
501 
502 if (self%nlocs /= other%nlocs) then
503  call abor1_ftn("ufo_geovals_diff: nlocs different between lhs and rhs")
504 endif
505 
506 do jv=1,self%nvar
507  iv = ufo_vars_getindex(other%variables, self%variables(jv))
508  if (iv .ne. -1) then !Only subtract if exists in RHS
509  if (self%geovals(jv)%nval /= other%geovals(iv)%nval) then
510  write(err_msg,*) 'ufo_geovals_diff: nvals for var ', trim(self%variables(jv)), ' are different in lhs and rhs'
511  call abor1_ftn(trim(err_msg))
512  endif
513  do jo=1,self%nlocs
514  do jz = 1, self%geovals(jv)%nval
515  self%geovals(jv)%vals(jz,jo) = self%geovals(jv)%vals(jz,jo) - other%geovals(iv)%vals(jz,jo)
516  enddo
517  enddo
518  endif
519 enddo
520 
521 end subroutine ufo_geovals_diff
522 
523 ! ------------------------------------------------------------------------------
524 !> Schur product of two GeoVaLs objects
525 
526 subroutine ufo_geovals_schurmult(self, other)
527 implicit none
528 type(ufo_geovals), intent(inout) :: self
529 type(ufo_geovals), intent(in) :: other
530 integer :: jv, jo, jz
531 integer :: iv
532 character(max_string) :: err_msg
533 
534 if (.not. self%linit) then
535  call abor1_ftn("ufo_geovals_schurmult: geovals not allocated")
536 endif
537 if (.not. other%linit) then
538  call abor1_ftn("ufo_geovals_schurmult: geovals not allocated")
539 endif
540 
541 if (self%nlocs /= other%nlocs) then
542  call abor1_ftn("ufo_geovals_schurmult: nlocs different between lhs and rhs")
543 endif
544 
545 do jv=1,self%nvar
546  iv = ufo_vars_getindex(other%variables, self%variables(jv))
547  if (iv .ne. -1) then !Only mult if exists in RHS
548  if (self%geovals(jv)%nval /= other%geovals(iv)%nval) then
549  write(err_msg,*) 'ufo_geovals_schurmult: nvals for var ', trim(self%variables(jv)), ' are different in lhs and rhs'
550  call abor1_ftn(trim(err_msg))
551  endif
552  do jo=1,self%nlocs
553  do jz = 1, self%geovals(jv)%nval
554  self%geovals(jv)%vals(jz,jo) = self%geovals(jv)%vals(jz,jo) * other%geovals(iv)%vals(jz,jo)
555  enddo
556  enddo
557  endif
558 enddo
559 
560 end subroutine ufo_geovals_schurmult
561 
562 ! ------------------------------------------------------------------------------
563 !> Copy one GeoVaLs object into another
564 !!
565 
566 subroutine ufo_geovals_copy(self, other)
567 implicit none
568 type(ufo_geovals), intent(in) :: self
569 type(ufo_geovals), intent(inout) :: other
570 integer :: jv
571 
572 if (.not. self%linit) then
573  call abor1_ftn("ufo_geovals_copy: geovals not defined")
574 endif
575 
576 call ufo_geovals_delete(other)
577 
578 other%nlocs = self%nlocs
579 other%nvar = self%nvar
580 allocate(other%variables(other%nvar))
581 other%variables(:) = self%variables(:)
582 
583 allocate(other%geovals(other%nvar))
584 do jv = 1, other%nvar
585  other%geovals(jv)%nval = self%geovals(jv)%nval
586  other%geovals(jv)%nlocs = self%geovals(jv)%nlocs
587  allocate(other%geovals(jv)%vals(other%geovals(jv)%nval, other%geovals(jv)%nlocs))
588  other%geovals(jv)%vals(:,:) = self%geovals(jv)%vals(:,:)
589 enddo
590 
591 other%missing_value = self%missing_value
592 other%linit = .true.
593 
594 end subroutine ufo_geovals_copy
595 
596 ! ------------------------------------------------------------------------------
597 !> Copy one location from GeoVaLs into a new object
598 !!
599 
600 subroutine ufo_geovals_copy_one(self, other, loc_index)
601 implicit none
602 type(ufo_geovals), intent(inout) :: self !> GeoVaLs for one location
603 type(ufo_geovals), intent(in) :: other !> GeoVaLs for many location
604 integer, intent(in) :: loc_index !> Index of the location in the "other" geoval
605 integer :: jv
606 
607 if (.not. other%linit) then
608  call abor1_ftn("ufo_geovals_copy_one: geovals not defined")
609 endif
610 
611 call ufo_geovals_delete(self)
612 
613 self%nlocs = 1
614 self%nvar = other%nvar
615 allocate(self%variables(self%nvar))
616 self%variables(:) = other%variables(:)
617 
618 allocate(self%geovals(self%nvar))
619 do jv = 1, self%nvar
620  self%geovals(jv)%nval = other%geovals(jv)%nval
621  self%geovals(jv)%nlocs = 1
622  allocate(self%geovals(jv)%vals(self%geovals(jv)%nval, self%geovals(jv)%nlocs))
623  self%geovals(jv)%vals(:,self%nlocs) = other%geovals(jv)%vals(:,loc_index)
624 enddo
625 
626 self%missing_value = other%missing_value
627 self%linit = .true.
628 
629 end subroutine ufo_geovals_copy_one
630 
631 ! ------------------------------------------------------------------------------
632 !> Initialize a GeoVaLs object based on an analytic state
633 !!
634 !! \details **ufo_geovals_analytic_init_c()** takes an existing ufo::GeoVaLs object
635 !! and fills in values based on one of several analytic solutions. This initialization
636 !! is intended to be used with the **TestStateInterpolation()** test; see there for
637 !! further information.
638 !!
639 !! Currently implemented options for analytic_init include:
640 !! * dcmip-test-1-1: 3D deformational flow
641 !! * dcmip-test-1-2: 3D Hadley-like meridional circulation
642 !! * dcmip-test-3-1: Non-orographic gravity waves on a small planet
643 !! * dcmip-test-4-0: Baroclinic instability
644 !!
645 !! \warning Currently only temperature is implemented. For variables other than
646 !! temperature, the input GeoVaLs object is not changed. This effectively
647 !! disables the interpolation test for that variable by setting the normalized
648 !! error to zero.
649 !!
650 !! \warning Currently there is no conversion between temperature and virtual
651 !! temperature
652 !!
653 !! \date May, 2018: Created by M. Miesch (JCSDA)
654 !! \date June, 2018: Added dcmip-test-4.0 (M. Miesch, JCSDA)
655 !!
656 !! \sa test::TestStateInterpolation()
657 !!
658 
659 subroutine ufo_geovals_analytic_init(self, locs, ic)
660 use dcmip_initial_conditions_test_1_2_3, only : test1_advection_deformation, &
661  test1_advection_hadley, test3_gravity_wave
662 use dcmip_initial_conditions_test_4, only : test4_baroclinic_wave
664 use ufo_utils_mod, only: cmp_strings
665 
666 implicit none
667 type(ufo_geovals), intent(inout) :: self
668 type(ufo_locations), intent(in) :: locs
669 character(*), intent(in) :: ic
670 
671 real(kind_real) :: pi = acos(-1.0_kind_real)
672 real(kind_real) :: deg_to_rad,rlat, rlon
673 real(kind_real) :: p0, kz, u0, v0, w0, t0, phis0, ps0, rho0, hum0
674 real(kind_real) :: q1, q2, q3, q4
675 real(kind_real), allocatable, dimension(:) :: lons, lats
676 integer :: nlocs, ivar, iloc, ival
677 
678 if (.not. self%linit) then
679  call abor1_ftn("ufo_geovals_analytic_init: geovals not defined")
680 endif
681 
682 ! The last variable should be the ln pressure coordinate. That's
683 ! where we get the height information for the analytic init
684 if (self%variables(self%nvar) /= var_prs .and. &
685  self%variables(self%nvar) /= var_prsi) then
686  call abor1_ftn("ufo_geovals_analytic_init: pressure coordinate not defined")
687 endif
688 
689 deg_to_rad = pi/180.0_kind_real
690 
691 nlocs = locs%nlocs()
692 allocate(lons(nlocs), lats(nlocs))
693 call locs%get_lons(lons)
694 call locs%get_lats(lats)
695 
696 do ivar = 1, self%nvar-1
697 
698  do iloc = 1, self%geovals(ivar)%nlocs
699 
700  ! convert lat and lon to radians
701  rlat = deg_to_rad * lats(iloc)
702  rlon = deg_to_rad*modulo(lons(iloc)+180.0_kind_real,360.0_kind_real) - pi
703 
704  do ival = 1, self%geovals(ivar)%nval
705 
706  ! obtain height from the existing GeoVaLs object, which should be an
707  ! output of the State::getValues() method
708  ! should be delivered in units of Pa
709  p0 = self%geovals(self%nvar)%vals(ival,iloc)
710 
711  init_option: select case (trim(ic))
712 
713  case ("invent_state")
714 
715  t0 = cos(deg_to_rad * lons(iloc) ) * cos(rlat)
716 
717  case ("dcmip-test-1-1")
718 
719  call test1_advection_deformation(rlon,rlat,p0,kz,0,u0,v0,w0,&
720  t0,phis0,ps0,rho0,hum0,q1,q2,q3,q4)
721 
722  case ("dcmip-test-1-2")
723 
724  call test1_advection_hadley(rlon,rlat,p0,kz,0,u0,v0,w0,&
725  t0,phis0,ps0,rho0,hum0,q1)
726 
727  case ("dcmip-test-3-1")
728 
729  call test3_gravity_wave(rlon,rlat,p0,kz,0,u0,v0,w0,&
730  t0,phis0,ps0,rho0,hum0)
731 
732  case ("dcmip-test-4-0")
733 
734  call test4_baroclinic_wave(0,1.0_kind_real,rlon,rlat,p0,kz,0,u0,v0,w0,&
735  t0,phis0,ps0,rho0,hum0,q1,q2)
736 
737  case default
738 
739  call abor1_ftn("ufo_geovals_analytic_init: invalid analytic_init")
740 
741  end select init_option
742 
743  ! currently only temperature is implemented
744  if (cmp_strings(self%variables(ivar), var_tv)) then
745  ! Warning: we may need a conversion from temperature to
746  ! virtual temperture here
747  self%geovals(ivar)%vals(ival,iloc) = t0
748  endif
749 
750  enddo
751  enddo
752 enddo
753 
754 deallocate(lons, lats)
755 
756 end subroutine ufo_geovals_analytic_init
757 
758 ! ------------------------------------------------------------------------------
759 !> Normalization of one GeoVaLs object by another
760 !!
761 !! \details This is a normalization operator that first computes the normalization
762 !! factor for each variable based on the rms amplitude of that variable across
763 !! all locations in the reference GeoVaLs object (other). Then each element of
764 !! the input GeoVals object (self) is divided by these normalization factors.
765 !! The operation is done in place. So, after execution, the input GeoVaLs
766 !! object will be nondimensional.
767 !!
768 !! \warning If the reference variable is identially zero across all
769 !! locations, then the result of this operatution is set to zero for that
770 !! variable. This is to used to bypass variables that do not have a reference
771 !! value in the State interpolation test.
772 !!
773 
774 subroutine ufo_geovals_normalize(self, other)
775 implicit none
776 type(ufo_geovals), intent(inout) :: self !> Input GeoVaLs object (LHS)
777 type(ufo_geovals), intent(in) :: other !> Reference GeoVaLs object (RHS)
778 integer :: jv, jo, jz
779 real(kind_real) :: over_nloc, vrms, norm
780 
781 if (.not. self%linit) then
782  call abor1_ftn("ufo_geovals_normalize: geovals not allocated")
783 endif
784 if (.not. other%linit) then
785  call abor1_ftn("ufo_geovals_normalize: geovals not allocated")
786 endif
787 if (self%nvar /= other%nvar) then
788  call abor1_ftn("ufo_geovals_normalize: reference geovals object must have the same variables as the original")
789 endif
790 
791 
792 do jv=1,self%nvar
793 
794  !> Compute normalization factors for the errors based on the rms amplitude of
795  !! each variable across all of the selected locations. Use the "other" GeoVaLs
796  !! object as a reference, since this may be the exact analytic answer
797 
798  over_nloc = 1.0_kind_real / &
799  (real(other%nlocs,kind_real)*real(other%geovals(jv)%nval,kind_real))
800 
801  vrms = 0.0_kind_real
802  do jo = 1, other%nlocs
803  do jz = 1, other%geovals(jv)%nval
804  vrms = vrms + other%geovals(jv)%vals(jz,jo)**2
805  enddo
806  enddo
807 
808  if (vrms > 0.0_kind_real) then
809  norm = 1.0_kind_real / sqrt(vrms*over_nloc)
810  else
811  norm = 0.0_kind_real
812  endif
813 
814  ! Now loop through the LHS locations to compute the normalized value
815  do jo=1,self%nlocs
816  do jz = 1, self%geovals(jv)%nval
817  self%geovals(jv)%vals(jz,jo) = norm*self%geovals(jv)%vals(jz,jo)
818  enddo
819  enddo
820 enddo
821 
822 end subroutine ufo_geovals_normalize
823 
824 ! ------------------------------------------------------------------------------
825 
826 subroutine ufo_geovals_reset_sec_arg(self, other, nlocs)
827 implicit none
828 type(ufo_geovals), intent(in) :: self
829 type(ufo_geovals), intent(inout) :: other
830 integer, intent(in) :: nlocs
831 integer :: ivar
832 
833 if (other%linit) call abor1_ftn("ufo_geovals_reset_sec_arg: other already have data")
834 
835 other%nlocs = nlocs
836 other%nvar = self%nvar
837 other%missing_value = self%missing_value
838 allocate(other%variables(self%nvar))
839 allocate(other%geovals(self%nvar))
840 do ivar = 1, self%nvar
841  other%variables(ivar) = self%variables(ivar)
842  other%geovals(ivar)%nlocs = nlocs
843  other%geovals(ivar)%nval = self%geovals(ivar)%nval
844  allocate(other%geovals(ivar)%vals(self%geovals(ivar)%nval, nlocs))
845  other%geovals(ivar)%vals(:,:) = 0.0
846 enddo
847 other%linit = .false.
848 
849 end subroutine ufo_geovals_reset_sec_arg
850 ! ------------------------------------------------------------------------------
851 
852 subroutine ufo_geovals_split(self, other1, other2)
853 implicit none
854 type(ufo_geovals), intent(in) :: self
855 type(ufo_geovals), intent(inout) :: other1
856 type(ufo_geovals), intent(inout) :: other2
857 
858 integer :: ivar, iobs
859 
860 if (.not. self%linit) &
861  call abor1_ftn("ufo_geovals_split: geovals self is not allocated or has no data")
862 
863 if (other1%linit .or. other2%linit) &
864  call abor1_ftn("ufo_geovals_split: geovals other1 or other2 already have data")
865 
866 call ufo_geovals_delete(other1)
867 call ufo_geovals_delete(other2)
868 call ufo_geovals_reset_sec_arg(self, other1, self%nlocs/2)
869 call ufo_geovals_reset_sec_arg(self, other2, self%nlocs - self%nlocs/2)
870 
871 do ivar = 1, self%nvar
872  do iobs = 1, self%nlocs/2
873  other1%geovals(ivar)%vals(:,iobs) = self%geovals(ivar)%vals(:,iobs)
874  enddo
875  do iobs = self%nlocs/2 + 1, self%nlocs
876  other2%geovals(ivar)%vals(:,iobs - self%nlocs/2) = self%geovals(ivar)%vals(:,iobs)
877  enddo
878 enddo
879 other1%linit = .true.
880 other2%linit = .true.
881 
882 end subroutine ufo_geovals_split
883 ! ------------------------------------------------------------------------------
884 
885 subroutine ufo_geovals_merge(self, other1, other2)
886 implicit none
887 type(ufo_geovals), intent(inout) :: self
888 type(ufo_geovals), intent(in) :: other1
889 type(ufo_geovals), intent(in) :: other2
890 
891 integer :: ivar, iobs
892 
893 if ((.not. other1%linit) .or. (.not. other2%linit)) &
894  call abor1_ftn("ufo_geovals_merge: geovals other1 or other2 is not allocated or has no data")
895 
896 call ufo_geovals_delete(self)
897 call ufo_geovals_reset_sec_arg(other1, self, other1%nlocs + other2%nlocs)
898 
899 do ivar = 1, self%nvar
900  do iobs = 1, other1%nlocs
901  self%geovals(ivar)%vals(:,iobs) = other1%geovals(ivar)%vals(:,iobs)
902  enddo
903  do iobs = other1%nlocs + 1, self%nlocs
904  self%geovals(ivar)%vals(:,iobs) = &
905  other2%geovals(ivar)%vals(:,iobs - other1%nlocs)
906  enddo
907 enddo
908 self%linit = .true.
909 
910 end subroutine ufo_geovals_merge
911 ! ------------------------------------------------------------------------------
912 
913 subroutine ufo_geovals_minmaxavg(self, kobs, kvar, pmin, pmax, prms)
914 implicit none
915 integer, intent(inout) :: kobs
916 integer, intent(in) :: kvar
917 real(kind_real), intent(inout) :: pmin, pmax, prms
918 type(ufo_geovals), intent(in) :: self
919 integer :: jo, jz, jv
920 
921 jv = kvar+1
922 kobs = 0
923 pmin = huge(pmin)
924 pmax = -huge(pmax)
925 prms = 0.0_kind_real
926 do jo = 1, self%nlocs
927  do jz = 1, self%geovals(jv)%nval
928  if (self%geovals(jv)%vals(jz,jo) .ne. self%missing_value) then
929  kobs = kobs + 1
930  if (self%geovals(jv)%vals(jz,jo) < pmin) pmin = self%geovals(jv)%vals(jz,jo)
931  if (self%geovals(jv)%vals(jz,jo) > pmax) pmax = self%geovals(jv)%vals(jz,jo)
932  prms = prms + self%geovals(jv)%vals(jz,jo) * self%geovals(jv)%vals(jz,jo)
933  endif
934  enddo
935 enddo
936 if (kobs > 0) prms = sqrt(prms/real(kobs,kind_real))
937 
938 end subroutine ufo_geovals_minmaxavg
939 
940 ! ------------------------------------------------------------------------------
941 !> Location where the summed geovals value is maximum
942 !!
943 !! \details This routine computes the rms value over the vertical profile for
944 !! each location and observation then returns the location number and the
945 !! variable number where this rms value is maximum. Intended for use with
946 !! the State interpotation test in which the input GeoVaLs object is a
947 !! nondimensional, positive-definite error measurement.
948 
949 subroutine ufo_geovals_maxloc(self, mxval, iobs, ivar)
950 implicit none
951 real(kind_real), intent(inout) :: mxval
952 integer, intent(inout) :: iobs, ivar
953 
954 type(ufo_geovals), intent(in) :: self
955 real(kind_real) :: vrms
956 integer :: jv, jo, jz
957 
958 if (.not. self%linit) then
959  call abor1_ftn("ufo_geovals_maxloc: geovals not allocated")
960 endif
961 
962 mxval = 0.0_kind_real
963 iobs = 1
964 ivar = 1
965 
966 do jv = 1,self%nvar
967  do jo = 1, self%nlocs
968 
969  vrms = 0.0_kind_real
970  do jz = 1, self%geovals(jv)%nval
971  vrms = vrms + self%geovals(jv)%vals(jz,jo)**2
972  enddo
973 
974  if ( self%geovals(jv)%nval > 0 ) then
975  vrms = sqrt(vrms/real(self%geovals(jv)%nval,kind_real))
976  end if
977 
978  if (vrms > mxval) then
979  mxval = vrms
980  iobs = jo
981  ivar = jv
982  endif
983 
984  enddo
985 enddo
986 
987 end subroutine ufo_geovals_maxloc
988 
989 ! ------------------------------------------------------------------------------
990 
991 subroutine ufo_geovals_read_netcdf(self, filename, loc_multiplier, c_obspace, vars)
992 use netcdf
993 use oops_variables_mod
994 implicit none
995 type(ufo_geovals), intent(inout) :: self
996 character(max_string), intent(in) :: filename
997 integer, intent(in) :: loc_multiplier
998 type(c_ptr), intent(in) :: c_obspace
999 type(oops_variables), intent(in) :: vars
1000 
1001 integer :: nlocs, gv_all_nlocs, nlocs_var
1002 integer :: nval
1003 integer :: obs_nlocs
1004 integer :: obs_all_nlocs
1005 integer :: iloc
1006 integer :: jloc, jloc_start, jloc_end
1007 integer :: iloc_new
1008 
1009 integer :: ncid, dimid, varid, vartype, ndims
1010 integer, dimension(3) :: dimids
1011 integer :: ivar
1012 integer :: ierr
1013 
1014 character(max_string) :: err_msg
1015 character(len=30) :: obs_nlocs_str
1016 character(len=30) :: geo_nlocs_str
1017 
1018 integer(c_size_t), allocatable, dimension(:) :: dist_indx
1019 integer(c_size_t), allocatable, dimension(:) :: obs_dist_indx
1020 
1021 real, allocatable :: field2d(:,:), field1d(:)
1022 
1023 ! open netcdf file
1024 call check('nf90_open', nf90_open(trim(filename),nf90_nowrite,ncid))
1025 
1026 ! find how many locs are in the file
1027 ierr = nf90_inq_dimid(ncid, "nlocs", dimid)
1028 if(ierr /= nf90_noerr) then
1029  write(err_msg,*) "Error: Dimension nlocs not found in ", trim(filename)
1030  call abor1_ftn(err_msg)
1031 endif
1032 call check('nf90_inq_dimid', nf90_inq_dimid(ncid, "nlocs", dimid))
1033 call check('nf90_inquire_dimension', nf90_inquire_dimension(ncid, dimid, len = gv_all_nlocs))
1034 
1035 !> round-robin distribute the observations to PEs
1036 !> Calculate how many obs. on each PE
1037 obs_all_nlocs = obsspace_get_gnlocs(c_obspace)
1038 obs_nlocs = obsspace_get_nlocs(c_obspace)
1039 allocate(obs_dist_indx(obs_nlocs))
1040 call obsspace_get_index(c_obspace, obs_dist_indx)
1041 
1042 ! loc_multiplier specifies how many locations in the geovals file per
1043 ! single location in the obs file. There needs to be at least
1044 ! loc_multiplier * obs_all_nlocs locations in the geovals file.
1045 
1046 if (gv_all_nlocs .lt. (loc_multiplier * obs_all_nlocs)) then
1047  write(obs_nlocs_str, *) loc_multiplier * obs_all_nlocs
1048  write(geo_nlocs_str, *) gv_all_nlocs
1049  write(err_msg,'(7a)') &
1050  "Error: Number of locations in the geovals file (", &
1051  trim(adjustl(geo_nlocs_str)), ") must be greater than or equal to ", &
1052  "the product of loc_multiplier and number of locations in the ", &
1053  "obs file (", trim(adjustl(obs_nlocs_str)), ")"
1054  call abor1_ftn(err_msg)
1055 endif
1056 
1057 ! We have enough locations in the geovals file to cover the span of the
1058 ! number of locations in the obs file. Generate the dist_indx according
1059 ! to the loc_multiplier and obs_nlocs values.
1060 if (loc_multiplier >= 0) then
1061  nlocs = loc_multiplier * obs_nlocs
1062  allocate(dist_indx(nlocs))
1063  iloc_new = 1
1064  do iloc = 1,obs_nlocs
1065  jloc_start = ((obs_dist_indx(iloc) - 1) * loc_multiplier) + 1
1066  jloc_end = obs_dist_indx(iloc) * loc_multiplier
1067  do jloc = jloc_start, jloc_end
1068  dist_indx(iloc_new) = jloc
1069  iloc_new = iloc_new + 1
1070  enddo
1071  enddo
1072 else
1073  nlocs = - loc_multiplier * obs_nlocs
1074  allocate(dist_indx(nlocs))
1075  iloc_new = 1
1076  do jloc = 1, - loc_multiplier
1077  do iloc = 1, obs_nlocs
1078  dist_indx(iloc_new) = obs_dist_indx(iloc) + (jloc - 1) * obs_all_nlocs
1079  iloc_new = iloc_new + 1
1080  enddo
1081  enddo
1082 end if
1083 
1084 ! allocate geovals structure
1085 call ufo_geovals_partial_setup(self, vars, nlocs)
1086 
1087 do ivar = 1, self%nvar
1088 
1089  ierr = nf90_inq_varid(ncid, self%variables(ivar), varid)
1090  if(ierr /= nf90_noerr) then
1091  write(err_msg,*) "Error: Variable ", trim(self%variables(ivar)), " not found in ", trim(filename)
1092  call abor1_ftn(err_msg)
1093  endif
1094 
1095  call check('nf90_inquire_variable', nf90_inquire_variable(ncid, varid, xtype = vartype, &
1096  ndims = ndims, dimids = dimids))
1097  !> read 1d variable
1098  if (ndims == 1) then
1099  call check('nf90_inquire_dimension', nf90_inquire_dimension(ncid, dimids(1), len = nlocs_var))
1100  if (nlocs_var /= gv_all_nlocs) then
1101  call abor1_ftn('ufo_geovals_read_netcdf: var dim /= gv_all_nlocs')
1102  endif
1103  nval = 1
1104  !> allocate geoval for this variable
1105  self%geovals(ivar)%nval = nval
1106  allocate(self%geovals(ivar)%vals(nval,nlocs))
1107 
1108  allocate(field1d(nlocs_var))
1109  call check('nf90_get_var', nf90_get_var(ncid, varid, field1d))
1110  self%geovals(ivar)%vals(1,:) = field1d(dist_indx)
1111  deallocate(field1d)
1112  !> read 2d variable
1113  elseif (ndims == 2) then
1114  call check('nf90_inquire_dimension', nf90_inquire_dimension(ncid, dimids(1), len = nval))
1115  call check('nf90_inquire_dimension', nf90_inquire_dimension(ncid, dimids(2), len = nlocs_var))
1116  if (nlocs_var /= gv_all_nlocs) then
1117  call abor1_ftn('ufo_geovals_read_netcdf: var dim /= gv_all_nlocs')
1118  endif
1119  !> allocate geoval for this variable
1120  self%geovals(ivar)%nval = nval
1121  allocate(self%geovals(ivar)%vals(nval,nlocs))
1122  allocate(field2d(nval, nlocs_var))
1123  call check('nf90_get_var', nf90_get_var(ncid, varid, field2d))
1124  self%geovals(ivar)%vals(:,:) = field2d(:,dist_indx)
1125  deallocate(field2d)
1126  !> only 1d & 2d vars
1127  else
1128  call abor1_ftn('ufo_geovals_read_netcdf: can only read 1d and 2d fields')
1129  endif
1130 
1131  ! set the missing value equal to IODA missing_value
1132  where (self%geovals(ivar)%vals > 1.0e08) self%geovals(ivar)%vals = self%missing_value
1133 
1134 enddo
1135 
1136 if (allocated(dist_indx)) deallocate(dist_indx)
1137 if (allocated(obs_dist_indx)) deallocate(obs_dist_indx)
1138 
1139 self%linit = .true.
1140 
1141 call check('nf90_close', nf90_close(ncid))
1142 
1143 end subroutine ufo_geovals_read_netcdf
1144 
1145 ! ------------------------------------------------------------------------------
1146 subroutine ufo_geovals_write_netcdf(self, filename)
1147 use netcdf
1148 implicit none
1149 type(ufo_geovals), intent(inout) :: self
1150 character(max_string), intent(in) :: filename
1151 
1152 integer :: i
1153 integer :: ncid, dimid_nlocs, dimid_nval, dims(2)
1154 integer, allocatable :: ncid_var(:)
1155 
1156 allocate(ncid_var(self%nvar))
1157 
1158 call check('nf90_create', nf90_create(trim(filename),nf90_hdf5,ncid))
1159 call check('nf90_def_dim', nf90_def_dim(ncid,'nlocs',self%nlocs, dimid_nlocs))
1160 dims(2) = dimid_nlocs
1161 
1162 do i = 1, self%nvar
1163  call check('nf90_def_dim', &
1164  nf90_def_dim(ncid,trim(self%variables(i))//"_nval",self%geovals(i)%nval, dimid_nval))
1165  dims(1) = dimid_nval
1166  call check('nf90_def_var', &
1167  nf90_def_var(ncid,trim(self%variables(i)),nf90_float,dims,ncid_var(i)))
1168 enddo
1169 
1170 call check('nf90_enddef', nf90_enddef(ncid))
1171 
1172 do i = 1, self%nvar
1173  call check('nf90_put_var', nf90_put_var(ncid,ncid_var(i),self%geovals(i)%vals(:,:)))
1174 enddo
1175 
1176 call check('nf90_close', nf90_close(ncid))
1177 deallocate(ncid_var)
1178 
1179 end subroutine ufo_geovals_write_netcdf
1180 
1181 ! ------------------------------------------------------------------------------
1182 subroutine check(action, status)
1183 use netcdf, only: nf90_noerr, nf90_strerror
1184 implicit none
1185 
1186 integer, intent (in) :: status
1187 character (len=*), intent (in) :: action
1188 character(max_string) :: err_msg
1189 
1190 if(status /= nf90_noerr) then
1191  write(err_msg,*) "During action: ", trim(action), ", received error: ", trim(nf90_strerror(status))
1192  call abor1_ftn(err_msg)
1193 end if
1194 
1195 end subroutine check
1196 
1197 ! ------------------------------------------------------------------------------
1198 
1199 subroutine ufo_geovals_print(self, iobs)
1200 implicit none
1201 type(ufo_geovals), intent(in) :: self
1202 integer, intent(in) :: iobs
1203 
1204 type(ufo_geoval), pointer :: geoval
1205 character(MAXVARLEN) :: varname
1206 integer :: ivar
1207 
1208 do ivar = 1, self%nvar
1209  varname = self%variables(ivar)
1210  call ufo_geovals_get_var(self, varname, geoval)
1211  if (associated(geoval)) then
1212  print *, 'geoval test: ', trim(varname), geoval%nval, geoval%vals(:,iobs)
1213  else
1214  print *, 'geoval test: ', trim(varname), ' doesnt exist'
1215  endif
1216 enddo
1217 
1218 end subroutine ufo_geovals_print
1219 
1220 ! ------------------------------------------------------------------------------
1221 
1222 end module ufo_geovals_mod
integer, parameter max_string
subroutine, public ufo_geovals_reorderzdir(self, varname, zdir)
subroutine check(action, status)
subroutine, public ufo_geovals_maxloc(self, mxval, iobs, ivar)
Location where the summed geovals value is maximum.
subroutine, public ufo_geovals_default_constr(self)
subroutine, public ufo_geovals_allocate(self, vars, nlevels)
Deprecated. Rely on ufo_geovals_setup to allocate GeoVaLs instead. Allocates GeoVaLs for vars variabl...
subroutine, public ufo_geovals_random(self)
subroutine, public ufo_geovals_setup(self, vars, nlocs, nvars, nvals)
Initializes and allocates self GeoVaLs with nlocs number of locations for vars variables....
subroutine, public ufo_geovals_write_netcdf(self, filename)
subroutine, public ufo_geovals_normalize(self, other)
Normalization of one GeoVaLs object by another.
subroutine, private ufo_geovals_reset_sec_arg(self, other, nlocs)
subroutine, public ufo_geovals_split(self, other1, other2)
subroutine, public ufo_geovals_read_netcdf(self, filename, loc_multiplier, c_obspace, vars)
subroutine, public ufo_geovals_rms(self, vrms)
subroutine, public ufo_geovals_copy(self, other)
Copy one GeoVaLs object into another.
subroutine, public ufo_geovals_delete(self)
subroutine, public ufo_geovals_analytic_init(self, locs, ic)
Initialize a GeoVaLs object based on an analytic state.
subroutine, public ufo_geovals_diff(self, other)
Difference between two GeoVaLs objects.
subroutine, public ufo_geovals_profmult(self, nlocs, values)
subroutine, public ufo_geovals_copy_one(self, other, loc_index)
Copy one location from GeoVaLs into a new object.
subroutine, public ufo_geovals_minmaxavg(self, kobs, kvar, pmin, pmax, prms)
subroutine, public ufo_geovals_scalmult(self, zz)
subroutine, public ufo_geovals_assign(self, rhs)
subroutine, public ufo_geovals_add(self, other)
Sum of two GeoVaLs objects.
subroutine, public ufo_geovals_print(self, iobs)
subroutine, public ufo_geovals_get_var(self, varname, geoval)
subroutine, public ufo_geovals_partial_setup(self, vars, nlocs)
Deprecated, use ufo_geovals_setup instead. Partially initializes self GeoVaLs with nlocs number of lo...
subroutine, public ufo_geovals_zero(self)
subroutine, public ufo_geovals_abs(self)
subroutine, public ufo_geovals_schurmult(self, other)
Schur product of two GeoVaLs objects.
subroutine, public ufo_geovals_merge(self, other1, other2)
Fortran interface to ufo::Locations.
integer function nlocs(this)
Return the number of observational locations in this Locations object.
Fortran module with various useful routines.
logical function, public cmp_strings(str1, str2)
character(len=maxvarlen), parameter, public var_prsi
integer function, public ufo_vars_getindex(vars, varname)
character(len=maxvarlen), parameter, public var_prs
character(len=maxvarlen), parameter, public var_tv
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators