UFO
ufo_adt_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 !> Fortran adt module for observation operator
7 
8 module ufo_adt_mod
9 
10  use fckit_configuration_module, only: fckit_configuration
11  use iso_c_binding
12  use kinds
13 
15  use ufo_basis_mod, only: ufo_basis
16  use ufo_vars_mod
17  use obsspace_mod
18  use missing_values_mod
19 
20  implicit none
21  private
22 
23  integer, parameter :: max_string=800
24 
25 !> Fortran derived type for the observation type
26  type, extends(ufo_basis), public :: ufo_adt
27  private
28  contains
29  procedure :: setup => ufo_adt_setup
30  procedure :: delete => ufo_adt_delete
31  procedure :: simobs => ufo_adt_simobs
32  end type ufo_adt
33 
34 contains
35 
36 ! ------------------------------------------------------------------------------
37 subroutine ufo_adt_setup(self, f_conf)
38 implicit none
39 class(ufo_adt), intent(inout) :: self
40 type(fckit_configuration), intent(in) :: f_conf
41 
42 end subroutine ufo_adt_setup
43 
44 ! ------------------------------------------------------------------------------
45 subroutine ufo_adt_delete(self)
46 implicit none
47 class(ufo_adt), intent(inout) :: self
48 
49 end subroutine ufo_adt_delete
50 
51 ! ------------------------------------------------------------------------------
52 subroutine ufo_adt_simobs(self, geovals, hofx, obss)
53 use fckit_mpi_module, only: fckit_mpi_comm, fckit_mpi_sum
54 implicit none
55  class(ufo_adt), intent(in) :: self
56  type(ufo_geovals), intent(in) :: geovals
57  real(c_double), intent(inout) :: hofx(:)
58  type(c_ptr), value, intent(in) :: obss
59 
60  character(len=*), parameter :: myname_="ufo_adt_simobs"
61  character(max_string) :: err_msg
62  type(ufo_geoval), pointer :: geoval_adt
63  real(kind_real), allocatable :: obs_adt(:)
64  integer :: obss_nlocs
65  integer :: iobs, cnt, cnt_glb
66  real(kind_real) :: offset_hofx, pe_offset_hofx
67  real(kind_real) :: offset_obs, pe_offset_obs
68  type(fckit_mpi_comm) :: f_comm
69  real(c_double) :: missing
70 
71  call obsspace_get_comm(obss, f_comm)
72 
73  ! Set missing flag
74  missing = missing_value(missing)
75 
76  ! check if nlocs is consistent in geovals & hofx
77  obss_nlocs = obsspace_get_nlocs(obss)
78 
79  !nlocs = size(hofx,1)
80  if (geovals%nlocs /= size(hofx,1)) then
81  write(err_msg,*) myname_, ' error: nlocs inconsistent!'
82  call abor1_ftn(err_msg)
83  endif
84 
85  ! check if adt variable is in geovals and get it
86  call ufo_geovals_get_var(geovals, var_abs_topo, geoval_adt)
87 
88  ! Read in obs data
89  allocate(obs_adt(obss_nlocs))
90 
91  call obsspace_get_db(obss, "ObsValue", "obs_absolute_dynamic_topography", obs_adt)
92 
93  ! Local offset
94  pe_offset_hofx = 0.0
95  pe_offset_obs = 0.0
96  cnt = 0
97  do iobs = 1, obss_nlocs
98  if (hofx(iobs)/=missing) then
99  pe_offset_hofx = pe_offset_hofx + geoval_adt%vals(1,iobs)
100  pe_offset_obs = pe_offset_obs + obs_adt(iobs)
101  cnt = cnt + 1
102  end if
103  end do
104 
105  ! Global offsets
106  call f_comm%allreduce(pe_offset_hofx, offset_hofx, fckit_mpi_sum())
107  call f_comm%allreduce(pe_offset_obs, offset_obs, fckit_mpi_sum())
108  call f_comm%allreduce(cnt, cnt_glb, fckit_mpi_sum())
109  offset_hofx = offset_hofx/cnt_glb
110  offset_obs = offset_obs/cnt_glb
111 
112  ! Adjust simulated obs to obs offset
113  do iobs = 1, obss_nlocs
114  hofx(iobs) = geoval_adt%vals(1,iobs) + (offset_obs-offset_hofx)
115  enddo
116 
117  deallocate(obs_adt)
118 
119  end subroutine ufo_adt_simobs
120 
121 end module ufo_adt_mod
Fortran adt module for observation operator.
Definition: ufo_adt_mod.F90:8
subroutine ufo_adt_setup(self, f_conf)
Definition: ufo_adt_mod.F90:38
subroutine ufo_adt_delete(self)
Definition: ufo_adt_mod.F90:46
subroutine ufo_adt_simobs(self, geovals, hofx, obss)
Definition: ufo_adt_mod.F90:53
integer, parameter max_string
Definition: ufo_adt_mod.F90:23
subroutine, public ufo_geovals_get_var(self, varname, geoval)
character(len=maxvarlen), public var_abs_topo
Fortran derived type for the observation type.
Definition: ufo_adt_mod.F90:26
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators