UFO
CloudDetectMinResidualIR.cc
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2019 UCAR
3  *
4  * This software is licensed under the terms of the Apache Licence Version 2.0
5  * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6  */
7 
9 
10 #include <cmath>
11 
12 #include <algorithm>
13 #include <iomanip>
14 #include <iostream>
15 #include <set>
16 #include <string>
17 #include <vector>
18 
19 #include "ioda/ObsDataVector.h"
20 #include "oops/util/IntSetParser.h"
21 #include "oops/util/missingValues.h"
22 #include "ufo/filters/Variable.h"
23 #include "ufo/utils/Constants.h"
24 
25 namespace ufo {
26 
27 static ObsFunctionMaker<CloudDetectMinResidualIR>
28  makerCloudDetectMinResidualIR_("CloudDetectMinResidualIR");
29 
30 // -----------------------------------------------------------------------------
31 
32 CloudDetectMinResidualIR::CloudDetectMinResidualIR(const eckit::LocalConfiguration & conf)
33  : invars_() {
34  // Check options
35  options_.deserialize(conf);
36 
37  // Get channels from options
38  std::set<int> channelset = oops::parseIntSet(options_.channelList);
39  std::copy(channelset.begin(), channelset.end(), std::back_inserter(channels_));
40  ASSERT(channels_.size() > 0);
41 
42  // Get test groups from options
43  const std::string &flaggrp = options_.testQCflag.value();
44  const std::string &errgrp = options_.testObserr.value();
45  const std::string &biasgrp = options_.testBias.value();
46  const std::string &hofxgrp = options_.testHofX.value();
47 
48  // Include required variables from ObsDiag
49  invars_ += Variable("brightness_temperature_jacobian_surface_temperature@ObsDiag", channels_);
50  invars_ += Variable("brightness_temperature_jacobian_air_temperature@ObsDiag", channels_);
51  invars_ += Variable("transmittances_of_atmosphere_layer@ObsDiag", channels_);
52  invars_ += Variable("pressure_level_at_peak_of_weightingfunction@ObsDiag", channels_);
53 
54  // Include list of required data from ObsSpace
55  invars_ += Variable("brightness_temperature@"+flaggrp, channels_);
56  invars_ += Variable("brightness_temperature@"+errgrp, channels_);
57  invars_ += Variable("brightness_temperature@"+biasgrp, channels_);
58  invars_ += Variable("brightness_temperature@"+hofxgrp, channels_);
59  invars_ += Variable("brightness_temperature@ObsValue", channels_);
60  invars_ += Variable("brightness_temperature@ObsError", channels_);
61 
62  // Include list of required data from GeoVaLs
63  invars_ += Variable("water_area_fraction@GeoVaLs");
64  invars_ += Variable("land_area_fraction@GeoVaLs");
65  invars_ += Variable("ice_area_fraction@GeoVaLs");
66  invars_ += Variable("surface_snow_area_fraction@GeoVaLs");
67  invars_ += Variable("average_surface_temperature_within_field_of_view@GeoVaLs");
68  invars_ += Variable("air_pressure@GeoVaLs");
69  invars_ += Variable("air_temperature@GeoVaLs");
70  invars_ += Variable("tropopause_pressure@GeoVaLs");
71 }
72 
73 // -----------------------------------------------------------------------------
74 
76 
77 // -----------------------------------------------------------------------------
78 
80  ioda::ObsDataVector<float> & out) const {
81  // Get channel use flags from options
82  std::vector<int> use_flag = options_.useflagChannel.value();
83  std::vector<int> use_flag_clddet = options_.useflagCloudDetect.value();
84 
85  // Get tuning parameters for surface sensitivity over sea/land/oce/snow/mixed from options
86  std::vector<float> dtempf_in = options_.obserrScaleFactorTsfc.value();
87 
88  // Get dimensions
89  size_t nlocs = in.nlocs();
90  size_t nchans = channels_.size();
91  size_t nlevs = in.nlevs(Variable("air_pressure@GeoVaLs"));
92 
93  // Get test groups from options
94  const std::string &flaggrp = options_.testQCflag.value();
95  const std::string &errgrp = options_.testObserr.value();
96  const std::string &biasgrp = options_.testBias.value();
97  const std::string &hofxgrp = options_.testHofX.value();
98 
99  // Get variables from ObsDiag
100  // Load surface temperature jacobian
101  std::vector<std::vector<float>> dbtdts(nchans, std::vector<float>(nlocs));
102  for (size_t ichan = 0; ichan < nchans; ++ichan) {
103  in.get(Variable("brightness_temperature_jacobian_surface_temperature@ObsDiag",
104  channels_)[ichan], dbtdts[ichan]);
105  }
106 
107  // Get temperature jacobian
108  std::vector<std::vector<std::vector<float>>>
109  dbtdt(nchans, std::vector<std::vector<float>>(nlevs, std::vector<float>(nlocs)));
110  for (size_t ichan = 0; ichan < nchans; ++ichan) {
111  for (size_t ilev = 0; ilev < nlevs; ++ilev) {
112  int level = nlevs - ilev;
113  in.get(Variable("brightness_temperature_jacobian_air_temperature@ObsDiag",
114  channels_)[ichan], level, dbtdt[ichan][ilev]);
115  }
116  }
117 
118  // Get layer-to-space transmittance
119  std::vector<std::vector<std::vector<float>>>
120  tao(nchans, std::vector<std::vector<float>>(nlevs, std::vector<float>(nlocs)));
121  for (size_t ichan = 0; ichan < nchans; ++ichan) {
122  for (size_t ilev = 0; ilev < nlevs; ++ilev) {
123  int level = nlevs - ilev;
124  in.get(Variable("transmittances_of_atmosphere_layer@ObsDiag",
125  channels_)[ichan], level, tao[ichan][ilev]);
126  }
127  }
128 
129  // Get pressure level at the peak of the weighting function
130  std::vector<float> values(nlocs, 0.0);
131  std::vector<std::vector<float>> wfunc_pmaxlev(nchans, std::vector<float>(nlocs));
132  for (size_t ichan = 0; ichan < nchans; ++ichan) {
133  in.get(Variable("pressure_level_at_peak_of_weightingfunction@ObsDiag",
134  channels_)[ichan], values);
135  for (size_t iloc = 0; iloc < nlocs; ++iloc) {
136  wfunc_pmaxlev[ichan][iloc] = nlevs - values[iloc] + 1;
137  }
138  }
139 
140  // Get variables from ObsSpace
141  // Get effective observation error and convert it to inverse of the error variance
142  const float missing = util::missingValue(missing);
143  std::vector<int> qcflag(nlocs, 0);
144  std::vector<std::vector<float>> varinv_use(nchans, std::vector<float>(nlocs, 0.0));
145  for (size_t ichan = 0; ichan < nchans; ++ichan) {
146  in.get(Variable("brightness_temperature@"+errgrp, channels_)[ichan], values);
147  in.get(Variable("brightness_temperature@"+flaggrp, channels_)[ichan], qcflag);
148  for (size_t iloc = 0; iloc < nlocs; ++iloc) {
149  if (flaggrp == "PreQC") values[iloc] == missing ? qcflag[iloc] = 100 : qcflag[iloc] = 0;
150  (qcflag[iloc] == 0) ? (values[iloc] = 1.0 / pow(values[iloc], 2)) : (values[iloc] = 0.0);
151  if (use_flag_clddet[ichan] > 0) varinv_use[ichan][iloc] = values[iloc];
152  }
153  }
154 
155  // Get bias corrected innovation (tbobs - hofx - bias)
156  std::vector<std::vector<float>> innovation(nchans, std::vector<float>(nlocs));
157  for (size_t ichan = 0; ichan < nchans; ++ichan) {
158  in.get(Variable("brightness_temperature@ObsValue", channels_)[ichan], innovation[ichan]);
159  in.get(Variable("brightness_temperature@"+hofxgrp, channels_)[ichan], values);
160  for (size_t iloc = 0; iloc < nlocs; ++iloc) {
161  innovation[ichan][iloc] = innovation[ichan][iloc] - values[iloc];
162  }
163  in.get(Variable("brightness_temperature@"+biasgrp, channels_)[ichan], values);
164  for (size_t iloc = 0; iloc < nlocs; ++iloc) {
165  innovation[ichan][iloc] = innovation[ichan][iloc] - values[iloc];
166  }
167  }
168 
169  // Get original observation error (uninflated) from ObsSpaec
170  std::vector<std::vector<float>> obserr(nchans, std::vector<float>(nlocs));
171  for (size_t ichan = 0; ichan < nchans; ++ichan) {
172  in.get(Variable("brightness_temperature@ObsError", channels_)[ichan], obserr[ichan]);
173  }
174 
175  // Get variables from GeoVaLS
176  // Get tropopause pressure [Pa]
177  std::vector<float> tropprs(nlocs);
178  in.get(Variable("tropopause_pressure@GeoVaLs"), tropprs);
179 
180  // Get average surface temperature within FOV
181  std::vector<float> tsavg(nlocs);
182  in.get(Variable("average_surface_temperature_within_field_of_view@GeoVaLs"), tsavg);
183 
184  // Get area fraction of each surface type
185  std::vector<float> water_frac(nlocs);
186  std::vector<float> land_frac(nlocs);
187  std::vector<float> ice_frac(nlocs);
188  std::vector<float> snow_frac(nlocs);
189  in.get(Variable("water_area_fraction@GeoVaLs"), water_frac);
190  in.get(Variable("land_area_fraction@GeoVaLs"), land_frac);
191  in.get(Variable("ice_area_fraction@GeoVaLs"), ice_frac);
192  in.get(Variable("surface_snow_area_fraction@GeoVaLs"), snow_frac);
193 
194  // Determine dominant surface type in each FOV
195  std::vector<bool> land(nlocs, false);
196  std::vector<bool> sea(nlocs, false);
197  std::vector<bool> ice(nlocs, false);
198  std::vector<bool> snow(nlocs, false);
199  std::vector<bool> mixed(nlocs, false);
200  for (size_t iloc = 0; iloc < nlocs; ++iloc) {
201  sea[iloc] = water_frac[iloc] >= 0.99;
202  land[iloc] = land_frac[iloc] >= 0.99;
203  ice[iloc] = ice_frac[iloc] >= 0.99;
204  snow[iloc] = snow_frac[iloc] >= 0.99;
205  mixed[iloc] = (!sea[iloc] && !land[iloc] && !ice[iloc] && !snow[iloc]);
206  }
207 
208  // Setup weight given to each surface type
209  std::vector<float> dtempf(nlocs);
210  for (size_t iloc = 0; iloc < nlocs; ++iloc) {
211  if (sea[iloc]) {
212  dtempf[iloc] = dtempf_in[0];
213  } else if (land[iloc]) {
214  dtempf[iloc] = dtempf_in[1];
215  } else if (ice[iloc]) {
216  dtempf[iloc] = dtempf_in[2];
217  } else if (snow[iloc]) {
218  dtempf[iloc] = dtempf_in[3];
219  } else {
220  dtempf[iloc] = dtempf_in[4];
221  }
222  }
223 
224  // Get air pressure [Pa]
225  std::vector<std::vector<float>> prsl(nlevs, std::vector<float>(nlocs));
226  for (size_t ilev = 0; ilev < nlevs; ++ilev) {
227  size_t level = nlevs - ilev;
228  in.get(Variable("air_pressure@GeoVaLs"), level, prsl[ilev]);
229  }
230 
231  // Get air temperature
232  std::vector<std::vector<float>> tair(nlevs, std::vector<float>(nlocs));
233  for (size_t ilev = 0; ilev < nlevs; ++ilev) {
234  size_t level = nlevs - ilev;
235  in.get(Variable("air_temperature@GeoVaLs"), level, tair[ilev]);
236  }
237 
238  // Minimum Residual Method (MRM) for Cloud Detection:
239  // Determine model level index of the cloud top (lcloud)
240  // Find pressure of the cloud top (cldprs)
241  // Estimate cloud fraction (cldfrac)
242  // output: out = 0 clear channel
243  // out = 1 cloudy channel
244  // out = 2 clear channel but too sensitive to surface condition
245 
246  // Set vectors to hold cloud infomation from cloud detection (cab be part of the output)
247  // std::vector<float> cloud_prsl(nlocs);
248  // std::vector<float> cloud_frac(nlocs);
249  // std::vector<int> cloud_lev(nlocs);
250 
251  // Loop through locations
252  const float btmax = 550.f, btmin = 50.f;
253  for (size_t iloc=0; iloc < nlocs; ++iloc) {
254  // Initialize at each location
255  // cloud_lev[iloc] = 0;
256  // cloud_prsl[iloc] = 0.0;
257  // cloud_frac[iloc] = 0.0;
258  float sum = 0.0, sum2 = 0.0, sum3 = 0.0;
259  float tmp = 0.0;
260  float cloudp = 0.0;
261  std::vector<float> dbt(nchans);
262  for (size_t ichan=0; ichan < nchans; ++ichan) {
263  if (varinv_use[ichan][iloc] > 0) {
264  sum3 = sum3 + innovation[ichan][iloc] * innovation[ichan][iloc] * varinv_use[ichan][iloc];
265  }
266  }
267  sum3 = 0.75 * sum3;
268 
269  // Set initial cloud condition
270  int lcloud = 0;
271  float cldfrac = 0.0;
272  float cldprs = prsl[0][iloc] * 0.01; // convert from [Pa] to [hPa]
273 
274  // Loop through vertical layer from surface to model top
275  for (size_t k = 0 ; k < nlevs ; ++k) {
276  // Perform cloud detection within troposphere
277  if (prsl[k][iloc] * 0.01 > tropprs[iloc] * 0.01) {
278  for (size_t ichan = 0; ichan < nchans; ++ichan) {
279  dbt[ichan] = (tair[k][iloc] - tsavg[iloc]) * dbtdts[ichan][iloc];
280  }
281  for (size_t kk = 0; kk < k; ++kk) {
282  for (size_t ichan = 0; ichan < nchans; ++ichan) {
283  dbt[ichan] = dbt[ichan] + (tair[k][iloc] - tair[kk][iloc]) * dbtdt[ichan][kk][iloc];
284  }
285  }
286  sum = 0.0;
287  sum2 = 0.0;
288  for (size_t ichan = 0; ichan < nchans; ++ichan) {
289  if (varinv_use[ichan][iloc] > 0.0) {
290  sum = sum + innovation[ichan][iloc] * dbt[ichan] * varinv_use[ichan][iloc];
291  sum2 = sum2 + dbt[ichan] * dbt[ichan] * varinv_use[ichan][iloc];
292  }
293  }
294  cloudp = std::min(std::max((sum/sum2), 0.f), 1.f);
295  sum = 0.0;
296  for (size_t ichan = 0; ichan < nchans; ++ichan) {
297  if (varinv_use[ichan][iloc] > 0.0) {
298  tmp = innovation[ichan][iloc] - cloudp * dbt[ichan];
299  sum = sum + tmp * tmp * varinv_use[ichan][iloc];
300  }
301  }
302  if (sum < sum3) {
303  sum3 = sum;
304  lcloud = k + 1; // array index + 1 -> model coordinate index
305  cldfrac = cloudp;
306  cldprs = prsl[k][iloc] * 0.01;
307  }
308  }
309  // end of vertical loop
310  }
311  // If more than 2% of the transmittance comes from the cloud layer,
312  // reject the channel (marked as cloudy channel)
313  float tao_cld = -999.0;
314  for (size_t ichan = 0; ichan < nchans; ++ichan) out[ichan][iloc] = 0;
315  if (lcloud > 0) {
316  for (size_t ichan = 0; ichan < nchans; ++ichan) {
317  // Get cloud top transmittance
318  tao_cld = tao[ichan][lcloud-1][iloc];
319  // Passive channels
320  if (use_flag[ichan] < 0 && lcloud >= wfunc_pmaxlev[ichan][iloc]) out[ichan][iloc] = 1;
321  // Active channels
322  if (out[ichan][iloc] < 1 && tao_cld > 0.02) out[ichan][iloc] = 1;
323  }
324  // cloud infomation output at model level
325  // cloud_lev[iloc] = lcloud;
326  // cloud_prsl[iloc] = cldprs;
327  // cloud_frac[iloc] = cldfrac;
328  } else {
329  // If no clouds is detected, do sensivity to surface temperature check
330  // Initialize at each location
331  float sum = 0.0, sum2 = 0.0;
332  float dts = 0.0;
333  float delta = 0.0;
334  const float dts_threshold = 3.0;
335  for (size_t ichan=0; ichan < nchans; ++ichan) {
336  delta = 0.0;
337  sum = sum + innovation[ichan][iloc] * dbtdts[ichan][iloc] * varinv_use[ichan][iloc];
338  sum2 = sum2 + dbtdts[ichan][iloc] * dbtdts[ichan][iloc] * varinv_use[ichan][iloc];
339  }
340  dts = std::fabs(sum / sum2);
341  if (std::abs(dts) > 1.0) {
342  if (sea[iloc] == false) {
343  dts = std::min(dtempf[iloc], dts);
344  } else {
345  dts = std::min(dts_threshold, dts);
346  }
347  for (size_t ichan=0; ichan < nchans; ++ichan) {
348  delta = std::max(0.05 * obserr[ichan][iloc], 0.02);
349  if (std::abs(dts * dbtdts[ichan][iloc]) > delta) out[ichan][iloc] = 2;
350  }
351  }
352  }
353  // end of location loop
354  }
355 }
356 
357 // -----------------------------------------------------------------------------
358 
360  return invars_;
361 }
362 
363 // -----------------------------------------------------------------------------
364 
365 } // namespace ufo
ufo::ObsFilterData::nlocs
size_t nlocs() const
Returns number of locations.
Definition: ObsFilterData.cc:66
ufo::Variables
Definition: src/ufo/filters/Variables.h:24
CloudDetectMinResidualIR.h
ufo::CloudDetectMinResidualIR::options_
CloudDetectMinResidualIRParameters options_
Definition: CloudDetectMinResidualIR.h:75
ufo
Definition: RunCRTM.h:27
ufo::CloudDetectMinResidualIR::CloudDetectMinResidualIR
CloudDetectMinResidualIR(const eckit::LocalConfiguration &)
Definition: CloudDetectMinResidualIR.cc:32
ufo::ObsFilterData::nlevs
size_t nlevs(const Variable &) const
Returns number of levels for specified variable if 3D GeoVaLs or ObsDiags.
Definition: ObsFilterData.cc:324
ufo::CloudDetectMinResidualIRParameters::testQCflag
oops::Parameter< std::string > testQCflag
Name of the data group to which the QC flag is applied (default is QCflagsData)
Definition: CloudDetectMinResidualIR.h:55
ufo::QCflags::missing
constexpr int missing
Definition: QCflags.h:15
ufo::CloudDetectMinResidualIRParameters::testHofX
oops::Parameter< std::string > testHofX
Name of the HofX group used to replace the default group (default is HofX)
Definition: CloudDetectMinResidualIR.h:49
ufo::CloudDetectMinResidualIRParameters::testBias
oops::Parameter< std::string > testBias
Name of the bias correction group used to replace the default group (default is ObsBias)
Definition: CloudDetectMinResidualIR.h:52
ufo::abs
util::Duration abs(const util::Duration &duration)
Definition: TrackCheckUtils.h:31
ufo::CloudDetectMinResidualIRParameters::useflagChannel
oops::RequiredParameter< std::vector< int > > useflagChannel
Useflag (-1: not used; 0: monitoring; 1: used) for each channel in channelList.
Definition: CloudDetectMinResidualIR.h:36
ufo_constants_mod::sea
real(kind_real), parameter, public sea
Definition: ufo_constants_mod.F90:34
ioda::ObsDataVector
Definition: BackgroundCheck.h:26
ufo::CloudDetectMinResidualIR::~CloudDetectMinResidualIR
~CloudDetectMinResidualIR()
Definition: CloudDetectMinResidualIR.cc:75
ufo::makerCloudDetectMinResidualIR_
static ObsFunctionMaker< CloudDetectMinResidualIR > makerCloudDetectMinResidualIR_("CloudDetectMinResidualIR")
ufo::CloudDetectMinResidualIR::invars_
ufo::Variables invars_
Definition: CloudDetectMinResidualIR.h:76
ufo::CloudDetectMinResidualIRParameters::obserrScaleFactorTsfc
oops::RequiredParameter< std::vector< float > > obserrScaleFactorTsfc
Definition: CloudDetectMinResidualIR.h:43
ufo::CloudDetectMinResidualIRParameters::testObserr
oops::Parameter< std::string > testObserr
Name of the data group to which the observation error is applied (default: ObsErrorData)
Definition: CloudDetectMinResidualIR.h:46
ufo::CloudDetectMinResidualIR::channels_
std::vector< int > channels_
Definition: CloudDetectMinResidualIR.h:77
Constants.h
ufo::ObsFilterData::get
void get(const Variable &, std::vector< float > &) const
Gets requested data from ObsFilterData.
Definition: ObsFilterData.cc:130
ufo::CloudDetectMinResidualIRParameters::useflagCloudDetect
oops::RequiredParameter< std::vector< int > > useflagCloudDetect
Useflag (-1: not used; 1: used) indicating if the channel is used for cloud detection.
Definition: CloudDetectMinResidualIR.h:39
ufo::Variable
Definition: Variable.h:23
ufo::ObsFilterData
ObsFilterData provides access to all data related to an ObsFilter.
Definition: src/ufo/filters/ObsFilterData.h:40
ufo::CloudDetectMinResidualIR::compute
void compute(const ObsFilterData &, ioda::ObsDataVector< float > &) const
compute the result of the function
Definition: CloudDetectMinResidualIR.cc:79
conf
Definition: conf.py:1
ufo::CloudDetectMinResidualIR::requiredVariables
const ufo::Variables & requiredVariables() const
geovals required to compute the function
Definition: CloudDetectMinResidualIR.cc:359
ufo::CloudDetectMinResidualIRParameters::channelList
oops::RequiredParameter< std::string > channelList
List of channels available for assimilation.
Definition: CloudDetectMinResidualIR.h:33
Variable.h