Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

avoid duplicates in HLTPMMassFilter plugin #40918

Merged
merged 1 commit into from Mar 2, 2023

Conversation

missirol
Copy link
Contributor

@missirol missirol commented Mar 1, 2023

PR description:

This PR modifies the plugin HLTPMMassFilter in order to (1) consider a given pair of candidates only once when counting how many pairs pass the selection requirements, and (2) avoid duplicates in the candidates added to the TriggerFilterObjectWithRefs output product of this plugin.

Discussed in CMSHLT-2635.

PR description:

None.

If this PR is a backport, please specify the original PR and why you need to backport that PR. If this PR will be backported, please specify to which release cycle the backport is meant for:

CMSSW_13_0_X

@cmsbuild
Copy link
Contributor

cmsbuild commented Mar 1, 2023

+code-checks

Logs: https://cmssdt.cern.ch/SDT/code-checks/cms-sw-PR-40918/34403

  • This PR adds an extra 16KB to repository

@cmsbuild
Copy link
Contributor

cmsbuild commented Mar 1, 2023

A new Pull Request was created by @missirol (Marino Missiroli) for master.

It involves the following packages:

  • HLTrigger/Egamma (hlt)

@cmsbuild, @missirol, @Martin-Grunewald can you please review it and eventually sign? Thanks.
@Martin-Grunewald, @silviodonato this is something you requested to watch as well.
@perrotta, @dpiparo, @rappoccio you are the release manager for this.

cms-bot commands are listed here

