UFO
ufo_avgkernel_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2020 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 averaging kernel observation operator
7 
9 
10  use oops_variables_mod
11  use ufo_vars_mod
12  use kinds
13  use iso_c_binding
14 
15  implicit none
16  private
17  integer, parameter :: max_string=800
18 
19 !> Fortran derived type for the observation type
20  type, public :: ufo_avgkernel
21  private
22  type(oops_variables), public :: obsvars
23  type(oops_variables), public :: geovars
24  integer :: nlayers_kernel
25  real(kind_real), allocatable, dimension(:) :: ak_kernel, bk_kernel
26  character(kind=c_char,len=:), allocatable :: obskernelvar, tracervars(:)
27  logical :: troposphere, totalcolumn
28  real(kind_real) :: convert_factor_model, convert_factor_hofx
29  contains
30  procedure :: setup => ufo_avgkernel_setup
31  procedure :: simobs => ufo_avgkernel_simobs
32  final :: destructor
33  end type ufo_avgkernel
34 
35 contains
36 
37 ! ------------------------------------------------------------------------------
38 subroutine ufo_avgkernel_setup(self, f_conf)
39  use fckit_configuration_module, only: fckit_configuration
40  use ufo_constants_mod, only: one
41  implicit none
42  class(ufo_avgkernel), intent(inout) :: self
43  type(fckit_configuration), intent(in) :: f_conf
44  type(fckit_configuration) :: f_confOpts
45  integer :: nlevs_yaml
46  integer :: ivar, nvars
47  character(len=max_string) :: err_msg
48 
49  ! get configuration for the averaging kernel operator
50  call f_conf%get_or_die("obs options",f_confopts)
51 
52  call f_confopts%get_or_die("nlayers_kernel", self%nlayers_kernel)
53  nlevs_yaml = f_confopts%get_size("ak")
54  if (nlevs_yaml /= self%nlayers_kernel+1) then
55  write(err_msg, *) "ufo_avgkernel_setup error: YAML nlayers_kernel != size of ak array"
56  call abor1_ftn(err_msg)
57  end if
58  nlevs_yaml = f_confopts%get_size("bk")
59  if (nlevs_yaml /= self%nlayers_kernel+1) then
60  write(err_msg, *) "ufo_avgkernel_setup error: YAML nlayers_kernel != size of bk array"
61  call abor1_ftn(err_msg)
62  end if
63 
64  ! allocate and read ak/bk for the averaging kernel
65  allocate(self%ak_kernel(nlevs_yaml))
66  allocate(self%bk_kernel(nlevs_yaml))
67  call f_confopts%get_or_die("ak", self%ak_kernel)
68  call f_confopts%get_or_die("bk", self%bk_kernel)
69 
70  ! get variable name from IODA for observation averaging kernel
71  if (.not. f_confopts%get("AvgKernelVar", self%obskernelvar)) then
72  self%obskernelvar = "averaging_kernel" ! default option
73  end if
74 
75  ! get name of geoval/tracer to use from the model
76  nvars = self%obsvars%nvars()
77  call f_confopts%get_or_die("tracer variables", self%tracervars)
78 
79  ! determine if this is a total column or troposphere calculation
80  ! support stratosphere, etc. later?
81  if (.not. f_confopts%get("tropospheric column", self%troposphere)) self%troposphere = .false.
82  if (.not. f_confopts%get("total column", self%totalcolumn)) self%totalcolumn = .false.
83 
84  ! do we need a conversion factor, say between ppmv and unity?
85  if (.not. f_confopts%get("model units coeff", self%convert_factor_model)) self%convert_factor_model = one
86  if (.not. f_confopts%get("hofx units coeff", self%convert_factor_hofx)) self%convert_factor_hofx = one
87 
88  ! add variables to geovars that are needed
89  ! specified tracers
90  do ivar = 1, nvars
91  call self%geovars%push_back(self%tracervars(ivar))
92  end do
93  ! surface pressure
94  call self%geovars%push_back(var_ps)
95  ! column pressure both layer and interface
96  call self%geovars%push_back(var_prs)
97  call self%geovars%push_back(var_prsi)
98  ! need air temperature for number density conversion
99  call self%geovars%push_back(var_ts)
100  ! need geopotential height on levels to convert units
101  call self%geovars%push_back(var_zi)
102 
103 end subroutine ufo_avgkernel_setup
104 
105 ! ------------------------------------------------------------------------------
106 subroutine destructor(self)
107  implicit none
108  type(ufo_avgkernel), intent(inout) :: self
109 
110 end subroutine destructor
111 
112 ! ------------------------------------------------------------------------------
113 ! averaging kernel observation operator
114 subroutine ufo_avgkernel_simobs(self, geovals_in, obss, nvars, nlocs, hofx)
115  use kinds
118  use ufo_constants_mod, only: zero, half
120  use iso_c_binding
121  use obsspace_mod
122  implicit none
123  class(ufo_avgkernel), intent(in) :: self
124  integer, intent(in) :: nvars, nlocs
125  type(ufo_geovals), intent(in) :: geovals_in
126  real(c_double), intent(inout) :: hofx(nvars, nlocs)
127  type(c_ptr), value, intent(in) :: obss
128 
129  ! Local variables
130  type(ufo_geoval), pointer :: prsi, prsl, psfc, temp, phii, tracer
131  integer :: ivar, iobs, ilev
132  character(len=MAXVARLEN) :: geovar, varstring
133  character(len=4) :: levstr
134  real(kind_real), allocatable, dimension(:,:) :: avgkernel_obs, prsl_obs, prsi_obs
135  real(kind_real), allocatable, dimension(:) :: airmass_tot, airmass_trop
136  integer, allocatable, dimension(:) :: troplev_obs
137  real(kind_real) :: hofx_tmp
138  type(ufo_geovals) :: geovals
139 
140  ! get geovals of atmospheric pressure
141  call ufo_geovals_copy(geovals_in, geovals) ! dont want to change geovals_in
142  call ufo_geovals_reorderzdir(geovals, self%geovars%variable(nvars+2), "bottom2top")
143  call ufo_geovals_get_var(geovals, self%geovars%variable(nvars+1), psfc)
144  call ufo_geovals_get_var(geovals, self%geovars%variable(nvars+2), prsl)
145  call ufo_geovals_get_var(geovals, self%geovars%variable(nvars+3), prsi)
146  call ufo_geovals_get_var(geovals, self%geovars%variable(nvars+4), temp)
147  call ufo_geovals_get_var(geovals, self%geovars%variable(nvars+5), phii)
148 
149  ! grab necesary metadata from IODA
150  ! get observation averaging kernel
151  ! once 2D arrays are allowed, rewrite/simplify this part
152  allocate(avgkernel_obs(self%nlayers_kernel, nlocs))
153  do ilev = 1, self%nlayers_kernel
154  write(levstr, fmt = "(I3)") ilev
155  levstr = adjustl(levstr)
156  varstring = trim(self%obskernelvar)//"_"//trim(levstr)
157  call obsspace_get_db(obss, "MetaData", trim(varstring), avgkernel_obs(ilev, :))
158  end do
159 
160  ! compute prsl_obs/prsi_obs from ak/bk/psfc
161  allocate(prsl_obs(self%nlayers_kernel, nlocs))
162  allocate(prsi_obs(self%nlayers_kernel+1, nlocs))
163  ! prsi_obs calculation
164  do ilev = 1, self%nlayers_kernel+1
165  prsi_obs(ilev,:) = self%ak_kernel(ilev) + self%bk_kernel(ilev) * psfc%vals(1,:)
166  end do
167  ! using simple averaging for now for prsl, can use more complex way later
168  do ilev = 1, self%nlayers_kernel
169  prsl_obs(ilev,:) = (prsi_obs(ilev,:) + prsi_obs(ilev+1,:)) * half
170  end do
171 
172  if (self%troposphere) then
173  allocate(troplev_obs(nlocs))
174  allocate(airmass_trop(nlocs))
175  allocate(airmass_tot(nlocs))
176  call obsspace_get_db(obss, "MetaData", "troposphere_layer_index", troplev_obs)
177  call obsspace_get_db(obss, "MetaData", "air_mass_factor_troposphere", airmass_trop)
178  call obsspace_get_db(obss, "MetaData", "air_mass_factor_total", airmass_tot)
179  end if
180 
181  ! loop through all variables
182  do ivar = 1, nvars
183  geovar = self%tracervars(ivar)
184  call ufo_geovals_get_var(geovals, geovar, tracer)
185  do iobs = 1, nlocs
186  if (avgkernel_obs(1,iobs) > -1.0e9_kind_real) then ! take care of missing obs
187  if (self%troposphere) then
188  call simulate_column_ob(self%nlayers_kernel, tracer%nval, avgkernel_obs(:,iobs), &
189  prsl_obs(:,iobs), prsl%vals(:,iobs), temp%vals(:,iobs),&
190  phii%vals(:,iobs), tracer%vals(:,iobs)*self%convert_factor_model, &
191  hofx_tmp, troplev_obs(iobs), airmass_tot(iobs), airmass_trop(iobs))
192  hofx(ivar,iobs) = hofx_tmp * self%convert_factor_hofx
193  end if
194  else
195  hofx(ivar,iobs) = zero ! default if we are unable to compute averaging kernel
196  end if
197  end do
198  end do
199 
200 
201 end subroutine ufo_avgkernel_simobs
202 
203 
204 ! ------------------------------------------------------------------------------
205 
206 end module ufo_avgkernel_mod
ufo_avgkernel_mod::max_string
integer, parameter max_string
Definition: ufo_avgkernel_mod.F90:17
satcolumn_mod
Definition: ufo_satcolumn_mod.F90:1
ufo_vars_mod::var_zi
character(len=maxvarlen), parameter, public var_zi
Definition: ufo_variables_mod.F90:30
ufo_atmvertinterp_tlad_mod::destructor
subroutine destructor(self)
Definition: ufo_atmvertinterp_tlad_mod.F90:203
ufo_constants_mod::half
real(kind_real), parameter, public half
Definition: ufo_constants_mod.F90:26
ufo_vars_mod::var_ps
character(len=maxvarlen), parameter, public var_ps
Definition: ufo_variables_mod.F90:28
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_vars_mod::var_prsi
character(len=maxvarlen), parameter, public var_prsi
Definition: ufo_variables_mod.F90:26
ufo_constants_mod
Definition: ufo_constants_mod.F90:2
satcolumn_mod::simulate_column_ob
subroutine simulate_column_ob(nlayers_obs, nlayers_model, avgkernel_obs, prsl_obs, prsl_model, temp_model, zi_model, profile_model, hofx, troplev_obs, airmass_tot, airmass_trop)
Definition: ufo_satcolumn_mod.F90:8
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_avgkernel_mod::ufo_avgkernel_simobs
subroutine ufo_avgkernel_simobs(self, geovals_in, obss, nvars, nlocs, hofx)
Definition: ufo_avgkernel_mod.F90:115
ufo_geovals_mod::ufo_geovals_get_var
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Definition: ufo_geovals_mod.F90:128
ufo_vars_mod::var_ts
character(len=maxvarlen), parameter, public var_ts
Definition: ufo_variables_mod.F90:19
ufo_avgkernel_mod::ufo_avgkernel_setup
subroutine ufo_avgkernel_setup(self, f_conf)
Definition: ufo_avgkernel_mod.F90:39
ufo_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
ufo_avgkernel_mod::ufo_avgkernel
Fortran derived type for the observation type.
Definition: ufo_avgkernel_mod.F90:20
ufo_constants_mod::zero
real(kind_real), parameter, public zero
Definition: ufo_constants_mod.F90:24
ufo_geovals_mod::ufo_geoval
type to hold interpolated field for one variable, one observation
Definition: ufo_geovals_mod.F90:40
ufo_vars_mod::var_prs
character(len=maxvarlen), parameter, public var_prs
Definition: ufo_variables_mod.F90:25
ufo_geovals_mod::ufo_geovals_reorderzdir
subroutine, public ufo_geovals_reorderzdir(self, varname, zdir)
Definition: ufo_geovals_mod.F90:334
ufo_geovals_mod::ufo_geovals_copy
subroutine, public ufo_geovals_copy(self, other)
Copy one GeoVaLs object into another.
Definition: ufo_geovals_mod.F90:506
ufo_avgkernel_mod
Fortran module for averaging kernel observation operator.
Definition: ufo_avgkernel_mod.F90:8
ufo_constants_mod::one
real(kind_real), parameter, public one
Definition: ufo_constants_mod.F90:27