UFO
ufo_seaicethickness_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 module for seaicethickness 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  use oops_variables_mod
20 
21  implicit none
22  private
23 
24  integer, parameter :: max_string=800
25 
26  !> Fortran derived type for the tl/ad observation operator
27  type, extends(ufo_basis_tlad), public :: ufo_seaicethickness_tlad
28  private
29  type(oops_variables), public :: obsvars
30  character(max_string) :: thickness_sim_option
31  type(ufo_geoval) :: icethick !< ice thickness (traj)
32  type(ufo_geoval) :: icefrac !< ice fraction (traj)
33  type(ufo_geoval) :: snowthick!< snow thickness(traj)
34  real(kind=kind_real) :: rho_ice = 905.0 !< [kg/m3]
35  real(kind=kind_real) :: rho_snow = 330.0 !< [kg/m3]
36  real(kind=kind_real) :: rho_water= 1000.0!< [kg/m3]
37 
38  contains
39  procedure :: setup => ufo_seaicethickness_tlad_setup
40  procedure :: delete => ufo_seaicethickness_tlad_delete
41  procedure :: settraj => ufo_seaicethickness_tlad_settraj
42  procedure :: simobs_tl => ufo_seaicethickness_simobs_tl
43  procedure :: simobs_ad => ufo_seaicethickness_simobs_ad
45 
46 contains
47 
48 ! ------------------------------------------------------------------------------
49 subroutine ufo_seaicethickness_tlad_setup(self, f_conf)
50 implicit none
51 class(ufo_seaicethickness_tlad), intent(inout) :: self
52 type(fckit_configuration), intent(in) :: f_conf
53 real(kind=kind_real) :: rho_ice, rho_snow, rho_water
54 integer :: ivar, nvars
55 character(max_string) :: err_msg
56 
57 nvars = self%obsvars%nvars()
58 if (nvars /= 1) then
59  write(err_msg,*) 'ufo_seaicethickness_tlad_setup error: only variables size 1 supported!'
60  call abor1_ftn(err_msg)
61 endif
62 
63 ! Set thickness-simulate option from ymal file
64 !self%thickness_sim_option = self%obsvars%variable(1)
65 
66 end subroutine ufo_seaicethickness_tlad_setup
67 
68 ! ------------------------------------------------------------------------------
70 implicit none
71 class(ufo_seaicethickness_tlad), intent(inout) :: self
72 
73 self%ltraj = .false.
74 
76 
77 ! ------------------------------------------------------------------------------
78 subroutine ufo_seaicethickness_tlad_settraj(self, geovals, obss)
79 implicit none
80 class(ufo_seaicethickness_tlad), intent(inout) :: self
81 type(ufo_geovals), intent(in) :: geovals
82 type(c_ptr), value, intent(in) :: obss
83 
84 character(len=*), parameter :: myname_="ufo_seaicethick_tlad_settraj"
85 
86 type(ufo_geoval), pointer :: icethick, icefrac, snowthick
87 
88 ! check if sea ice thickness variables is in geovals and get it
89 call ufo_geovals_get_var(geovals, var_seaicethick, icethick)
90 
91 ! check if sea ice fraction variables is in geovals and get it
92 call ufo_geovals_get_var(geovals, var_seaicefrac, icefrac)
93 
94 if (trim(self%obsvars%variable(1)) == "sea_ice_freeboard") then
95  call ufo_geovals_get_var(geovals, var_seaicesnowthick, snowthick)
96  self%snowthick= snowthick
97 endif
98 
99 self%icethick = icethick
100 self%icefrac = icefrac
101 self%ltraj = .true.
102 
104 
105 ! ------------------------------------------------------------------------------
106 subroutine ufo_seaicethickness_simobs_tl(self, geovals, hofx, obss)
107 implicit none
108 class(ufo_seaicethickness_tlad), intent(in) :: self
109 type(ufo_geovals), intent(in) :: geovals
110 real(c_double), intent(inout) :: hofx(:)
111 type(c_ptr), value, intent(in) :: obss
112 
113 character(len=*), parameter :: myname_="ufo_seaicethick_simobs_tl"
114 character(max_string) :: err_msg
115 
116 integer :: iobs, icat, ncat
117 type(ufo_geoval), pointer :: icethick_d, icefrac_d, snowthick
118 real(kind=kind_real) :: rho_wiw, rho_wsw
119 
120 ! check if trajectory was set
121 if (.not. self%ltraj) then
122  write(err_msg,*) myname_, ' trajectory wasnt set!'
123  call abor1_ftn(err_msg)
124 endif
125 
126 ! check if nlocs is consistent in geovals & hofx
127 if (geovals%nlocs /= size(hofx,1)) then
128  write(err_msg,*) myname_, ' error: nlocs inconsistent!'
129  call abor1_ftn(err_msg)
130 endif
131 
132 ! check if sea ice fraction variable is in geovals and get it
133 call ufo_geovals_get_var(geovals, var_seaicefrac, icefrac_d)
134 
135 ! check if sea ice thickness variable is in geovals and get it
136 call ufo_geovals_get_var(geovals, var_seaicethick, icethick_d)
137 
138 if (trim(self%obsvars%variable(1)) == "sea_ice_freeboard") then
139  rho_wiw = (self%rho_water-self%rho_ice)/self%rho_water
140  rho_wsw = (self%rho_water-self%rho_snow)/self%rho_water
141 endif
142 
143 ! sea ice thickness obs operator
144 ncat = icefrac_d%nval
145 hofx = 0.0
146 
147 select case (trim(self%obsvars%variable(1)))
148 case ("sea_ice_freeboard")
149  do iobs = 1, size(hofx,1)
150  do icat = 1, ncat
151  hofx(iobs) = hofx(iobs) + &
152  rho_wiw * self%icefrac%vals(icat,iobs) * icethick_d%vals(icat,iobs) + &
153  rho_wiw * icefrac_d%vals(icat,iobs) * self%icethick%vals(icat,iobs) + &
154  rho_wsw * icefrac_d%vals(icat,iobs) * self%snowthick%vals(icat,iobs)
155  enddo
156  enddo
157 case ("sea_ice_thickness")
158  do iobs = 1, size(hofx,1)
159  do icat = 1, ncat
160  hofx(iobs) = hofx(iobs) + &
161  self%icefrac%vals(icat,iobs) * icethick_d%vals(icat,iobs) + &
162  icefrac_d%vals(icat,iobs) * self%icethick%vals(icat,iobs)
163  enddo
164  enddo
165 case default
166  write(err_msg,*) myname_, ' error: no match seaice thickness_option!'
167  call abor1_ftn(err_msg)
168 end select
169 
170 end subroutine ufo_seaicethickness_simobs_tl
171 
172 ! ------------------------------------------------------------------------------
173 subroutine ufo_seaicethickness_simobs_ad(self, geovals, hofx, obss)
174 implicit none
175 class(ufo_seaicethickness_tlad), intent(in) :: self
176 type(ufo_geovals), intent(inout) :: geovals
177 real(c_double), intent(in) :: hofx(:)
178 type(c_ptr), value, intent(in) :: obss
179 
180 character(len=*), parameter :: myname_="ufo_seaicethick_simobs_ad"
181 character(max_string) :: err_msg
182 
183 integer :: iobs, icat, ncat
184 type(ufo_geoval), pointer :: icefrac_d, icethick_d
185 real(c_double) :: missing
186 real(kind=kind_real) :: rho_wiw, rho_wsw
187 
188 !> Set missing value
189 missing = missing_value(missing)
190 
191 ! check if trajectory was set
192 if (.not. self%ltraj) then
193  write(err_msg,*) myname_, ' trajectory wasnt set!'
194  call abor1_ftn(err_msg)
195 endif
196 
197 ! check if nlocs is consistent in geovals & hofx
198 if (geovals%nlocs /= size(hofx,1)) then
199  write(err_msg,*) myname_, ' error: nlocs inconsistent!'
200  call abor1_ftn(err_msg)
201 endif
202 
203 if (trim(self%obsvars%variable(1)) == "sea_ice_freeboard") then
204  rho_wiw = (self%rho_water-self%rho_ice)/self%rho_water
205  rho_wsw = (self%rho_water-self%rho_snow)/self%rho_water
206 endif
207 
208 if (.not. geovals%linit ) geovals%linit=.true.
209 
210 ! Get sea-ice fraction & thickness geovals
211 call ufo_geovals_get_var(geovals, var_seaicefrac, icefrac_d)
212 call ufo_geovals_get_var(geovals, var_seaicethick, icethick_d)
213 
214 ncat = self%icethick%nval
215 if (.not.(allocated(icefrac_d%vals) .or. .not. allocated(icethick_d%vals))) then
216  if (ncat < 1) then
217  write(err_msg,*) myname_, ' unknown number of categories'
218  call abor1_ftn(err_msg)
219  endif
220  if (.not. allocated(icefrac_d%vals)) allocate(icefrac_d%vals(ncat,size(hofx,1)))
221  if (.not. allocated(icethick_d%vals)) allocate(icethick_d%vals(ncat, size(hofx,1)))
222 end if
223 
224 ! backward sea ice thickness obs operator
225 
226 if (.not. allocated(icefrac_d%vals)) allocate(icefrac_d%vals(ncat,size(hofx,1)))
227 if (.not. allocated(icethick_d%vals)) allocate(icethick_d%vals(ncat, size(hofx,1)))
228 
229 icethick_d%vals = 0.0
230 icefrac_d%vals = 0.0
231 
232 select case (trim(self%obsvars%variable(1)))
233 case ("sea_ice_freeboard")
234  do iobs = 1, size(hofx,1)
235  if (hofx(iobs) /= missing) then
236  do icat = 1, ncat
237  icefrac_d%vals(icat,iobs) = icefrac_d%vals(icat,iobs)&
238  + rho_wiw*self%icethick%vals(icat,iobs) * hofx(iobs)&
239  + rho_wsw*self%snowthick%vals(icat,iobs) * hofx(iobs)
240  icethick_d%vals(icat,iobs) = icethick_d%vals(icat,iobs)&
241  + rho_wiw*self%icefrac%vals(icat,iobs) * hofx(iobs)
242  end do
243  end if
244  enddo
245 case ("sea_ice_thickness")
246  do iobs = 1, size(hofx,1)
247  if (hofx(iobs) /= missing) then
248  do icat = 1, ncat
249  icefrac_d%vals(icat,iobs) = icefrac_d%vals(icat,iobs) + self%icethick%vals(icat,iobs) * hofx(iobs)
250  icethick_d%vals(icat,iobs) = icethick_d%vals(icat,iobs) + self%icefrac%vals(icat,iobs) * hofx(iobs)
251  end do
252  end if
253  enddo
254 case default
255  write(err_msg,*) myname_, ' error: no match seaice thickness_option!'
256  call abor1_ftn(err_msg)
257 end select
258 
259 end subroutine ufo_seaicethickness_simobs_ad
260 
261 ! ------------------------------------------------------------------------------
262 
ufo_avgkernel_mod::max_string
integer, parameter max_string
Definition: ufo_avgkernel_mod.F90:17
ufo_basis_tlad_mod
Definition: ufo_basis_tlad_mod.F90:6
ufo_seaicethickness_tlad_mod
Fortran module for seaicethickness tl/ad observation operator.
Definition: ufo_seaicethickness_tlad_mod.F90:8
ufo_vars_mod::var_seaicethick
character(len=maxvarlen), public var_seaicethick
Definition: ufo_variables_mod.F90:138
ufo_vars_mod::var_seaicefrac
character(len=maxvarlen), public var_seaicefrac
Definition: ufo_variables_mod.F90:137
ufo_seaicethickness_tlad_mod::ufo_seaicethickness_tlad_setup
subroutine ufo_seaicethickness_tlad_setup(self, f_conf)
Definition: ufo_seaicethickness_tlad_mod.F90:50
ufo_seaicethickness_tlad_mod::ufo_seaicethickness_tlad_delete
subroutine ufo_seaicethickness_tlad_delete(self)
Definition: ufo_seaicethickness_tlad_mod.F90:70
ufo_basis_tlad_mod::ufo_basis_tlad
Definition: ufo_basis_tlad_mod.F90:12
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_seaicethickness_tlad_mod::ufo_seaicethickness_simobs_ad
subroutine ufo_seaicethickness_simobs_ad(self, geovals, hofx, obss)
Definition: ufo_seaicethickness_tlad_mod.F90:174
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_geovals_mod::ufo_geovals_get_var
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Definition: ufo_geovals_mod.F90:128
ufo_seaicethickness_tlad_mod::ufo_seaicethickness_simobs_tl
subroutine ufo_seaicethickness_simobs_tl(self, geovals, hofx, obss)
Definition: ufo_seaicethickness_tlad_mod.F90:107
ufo_seaicethickness_tlad_mod::ufo_seaicethickness_tlad
Fortran derived type for the tl/ad observation operator.
Definition: ufo_seaicethickness_tlad_mod.F90:27
ufo_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
ufo_geovals_mod::ufo_geoval
type to hold interpolated field for one variable, one observation
Definition: ufo_geovals_mod.F90:40
ufo_seaicethickness_tlad_mod::ufo_seaicethickness_tlad_settraj
subroutine ufo_seaicethickness_tlad_settraj(self, geovals, obss)
Definition: ufo_seaicethickness_tlad_mod.F90:79
ufo_vars_mod::var_seaicesnowthick
character(len=maxvarlen), public var_seaicesnowthick
Definition: ufo_variables_mod.F90:139