FV3 Bundle
ufo_adt_tlad_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017 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 module to handle adt observations
7 
9 
10 use ufo_vars_mod
11 use ioda_locs_mod
13 use kinds
14 use iso_c_binding
15 use obsspace_mod
16 
17 implicit none
18 public :: ufo_adt_tlad
19 public :: ufo_adt_tlad_delete
20 public :: ufo_adt_tlad_settraj
21 public :: ufo_adt_simobs_tl
22 public :: ufo_adt_simobs_ad
23 private
24 integer, parameter :: max_string=800
25 
26 !> Fortran derived type for adt observation operator
27 type :: ufo_adt_tlad
28  type(ufo_geoval) :: geoval_adt !< adt (traj)
29  logical :: ltraj = .false. !< trajectory set?
30 end type ufo_adt_tlad
31 
32 
33 ! ------------------------------------------------------------------------------
34 
35 contains
36 
37 ! ------------------------------------------------------------------------------
38 
39 subroutine ufo_adt_tlad_delete(self)
40 implicit none
41 type(ufo_adt_tlad), intent(inout) :: self
42 
43 self%ltraj = .false.
44 
45 end subroutine ufo_adt_tlad_delete
46 
47 ! ------------------------------------------------------------------------------
48 
49 subroutine ufo_adt_tlad_settraj(self, geovals)
50 implicit none
51 type(ufo_adt_tlad), intent(inout) :: self
52 type(ufo_geovals), intent(in) :: geovals
53 
54 character(len=*), parameter :: myname_="ufo_adt_tlad_settraj"
55 character(max_string) :: err_msg
56 
57 type(ufo_geoval), pointer :: geoval_adt
58 
59 ! check if adt variables is in geovals and get it
60 call ufo_geovals_get_var(geovals, var_abs_topo, geoval_adt)
61 
62 self%geoval_adt = geoval_adt
63 self%ltraj = .true.
64 
65 end subroutine ufo_adt_tlad_settraj
66 
67 
68 ! ------------------------------------------------------------------------------
69 
70 subroutine ufo_adt_simobs_tl(self, geovals, hofx)
71 implicit none
72 type(ufo_adt_tlad), intent(in) :: self
73 type(ufo_geovals), intent(in) :: geovals
74 real(c_double), intent(inout) :: hofx(:)
75 
76 character(len=*), parameter :: myname_="ufo_adt_simobs_tl"
77 character(max_string) :: err_msg
78 
79 integer :: iobs, nobs
80 type(ufo_geoval), pointer :: geoval_adt
81 real(kind_real) :: offset_obs, offset_hofx
82 
83 ! check if trajectory was set
84 if (.not. self%ltraj) then
85  write(err_msg,*) myname_, ' trajectory wasnt set!'
86  call abor1_ftn(err_msg)
87 endif
88 
89 ! check if nobs is consistent in geovals & hofx
90 nobs = size(hofx,1)
91 if (geovals%nobs /= nobs) then
92  write(err_msg,*) myname_, ' error: nobs inconsistent!'
93  call abor1_ftn(err_msg)
94 endif
95 
96 ! check if adt variable is in geovals and get it
97 call ufo_geovals_get_var(geovals, var_abs_topo, geoval_adt)
98 
99 ! Compute offset
100 offset_hofx=sum(geoval_adt%vals(1,:))/nobs
101 
102 ! adt obs operator
103 hofx = 0.0
104 do iobs = 1, nobs
105  hofx(iobs) = geoval_adt%vals(1,iobs) - offset_hofx
106 enddo
107 
108 end subroutine ufo_adt_simobs_tl
109 
110 ! ------------------------------------------------------------------------------
111 
112 subroutine ufo_adt_simobs_ad(self, geovals, hofx)
113 implicit none
114 type(ufo_adt_tlad), intent(in) :: self
115 type(ufo_geovals), intent(inout) :: geovals
116 real(c_double), intent(inout) :: hofx(:)
117 
118 character(len=*), parameter :: myname_="ufo_adt_simobs_ad"
119 character(max_string) :: err_msg
120 
121 integer :: iobs, nobs
122 type(ufo_geoval), pointer :: geoval_adt
123 real(kind_real) :: offset_hofx
124 
125 ! check if trajectory was set
126 if (.not. self%ltraj) then
127  write(err_msg,*) myname_, ' trajectory wasnt set!'
128  call abor1_ftn(err_msg)
129 endif
130 
131 ! check if nobs is consistent in geovals & hofx
132 nobs = size(hofx,1)
133 if (geovals%nobs /= nobs) then
134  write(err_msg,*) myname_, ' error: nobs inconsistent!'
135  call abor1_ftn(err_msg)
136 endif
137 
138 if (.not. geovals%linit ) geovals%linit=.true.
139 
140 ! check if adt variable is in geovals and get it
141 call ufo_geovals_get_var(geovals, var_abs_topo, geoval_adt)
142 
143 ! backward adt obs operator
144 
145 ! Compute offset
146 offset_hofx=sum(hofx)/nobs
147 
148 if (.not. allocated(geoval_adt%vals)) allocate(geoval_adt%vals(1,nobs))
149 geoval_adt%vals = 0.0
150 do iobs = 1, nobs
151  geoval_adt%vals(1,iobs) = geoval_adt%vals(1,iobs) + hofx(iobs) - offset_hofx
152 enddo
153 
154 end subroutine ufo_adt_simobs_ad
155 
156 end module ufo_adt_tlad_mod
subroutine, public ufo_geovals_get_var(self, varname, geoval, status)
integer, parameter max_string
subroutine, public ufo_adt_simobs_ad(self, geovals, hofx)
Fortran derived type for adt observation operator.
subroutine, public ufo_adt_tlad_settraj(self, geovals)
type to hold interpolated fields required by the obs operators
character(len=maxvarlen), public var_abs_topo
Fortran module to handle adt observations.
subroutine, public ufo_adt_simobs_tl(self, geovals, hofx)
Fortran module handling observation locations.
subroutine, public ufo_adt_tlad_delete(self)
Fortran interface to ObsSpace.
Definition: obsspace_mod.F90:9
type to hold interpolated field for one variable, one observation