if (mass >= lowerMassCut_ && mass <= upperMassCut_) {
if (isGoodPair(pEleCh1[jj], pEleCh1[ii])) {
n++;
refele = electrons[ii];
filterproduct.addObject(TriggerElectron, refele);
refele = electrons[jj];
filterproduct.addObject(TriggerElectron, refele);
for (auto const idx : {jj, ii}) {
if (save_cand[idx])
filterproduct.addObject(TriggerElectron, electrons[idx]);
save_cand[idx] = false;
}
Copy link
Contributor Author

@missirol missirol Mar 1, 2023

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

In the isElectron1_ && isElectron2_ case, there is no change in how candidate pairs are counted (as the nested loop was already considering only unique pairs). The only change is that the same TriggerElectron is not added more than once to the output product.

Comment on lines +134 to +143
std::vector<bool> save_cand(scs.size(), true);
for (unsigned int jj = 0; jj < scs.size(); jj++) {
TLorentzVector p1Ele = pEleCh1.at(jj);
for (unsigned int ii = 0; ii < scs.size(); ii++) {
TLorentzVector p2Ele = pEleCh2.at(ii);

if (fabs(p1Ele.E() - p2Ele.E()) < 0.00001)
continue;

TLorentzVector pTot = p1Ele + p2Ele;
double mass = pTot.M();

if (mass >= lowerMassCut_ && mass <= upperMassCut_) {
for (unsigned int ii = jj + 1; ii < scs.size(); ii++) {
if (isGoodPair(pEleCh1[jj], pEleCh2[ii]) or isGoodPair(pEleCh2[jj], pEleCh1[ii])) {
n++;
refsc = scs[ii];
filterproduct.addObject(TriggerCluster, refsc);
refsc = scs[jj];
filterproduct.addObject(TriggerCluster, refsc);
for (auto const idx : {jj, ii}) {
if (save_cand[idx])
filterproduct.addObject(TriggerCluster, scs[idx]);
save_cand[idx] = false;
}
Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

In the non-isElectron1_ && isElectron2_ case, there is a change in how candidate pairs are counted, as (1) the case jj == ii is explicitly skipped (supposedly it was effectively skipped pre-PR too, because of the requirement on delta-Energy), and (2) a given pair of TriggerClusters is only counted at most once: a given pair is considered 'good' if it passes the selection under the "plus minus" or "minus plus" hypothesis, while pre-PR the pair was counted twice when both hypotheses were passing the selection (because both indices of the nested loop ranged from 0 to size-1). Here too, a given TriggerCluster is added at most once to the output product, unlike in the pre-PR implementation.

@missirol
Copy link
Contributor Author

missirol commented Mar 1, 2023

enable gpu

Nothing to do with GPUs, but to test this change with the HLT GPU reconstruction on. (potentially different inputs to the filter in question.)

@missirol
Copy link
Contributor Author

missirol commented Mar 1, 2023

please test

@missirol
Copy link
Contributor Author

missirol commented Mar 1, 2023

type bugfix

@cmsbuild
Copy link
Contributor

cmsbuild commented Mar 1, 2023

+1

Summary: https://cmssdt.cern.ch/SDT/jenkins-artifacts/pull-request-integration/PR-7258b4/30977/summary.html
COMMIT: 7732627
CMSSW: CMSSW_13_1_X_2023-03-01-1100/el8_amd64_gcc11
Additional Tests: GPU
User test area: For local testing, you can use /cvmfs/cms-ci.cern.ch/week0/cms-sw/cmssw/40918/30977/install.sh to create a dev area with all the needed externals and cmssw changes.

Comparison Summary

Summary:

  • You potentially removed 27 lines from the logs
  • Reco comparison results: 4 differences found in the comparisons
  • DQMHistoTests: Total files compared: 49
  • DQMHistoTests: Total histograms compared: 3529699
  • DQMHistoTests: Total failures: 4
  • DQMHistoTests: Total nulls: 0
  • DQMHistoTests: Total successes: 3529673
  • DQMHistoTests: Total skipped: 22
  • DQMHistoTests: Total Missing objects: 0
  • DQMHistoSizes: Histogram memory added: 0.0 KiB( 48 files compared)
  • Checked 213 log files, 164 edm output root files, 49 DQM output files
  • TriggerResults: no differences found

GPU Comparison Summary

Summary:

  • No significant changes to the logs found
  • Reco comparison results: 8 differences found in the comparisons
  • DQMHistoTests: Total files compared: 4
  • DQMHistoTests: Total histograms compared: 19862
  • DQMHistoTests: Total failures: 457
  • DQMHistoTests: Total nulls: 0
  • DQMHistoTests: Total successes: 19405
  • DQMHistoTests: Total skipped: 0
  • DQMHistoTests: Total Missing objects: 0
  • DQMHistoSizes: Histogram memory added: 0.0 KiB( 3 files compared)
  • Checked 12 log files, 9 edm output root files, 4 DQM output files
  • TriggerResults: no differences found

@missirol
Copy link
Contributor Author

missirol commented Mar 1, 2023

The differences in DQM outputs for wf 11634.7 (EgammaV/ConversionValidator) are most likely unrelated to this PR. They were recently seen also elsewhere (#40873 (comment)).

@missirol
Copy link
Contributor Author

missirol commented Mar 1, 2023

assign egamma-pog

Since this HLT filter is mostly related to E/gamma, I would kindly ask the POG to review this PR.

It is not really urgent, but if this PR proceeds quickly, its backport might be included in CMSSW_13_0_0.

@cmsbuild
Copy link
Contributor

cmsbuild commented Mar 1, 2023

New categories assigned: egamma-pog

@a-kapoor,@swagata87 you have been requested to review this Pull request/Issue and eventually sign? Thanks


if (mass >= lowerMassCut_ && mass <= upperMassCut_) {
for (unsigned int ii = jj + 1; ii < scs.size(); ii++) {
if (isGoodPair(pEleCh1[jj], pEleCh2[ii]) or isGoodPair(pEleCh2[jj], pEleCh1[ii])) {
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

if (isGoodPair(pEleCh1[jj], pEleCh2[ii]) ) { will suffice, and the or with isGoodPair(pEleCh2[jj], pEleCh1[ii]) can be dropped; because the way isGoodPair function is written, it does not differentiate between electron and positron

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

or maybe keep it to be safe, running on 1K Zee events I do not see any change in results with and w/o the or, but I'm not sure if low-pT would be affected by a wrong charge hypothesis. The way it is now seems safer approach.

@swagata87
Copy link
Contributor

+1

  • Fixes a double-counting issue in the filter.
  • Double-checked by running HLT_Ele28_HighEta_SC20_Mass55_v, the only non-parking HLT path to my knowledge that uses this filter
  • No difference in result before and after this PR, running on 1K DY events
  • Note that the momentum vectors used in this filter code to get the invariant mass is approximate, because it uses beamspot as a proxy of primary vertex. Given this high uncertainty on Z position of vertex, it's unclear if different charge hypothesis are really needed. But the way it is done now seems the safest approach.
  • The mass cut at HLT level should be loose enough to account for the unavoidable momentum approximations.

@missirol
Copy link
Contributor Author

missirol commented Mar 2, 2023

Thanks for the checks and comments, @swagata87 .

@missirol
Copy link
Contributor Author

missirol commented Mar 2, 2023

+hlt

@cmsbuild
Copy link
Contributor

cmsbuild commented Mar 2, 2023

This pull request is fully signed and it will be integrated in one of the next master IBs (tests are also fine). This pull request will now be reviewed by the release team before it's merged. @perrotta, @dpiparo, @rappoccio (and backports should be raised in the release meeting by the corresponding L2)

@rappoccio
Copy link
Contributor

+1

@cmsbuild cmsbuild merged commit aec9b0a into cms-sw:master Mar 2, 2023
@missirol missirol deleted the devel_fixDuplHLTPMMassFilter branch March 2, 2023 16:09
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

Successfully merging this pull request may close these issues.

None yet

4 participants