10 use fckit_configuration_module,
only: fckit_configuration
18 use missing_values_mod
29 real(c_double) :: r_miss_val
30 real (kind=kind_real),
allocatable :: jac(:,:)
44 type(fckit_configuration),
intent(in) :: f_conf
52 if (
allocated(self%jac))
deallocate(self%jac)
64 type(c_ptr),
value,
intent(in) :: obss
66 character(len=*),
parameter :: myname_=
"ufo_coolskin_tlad_settraj"
67 type(
ufo_geoval),
pointer :: S_ns,H_I,H_s,R_nl,Td,u
70 self%nlocs = obsspace_get_nlocs(obss)
83 self%r_miss_val = missing_value(self%r_miss_val)
86 allocate(self%jac(6,self%nlocs))
87 do iobs = 1, self%nlocs
105 real(c_double),
intent(inout) :: hofx(:)
106 type(c_ptr),
value,
intent(in) :: obss
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
123 if (.not. self%ltraj)
then
124 write(err_msg,*) myname_,
' trajectory wasnt set!'
125 call abor1_ftn(err_msg)
131 if (geovals%nlocs /= nobs)
then
132 write(err_msg,*) myname_,
' error: nobs inconsistent!'
133 call abor1_ftn(err_msg)
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)
156 real(c_double),
intent(in) :: hofx(:)
157 type(c_ptr),
value,
intent(in) :: obss
159 character(len=*),
parameter :: myname_=
"ufo_coolskin_simobs_ad"
160 character(max_string) :: err_msg
162 integer :: iobs, nobs
164 type(
ufo_geoval),
pointer :: S_ns, H_I, H_s, R_nl, Td, u
168 if (.not. self%ltraj)
then
169 write(err_msg,*) myname_,
' trajectory wasnt set!'
170 call abor1_ftn(err_msg)
175 if (geovals%nlocs /= nobs)
then
176 write(err_msg,*) myname_,
' error: nobs inconsistent!'
177 call abor1_ftn(err_msg)
180 if (.not. geovals%linit ) geovals%linit=.true.
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)
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