19 #include "ioda/ObsDataVector.h"
20 #include "oops/util/IntSetParser.h"
21 #include "oops/util/missingValues.h"
27 static ObsFunctionMaker<CloudDetectMinResidualIR>
39 std::copy(channelset.begin(), channelset.end(), std::back_inserter(
channels_));
67 invars_ +=
Variable(
"average_surface_temperature_within_field_of_view@GeoVaLs");
89 size_t nlocs = in.
nlocs();
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",
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]);
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]);
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",
135 for (
size_t iloc = 0; iloc < nlocs; ++iloc) {
136 wfunc_pmaxlev[ichan][iloc] = nlevs - values[iloc] + 1;
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) {
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];
156 std::vector<std::vector<float>> innovation(nchans, std::vector<float>(nlocs));
157 for (
size_t ichan = 0; ichan < nchans; ++ichan) {
160 for (
size_t iloc = 0; iloc < nlocs; ++iloc) {
161 innovation[ichan][iloc] = innovation[ichan][iloc] - values[iloc];
164 for (
size_t iloc = 0; iloc < nlocs; ++iloc) {
165 innovation[ichan][iloc] = innovation[ichan][iloc] - values[iloc];
170 std::vector<std::vector<float>> obserr(nchans, std::vector<float>(nlocs));
171 for (
size_t ichan = 0; ichan < nchans; ++ichan) {
177 std::vector<float> tropprs(nlocs);
178 in.
get(
Variable(
"tropopause_pressure@GeoVaLs"), tropprs);
181 std::vector<float> tsavg(nlocs);
182 in.
get(
Variable(
"average_surface_temperature_within_field_of_view@GeoVaLs"), tsavg);
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);
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]);
209 std::vector<float> dtempf(nlocs);
210 for (
size_t iloc = 0; iloc < nlocs; ++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];
220 dtempf[iloc] = dtempf_in[4];
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]);
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]);
252 const float btmax = 550.f, btmin = 50.f;
253 for (
size_t iloc=0; iloc < nlocs; ++iloc) {
258 float sum = 0.0, sum2 = 0.0, sum3 = 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];
272 float cldprs = prsl[0][iloc] * 0.01;
275 for (
size_t k = 0 ; k < nlevs ; ++k) {
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];
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];
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];
294 cloudp = std::min(std::max((sum/sum2), 0.f), 1.f);
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];
306 cldprs = prsl[k][iloc] * 0.01;
313 float tao_cld = -999.0;
314 for (
size_t ichan = 0; ichan < nchans; ++ichan) out[ichan][iloc] = 0;
316 for (
size_t ichan = 0; ichan < nchans; ++ichan) {
318 tao_cld = tao[ichan][lcloud-1][iloc];
320 if (use_flag[ichan] < 0 && lcloud >= wfunc_pmaxlev[ichan][iloc]) out[ichan][iloc] = 1;
322 if (out[ichan][iloc] < 1 && tao_cld > 0.02) out[ichan][iloc] = 1;
331 float sum = 0.0, sum2 = 0.0;
334 const float dts_threshold = 3.0;
335 for (
size_t ichan=0; ichan < nchans; ++ichan) {
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];
340 dts = std::fabs(sum / sum2);
342 if (
sea[iloc] ==
false) {
343 dts = std::min(dtempf[iloc], dts);
345 dts = std::min(dts_threshold, dts);
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;