10 use oops_variables_mod
20 type(oops_variables),
public :: obsvars
21 type(oops_variables),
public :: geovars
23 real(kind_real) :: magl
33 use fckit_configuration_module,
only: fckit_configuration
36 type(fckit_configuration),
intent(in) :: f_conf
37 integer :: nvars, ivar, fact10tmp
40 self%use_fact10 = .false.
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.
48 call self%geovars%push_back(
var_z)
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)
66 nvars = self%obsvars%nvars()
68 call self%geovars%push_back(self%obsvars%variable(ivar))
85 integer,
intent(in) :: nvars, nlocs
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, &
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
104 real(kind_real),
parameter :: ka = 2.4e-5_kind_real
131 allocate(obselev(nlocs))
132 call obsspace_get_db(obss,
"MetaData",
"station_elevation", obselev)
135 allocate(obshgt(nlocs))
136 call obsspace_get_db(obss,
"MetaData",
"height", obshgt)
140 z0 = roughlen%vals(1,iobs)
141 if (z0 < minroughlen) z0 = minroughlen
144 if (landmask%vals(1,iobs) < 0.01) zq0 = z0
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)
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)
161 rib = (
grav * phi%vals(1,iobs) / th1) * (thv1 - thvg) / v2
163 gzsoz0 = log(phi%vals(1,iobs)/z0)
164 gzzoz0 = log((obshgt(iobs)-obselev(iobs))/z0)
168 thg, phi%vals(1,iobs), obshgt(iobs)-obselev(iobs),&
169 psim, psih, psimz, psihz)
172 geovar = self%obsvars%variable(ivar)
176 select case(trim(geovar))
177 case(
"air_temperature",
"virtual_temperature")
179 psitz = gzzoz0 - psihz
180 if (trim(geovar) ==
"air_temperature")
then
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
190 hofx(ivar,iobs) = profile%vals(1,iobs) * rad10%vals(1,iobs)
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
197 case(
"specific_humidity")
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
206 deallocate(obshgt,obselev)