FV3 Bundle
ufo_seasurfacetemp_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 sea-surface temperature observations
7 
9 
10 use iso_c_binding
11 use ufo_vars_mod
13 use kinds
14 
15 implicit none
20 private
21 integer, parameter :: max_string=800
22 
23 !> Fortran derived type for sea-surface temperature observation operator
26 
27 
28 ! ------------------------------------------------------------------------------
29 
30 contains
31 
32 
33 ! ------------------------------------------------------------------------------
34 
35 subroutine ufo_seasurfacetemp_tlad_settraj(self, geovals)
36 implicit none
37 type(ufo_seasurfacetemp_tlad), intent(inout) :: self
38 type(ufo_geovals), intent(in) :: geovals
39 
40 character(len=*), parameter :: myname_="ufo_seasurfacetemp_tlad_settraj"
41 character(max_string) :: err_msg
42 
43 type(ufo_geoval), pointer :: geoval_sst
44 
45 ! since observation operator is linear, don't care about trajectory itself
46 
48 
49 
50 ! ------------------------------------------------------------------------------
51 
52 subroutine ufo_seasurfacetemp_simobs_tl(self, geovals, hofx)
53 implicit none
54 type(ufo_seasurfacetemp_tlad), intent(in) :: self
55 type(ufo_geovals), intent(in) :: geovals
56 real(c_double),intent(inout) :: hofx(:)
57 
58 character(len=*), parameter :: myname_="ufo_seasurfacetemp_simobs_tl"
59 character(max_string) :: err_msg
60 
61 integer :: iobs
62 type(ufo_geoval), pointer :: geoval_sst
63 
64 print *, myname_, ' nobs: ', geovals%nobs, size(hofx,1)
65 
66 ! check if nobs is consistent in geovals & hofx
67 if (geovals%nobs /= size(hofx,1)) then
68  write(err_msg,*) myname_, ' error: nobs inconsistent!'
69  call abor1_ftn(err_msg)
70 endif
71 
72 ! check if sst variables is in geovals and get it
73 call ufo_geovals_get_var(geovals, var_ocn_sst, geoval_sst)
74 
75 ! sst obs operator
76 do iobs = 1, size(hofx,1)
77  hofx(iobs) = geoval_sst%vals(1,iobs)
78 enddo
79 
80 end subroutine ufo_seasurfacetemp_simobs_tl
81 
82 ! ------------------------------------------------------------------------------
83 
84 subroutine ufo_seasurfacetemp_simobs_ad(self, geovals, hofx)
85 implicit none
86 type(ufo_seasurfacetemp_tlad), intent(in) :: self
87 type(ufo_geovals), intent(inout) :: geovals
88 real(c_double),intent(inout) :: hofx(:)
89 
90 character(len=*), parameter :: myname_="ufo_seasurfacetemp_simobs_ad"
91 character(max_string) :: err_msg
92 
93 integer :: iobs
94 type(ufo_geoval), pointer :: geoval_sst
95 
96 ! check if nobs is consistent in geovals & hofx
97 if (geovals%nobs /= size(hofx,1)) then
98  write(err_msg,*) myname_, ' error: nobs inconsistent!'
99  call abor1_ftn(err_msg)
100 endif
101 
102 if (.not. geovals%linit ) geovals%linit=.true.
103 
104 ! check if sst variables is in geovals and get it
105 call ufo_geovals_get_var(geovals, var_ocn_sst, geoval_sst)
106 
107 if (.not.(allocated(geoval_sst%vals))) then
108  geoval_sst%nval=1
109  allocate(geoval_sst%vals(1,size(hofx,1)))
110 end if
111 
112 ! backward sst obs operator
113 geoval_sst%vals=0.0
114 do iobs = 1, size(hofx,1)
115  geoval_sst%vals(1,iobs) = geoval_sst%vals(1,iobs) + hofx(iobs)
116 enddo
117 
118 end subroutine ufo_seasurfacetemp_simobs_ad
119 
character(len=maxvarlen), public var_ocn_sst
subroutine, public ufo_geovals_get_var(self, varname, geoval, status)
integer, parameter max_string
Fortran module to handle sea-surface temperature observations.
subroutine, public ufo_seasurfacetemp_simobs_tl(self, geovals, hofx)
subroutine, public ufo_seasurfacetemp_tlad_settraj(self, geovals)
type to hold interpolated fields required by the obs operators
Fortran derived type for sea-surface temperature observation operator.
subroutine, public ufo_seasurfacetemp_simobs_ad(self, geovals, hofx)
type to hold interpolated field for one variable, one observation