UFO
ufo_coolskin_tlad_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 coolskin module for tl/ad observation operator
7 
9 
10  use fckit_configuration_module, only: fckit_configuration
11  use iso_c_binding
12  use kinds
13 
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 tl/ad observation operator
26  type, extends(ufo_basis_tlad), public :: ufo_coolskin_tlad
27  private
28  integer :: nlocs !< Local number of obs
29  real(c_double) :: r_miss_val !< Missing value flag
30  real (kind=kind_real), allocatable :: jac(:,:) !< Jacobian [6*nobs]
31  contains
32  procedure :: setup => ufo_coolskin_tlad_setup
33  procedure :: delete => ufo_coolskin_tlad_delete
34  procedure :: settraj => ufo_coolskin_tlad_settraj
35  procedure :: simobs_tl => ufo_coolskin_simobs_tl
36  procedure :: simobs_ad => ufo_coolskin_simobs_ad
37  end type ufo_coolskin_tlad
38 
39 contains
40 
41 ! ------------------------------------------------------------------------------
42 subroutine ufo_coolskin_tlad_setup(self, f_conf)
43 class(ufo_coolskin_tlad), intent(inout) :: self
44 type(fckit_configuration), intent(in) :: f_conf
45 
46 end subroutine ufo_coolskin_tlad_setup
47 
48 ! ------------------------------------------------------------------------------
49 subroutine ufo_coolskin_tlad_delete(self)
50 class(ufo_coolskin_tlad), intent(inout) :: self
51 
52 if (allocated(self%jac)) deallocate(self%jac)
53 self%ltraj = .false.
54 
55 end subroutine ufo_coolskin_tlad_delete
56 
57 ! ------------------------------------------------------------------------------
58 subroutine ufo_coolskin_tlad_settraj(self, geovals, obss)
60 
61 implicit none
62 class(ufo_coolskin_tlad), intent(inout) :: self
63 type(ufo_geovals), intent(in) :: geovals
64 type(c_ptr), value, intent(in) :: obss
65 
66 character(len=*), parameter :: myname_="ufo_coolskin_tlad_settraj"
67 type(ufo_geoval), pointer :: S_ns,H_I,H_s,R_nl,Td,u
68 
69 integer :: iobs
70 self%nlocs = obsspace_get_nlocs(obss)
71 
72 ! Get trajectory geovals (state at which the Jacobian is computed)
73 call ufo_geovals_get_var(geovals, var_ocn_sst, td)
74 call ufo_geovals_get_var(geovals, var_sw_rad , r_nl )
75 call ufo_geovals_get_var(geovals, var_latent_heat , h_i )
76 call ufo_geovals_get_var(geovals, var_sens_heat , h_s )
77 call ufo_geovals_get_var(geovals, var_lw_rad , s_ns )
78 call ufo_geovals_get_var(geovals, var_sea_fric_vel , u )
79 
80 self%ltraj = .true.
81 
82 ! Set missing flag
83 self%r_miss_val = missing_value(self%r_miss_val)
84 
85 ! Initialize and compute traj dependent Jacobian
86 allocate(self%jac(6,self%nlocs))
87 do iobs = 1, self%nlocs
88  call ufo_coolskin_jac(self%jac(:,iobs),&
89  s_ns%vals(1,iobs),&
90  h_i%vals(1,iobs),&
91  h_s%vals(1,iobs),&
92  r_nl%vals(1,iobs),&
93  td%vals(1,iobs),&
94  u%vals(1,iobs))
95 enddo
96 
97 end subroutine ufo_coolskin_tlad_settraj
98 
99 ! ------------------------------------------------------------------------------
100 subroutine ufo_coolskin_simobs_tl(self, geovals, hofx, obss)
101 
102 implicit none
103 class(ufo_coolskin_tlad), intent(in) :: self
104 type(ufo_geovals), intent(in) :: geovals
105 real(c_double), intent(inout) :: hofx(:)
106 type(c_ptr), value, intent(in) :: obss
107 
108 character(len=*), parameter :: myname_="ufo_coolskin_simobs_tl"
109 character(max_string) :: err_msg
110 integer :: iobs, nobs
111 type(ufo_geoval), pointer :: S_ns,H_I,H_s,R_nl,Td,u
112 
113 
114 ! Get geovals increments
115 call ufo_geovals_get_var(geovals, var_ocn_sst, td)
116 call ufo_geovals_get_var(geovals, var_sw_rad , r_nl )
117 call ufo_geovals_get_var(geovals, var_latent_heat , h_i )
118 call ufo_geovals_get_var(geovals, var_sens_heat , h_s )
119 call ufo_geovals_get_var(geovals, var_lw_rad , s_ns )
120 call ufo_geovals_get_var(geovals, var_sea_fric_vel , u )
121 
122 ! check if trajectory was set
123 if (.not. self%ltraj) then
124  write(err_msg,*) myname_, ' trajectory wasnt set!'
125  call abor1_ftn(err_msg)
126 endif
127 
128 ! check if nobs is consistent in geovals & hofx
129 nobs = self%nlocs
130 
131 if (geovals%nlocs /= nobs) then
132  write(err_msg,*) myname_, ' error: nobs inconsistent!'
133  call abor1_ftn(err_msg)
134 endif
135 
136 
137 ! Perturbation coolskin obs operator
138 hofx = 0.0
139 do iobs = 1, self%nlocs
140  hofx(iobs) = self%jac(1,iobs)*s_ns%vals(1,iobs) + &
141  self%jac(2,iobs)*h_i%vals(1,iobs) + &
142  self%jac(3,iobs)*h_s%vals(1,iobs) + &
143  self%jac(4,iobs)*r_nl%vals(1,iobs) + &
144  self%jac(5,iobs)*td%vals(1,iobs) + &
145  self%jac(6,iobs)*u%vals(1,iobs)
146 enddo
147 
148 end subroutine ufo_coolskin_simobs_tl
149 
150 ! ------------------------------------------------------------------------------
151 subroutine ufo_coolskin_simobs_ad(self, geovals, hofx, obss)
152 
153 implicit none
154 class(ufo_coolskin_tlad), intent(in) :: self
155 type(ufo_geovals), intent(inout) :: geovals
156 real(c_double), intent(in) :: hofx(:)
157 type(c_ptr), value, intent(in) :: obss
158 
159 character(len=*), parameter :: myname_="ufo_coolskin_simobs_ad"
160 character(max_string) :: err_msg
161 
162 integer :: iobs, nobs
163 
164 type(ufo_geoval), pointer :: S_ns, H_I, H_s, R_nl, Td, u
165 
166 
167 ! check if trajectory was set
168 if (.not. self%ltraj) then
169  write(err_msg,*) myname_, ' trajectory wasnt set!'
170  call abor1_ftn(err_msg)
171 endif
172 
173 ! check if nobs is consistent in geovals & hofx
174 nobs = self%nlocs
175 if (geovals%nlocs /= nobs) then
176  write(err_msg,*) myname_, ' error: nobs inconsistent!'
177  call abor1_ftn(err_msg)
178 endif
179 
180 if (.not. geovals%linit ) geovals%linit=.true.
181 
182 ! Associate geovals pointers
183 call ufo_geovals_get_var(geovals, var_ocn_sst, td)
184 call ufo_geovals_get_var(geovals, var_sw_rad , r_nl )
185 call ufo_geovals_get_var(geovals, var_latent_heat , h_i )
186 call ufo_geovals_get_var(geovals, var_sens_heat , h_s )
187 call ufo_geovals_get_var(geovals, var_lw_rad , s_ns )
188 call ufo_geovals_get_var(geovals, var_sea_fric_vel , u )
189 
190 ! Apply adjoint obs operator
191 do iobs = 1, nobs
192  if (hofx(iobs)/=self%r_miss_val) then
193  s_ns%vals(1,iobs) = s_ns%vals(1,iobs) + self%jac(1,iobs)*hofx(iobs)
194  h_i%vals(1,iobs) = h_i%vals(1,iobs) + self%jac(2,iobs)*hofx(iobs)
195  h_s%vals(1,iobs) = h_s%vals(1,iobs) + self%jac(3,iobs)*hofx(iobs)
196  r_nl%vals(1,iobs) = r_nl%vals(1,iobs) + self%jac(4,iobs)*hofx(iobs)
197  td%vals(1,iobs) = td%vals(1,iobs) + self%jac(5,iobs)*hofx(iobs)
198  u%vals(1,iobs) = u%vals(1,iobs) + self%jac(6,iobs)*hofx(iobs)
199  end if
200 enddo
201 
202 end subroutine ufo_coolskin_simobs_ad
203 
204 !--------------------------------------------------
205 
206 end module ufo_coolskin_tlad_mod
subroutine, public ufo_coolskin_jac(jac, S_ns, H_I, H_s, R_nl, Tdc, u0)
Fortran coolskin module for tl/ad observation operator.
subroutine ufo_coolskin_tlad_setup(self, f_conf)
subroutine ufo_coolskin_tlad_delete(self)
subroutine ufo_coolskin_simobs_ad(self, geovals, hofx, obss)
subroutine ufo_coolskin_simobs_tl(self, geovals, hofx, obss)
integer, parameter max_string
subroutine ufo_coolskin_tlad_settraj(self, geovals, obss)
subroutine, public ufo_geovals_get_var(self, varname, geoval)
character(len=maxvarlen), public var_sens_heat
character(len=maxvarlen), public var_latent_heat
character(len=maxvarlen), public var_lw_rad
character(len=maxvarlen), public var_sw_rad
character(len=maxvarlen), public var_ocn_sst
character(len=maxvarlen), parameter, public var_sea_fric_vel
Fortran derived type for the tl/ad observation operator.
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators