71 : public edm::global::EDAnalyzer<edm::
RunCache<gxsec::
RunCache>, edm::LuminosityBlockCache<gxsec::LumiCache>> {
89 void combine(
double &,
double &,
double &,
const double &,
const double &,
const double &)
const;
138 xsecPreFilter_(-1, -1),
140 filterOnlyEffStat_(0, 0, 0, 0, 0., 0., 0., 0.),
141 hepMCFilterEffStat_(0, 0, 0, 0, 0., 0., 0., 0.) {
142 genFilterInfoToken_ = consumes<GenFilterInfo, edm::InLumi>(
edm::InputTag(
"genFilterEfficiencyProducer",
""));
143 hepMCFilterInfoToken_ = consumes<GenFilterInfo, edm::InLumi>(
edm::InputTag(
"generator",
""));
144 genLumiInfoToken_ = consumes<GenLumiInfoProduct, edm::InLumi>(
edm::InputTag(
"generator",
""));
145 lheRunInfoToken_ = consumes<LHERunInfoProduct, edm::InRun>(
edm::InputTag(
"externalLHEProducer",
""));
158 std::lock_guard
l{mutex_};
159 xsecBeforeMatching_.clear();
160 xsecAfterMatching_.clear();
161 jetMatchEffStat_.clear();
163 return std::make_shared<gxsec::RunCache>();
168 return std::shared_ptr<gxsec::LumiCache>();
175 iLumi.
getByToken(genLumiInfoToken_, genLumiInfo);
178 hepidwtup_ = genLumiInfo->getHEPIDWTUP();
180 std::vector<GenLumiInfoProduct::ProcessInfo>
const &theProcesses = genLumiInfo->getProcessInfos();
182 unsigned int theProcesses_size = theProcesses.size();
186 if (hepidwtup_ == -1) {
187 if (theProcesses_size != 1) {
188 edm::LogError(
"GenXSecAnalyzer::endLuminosityBlock") <<
"Pure parton shower has thisProcessInfos size!=1";
193 for (
unsigned int ip = 0; ip < theProcesses_size; ip++) {
194 if (theProcesses[ip].lheXSec().
value() < 0) {
196 <<
"cross section of process " << ip <<
" value = " << theProcesses[ip].lheXSec().value();
203 std::lock_guard
g{mutex_};
204 runC->product_.mergeProduct(*genLumiInfo);
207 iLumi.
getByToken(genFilterInfoToken_, genFilter);
210 std::lock_guard
g{mutex_};
211 filterOnlyEffStat_.mergeProduct(*genFilter);
212 runC->filterOnlyEffRun_.mergeProduct(*genFilter);
213 runC->thisRunWeight_ += genFilter->sumPassWeights();
217 iLumi.
getByToken(hepMCFilterInfoToken_, hepMCFilter);
220 std::lock_guard
g{mutex_};
221 hepMCFilterEffStat_.mergeProduct(*hepMCFilter);
222 runC->hepMCFilterEffRun_.mergeProduct(*hepMCFilter);
225 std::lock_guard
g{mutex_};
228 for (
unsigned int ip = 0; ip < theProcesses_size; ip++) {
229 int id = theProcesses[ip].process();
234 double passw = temp_killed.
sum();
235 double passw2 = temp_killed.
sum2();
236 double totalw = temp_selected.
sum();
237 double totalw2 = temp_selected.
sum2();
239 theProcesses[ip].nPassNeg(),
240 theProcesses[ip].nTotalPos(),
241 theProcesses[ip].nTotalNeg(),
248 jetMatchEffStat_[10000].mergeProduct(tempInfo);
249 double currentValue = theProcesses[ip].lheXSec().value();
250 double currentError = theProcesses[ip].lheXSec().error();
253 auto &thisRunWeightPre = runC->thisRunWeightPre_;
256 double previousValue = runC->previousLumiBlockLHEXSec_[
id].value();
258 if (currentValue != previousValue)
260 double xsec = y.
value();
262 combine(xsec, err, thisRunWeightPre, currentValue, currentError, totalw);
265 thisRunWeightPre += totalw;
271 y = theProcesses[ip].lheXSec();
272 thisRunWeightPre += totalw;
275 runC->previousLumiBlockLHEXSec_[
id] = theProcesses[ip].lheXSec();
288 for (
unsigned int iSize = 0; iSize < thisHeprup.
XSECUP.size(); iSize++) {
289 std::cout << std::setw(14) << std::fixed << thisHeprup.
XSECUP[iSize] << std::setw(14) << std::fixed
290 << thisHeprup.
XERRUP[iSize] << std::setw(14) << std::fixed << thisHeprup.
XMAXUP[iSize] << std::setw(14)
291 << std::fixed << thisHeprup.
LPRUP[iSize] << std::endl;
296 auto runC = runCache(iRun.
index());
297 std::lock_guard
l{mutex_};
302 std::vector<GenLumiInfoProduct::ProcessInfo> newInfos;
303 for (std::map<int, GenLumiInfoProduct::XSec>::const_iterator iter = runC->currentLumiBlockLHEXSec_.begin();
304 iter != runC->currentLumiBlockLHEXSec_.end();
307 temp.
setLheXSec(iter->second.value(), iter->second.error());
308 newInfos.push_back(temp);
310 runC->product_.setProcessInfo(newInfos);
314 combine(xsecPreFilter_, totalWeightPre_, thisRunXSecPre, runC->thisRunWeightPre_);
316 double thisHepFilterEff = 1;
317 double thisHepFilterErr = 0;
319 if (runC->hepMCFilterEffRun_.sumWeights2() > 0) {
320 thisHepFilterEff = runC->hepMCFilterEffRun_.filterEfficiency(hepidwtup_);
321 thisHepFilterErr = runC->hepMCFilterEffRun_.filterEfficiencyError(hepidwtup_);
322 if (thisHepFilterEff < 0) {
323 thisHepFilterEff = 1;
324 thisHepFilterErr = 0;
328 double thisGenFilterEff = 1;
329 double thisGenFilterErr = 0;
331 if (runC->filterOnlyEffRun_.sumWeights2() > 0) {
332 thisGenFilterEff = runC->filterOnlyEffRun_.filterEfficiency(hepidwtup_);
333 thisGenFilterErr = runC->filterOnlyEffRun_.filterEfficiencyError(hepidwtup_);
334 if (thisGenFilterEff < 0) {
335 thisGenFilterEff = 1;
336 thisGenFilterErr = 0;
339 double thisXsec = thisRunXSecPre.
value() > 0 ? thisHepFilterEff * thisGenFilterEff * thisRunXSecPre.
value() : 0;
341 thisRunXSecPre.
value() > 0
343 pow(thisHepFilterErr / thisHepFilterEff, 2) +
pow(thisGenFilterErr / thisGenFilterEff, 2))
346 combine(xsec_, totalWeight_, thisRunXSec, runC->thisRunWeight_);
352 const double ¤tValue,
353 const double ¤tError,
354 const double ¤tWeight)
const {
355 if (finalValue <= 0) {
356 finalValue = currentValue;
357 finalError = currentError;
358 finalWeight += currentWeight;
360 double wgt1 = (finalError <= 0 || currentError <= 0) ? finalWeight : 1 / (finalError * finalError);
361 double wgt2 = (finalError <= 0 || currentError <= 0) ? currentWeight : 1 / (currentError * currentError);
362 double xsec = (wgt1 * finalValue + wgt2 * currentValue) / (wgt1 + wgt2);
363 double err = (finalError <= 0 || currentError <= 0) ? 0 : 1.0 /
std::sqrt(wgt1 + wgt2);
366 finalWeight += currentWeight;
374 const double &thisw)
const {
377 double thisValue = thisRunXSec.
value();
378 double thisError = thisRunXSec.
error();
379 combine(value, error, totalw, thisValue, thisError, thisw);
386 double sigSelSum = 0.0;
387 double err2SelSum = 0.0;
389 double err2Sum = 0.0;
391 std::vector<GenLumiInfoProduct::XSec> tempVector_before;
392 std::vector<GenLumiInfoProduct::XSec> tempVector_after;
396 for (
unsigned int ip = 0; ip < vectorSize; ip++) {
402 sigSelSum += hepxsec_value;
403 err2SelSum += hepxsec_error * hepxsec_error;
415 switch (hepidwtup_) {
418 fracAcc = ntotal > 0 ? npass / ntotal : -1;
431 double sigmaFin = hepxsec_value * fracAcc;
436 switch (hepidwtup_) {
440 double effp = ntotal_pos > 0 ? (double)proc.
nPassPos() / ntotal_pos : 0;
441 double effp_err2 = ntotal_pos > 0 ? (1 - effp) * effp / ntotal_pos : 0;
444 double effn = ntotal_neg > 0 ? (double)proc.
nPassNeg() / ntotal_neg : 0;
445 double effn_err2 = ntotal_neg > 0 ? (1 - effn) * effn / ntotal_neg : 0;
448 ? (ntotal_pos * ntotal_pos * effp_err2 + ntotal_neg * ntotal_neg * effn_err2) / ntotal / ntotal
458 double numerator = (passw2 * failw * failw + failw2 * passw * passw);
460 efferr2 = denominator > 0 ? numerator / denominator : 0;
464 double delta2Veto = efferr2 / fracAcc / fracAcc;
468 double sigma2Sum, sigma2Err;
469 sigma2Sum = hepxsec_value * hepxsec_value;
470 sigma2Err = hepxsec_error * hepxsec_error;
472 double delta2Sum = delta2Veto + sigma2Err / sigma2Sum;
473 relErr = (delta2Sum > 0.0 ?
std::sqrt(delta2Sum) : 0.0);
474 double deltaFin = sigmaFin * relErr;
480 err2Sum += deltaFin * deltaFin;
485 double total_matcheff = jetMatchEffStat_[10000].filterEfficiency(hepidwtup_);
486 double total_matcherr = jetMatchEffStat_[10000].filterEfficiencyError(hepidwtup_);
488 double xsec_after = sigSelSum * total_matcheff;
489 double xsecerr_after = (total_matcheff > 0 && sigSelSum > 0)
490 ? xsec_after *
sqrt(err2SelSum / sigSelSum / sigSelSum +
491 total_matcherr * total_matcherr / total_matcheff / total_matcheff)
495 tempVector_after.push_back(result);
497 xsecBeforeMatching_ = tempVector_before;
498 xsecAfterMatching_ = tempVector_after;
505 <<
"------------------------------------"
507 <<
"GenXsecAnalyzer:"
509 <<
"------------------------------------";
511 if (jetMatchEffStat_.empty()) {
512 edm::LogPrint(
"GenXSecAnalyzer") <<
"------------------------------------"
514 <<
"Cross-section summary not available"
516 <<
"------------------------------------";
521 double final_fract_neg_w = 0;
522 double final_fract_neg_w_unc = 0;
526 if (nMCs_ == 1 && hepidwtup_ != -1) {
528 <<
"-----------------------------------------------------------------------------------------------------------"
529 "--------------------------------------------------------------- \n"
530 <<
"Overall cross-section summary \n"
531 <<
"-----------------------------------------------------------------------------------------------------------"
532 "---------------------------------------------------------------";
533 edm::LogPrint(
"GenXSecAnalyzer") <<
"Process\t\txsec_before [pb]\t\tpassed\tnposw\tnnegw\ttried\tnposw\tnnegw "
534 "\txsec_match [pb]\t\t\taccepted [%]\t event_eff [%]";
536 const unsigned sizeOfInfos = jetMatchEffStat_.size();
537 const unsigned last = sizeOfInfos - 1;
540 double jetmatch_eff = 0;
541 double jetmatch_err = 0;
542 double matching_eff = 1;
543 double matching_efferr = 1;
545 for (std::map<int, GenFilterInfo>::const_iterator iter = jetMatchEffStat_.begin(); iter != jetMatchEffStat_.end();
565 <<
"-------------------------------------------------------------------------------------------------------"
566 "------------------------------------------------------------------- ";
572 final_fract_neg_w_unc =
574 ? final_fract_neg_w * final_fract_neg_w / thisEventEffStat.
numEventsPassed() *
580 title[
i] = Form(
"%d", i);
583 edm::LogPrint(
"GenXSecAnalyzer") << title[
i] <<
"\t\t" << std::scientific << std::setprecision(3)
584 << xsecBeforeMatching_[
i].value() <<
" +/- " << xsecBeforeMatching_[
i].error()
591 << std::setprecision(3) << xsecAfterMatching_[
i].value() <<
" +/- "
592 << xsecAfterMatching_[
i].error() <<
"\t\t" << std::fixed << std::setprecision(1)
593 << (jetmatch_eff * 100) <<
" +/- " << (jetmatch_err * 100) <<
"\t" << std::fixed
603 <<
"-----------------------------------------------------------------------------------------------------------"
604 "---------------------------------------------------------------";
606 edm::LogPrint(
"GenXSecAnalyzer") <<
"Before matching: total cross section = " << std::scientific
607 << std::setprecision(3) << xsecBeforeMatching_[
last].value() <<
" +- "
608 << xsecBeforeMatching_[
last].error() <<
" pb";
610 edm::LogPrint(
"GenXSecAnalyzer") <<
"After matching: total cross section = " << std::scientific
611 << std::setprecision(3) << xsecAfterMatching_[
last].value() <<
" +- "
612 << xsecAfterMatching_[
last].error() <<
" pb";
614 edm::LogPrint(
"GenXSecAnalyzer") <<
"Matching efficiency = " << std::fixed << std::setprecision(1) << matching_eff
615 <<
" +/- " << matching_efferr <<
" [TO BE USED IN MCM]";
617 }
else if (hepidwtup_ == -1)
618 edm::LogPrint(
"GenXSecAnalyzer") <<
"Before Filter: total cross section = " << std::scientific
619 << std::setprecision(3) << xsecPreFilter_.value() <<
" +- "
620 << xsecPreFilter_.error() <<
" pb";
623 double hepMCFilter_eff = 1.0;
624 double hepMCFilter_err = 0.0;
625 if (hepMCFilterEffStat_.sumWeights2() > 0) {
626 hepMCFilter_eff = hepMCFilterEffStat_.filterEfficiency(-1);
627 hepMCFilter_err = hepMCFilterEffStat_.filterEfficiencyError(-1);
628 edm::LogPrint(
"GenXSecAnalyzer") <<
"HepMC filter efficiency (taking into account weights)= "
629 <<
"(" << hepMCFilterEffStat_.sumPassWeights() <<
")"
631 <<
"(" << hepMCFilterEffStat_.sumWeights() <<
")"
632 <<
" = " << std::scientific << std::setprecision(3) << hepMCFilter_eff <<
" +- "
635 double hepMCFilter_event_total =
636 hepMCFilterEffStat_.numTotalPositiveEvents() + hepMCFilterEffStat_.numTotalNegativeEvents();
637 double hepMCFilter_event_pass =
638 hepMCFilterEffStat_.numPassPositiveEvents() + hepMCFilterEffStat_.numPassNegativeEvents();
639 double hepMCFilter_event_eff = hepMCFilter_event_total > 0 ? hepMCFilter_event_pass / hepMCFilter_event_total : 0;
640 double hepMCFilter_event_err =
641 hepMCFilter_event_total > 0
642 ?
sqrt((1 - hepMCFilter_event_eff) * hepMCFilter_event_eff / hepMCFilter_event_total)
644 edm::LogPrint(
"GenXSecAnalyzer") <<
"HepMC filter efficiency (event-level)= "
645 <<
"(" << hepMCFilter_event_pass <<
")"
647 <<
"(" << hepMCFilter_event_total <<
")"
648 <<
" = " << std::scientific << std::setprecision(3) << hepMCFilter_event_eff
649 <<
" +- " << hepMCFilter_event_err;
653 if (filterOnlyEffStat_.sumWeights2() > 0) {
654 double filterOnly_eff = filterOnlyEffStat_.filterEfficiency(-1);
655 double filterOnly_err = filterOnlyEffStat_.filterEfficiencyError(-1);
657 edm::LogPrint(
"GenXSecAnalyzer") <<
"Filter efficiency (taking into account weights)= "
658 <<
"(" << filterOnlyEffStat_.sumPassWeights() <<
")"
660 <<
"(" << filterOnlyEffStat_.sumWeights() <<
")"
661 <<
" = " << std::scientific << std::setprecision(3) << filterOnly_eff <<
" +- "
664 double filterOnly_event_total =
665 filterOnlyEffStat_.numTotalPositiveEvents() + filterOnlyEffStat_.numTotalNegativeEvents();
666 double filterOnly_event_pass =
667 filterOnlyEffStat_.numPassPositiveEvents() + filterOnlyEffStat_.numPassNegativeEvents();
668 double filterOnly_event_eff = filterOnly_event_total > 0 ? filterOnly_event_pass / filterOnly_event_total : 0;
669 double filterOnly_event_err = filterOnly_event_total > 0
670 ?
sqrt((1 - filterOnly_event_eff) * filterOnly_event_eff / filterOnly_event_total)
672 edm::LogPrint(
"GenXSecAnalyzer") <<
"Filter efficiency (event-level)= "
673 <<
"(" << filterOnly_event_pass <<
")"
675 <<
"(" << filterOnly_event_total <<
")"
676 <<
" = " << std::scientific << std::setprecision(3) << filterOnly_event_eff
677 <<
" +- " << filterOnly_event_err <<
" [TO BE USED IN MCM]";
681 filterOnly_event_pass > 0 ? filterOnlyEffStat_.numPassNegativeEvents() / (filterOnly_event_pass) : 0;
682 final_fract_neg_w_unc =
683 filterOnlyEffStat_.numPassNegativeEvents() > 0
684 ? final_fract_neg_w * final_fract_neg_w / filterOnly_event_pass *
685 sqrt(filterOnlyEffStat_.numPassPositiveEvents() * filterOnlyEffStat_.numPassPositiveEvents() /
686 filterOnlyEffStat_.numPassNegativeEvents() +
687 filterOnlyEffStat_.numPassPositiveEvents())
691 edm::LogPrint(
"GenXSecAnalyzer") <<
"\nAfter filter: final cross section = " << std::scientific
692 << std::setprecision(3) << xsec_.value() <<
" +- " << xsec_.error() <<
" pb";
694 edm::LogPrint(
"GenXSecAnalyzer") <<
"After filter: final fraction of events with negative weights = "
695 << std::scientific << std::setprecision(3) << final_fract_neg_w <<
" +- "
696 << final_fract_neg_w_unc;
699 double lumi_1M_evts =
700 xsec_.value() > 0 ? 1e6 * (1 - 2 * final_fract_neg_w) * (1 - 2 * final_fract_neg_w) / xsec_.value() / 1e3 : 0;
701 double lumi_1M_evts_unc =
702 xsec_.value() > 0 ? (1 - 2 * final_fract_neg_w) * lumi_1M_evts *
703 sqrt(1
e-6 + 16 *
pow(final_fract_neg_w_unc, 2) /
pow(1 - 2 * final_fract_neg_w, 2) +
704 pow(xsec_.error() / xsec_.value(), 2))
706 edm::LogPrint(
"GenXSecAnalyzer") <<
"After filter: final equivalent lumi for 1M events (1/fb) = " << std::scientific
707 << std::setprecision(3) << lumi_1M_evts <<
" +- " << lumi_1M_evts_unc;
std::shared_ptr< gxsec::LumiCache > globalBeginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) const final
const std::vector< ProcessInfo > & getProcessInfos() const
unsigned int nTotalNeg() const
edm::EDGetTokenT< GenFilterInfo > hepMCFilterInfoToken_
uint16_t *__restrict__ id
unsigned int numTotalPositiveEvents() const
GenLumiInfoProduct product_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
double filterEfficiency(int idwtup=+3) const
edm::EDGetTokenT< GenLumiInfoProduct > genLumiInfoToken_
unsigned int nPassPos() const
GenFilterInfo filterOnlyEffRun_
Log< level::Error, false > LogError
void globalEndRun(edm::Run const &, edm::EventSetup const &) const final
GenLumiInfoProduct::XSec compute(const GenLumiInfoProduct &) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
std::atomic< int > hepidwtup_
void globalEndLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) const final
unsigned int nPassNeg() const
bool mergeProduct(GenFilterInfo const &other)
void analyze(edm::StreamID, const edm::Event &, const edm::EventSetup &) const final
XSec const & lheXSec() const
unsigned int numEventsPassed() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
static type combine(const A &_1, const B &_2)
#define CMS_THREAD_GUARD(_var_)
unsigned int numTotalNegativeEvents() const
Log< level::Warning, true > LogPrint
std::vector< double > XERRUP
unsigned int numEventsTotal() const
std::vector< double > XMAXUP
void combine(GenLumiInfoProduct::XSec &, double &, const GenLumiInfoProduct::XSec &, const double &) const
double filterEfficiencyError(int idwtup=+3) const
FinalStat const & killed() const
unsigned int numPassPositiveEvents() const
void setLheXSec(double value, double err)
edm::EDGetTokenT< GenFilterInfo > genFilterInfoToken_
~GenXSecAnalyzer() override
edm::EDGetTokenT< LHERunInfoProduct > lheRunInfoToken_
std::vector< double > XSECUP
unsigned int nTotalPos() const
std::map< int, GenLumiInfoProduct::XSec > currentLumiBlockLHEXSec_
std::shared_ptr< gxsec::RunCache > globalBeginRun(edm::Run const &, edm::EventSetup const &) const final
GenFilterInfo hepMCFilterEffRun_
Power< A, B >::type pow(const A &a, const B &b)
unsigned int numPassNegativeEvents() const
FinalStat const & selected() const
std::map< int, GenLumiInfoProduct::XSec > previousLumiBlockLHEXSec_
Run const & getRun() const