Belle II Software  release-06-02-00
FilterCalculator.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 #include <hlt/softwaretrigger/calculations/FilterCalculator.h>
9 #include <hlt/softwaretrigger/calculations/utilities.h>
10 
11 #include <TLorentzVector.h>
12 
13 #include <mdst/dataobjects/TrackFitResult.h>
14 #include <mdst/dataobjects/HitPatternCDC.h>
15 
16 #include <numeric>
17 #include <stdexcept>
18 
19 using namespace Belle2;
20 using namespace SoftwareTrigger;
21 
25  double pT = NAN;
27  const Track* track = nullptr;
29  double clusterEnergySumCMS = 0;
31  double clusterEnergySumLab = 0;
33  double pCMS = NAN;
35  double pLab = NAN;
37  TLorentzVector p4CMS;
39  TLorentzVector p4Lab;
40 };
41 
45  const ECLCluster* cluster = nullptr;
47  double energyCMS = NAN;
49  double energyLab = NAN;
51  TLorentzVector p4CMS;
53  TLorentzVector p4Lab;
55  bool isTrack = false;
57  double clusterTime = NAN;
58 };
59 
61 const double flatBoundaries[10] = {0., 19., 22., 25., 30., 35., 45., 60., 90., 180.};
62 
63 
64 FilterCalculator::FilterCalculator() : m_bitsNN("CDCTriggerNNBits")
65 {
66 
67 }
68 
70 {
71  m_tracks.isRequired();
72  m_eclClusters.isRequired();
73  m_l1Trigger.isOptional();
74  m_bitsNN.isOptional();
75 }
76 
77 void FilterCalculator::doCalculation(SoftwareTriggerObject& calculationResult)
78 {
79  calculationResult["nTrkLoose"] = 0;
80  calculationResult["nTrkTight"] = 0;
81  calculationResult["ee2leg"] = 0;
82  calculationResult["nEmedium"] = 0;
83  calculationResult["nElow"] = 0;
84  calculationResult["nEhigh"] = 0;
85  calculationResult["nE180Lab"] = 0;
86  calculationResult["nE300Lab"] = 0;
87  calculationResult["nE500Lab"] = 0;
88  calculationResult["nE2000CMS"] = 0;
89  calculationResult["nE4000CMS"] = 0;
90  calculationResult["nE250Lab"] = 0;
91  calculationResult["nMaxEPhotonAcc"] = 0;
92  calculationResult["dphiCmsClust"] = NAN;
94  calculationResult["netChargeLoose"] = 0;
95  calculationResult["maximumPCMS"] = NAN;
96  calculationResult["eexx"] = 0;
97  calculationResult["ee1leg1trk"] = 0;
98  calculationResult["nEhighLowAng"] = 0;
99  calculationResult["nEsingleClust"] = 0;
100  calculationResult["nEsinglePhotonBarrel"] = 0;
101  calculationResult["nEsinglePhotonExtendedBarrel"] = 0;
102  calculationResult["nEsinglePhotonEndcap"] = 0;
103  calculationResult["nEsingleElectronBarrel"] = 0;
104  calculationResult["nEsingleElectronExtendedBarrel"] = 0;
105  calculationResult["nReducedEsinglePhotonReducedBarrel"] = 0;
106  calculationResult["nVetoClust"] = 0;
107  calculationResult["chrgClust2GeV"] = 0;
108  calculationResult["neutClust045GeVAcc"] = 0;
109  calculationResult["neutClust045GeVBarrel"] = 0;
110  calculationResult["singleTagLowMass"] = 0;
111  calculationResult["singleTagHighMass"] = 0;
112  calculationResult["n2GeVNeutBarrel"] = 0;
113  calculationResult["n2GeVNeutEndcap"] = 0;
114  calculationResult["n2GeVChrg"] = 0;
115  calculationResult["n2GeVPhotonBarrel"] = 0;
116  calculationResult["n2GeVPhotonEndcap"] = 0;
117  calculationResult["ee1leg"] = 0;
118  calculationResult["ee1leg1clst"] = 0;
119  calculationResult["ee1leg1e"] = 0;
120  calculationResult["ee2clst"] = 0;
121  calculationResult["gg2clst"] = 0;
122  calculationResult["eeee"] = 0;
123  calculationResult["eemm"] = 0;
124  calculationResult["eexxSelect"] = 0;
125  calculationResult["radBhabha"] = 0;
126  calculationResult["eeBrem"] = 0;
127  calculationResult["isrRadBhabha"] = 0;
128  calculationResult["muonPairECL"] = 0;
129  calculationResult["selectee1leg1trk"] = 0;
130  calculationResult["selectee1leg1clst"] = 0;
131  calculationResult["selectee"] = 0;
132  calculationResult["ggBarrelVL"] = 0;
133  calculationResult["ggBarrelLoose"] = 0;
134  calculationResult["ggBarrelTight"] = 0;
135  calculationResult["ggEndcapVL"] = 0;
136  calculationResult["ggEndcapLoose"] = 0;
137  calculationResult["ggEndcapTight"] = 0;
138  calculationResult["muonPairV"] = 0;
139  calculationResult["selectmumu"] = 0;
140  calculationResult["singleMuon"] = 0;
141  calculationResult["cosmic"] = 0;
142  calculationResult["eeFlat0"] = 0;
143  calculationResult["eeFlat1"] = 0;
144  calculationResult["eeFlat2"] = 0;
145  calculationResult["eeFlat3"] = 0;
146  calculationResult["eeFlat4"] = 0;
147  calculationResult["eeFlat5"] = 0;
148  calculationResult["eeFlat6"] = 0;
149  calculationResult["eeFlat7"] = 0;
150  calculationResult["eeFlat8"] = 0;
151  calculationResult["eeOneClust"] = 0;
152  // Passed on L1 information
153  if (m_l1Trigger.isValid()) {
154  calculationResult["l1_trigger_random"] = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_RAND;
155  calculationResult["l1_trigger_delayed_bhabha"] = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_DPHY;
156  calculationResult["l1_trigger_poisson"] = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_POIS;
157  bool bha3d;
158  try {
159  bha3d = m_l1Trigger->testPsnm("bha3d");
160  } catch (const std::exception&) {
161  bha3d = false;
162  }
163  bool bhapurPsnm;
164  try {
165  bhapurPsnm = m_l1Trigger->testPsnm("bhapur");
166  } catch (const std::exception&) {
167  bhapurPsnm = false;
168  }
169  bool bhapurFtdl;
170  try {
171  bhapurFtdl = m_l1Trigger->testFtdl("bhapur");
172  } catch (const std::exception&) {
173  bhapurFtdl = false;
174  }
175  bool lml1;
176  try {
177  lml1 = m_l1Trigger->testPsnm("lml1");
178  } catch (const std::exception&) {
179  lml1 = false;
180  }
181  calculationResult["bha3d"] = bha3d;
182  calculationResult["bhapur"] = bhapurPsnm;
183  calculationResult["bhapur_lml1"] = lml1 and bhapurFtdl;
184  } else {
185  calculationResult["l1_trigger_random"] = 1; // save every event if no L1 trigger info
186  calculationResult["l1_trigger_delayed_bhabha"] = 0;
187  calculationResult["l1_trigger_poisson"] = 0;
188  calculationResult["bha3d"] = 0;
189  calculationResult["bhapur"] = 0;
190  calculationResult["bhapur_lml1"] = 0;
191  }
192 
193  // Every 256th event has CDC NN trigger information
194  calculationResult["l1_trg_NN_info"] = 0;
195  if (m_bitsNN.isValid() and m_bitsNN.getEntries() > 0) {calculationResult["l1_trg_NN_info"] = 1;}
196 
197  calculationResult["true"] = 1;
198  calculationResult["false"] = 0;
199 
200  // Some utilities
201  // TODO: into constructor
202  ClusterUtils cUtil;
203  const TVector3 clustervertex = cUtil.GetIPPosition();
204  PCmsLabTransform boostrotate;
205  TLorentzVector p4ofCOM;
206  p4ofCOM.SetPxPyPzE(0, 0, 0, boostrotate.getCMSEnergy());
207 
208  // Pointers to the two tracks with the maximum pt
209  std::map<short, std::optional<MaximumPtTrack>> maximumPtTracks = {
210  { -1, {}}, {1, {}}
211  };
212 
213  // Pointer to the two tracks with maximum pt without a cut applied on z0 (used for cosmic trigger)
214  std::map<short, std::optional<MaximumPtTrack>> maximumPtTracksWithoutZCut = {
215  { -1, {}}, {1, {}}
216  };
217 
218  // --- Track variables -- //
219  for (const Track& track : m_tracks) {
220  const TrackFitResult* trackFitResult = track.getTrackFitResultWithClosestMass(Const::pion);
221  if (not trackFitResult) {
222  // Hu? No track fit result for pion? Ok, lets skip this track
223  continue;
224  }
225 
226  // TODO: Temporary cut on number of CDC hits
227  if (trackFitResult->getHitPatternCDC().getNHits() == 0) {
228  continue;
229  }
230 
231  // Count loose and tight tracks
232  const double z0 = trackFitResult->getZ0();
233  if (std::abs(z0) < m_tightTrkZ0) {
234  calculationResult["nTrkTight"] += 1;
235  }
236 
237  // From here on use only tracks with defined charge
238  const short charge = trackFitResult->getChargeSign();
239  if (charge == 0) {
240  continue;
241  }
242 
243  const TLorentzVector& momentumLab = trackFitResult->get4Momentum();
244  const TLorentzVector momentumCMS = boostrotate.rotateLabToCms() * momentumLab;
245  double pCMS = momentumCMS.Rho();
246 
247  // Find the maximum pt negative [0] and positive [1] tracks without z0 cut
248  const double pT = trackFitResult->getTransverseMomentum();
249  const auto& currentMaximum = maximumPtTracksWithoutZCut.at(charge);
250  if (not currentMaximum or pT > currentMaximum->pT) {
251  MaximumPtTrack newMaximum;
252  newMaximum.pT = pT;
253  newMaximum.track = &track;
254  newMaximum.pCMS = pCMS;
255  newMaximum.pLab = momentumLab.Rho();
256  newMaximum.p4CMS = momentumCMS;
257  newMaximum.p4Lab = momentumLab;
258  maximumPtTracksWithoutZCut[charge] = newMaximum;
259  }
260 
261  // Loose tracks
262  if (std::abs(z0) < m_looseTrkZ0) {
263  calculationResult["nTrkLoose"] += 1;
264  calculationResult["netChargeLoose"] += charge;
265 
266  if (std::isnan(calculationResult["maximumPCMS"]) or pCMS > calculationResult["maximumPCMS"]) {
267  calculationResult["maximumPCMS"] = pCMS;
268  }
269 
270  // Find the maximum pt negative [0] and positive [1] tracks
271  const double pTLoose = trackFitResult->getTransverseMomentum();
272  const auto& currentMaximumLoose = maximumPtTracks.at(charge);
273  if (not currentMaximumLoose or pTLoose > currentMaximumLoose->pT) {
274  MaximumPtTrack newMaximum;
275  newMaximum.pT = pTLoose;
276  newMaximum.track = &track;
277  newMaximum.pCMS = pCMS;
278  newMaximum.pLab = momentumLab.Rho();
279  newMaximum.p4CMS = momentumCMS;
280  newMaximum.p4Lab = momentumLab;
281  maximumPtTracks[charge] = newMaximum;
282  }
283  }
284  }
285 
286  // -- Cluster variables -- //
287  std::vector<SelectedECLCluster> selectedClusters;
288  for (const ECLCluster& cluster : m_eclClusters) {
289  // Use the photon hypothesis, and only look at clusters with times < 200 ns
290  const double time = cluster.getTime();
291  if (not cluster.hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons) or
292  std::abs(time) > 200 or
293  cluster.getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons) < 0.1) {
294  continue;
295  }
296 
297  const double dt99 = cluster.getDeltaTime99();
298 
299  // Store infos if above the single photon threshold
300  SelectedECLCluster selectedCluster;
301  selectedCluster.p4Lab = cUtil.Get4MomentumFromCluster(&cluster, clustervertex,
303  selectedCluster.p4CMS = boostrotate.rotateLabToCms() * selectedCluster.p4Lab; // was clustp4COM
304  selectedCluster.cluster = &cluster;
305  selectedCluster.clusterTime = time / dt99; // was clustT
306  selectedCluster.energyLab = selectedCluster.p4Lab.E();
307  selectedCluster.energyCMS = selectedCluster.p4CMS.E(); // was first of ECOMPair
308  selectedCluster.isTrack = cluster.isTrack(); // was tempCharge
309 
310  selectedClusters.push_back(selectedCluster);
311 
312  if (selectedCluster.energyCMS > m_E2min) {
313  calculationResult["nElow"] += 1;
314  }
315  if (selectedCluster.energyCMS > m_E0min) {
316  calculationResult["nEmedium"] += 1;
317  }
318  if (selectedCluster.energyCMS > m_Ehigh) {
319  calculationResult["nEhigh"] += 1;
320  }
321 
322  if (selectedCluster.energyCMS > m_E0min and std::abs(selectedCluster.clusterTime) < 10) {
323  calculationResult["nVetoClust"] += 1;
324  }
325 
326 
327  //..Single cluster trigger objects use charge, Zernike moment, and thetaLab
328  const double thetaLab = selectedCluster.p4Lab.Theta() * TMath::RadToDeg();
329  const double zmva = cluster.getZernikeMVA();
330  const bool photon = zmva > 0.5 and not selectedCluster.isTrack;
331  const bool electron = zmva > 0.5 and selectedCluster.isTrack;
332 
333  //..For 1 track radiative Bhabha control sample
334  if (selectedCluster.energyCMS > 2. and selectedCluster.isTrack) {
335  calculationResult["chrgClust2GeV"] += 1;
336  }
337  if (selectedCluster.energyCMS > 0.45 and not selectedCluster.isTrack) {
338  const bool isInAcceptance = 17. < thetaLab and thetaLab < 150.;
339  if (isInAcceptance) {calculationResult["neutClust045GeVAcc"] += 1;}
340  const bool isInBarrel = 30. < thetaLab and thetaLab < 130.;
341  if (isInBarrel) {calculationResult["neutClust045GeVBarrel"] += 1;}
342  }
343 
344 
345  // improved 3 cluster (4 cluster) trigger
346  const bool notInHighBackgroundEndcapRegion = 18.5 < thetaLab and thetaLab < 139.3;
347  if (selectedCluster.energyLab > m_EminLab and notInHighBackgroundEndcapRegion) {
348  calculationResult["nE180Lab"] += 1;
349  }
350 
351  if (selectedCluster.energyLab > m_EminLab4Cluster and notInHighBackgroundEndcapRegion) {
352  calculationResult["nE300Lab"] += 1;
353  }
354 
355  if (selectedCluster.energyLab > m_EminLab3Cluster and notInHighBackgroundEndcapRegion) {
356  calculationResult["nE500Lab"] += 1;
357  }
358 
359  if (selectedCluster.energyCMS > m_Ehigh and notInHighBackgroundEndcapRegion) {
360  calculationResult["nE2000CMS"] += 1;
361  }
362 
363  //..For two-photon fusion ALP trigger
364  if (selectedCluster.energyLab > 0.25) {
365  calculationResult["nE250Lab"] += 1;
366  }
367  if (selectedCluster.energyCMS > 4.) {
368  calculationResult["nE4000CMS"] += 1;
369  }
370 
371  //..Single cluster triggers
372  if (selectedCluster.energyCMS > m_EsinglePhoton) {
373  calculationResult["nEsingleClust"] += 1;
374 
375  const bool barrelRegion = thetaLab > 45 and thetaLab < 115;
376  const bool extendedBarrelRegion = thetaLab > 30 and thetaLab < 130;
377  const bool endcapRegion = (thetaLab > 22 and thetaLab < 45) or (thetaLab > 115 and thetaLab < 145);
378 
379  if (photon and barrelRegion) {
380  calculationResult["nEsinglePhotonBarrel"] += 1;
381  }
382 
383  if (photon and extendedBarrelRegion) {
384  calculationResult["nEsinglePhotonExtendedBarrel"] += 1;
385  }
386 
387  if (electron and barrelRegion) {
388  calculationResult["nEsingleElectronBarrel"] += 1;
389  }
390 
391  if (electron and extendedBarrelRegion) {
392  calculationResult["nEsingleElectronExtendedBarrel"] += 1;
393  }
394 
395  if (photon and endcapRegion) {
396  calculationResult["nEsinglePhotonEndcap"] += 1;
397  }
398  }
399 
400  if (selectedCluster.energyCMS > m_reducedEsinglePhoton) {
401  const bool reducedBarrelRegion = thetaLab > 44 and thetaLab < 98;
402 
403  if (photon and reducedBarrelRegion) {
404  calculationResult["nReducedEsinglePhotonReducedBarrel"] += 1;
405  }
406  }
407 
408  if (selectedCluster.energyCMS > m_Ehigh) {
409  // TODO: different definition!
410  const bool barrelRegion = thetaLab > 32 and thetaLab < 130;
411  const bool endcapRegion = (thetaLab > 22 and thetaLab < 32) or (thetaLab > 130 and thetaLab < 145);
412  const bool lowAngleRegion = thetaLab < 22 or thetaLab > 145;
413 
414  if (not selectedCluster.isTrack and barrelRegion) {
415  calculationResult["n2GeVNeutBarrel"] += 1;
416  }
417  if (not selectedCluster.isTrack and endcapRegion) {
418  calculationResult["n2GeVNeutEndcap"] += 1;
419  }
420  if (selectedCluster.isTrack and not lowAngleRegion) {
421  calculationResult["n2GeVChrg"] += 1;
422  }
423  if (lowAngleRegion) {
424  calculationResult["nEhighLowAng"] += 1;
425  }
426  if (photon and barrelRegion) {
427  calculationResult["n2GeVPhotonBarrel"] += 1;
428  }
429  if (photon and endcapRegion) {
430  calculationResult["n2GeVPhotonEndcap"] += 1;
431  }
432  }
433  }
434  std::sort(selectedClusters.begin(), selectedClusters.end(), [](const auto & lhs, const auto & rhs) {
435  return lhs.energyCMS > rhs.energyCMS;
436  });
437 
438  // -- Bhabha and two-photon lepton preparation -- //
439  for (short charge : { -1, 1}) {
440  auto& maximumPtTrack = maximumPtTracks.at(charge);
441  if (not maximumPtTrack) {
442  continue;
443  }
444 
445  // Prep max two tracks for Bhabha and two-photon lepton selections
446  // Track / cluster matching
447  for (auto& cluster : maximumPtTrack->track->getRelationsTo<ECLCluster>()) {
448  if (cluster.hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
449  maximumPtTrack->clusterEnergySumLab += cluster.getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
450  }
451  }
452  maximumPtTrack->clusterEnergySumCMS =
453  maximumPtTrack->clusterEnergySumLab * maximumPtTrack->pCMS / maximumPtTrack->pLab;
454 
455  // Single leg Bhabha selections
456  if (maximumPtTrack->clusterEnergySumCMS > 4.5) {
457  calculationResult["ee1leg"] = 1;
458  }
459 
460  // single muon
461  if (maximumPtTrack->pCMS > 3 and maximumPtTrack->clusterEnergySumLab > 0. and maximumPtTrack->clusterEnergySumLab < 1.) {
462  calculationResult["singleMuon"] = 1;
463  }
464  }
465 
466  // Bhabha selections using two tracks
467  if (maximumPtTracks.at(-1) and maximumPtTracks.at(1)) {
468  const MaximumPtTrack& negativeTrack = *maximumPtTracks.at(-1);
469  const MaximumPtTrack& positiveTrack = *maximumPtTracks.at(1);
470 
471  double dphi = std::abs(negativeTrack.p4CMS.Phi() - positiveTrack.p4CMS.Phi()) * TMath::RadToDeg();
472  if (dphi > 180) {
473  dphi = 360 - dphi;
474  }
475 
476  const double negativeClusterSum = negativeTrack.clusterEnergySumCMS;
477  const double negativeClusterSumLab = negativeTrack.clusterEnergySumLab;
478  const double negativeP = negativeTrack.pCMS;
479  const double positiveClusterSum = positiveTrack.clusterEnergySumCMS;
480  const double positiveClusterSumLab = positiveTrack.clusterEnergySumLab;
481  const double positiveP = positiveTrack.pCMS;
482 
483  // two leg Bhabha
484  const double thetaSum = (negativeTrack.p4CMS.Theta() + positiveTrack.p4CMS.Theta()) * TMath::RadToDeg();
485  const double dthetaSum = std::abs(thetaSum - 180);
486  const auto back2back = dphi > 175 and dthetaSum < 15;
487  if (back2back and negativeClusterSum > 3 and positiveClusterSum > 3 and
488  (negativeClusterSum > 4.5 or positiveClusterSum > 4.5)) {
489  calculationResult["ee2leg"] = 1;
490  }
491 
492  // one leg one track Bhabha
493  if (back2back and ((negativeClusterSum > 4.5 and positiveP > 3) or (positiveClusterSum > 4.5 and negativeP > 3))) {
494  calculationResult["ee1leg1trk"] = 1;
495  }
496 
497 
498  // one leg plus one electron
499  if ((negativeClusterSum > 4.5 and positiveClusterSum > 0.8 * positiveP) or
500  (positiveClusterSum > 4.5 and negativeClusterSum > 0.8 * negativeP)) {
501  calculationResult["ee1leg1e"] = 1;
502  }
503 
504  // eeee, eemm, mu mu, and radiative Bhabha selection
505  const TLorentzVector p4Miss = p4ofCOM - negativeTrack.p4CMS - positiveTrack.p4CMS;
506  const double pmissTheta = p4Miss.Theta() * TMath::RadToDeg();
507  const double pmissp = p4Miss.Rho();
508  const double relMissAngle0 = negativeTrack.p4CMS.Angle(p4Miss.Vect()) * TMath::RadToDeg();
509  const double relMissAngle1 = positiveTrack.p4CMS.Angle(p4Miss.Vect()) * TMath::RadToDeg();
510 
511  const bool electronEP = positiveClusterSum > 0.8 * positiveP or negativeClusterSum > 0.8 * negativeP;
512  const bool notMuonPair = negativeClusterSumLab > 1 or positiveClusterSumLab > 1;
513  const double highp = std::max(negativeP, positiveP);
514  const double lowp = std::min(negativeP, positiveP);
515  const bool lowEdep = negativeClusterSumLab < 0.5 and positiveClusterSumLab < 0.5;
516 
517 
518  if (calculationResult["maximumPCMS"] < 2 and dphi > 160 and (pmissTheta < 25. or pmissTheta > 155.)) {
519  calculationResult["eexxSelect"] = 1;
520  if (electronEP) {
521  calculationResult["eeee"] = 1;
522  } else {
523  calculationResult["eemm"] = 1;
524  }
525  }
526  if (1. < negativeP and 1. < positiveP and 2 == calculationResult["nTrkLoose"] and
527  calculationResult["nTrkTight"] >= 1 and dphi < 170. and
528  10. < pmissTheta and pmissTheta < 170. and 1. < pmissp and electronEP) {
529  calculationResult["radBhabha"] = 1;
530  }
531 
532  if (2. < negativeP and 2. < positiveP and 2 == calculationResult["nTrkLoose"] and
533  calculationResult["nTrkTight"] >= 1 and 175. < dphi and
534  (pmissTheta < 5. or pmissTheta > 175.) and electronEP) {
535  calculationResult["isrRadBhabha"] = 1;
536  }
537  if ((pmissTheta < 20. or pmissTheta > 160.) and
538  ((calculationResult["maximumPCMS"] < 1.2 and dphi > 150.) or
539  (calculationResult["maximumPCMS"] < 2. and 175. < dphi))) {
540  calculationResult["eexx"] = 1;
541  }
542  if (calculationResult["nTrkLoose"] == 2 and highp > 4.5 and notMuonPair and pmissp > 1. and
543  (relMissAngle0 < 5. or relMissAngle1 < 5.)) {
544  calculationResult["eeBrem"] = 1;
545  }
546 
547  if (calculationResult["nTrkLoose"] == 2 and highp > 4.5 and lowEdep and dphi > 175. and 175. < thetaSum and
548  thetaSum < 185.) {
549  calculationResult["muonPairV"] = 1;
550  }
551  if (3. < highp and 2.5 < lowp and 165. < dphi and
552  ((negativeClusterSumLab > 0. and negativeClusterSumLab < 1.) or
553  (positiveClusterSumLab > 0. and positiveClusterSumLab < 1.))) {
554  calculationResult["selectmumu"] = 1;
555  }
556  }
557 
558  if (selectedClusters.size() >= 2) {
559  const SelectedECLCluster& firstCluster = selectedClusters[0];
560  const SelectedECLCluster& secondCluster = selectedClusters[1];
561 
562  // Bhabha / gg vetoes using max two clusters
563  double dphi = std::abs(firstCluster.p4CMS.Phi() - secondCluster.p4CMS.Phi()) * TMath::RadToDeg();
564  if (dphi > 180) {
565  dphi = 360 - dphi;
566  }
567  double thetaSum = (firstCluster.p4CMS.Theta() + secondCluster.p4CMS.Theta()) * TMath::RadToDeg();
568  double dthetaSum = std::abs(thetaSum - 180);
569 
570  //..Quantities for two-photon fusion triggers
571  calculationResult["dphiCmsClust"] = dphi;
572  for (int ic = 0; ic < 2; ic++) {
573  const double thetaLab = selectedClusters[ic].p4Lab.Theta() * TMath::RadToDeg();
574  const bool isInAcceptance = 17. < thetaLab and thetaLab < 150.;
575  const ECLCluster* cluster = selectedClusters[ic].cluster;
576  const double zmva = cluster->getZernikeMVA();
577  const bool photon = zmva > 0.5 and not selectedClusters[ic].isTrack;
578  if (isInAcceptance and photon) {calculationResult["nMaxEPhotonAcc"] += 1;}
579  }
580 
581  const double firstEnergy = firstCluster.p4CMS.E();
582  const double secondEnergy = secondCluster.p4CMS.E();
583 
584  const bool highEnergetic = firstEnergy > 3 and secondEnergy > 3 and (firstEnergy > 4.5 or secondEnergy > 4.5);
585 
586  if (dphi > 160 and dphi < 178 and dthetaSum < 15 and highEnergetic) {
587  calculationResult["ee2clst"] = 1;
588  }
589 
590  if (dphi > 178 and dthetaSum < 15 and highEnergetic) {
591  calculationResult["gg2clst"] = 1;
592  }
593 
594  if ((calculationResult["ee2clst"] == 1 or calculationResult["gg2clst"] == 1) and
595  calculationResult["ee1leg"] == 1) {
596  calculationResult["ee1leg1clst"] = 1;
597  }
598 
599  const double Elab0 = firstCluster.p4Lab.E();
600  const double Elab1 = secondCluster.p4Lab.E();
601 
602  // gg and mumu accept triggers using ECL
603  if (firstEnergy > 2 and secondEnergy > 2) {
604  const double thetaLab0 = firstCluster.p4Lab.Theta() * TMath::RadToDeg();
605  const double thetaLab1 = secondCluster.p4Lab.Theta() * TMath::RadToDeg();
606 
607  const bool barrel0 = 32. < thetaLab0 and thetaLab0 < 130.;
608  const bool barrel1 = 32. < thetaLab1 and thetaLab1 < 130.;
609  const bool oneClustersAbove4 = firstEnergy > 4 or secondEnergy > 4;
610  const bool oneIsNeutral = not firstCluster.isTrack or not secondCluster.isTrack;
611  const bool bothAreNeutral = not firstCluster.isTrack and not secondCluster.isTrack;
612  const bool oneIsBarrel = barrel0 or barrel1;
613  const bool dphiCutExtraLoose = dphi > 175;
614  const bool dphiCutLoose = dphi > 177;
615  const bool dphiCutTight = dphi > 177.5;
616 
617  if (dphiCutExtraLoose and oneIsNeutral and oneIsBarrel) {
618  calculationResult["ggBarrelVL"] = 1;
619  }
620  if (oneClustersAbove4 and dphiCutLoose and oneIsNeutral and oneIsBarrel) {
621  calculationResult["ggBarrelLoose"] = 1;
622  }
623  if (oneClustersAbove4 and dphiCutTight and bothAreNeutral and oneIsBarrel) {
624  calculationResult["ggBarrelTight"] = 1;
625  }
626  if (dphiCutExtraLoose and oneIsNeutral and not oneIsBarrel) {
627  calculationResult["ggEndcapVL"] = 1;
628  }
629  if (oneClustersAbove4 and dphiCutLoose and oneIsNeutral and not oneIsBarrel) {
630  calculationResult["ggEndcapLoose"] = 1;
631  }
632  if (oneClustersAbove4 and dphiCutTight and bothAreNeutral and not oneIsBarrel) {
633  calculationResult["ggEndcapTight"] = 1;
634  }
635  }
636 
637  const double minEnergy = std::min(Elab0, Elab1);
638  const double maxEnergy = std::max(Elab0, Elab1);
639  if (dphi > 155 and thetaSum > 165 and thetaSum < 195 and minEnergy > 0.15 and minEnergy < 0.5 and
640  maxEnergy > 0.15 and maxEnergy < 0.5) {
641  calculationResult["muonPairECL"] = 1;
642  }
643 
644  }
645 
646  // Bhabha accept triggers.
647  // Use theta_lab of negative track if available; otherwise cluster.
648  double thetaFlatten = 0;
649 
650  // One leg, one cluster.
651  for (short charge : { -1, 1}) {
652  const auto& maximumPtTrack = maximumPtTracks.at(charge);
653  if (not maximumPtTrack) {
654  continue;
655  }
656 
657  if (maximumPtTrack->clusterEnergySumCMS > 1.5) {
658  double invMass = 0.;
659  double tempFlatten = 0.;
660  for (const SelectedECLCluster& cluster : selectedClusters) {
661  double tempInvMass = (maximumPtTrack->p4Lab + cluster.p4Lab).M();
662  if (tempInvMass > invMass) {
663  invMass = tempInvMass;
664  if (charge == 1) {
665  tempFlatten = cluster.p4Lab.Theta() * TMath::RadToDeg();
666  }
667  }
668  }
669  if (charge == -1) {
670  tempFlatten = maximumPtTrack->p4Lab.Theta() * TMath::RadToDeg();
671  }
672  if (invMass > 5.29) {
673  calculationResult["selectee1leg1clst"] = 1;
674  thetaFlatten = tempFlatten;
675  }
676  }
677  }
678 
679  // Two legs
680  if (maximumPtTracks.at(-1) and maximumPtTracks.at(1)) {
681  const MaximumPtTrack& negativeTrack = *maximumPtTracks.at(-1);
682  const MaximumPtTrack& positiveTrack = *maximumPtTracks.at(1);
683  const double invMass = (negativeTrack.p4Lab + positiveTrack.p4Lab).M();
684  if (invMass > 5.29 and (negativeTrack.clusterEnergySumCMS > 1.5 or positiveTrack.clusterEnergySumCMS > 1.5)) {
685  calculationResult["selectee1leg1trk"] = 1;
686  }
687 
688  thetaFlatten = negativeTrack.p4Lab.Theta() * TMath::RadToDeg();
689 
690  // Two tracks but (exactly) 1 high energy cluster
691  const bool missNegClust = positiveTrack.clusterEnergySumCMS > 4.5 and negativeTrack.clusterEnergySumCMS < 1.;
692  const bool missPosClust = negativeTrack.clusterEnergySumCMS > 4.5 and positiveTrack.clusterEnergySumCMS < 1.;
693  if ((invMass > 9.) and (missNegClust or missPosClust)) {
694  calculationResult["eeOneClust"] = 1;
695  }
696  }
697 
698  if (calculationResult["selectee1leg1trk"] == 1 or calculationResult["selectee1leg1clst"] == 1) {
699  for (int iflat = 0; iflat < 9; iflat++) {
700  const std::string& eeFlatName = "eeFlat" + std::to_string(iflat);
701  calculationResult[eeFlatName] =
702  thetaFlatten > flatBoundaries[iflat] and thetaFlatten < flatBoundaries[iflat + 1];
703  if (calculationResult[eeFlatName]) {
704  calculationResult["selectee"] = 1;
705  }
706  }
707  }
708 
709  // Single tag pi0 / eta dedicated lines
710  if (calculationResult["nTrkLoose"] == 1 and calculationResult["maximumPCMS"] > 0.8 and selectedClusters.size() >= 2) {
711 
712  decltype(selectedClusters) selectedSingleTagClusters(selectedClusters.size());
713  auto lastItem = std::copy_if(selectedClusters.begin(), selectedClusters.end(), selectedSingleTagClusters.begin(),
714  [](auto & cluster) {
715  const bool isNeutralCluster = not cluster.isTrack;
716  const bool hasEnoughEnergy = cluster.energyLab > 0.1;
717  const double clusterThetaLab = cluster.p4Lab.Theta() * TMath::RadToDeg();
718  const bool isInAcceptance = 17 < clusterThetaLab and clusterThetaLab < 150.;
719  return isNeutralCluster and hasEnoughEnergy and isInAcceptance;
720  });
721  selectedSingleTagClusters.resize(std::distance(selectedSingleTagClusters.begin(), lastItem));
722 
723  if (selectedSingleTagClusters.size() >= 2) { // One track and at least two clusters are found
724 
725  const auto& track = maximumPtTracks.at(-1) ? *maximumPtTracks.at(-1) : *maximumPtTracks.at(1);
726  // in real signal, the pi0 decay daughters are always the two most-energetic neutral clusters.
727  const auto& firstCluster = selectedSingleTagClusters[0];
728  const auto& secondCluster = selectedSingleTagClusters[1];
729 
730  const TLorentzVector trackP4CMS = track.p4CMS;
731  const TLorentzVector pi0P4CMS = firstCluster.p4CMS + secondCluster.p4CMS;
732 
733  const bool passPi0ECMS = pi0P4CMS.E() > 1. and pi0P4CMS.E() < 0.525 * p4ofCOM.M();
734  const double thetaSumCMS = (pi0P4CMS.Theta() + trackP4CMS.Theta()) * TMath::RadToDeg();
735  const bool passThetaSum = thetaSumCMS < 170. or thetaSumCMS > 190.;
736 
737  double dphiCMS = std::abs(trackP4CMS.Phi() - pi0P4CMS.Phi()) * TMath::RadToDeg();
738  if (dphiCMS > 180) {
739  dphiCMS = 360 - dphiCMS;
740  }
741  const bool passdPhi = dphiCMS > 160.;
742 
743  if (passPi0ECMS and passThetaSum and passdPhi and pi0P4CMS.M() < 0.7) {
744  calculationResult["singleTagLowMass"] = 1;
745  } else if (passPi0ECMS and passThetaSum and passdPhi and pi0P4CMS.M() > 0.7) {
746  calculationResult["singleTagHighMass"] = 1;
747  }
748  }
749  }
750 
751  //..Cosmic selection. Need exactly two tracks.
752  if (m_tracks.getEntries() == 2) {
753 
754  const auto negTrack = maximumPtTracksWithoutZCut.at(-1);
755  const auto posTrack = maximumPtTracksWithoutZCut.at(1);
756 
757  if (negTrack and posTrack) {
758 
759  const double maxNegpT = negTrack->pT;
760  const double maxPospT = posTrack->pT;
761 
762  auto accumulatePhotonEnergy = [](double result, const auto & cluster) {
763  return result + (cluster.hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons) ? cluster.getEnergy(
765  };
766 
767  const auto& clustersOfNegTrack = negTrack->track->getRelationsTo<ECLCluster>();
768  const auto& clustersOfPosTrack = posTrack->track->getRelationsTo<ECLCluster>();
769 
770  const double maxClusterENeg = std::accumulate(clustersOfNegTrack.begin(), clustersOfNegTrack.end(), 0.0, accumulatePhotonEnergy);
771  const double maxClusterEPos = std::accumulate(clustersOfPosTrack.begin(), clustersOfPosTrack.end(), 0.0, accumulatePhotonEnergy);
772 
773  const TLorentzVector& momentumLabNeg(negTrack->p4Lab);
774  const TLorentzVector& momentumLabPos(posTrack->p4Lab);
775 
776  const double& z0Neg = negTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getZ0();
777  const double& d0Neg = negTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getD0();
778  const double& z0Pos = posTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getZ0();
779  const double& d0Pos = posTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getD0();
780 
781  // Select cosmic using these tracks
782  const bool goodMagneticRegion = (z0Neg<m_goodMagneticRegionZ0 or abs(d0Neg)>m_goodMagneticRegionD0)
783  and (z0Pos<m_goodMagneticRegionZ0 or abs(d0Pos)>m_goodMagneticRegionD0);
784  if (maxNegpT > m_cosmicMinPt and maxPospT > m_cosmicMinPt and maxClusterENeg < m_cosmicMaxClusterEnergy
785  and maxClusterEPos < m_cosmicMaxClusterEnergy and goodMagneticRegion) {
786  double dphiLab = std::abs(momentumLabNeg.Phi() - momentumLabPos.Phi()) * TMath::RadToDeg();
787  if (dphiLab > 180) {
788  dphiLab = 360 - dphiLab;
789  }
790 
791  const double thetaSumLab = (momentumLabNeg.Theta() + momentumLabPos.Theta()) * TMath::RadToDeg();
792 
793  constexpr double phiBackToBackTolerance = 2.;
794  constexpr double thetaBackToBackTolerance = 2.;
795  if ((180 - dphiLab) < phiBackToBackTolerance and std::abs(180 - thetaSumLab) < thetaBackToBackTolerance) {
796  calculationResult["cosmic"] = 1;
797  }
798  }
799  }
800  }
801 
802 }
Class to provide momentum-related information from ECLClusters.
Definition: ClusterUtils.h:34
const TLorentzVector Get4MomentumFromCluster(const ECLCluster *cluster, ECLCluster::EHypothesisBit hypo)
Returns four momentum vector.
Definition: ClusterUtils.cc:24
const TVector3 GetIPPosition()
Returns default IP position from beam parameters.
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
ECL cluster data.
Definition: ECLCluster.h:27
@ c_nPhotons
CR is split into n photons (N1)
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
Class to hold Lorentz transformations from/to CMS and boost vector.
double getCMSEnergy() const
Returns CMS energy of e+e- (aka.
const TLorentzRotation rotateLabToCms() const
Returns Lorentz transformation from Lab to CMS.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
double m_EminLab
which lab energy defines nE180Lab
double m_reducedEsinglePhoton
which CMS energy defines nReducedEsingle clusters
StoreObjPtr< TRGSummary > m_l1Trigger
Store Object with the trigger result.
double m_E0min
which CMS energy defines nEmedium
void requireStoreArrays() override
Require the particle list. We do not need more here.
void doCalculation(SoftwareTriggerObject &calculationResult) override
Actually write out the variables into the map.
double m_goodMagneticRegionZ0
maximum z0 for well understood magnetic field (cm)
double m_cosmicMinPt
which LAB pt defines a cosmic
double m_EsinglePhoton
which CMS energy defines nEsingleClust
double m_E2min
which CMS energy defines nElow
double m_tightTrkZ0
which Z0 defines a tight track
StoreArray< Track > m_tracks
Store Array of the tracks to be used.
double m_EminLab3Cluster
which lab energy defines nE500Lab
double m_Ehigh
which CMS energy defines nEhigh
StoreArray< ECLCluster > m_eclClusters
Store Array of the ecl clusters to be used.
double m_EminLab4Cluster
which lab energy defines nE300Lab
StoreArray< CDCTriggerUnpacker::NNBitStream > m_bitsNN
Store Object with the trigger NN bits.
FilterCalculator()
Set the default names for the store object particle lists.
double m_cosmicMaxClusterEnergy
which LAB cluster energy vetoes a cosmic candidate
double m_goodMagneticRegionD0
minimum d0 for well understood magnetic field, if z0 is large (cm)
double m_looseTrkZ0
which Z0 defines a loose track
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
TLorentzVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
double getTransverseMomentum() const
Getter for the absolute value of the transverse momentum at the perigee.
double getZ0() const
Getter for z0.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
Abstract base class for different kinds of events.
Temporary data structure holding the track(s) with the maximum pT.
double pLab
the momentum magnitude in lab system
double pT
the pT of the track
double clusterEnergySumLab
the sum of related cluster energies in lab system
double clusterEnergySumCMS
the sum of related cluster energies in CMS system
double pCMS
the momentum magnitude in CMS system
TLorentzVector p4CMS
the 4 momentum in CMS system
TLorentzVector p4Lab
the 4 momentum in lab system
const Track * track
the track
Temporary data structure holding the ECL clusters used for this analysis.
double energyLab
the energy in Lab system
double clusterTime
the time of the cluster
TLorentzVector p4CMS
the 4 momentum in CMS system
TLorentzVector p4Lab
the 4 momentum in lab system
const ECLCluster * cluster
The ECL cluster.
bool isTrack
is this ECL cluster likely from a track (or a photon) = is it charged?
double energyCMS
the energy in CMS system