UFO
ufo_atmsfcinterp_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 atmsfcinterp observation operator
7 
9 
10  use oops_variables_mod
11  use ufo_vars_mod
12  use kinds
13 
14  implicit none
15  private
16 
17  !> Fortran derived type for the observation type
18  type, public :: ufo_atmsfcinterp
19  private
20  type(oops_variables), public :: obsvars
21  type(oops_variables), public :: geovars
22  logical :: use_fact10
23  real(kind_real) :: magl
24  contains
25  procedure :: setup => atmsfcinterp_setup_
26  procedure :: simobs => atmsfcinterp_simobs_
27  end type ufo_atmsfcinterp
28 
29 contains
30 
31 ! ------------------------------------------------------------------------------
32 subroutine atmsfcinterp_setup_(self, f_conf)
33  use fckit_configuration_module, only: fckit_configuration
34  implicit none
35  class(ufo_atmsfcinterp), intent(inout) :: self
36  type(fckit_configuration), intent(in) :: f_conf
37  integer :: nvars, ivar, fact10tmp
38 
39  ! check for if we need to look for wind reduction factor
40  self%use_fact10 = .false.
41  fact10tmp = 0
42  if (f_conf%has("use_fact10")) call f_conf%get_or_die("use_fact10",fact10tmp)
43  if (fact10tmp /= 0) then
44  self%use_fact10 = .true.
45  end if
46 
47  !> add geopotential height
48  call self%geovars%push_back(var_z)
49  !> need skin temperature for near-surface interpolations
50  call self%geovars%push_back(var_sfc_t)
51  !> need surface geopotential height to get difference from phi
52  call self%geovars%push_back(var_sfc_z)
53  !> need surface roughness
54  call self%geovars%push_back(var_sfc_rough)
55  !> need surface and atmospheric pressure for potential temperature
56  call self%geovars%push_back(var_ps)
57  call self%geovars%push_back(var_prs)
58  call self%geovars%push_back(var_prsi)
59  call self%geovars%push_back(var_ts)
60  call self%geovars%push_back(var_tv)
61  call self%geovars%push_back(var_q)
62  call self%geovars%push_back(var_u)
63  call self%geovars%push_back(var_v)
64  call self%geovars%push_back(var_sfc_lfrac)
65  if (self%use_fact10) call self%geovars%push_back(var_sfc_fact10)
66  nvars = self%obsvars%nvars()
67  do ivar = 1, nvars
68  call self%geovars%push_back(self%obsvars%variable(ivar))
69  enddo
70 
71 end subroutine atmsfcinterp_setup_
72 
73 ! ------------------------------------------------------------------------------
74 
75 subroutine atmsfcinterp_simobs_(self, geovals, obss, nvars, nlocs, hofx)
81  use obsspace_mod
82  use iso_c_binding
83  implicit none
84  class(ufo_atmsfcinterp), intent(in) :: self
85  integer, intent(in) :: nvars, nlocs
86  type(ufo_geovals), intent(in) :: geovals
87  real(c_double), intent(inout) :: hofx(nvars, nlocs)
88  type(c_ptr), value, intent(in) :: obss
89  type(ufo_geoval), pointer :: phi, hgt, tsfc, roughlen, psfc, prs, prsi, &
90  tsen, tv, q, u, v, landmask, &
91  profile, rad10
92  integer :: ivar, iobs
93  real(kind_real), allocatable :: obselev(:), obshgt(:)
94  real(kind_real), parameter :: minroughlen = 1.0e-4_kind_real
95  character(len=MAXVARLEN) :: geovar
96  real(kind_real) :: thv1, thv2, th1, thg, thvg, rib, V2
97  real(kind_real) :: gzsoz0, gzzoz0
98  real(kind_real) :: redfac, psim, psimz, psih, psihz
99  real(kind_real) :: ttmp1, ttmpg, eg, qg
100  real(kind_real) :: z0, zq0, tvsfc
101  real(kind_real), parameter :: fv = rv/rd - 1.0_kind_real
102  real(kind_real) :: psit, psitz, ust, psiq, psiqz
103  real(kind_real), parameter :: zint0 = 0.01_kind_real ! default roughness over land
104  real(kind_real), parameter :: ka = 2.4e-5_kind_real
105 
106  ! to compute the value near the surface we are going to use
107  ! similarity theory which requires a number of near surface parameters
108  ! that need to be grabbed from GeoVaLs regardless of the observation type
109 
110  ! for low altitudes, we can assume: phi = z
111  ! grabbing geopotential height profile and surface elevation,
112  ! surface temperature, surface roughness length
113  ! surface and profile air pressure
114  ! profiles of tsen, q, u, v, and surface land fraction
115  call ufo_geovals_get_var(geovals, var_z, phi)
116  call ufo_geovals_get_var(geovals, var_sfc_z, hgt)
117  call ufo_geovals_get_var(geovals, var_sfc_t, tsfc)
118  call ufo_geovals_get_var(geovals, var_sfc_rough, roughlen)
119  call ufo_geovals_get_var(geovals, var_ps, psfc)
120  call ufo_geovals_get_var(geovals, var_prs, prs)
121  call ufo_geovals_get_var(geovals, var_prsi, prsi)
122  call ufo_geovals_get_var(geovals, var_ts, tsen)
123  call ufo_geovals_get_var(geovals, var_tv, tv)
124  call ufo_geovals_get_var(geovals, var_q, q)
125  call ufo_geovals_get_var(geovals, var_u, u)
126  call ufo_geovals_get_var(geovals, var_v, v)
127  call ufo_geovals_get_var(geovals, var_sfc_lfrac, landmask)
128  if (self%use_fact10) call ufo_geovals_get_var(geovals, var_sfc_fact10, rad10)
129 
130  ! get station elevation from obs
131  allocate(obselev(nlocs))
132  call obsspace_get_db(obss, "MetaData", "station_elevation", obselev)
133 
134  ! get observation height (above sea level)
135  allocate(obshgt(nlocs))
136  call obsspace_get_db(obss, "MetaData", "height", obshgt)
137 
138  do iobs = 1, nlocs
139  ! minimum roughness length
140  z0 = roughlen%vals(1,iobs)
141  if (z0 < minroughlen) z0 = minroughlen
142  ! roughness length for over water
143  zq0 = zint0
144  if (landmask%vals(1,iobs) < 0.01) zq0 = z0
145 
146  ! get virtual temperature of the ground assuming saturation
147  call gsi_tp_to_qs(tsfc%vals(1,iobs), psfc%vals(1,iobs), eg, qg)
148  tvsfc = tsfc%vals(1,iobs) * (1.0_kind_real + fv * qg)
149 
150  ! get potential temperatures for calculating psi
151  call calc_theta(tv%vals(1,iobs), prs%vals(1,iobs), thv1)
152  call calc_theta(tv%vals(2,iobs), prs%vals(2,iobs), thv2)
153  call calc_theta(tsen%vals(1,iobs), prs%vals(1,iobs), th1)
154  call calc_theta(tsfc%vals(1,iobs), psfc%vals(1,iobs), thg)
155  call calc_theta(tvsfc, psfc%vals(1,iobs), thvg)
156 
157  ! calculate convective velocity
158  call calc_conv_vel_gsi(u%vals(1,iobs), v%vals(1,iobs), thvg, thv1, v2)
159 
160  ! calculate bulk richardson number
161  rib = (grav * phi%vals(1,iobs) / th1) * (thv1 - thvg) / v2
162 
163  gzsoz0 = log(phi%vals(1,iobs)/z0)
164  gzzoz0 = log((obshgt(iobs)-obselev(iobs))/z0)
165 
166  ! calculate parameters regardless of variable
167  call calc_psi_vars_gsi(rib, gzsoz0, gzzoz0, thv1, thv2, v2, th1,&
168  thg, phi%vals(1,iobs), obshgt(iobs)-obselev(iobs),&
169  psim, psih, psimz, psihz)
170  do ivar = 1, nvars
171  ! Get the name of input variable in geovals
172  geovar = self%obsvars%variable(ivar)
173  ! Get profile for this variable from geovals
174  call ufo_geovals_get_var(geovals, geovar, profile)
175 
176  select case(trim(geovar))
177  case("air_temperature", "virtual_temperature")
178  psit = gzsoz0 - psih
179  psitz = gzzoz0 - psihz
180  if (trim(geovar) == "air_temperature") then
181  ttmp1 = th1
182  ttmpg = thg
183  else
184  ttmp1 = thv1
185  ttmpg = thvg
186  end if
187  hofx(ivar,iobs) = (ttmpg + (ttmp1 - ttmpg)*psitz/psit)*(psfc%vals(1,iobs)/1.0e5_kind_real)**rd_over_cp
188  case("eastward_wind", "northward_wind")
189  if (self%use_fact10) then ! use provided fact10 from model
190  hofx(ivar,iobs) = profile%vals(1,iobs) * rad10%vals(1,iobs)
191  else ! compute wind reduction factor
192  call sfc_wind_fact_gsi(u%vals(1,iobs), v%vals(1,iobs), tsen%vals(1,iobs), q%vals(1,iobs),&
193  psfc%vals(1,iobs), prsi%vals(1,iobs), prsi%vals(2,iobs),&
194  tsfc%vals(1,iobs), z0, landmask%vals(1,iobs), redfac)
195  hofx(ivar,iobs) = profile%vals(1,iobs) * redfac
196  end if
197  case("specific_humidity")
198  ust = von_karman * sqrt(v2) / (gzsoz0 - psim)
199  psiq = log(von_karman*ust*phi%vals(1,iobs)/ka + phi%vals(1,iobs) / zq0) - psih
200  psiqz = log(von_karman*ust*(obshgt(iobs)-obselev(iobs))/ka + (obshgt(iobs)-obselev(iobs)) / zq0) - psihz
201  hofx(ivar,iobs) = qg + (q%vals(1,iobs) - qg)*psiqz/psiq
202  end select
203  end do
204  end do
205 
206  deallocate(obshgt,obselev)
207 
208 end subroutine atmsfcinterp_simobs_
209 
210 ! ------------------------------------------------------------------------------
211 
212 end module ufo_atmsfcinterp_mod
ufo_vars_mod::var_sfc_fact10
character(len=maxvarlen), parameter, public var_sfc_fact10
Definition: ufo_variables_mod.F90:73
thermo_utils_mod::gsi_tp_to_qs
subroutine gsi_tp_to_qs(t, p, es, qs)
Definition: thermo_utils.F90:35
ufo_constants_mod::grav
real(kind_real), parameter, public grav
Definition: ufo_constants_mod.F90:9
ufo_constants_mod::rv
real(kind_real), parameter, public rv
Definition: ufo_constants_mod.F90:13
ufo_vars_mod::var_v
character(len=maxvarlen), parameter, public var_v
Definition: ufo_variables_mod.F90:24
ufo_atmsfcinterp_mod
Fortran module for atmsfcinterp observation operator.
Definition: ufo_atmsfcinterp_mod.F90:8
atmsfc_mod::calc_conv_vel_gsi
subroutine calc_conv_vel_gsi(u1, v1, thvg, thv1, V2)
Definition: ufo_atmsfc_mod.F90:243
ufo_atmsfcinterp_mod::atmsfcinterp_setup_
subroutine atmsfcinterp_setup_(self, f_conf)
Definition: ufo_atmsfcinterp_mod.F90:33
ufo_vars_mod::var_sfc_t
character(len=maxvarlen), parameter, public var_sfc_t
Definition: ufo_variables_mod.F90:72
atmsfc_mod
Definition: ufo_atmsfc_mod.F90:1
atmsfc_mod::calc_psi_vars_gsi
subroutine calc_psi_vars_gsi(rib, gzsoz0, gzzoz0, thv1, thv2, V2, th1, thg, phi1, obshgt, psim, psih, psimz, psihz)
Definition: ufo_atmsfc_mod.F90:162
ufo_vars_mod::var_sfc_z
character(len=maxvarlen), parameter, public var_sfc_z
Definition: ufo_variables_mod.F90:31
ufo_vars_mod::var_ps
character(len=maxvarlen), parameter, public var_ps
Definition: ufo_variables_mod.F90:28
ufo_vars_mod::var_sfc_rough
character(len=maxvarlen), parameter, public var_sfc_rough
Definition: ufo_variables_mod.F90:71
ufo_atmsfcinterp_mod::atmsfcinterp_simobs_
subroutine atmsfcinterp_simobs_(self, geovals, obss, nvars, nlocs, hofx)
Definition: ufo_atmsfcinterp_mod.F90:76
ufo_vars_mod::var_tv
character(len=maxvarlen), parameter, public var_tv
Definition: ufo_variables_mod.F90:18
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
ufo_vars_mod::var_u
character(len=maxvarlen), parameter, public var_u
Definition: ufo_variables_mod.F90:23
ufo_vars_mod::var_sfc_lfrac
character(len=maxvarlen), parameter, public var_sfc_lfrac
Definition: ufo_variables_mod.F90:52
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_vars_mod::var_z
character(len=maxvarlen), parameter, public var_z
Definition: ufo_variables_mod.F90:29
thermo_utils_mod
Fortran module for thermodynamic computations and conversions for use in UFO.
Definition: thermo_utils.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_vars_mod::var_q
character(len=maxvarlen), parameter, public var_q
Definition: ufo_variables_mod.F90:22
ufo_vars_mod::var_ts
character(len=maxvarlen), parameter, public var_ts
Definition: ufo_variables_mod.F90:19
ufo_atmsfcinterp_mod::ufo_atmsfcinterp
Fortran derived type for the observation type.
Definition: ufo_atmsfcinterp_mod.F90:18
ufo_constants_mod::rd_over_cp
real(kind_real), parameter, public rd_over_cp
Definition: ufo_constants_mod.F90:19
ufo_constants_mod::von_karman
real(kind_real), parameter, public von_karman
Definition: ufo_constants_mod.F90:47
ufo_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
ufo_constants_mod::rd
real(kind_real), parameter, public rd
Definition: ufo_constants_mod.F90:12
thermo_utils_mod::calc_theta
subroutine calc_theta(t_in, p_in, t_out)
Definition: thermo_utils.F90:20
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
atmsfc_mod::sfc_wind_fact_gsi
subroutine sfc_wind_fact_gsi(u, v, tsen, q, psfc, prsi1, prsi2, skint, z0, lsmask, f10m)
Definition: ufo_atmsfc_mod.F90:7