Belle II Software release-09-00-02
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 <Math/Vector3D.h>
12#include <Math/Vector4D.h>
13#include <Math/VectorUtil.h>
14
15#include <mdst/dataobjects/TrackFitResult.h>
16#include <mdst/dataobjects/HitPatternCDC.h>
17
18#include <analysis/dataobjects/Particle.h>
19#include <analysis/VertexFitting/KFit/VertexFitKFit.h>
20
21#include <framework/logging/Logger.h>
22
23#include <numeric>
24#include <stdexcept>
25#include <optional>
26
27using namespace Belle2;
28using namespace SoftwareTrigger;
29
33 double pT = NAN;
35 const Track* track = nullptr;
41 double pCMS = NAN;
43 double pLab = NAN;
45 ROOT::Math::PxPyPzEVector p4CMS;
47 ROOT::Math::PxPyPzEVector p4Lab;
48};
49
53 const ECLCluster* cluster = nullptr;
55 double energyCMS = NAN;
57 double energyLab = NAN;
59 ROOT::Math::PxPyPzEVector p4CMS;
61 ROOT::Math::PxPyPzEVector p4Lab;
63 bool isTrack = false;
65 double clusterTime = NAN;
66};
67
69const double flatBoundaries[10] = {0., 19., 22., 25., 30., 35., 45., 60., 90., 180.};
70
71
72FilterCalculator::FilterCalculator() : m_bitsNN("CDCTriggerNNBits")
73{
74
75}
76
78{
79 m_tracks.isRequired();
80 m_eclClusters.isRequired();
81 m_l1Trigger.isOptional();
82 m_bitsNN.isOptional();
83}
84
85void FilterCalculator::doCalculation(SoftwareTriggerObject& calculationResult)
86{
87 calculationResult["nTrkLoose"] = 0;
88 calculationResult["nTrkTight"] = 0;
89 calculationResult["ee2leg"] = 0;
90 calculationResult["nEmedium"] = 0;
91 calculationResult["nElow"] = 0;
92 calculationResult["nEhigh"] = 0;
93 calculationResult["nE180Lab"] = 0;
94 calculationResult["nE300Lab"] = 0;
95 calculationResult["nE500Lab"] = 0;
96 calculationResult["nE2000CMS"] = 0;
97 calculationResult["nE4000CMS"] = 0;
98 calculationResult["nE250Lab"] = 0;
99 calculationResult["nMaxEPhotonAcc"] = 0;
100 calculationResult["dphiCmsClust"] = NAN;
102 calculationResult["netChargeLoose"] = 0;
103 calculationResult["maximumPCMS"] = NAN;
104 calculationResult["maximumPLab"] = NAN;
105 calculationResult["eexx"] = 0;
106 calculationResult["ee1leg1trk"] = 0;
107 calculationResult["nEhighLowAng"] = 0;
108 calculationResult["nEsingleClust"] = 0;
109 calculationResult["nEsinglePhotonBarrel"] = 0;
110 calculationResult["nEsinglePhotonExtendedBarrel"] = 0;
111 calculationResult["nEsinglePhotonEndcap"] = 0;
112 calculationResult["nEsingleElectronBarrel"] = 0;
113 calculationResult["nEsingleElectronExtendedBarrel"] = 0;
114 calculationResult["nReducedEsinglePhotonReducedBarrel"] = 0;
115 calculationResult["nVetoClust"] = 0;
116 calculationResult["chrgClust2GeV"] = 0;
117 calculationResult["neutClust045GeVAcc"] = 0;
118 calculationResult["neutClust045GeVBarrel"] = 0;
119 calculationResult["singleTagLowMass"] = 0;
120 calculationResult["singleTagHighMass"] = 0;
121 calculationResult["n2GeVNeutBarrel"] = 0;
122 calculationResult["n2GeVNeutEndcap"] = 0;
123 calculationResult["n2GeVChrg"] = 0;
124 calculationResult["n2GeVPhotonBarrel"] = 0;
125 calculationResult["n2GeVPhotonEndcap"] = 0;
126 calculationResult["ee1leg"] = 0;
127 calculationResult["ee1leg1clst"] = 0;
128 calculationResult["ee1leg1e"] = 0;
129 calculationResult["ee2clst"] = 0;
130 calculationResult["gg2clst"] = 0;
131 calculationResult["eeee"] = 0;
132 calculationResult["eemm"] = 0;
133 calculationResult["eexxSelect"] = 0;
134 calculationResult["radBhabha"] = 0;
135 calculationResult["eeBrem"] = 0;
136 calculationResult["isrRadBhabha"] = 0;
137 calculationResult["muonPairECL"] = 0;
138 calculationResult["ggHighPt"] = 0;
139 calculationResult["selectee1leg1trk"] = 0;
140 calculationResult["selectee1leg1clst"] = 0;
141 calculationResult["selectee"] = 0;
142 calculationResult["ggBarrelVL"] = 0;
143 calculationResult["ggBarrelLoose"] = 0;
144 calculationResult["ggBarrelTight"] = 0;
145 calculationResult["ggEndcapVL"] = 0;
146 calculationResult["ggEndcapLoose"] = 0;
147 calculationResult["ggEndcapTight"] = 0;
148 calculationResult["muonPairV"] = 0;
149 calculationResult["selectmumu"] = 0;
150 calculationResult["singleMuon"] = 0;
151 calculationResult["cosmic"] = 0;
152 calculationResult["displacedVertex"] = 0;
153 calculationResult["eeFlat0"] = 0;
154 calculationResult["eeFlat1"] = 0;
155 calculationResult["eeFlat2"] = 0;
156 calculationResult["eeFlat3"] = 0;
157 calculationResult["eeFlat4"] = 0;
158 calculationResult["eeFlat5"] = 0;
159 calculationResult["eeFlat6"] = 0;
160 calculationResult["eeFlat7"] = 0;
161 calculationResult["eeFlat8"] = 0;
162 calculationResult["eeOneClust"] = 0;
163
164 //..New track definitions
165 calculationResult["nTrkLooseB"] = 0;
166 calculationResult["nTrkTightB"] = 0;
167 calculationResult["maximumPCMSB"] = NAN;
168 calculationResult["netChargeLooseB"] = 0;
169 calculationResult["muonPairVB"] = 0;
170 calculationResult["eeBremB"] = 0;
171 calculationResult["singleTagLowMassB"] = 0;
172 calculationResult["singleTagHighMassB"] = 0;
173 calculationResult["radBhabhaB"] = 0;
175 //..Tracks without a requirement on nCDCHits
176 calculationResult["nTrackC"] = 0;
177 calculationResult["maximumPCMSC"] = NAN;
178 calculationResult["pCmsNegC"] = NAN;
179 calculationResult["clusterENegC"] = NAN;
180 calculationResult["pCmsPosC"] = NAN;
181 calculationResult["clusterEPosC"] = NAN;
182 calculationResult["dPhiCmsC"] = NAN;
185 // Passed on L1 information
186 if (m_l1Trigger.isValid()) {
187 calculationResult["l1_trigger_random"] = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_RAND;
188 calculationResult["l1_trigger_delayed_bhabha"] = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_DPHY;
189 calculationResult["l1_trigger_poisson"] = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_POIS;
190 bool bha3d;
191 try {
192 bha3d = m_l1Trigger->testPsnm("bha3d");
193 } catch (const std::exception&) {
194 bha3d = false;
195 }
196 bool bhapurPsnm;
197 try {
198 bhapurPsnm = m_l1Trigger->testPsnm("bhapur");
199 } catch (const std::exception&) {
200 bhapurPsnm = false;
201 }
202 bool bhapurFtdl;
203 try {
204 bhapurFtdl = m_l1Trigger->testFtdl("bhapur");
205 } catch (const std::exception&) {
206 bhapurFtdl = false;
207 }
208 bool lml1;
209 try {
210 lml1 = m_l1Trigger->testPsnm("lml1");
211 } catch (const std::exception&) {
212 lml1 = false;
213 }
214 bool l1_bit_f;
215 try {
216 l1_bit_f = m_l1Trigger->testPsnm("fpre");
217 } catch (const std::exception&) {
218 l1_bit_f = false;
219 }
220 calculationResult["bha3d"] = bha3d;
221 calculationResult["bhapur"] = bhapurPsnm;
222 calculationResult["bhapur_lml1"] = lml1 and bhapurFtdl;
223 calculationResult["l1_bit_f"] = l1_bit_f;
224 } else {
225 calculationResult["l1_trigger_random"] = 1; // save every event if no L1 trigger info
226 calculationResult["l1_trigger_delayed_bhabha"] = 0;
227 calculationResult["l1_trigger_poisson"] = 0;
228 calculationResult["bha3d"] = 0;
229 calculationResult["bhapur"] = 0;
230 calculationResult["bhapur_lml1"] = 0;
231 calculationResult["l1_bit_f"] = 0;
232 }
233
234 // Every 256th event has CDC NN trigger information
235 calculationResult["l1_trg_NN_info"] = 0;
236 if (m_bitsNN.isValid() and m_bitsNN.getEntries() > 0) {calculationResult["l1_trg_NN_info"] = 1;}
237
238 calculationResult["true"] = 1;
239 calculationResult["false"] = 0;
240
241 // Some utilities
242 // TODO: into constructor
243 ClusterUtils cUtil;
244 const ROOT::Math::XYZVector clustervertex = cUtil.GetIPPosition();
245 PCmsLabTransform boostrotate;
246 ROOT::Math::PxPyPzEVector p4ofCOM;
247 p4ofCOM.SetPxPyPzE(0, 0, 0, boostrotate.getCMSEnergy());
248
249 // Pointers to the two tracks with the maximum pt
250 std::map<short, std::optional<MaximumPtTrack>> maximumPtTracks = {
251 { -1, {}}, {1, {}}
252 };
253
254 // Pointer to the two tracks with maximum pt without a cut applied on z0 (used for cosmic trigger)
255 std::map<short, std::optional<MaximumPtTrack>> maximumPtTracksWithoutZCut = {
256 { -1, {}}, {1, {}}
257 };
258
259 //..Pointers to the two "C" tracks with maximum pCms (no requirement on nCDC hits)
260 std::map<short, std::optional<MaximumPtTrack>> maximumPCmsTracksC = {
261 { -1, {}}, {1, {}}
262 };
263
264 //------------------------------------------------------------
265 // --- Track variables -- //
266 for (const Track& track : m_tracks) {
267 const TrackFitResult* trackFitResult = track.getTrackFitResultWithClosestMass(Const::pion);
268 if (not trackFitResult) {
269 // Huh? No track fit result for pion? Ok, lets skip this track
270 continue;
271 }
272
273 //..Track properties. Only use tracks with defined charge.
274 const short charge = trackFitResult->getChargeSign();
275 if (charge == 0) {continue;}
276 const double z0 = trackFitResult->getZ0();
277 const ROOT::Math::PxPyPzEVector& momentumLab = trackFitResult->get4Momentum();
278 const ROOT::Math::PxPyPzEVector momentumCMS = boostrotate.rotateLabToCms() * momentumLab;
279 const double pCMS = momentumCMS.P();
280 const double pLab = momentumLab.P();
281 const double trackTime = track.getTrackTime();
282 const double pT = trackFitResult->getTransverseMomentum();
283 const int nCDCHits = (int)(0.5 + trackFitResult->getHitPatternCDC().getNHits());
284 double clusterELab = 0.;
285 for (auto& cluster : track.getRelationsTo<ECLCluster>()) {
286 if (cluster.hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
287 clusterELab += cluster.getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
288 }
289 }
290 const double clusterECMS = clusterELab * pCMS / pLab;
291
292
293 //------------------------------------------------------------
294 //..New list without a requirement on number of CDC hits
295 bool goodTrackCTime = true;
296 if (std::abs(trackTime) > 15.) {goodTrackCTime = false;}
297 if (std::abs(z0) < 1. and goodTrackCTime and pCMS > 0.2) {
298 calculationResult["nTrackC"] += 1;
299 if (std::isnan(calculationResult["maximumPCMSC"]) or pCMS > calculationResult["maximumPCMSC"]) {
300 calculationResult["maximumPCMSC"] = pCMS;
301 }
302 }
303
304 //..Store maximum pCms positive and negative tracks
305 const auto& currentMaximumC = maximumPCmsTracksC.at(charge);
306 if (not currentMaximumC or pCMS > currentMaximumC->pCMS) {
307 MaximumPtTrack newMaximum;
308 newMaximum.pT = pT;
309 newMaximum.track = &track;
310 newMaximum.pCMS = pCMS;
311 newMaximum.pLab = pLab;
312 newMaximum.p4CMS = momentumCMS;
313 newMaximum.p4Lab = momentumLab;
314 newMaximum.clusterEnergySumLab = clusterELab;
315 newMaximum.clusterEnergySumCMS = clusterECMS;
316 maximumPCmsTracksC[charge] = newMaximum;
317 }
318
319 //------------------------------------------------------------
320 //..Now original loose and tight tracks, which require nCDC hits
321 // TODO: Temporary cut on number of CDC hits
322 if (nCDCHits == 0) {
323 continue;
324 }
325
326 // Count tight tracks
327 if (std::abs(z0) < m_tightTrkZ0) {
328 calculationResult["nTrkTight"] += 1;
329 }
330
331
332 //------------------------------------------------------------
333 // Find the maximum pt negative [0] and positive [1] tracks without z0 cut
334 const auto& currentMaximum = maximumPtTracksWithoutZCut.at(charge);
335 if (not currentMaximum or pT > currentMaximum->pT) {
336 MaximumPtTrack newMaximum;
337 newMaximum.pT = pT;
338 newMaximum.track = &track;
339 newMaximum.pCMS = pCMS;
340 newMaximum.pLab = pLab;
341 newMaximum.p4CMS = momentumCMS;
342 newMaximum.p4Lab = momentumLab;
343 newMaximum.clusterEnergySumLab = clusterELab;
344 newMaximum.clusterEnergySumCMS = clusterECMS;
345 maximumPtTracksWithoutZCut[charge] = newMaximum;
346 }
347
348
349 //------------------------------------------------------------
350 // Loose tracks
351 if (std::abs(z0) < m_looseTrkZ0) {
352 calculationResult["nTrkLoose"] += 1;
353 calculationResult["netChargeLoose"] += charge;
354
355 if (std::isnan(calculationResult["maximumPCMS"]) or pCMS > calculationResult["maximumPCMS"]) {
356 calculationResult["maximumPCMS"] = pCMS;
357 }
358
359 if (std::isnan(calculationResult["maximumPLab"]) or pLab > calculationResult["maximumPLab"]) {
360 calculationResult["maximumPLab"] = pLab;
361 }
362
363 // Find the maximum pt negative [0] and positive [1] tracks
364 const double pTLoose = trackFitResult->getTransverseMomentum();
365 const auto& currentMaximumLoose = maximumPtTracks.at(charge);
366 if (not currentMaximumLoose or pTLoose > currentMaximumLoose->pT) {
367 MaximumPtTrack newMaximum;
368 newMaximum.pT = pTLoose;
369 newMaximum.track = &track;
370 newMaximum.pCMS = pCMS;
371 newMaximum.pLab = momentumLab.P();
372 newMaximum.p4CMS = momentumCMS;
373 newMaximum.p4Lab = momentumLab;
374 newMaximum.clusterEnergySumLab = clusterELab;
375 newMaximum.clusterEnergySumCMS = clusterECMS;
376 maximumPtTracks[charge] = newMaximum;
377 }
378 }
379
380 //..New tighter track definitions.
381 // We can use pCMS and charge from above, because every looseB track is also a loose track
382 if (nCDCHits >= 5) {
383 if (std::abs(z0) < 1.) {calculationResult["nTrkTightB"] += 1;}
384 if (std::abs(z0) < 5.) {
385 calculationResult["nTrkLooseB"] += 1;
386 calculationResult["netChargeLooseB"] += charge;
387 if (std::isnan(calculationResult["maximumPCMSB"]) or pCMS > calculationResult["maximumPCMSB"]) {
388 calculationResult["maximumPCMSB"] = pCMS;
389
390 }
391 }
392 }
393
394 } // End loop over tracks
395
396
397 //------------------------------------------------------------
398 //..Store the "C" track calculations
399 for (short charge : { -1, 1}) {
400 auto& maximumPcmsTrackC = maximumPCmsTracksC.at(charge);
401 if (not maximumPcmsTrackC) {
402 continue;
403 }
404
405 //..Calculations for each track
406 if (charge == -1) {
407 calculationResult["pCmsNegC"] = maximumPcmsTrackC->pCMS;
408 calculationResult["clusterENegC"] = maximumPcmsTrackC->clusterEnergySumLab;
409 } else {
410 calculationResult["pCmsPosC"] = maximumPcmsTrackC->pCMS;
411 calculationResult["clusterEPosC"] = maximumPcmsTrackC->clusterEnergySumLab;
412 }
413 }
414
415 //..Calculations using both highest-pCms C tracks
416 if (maximumPCmsTracksC.at(-1) and maximumPCmsTracksC.at(1)) {
417 const MaximumPtTrack& negativeTrack = *maximumPCmsTracksC.at(-1);
418 const MaximumPtTrack& positiveTrack = *maximumPCmsTracksC.at(1);
419
420 calculationResult["dPhiCmsC"] = std::abs(ROOT::Math::VectorUtil::DeltaPhi(negativeTrack.p4CMS,
421 positiveTrack.p4CMS)) * TMath::RadToDeg();
422 }
423
424
425 //------------------------------------------------------------
426 // -- Cluster variables -- //
427 std::vector<SelectedECLCluster> selectedClusters;
428 for (const ECLCluster& cluster : m_eclClusters) {
429 // Use the photon hypothesis, and only look at clusters with times < 200 ns
430 const double time = cluster.getTime();
431 if (not cluster.hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons) or
432 std::abs(time) > 200 or
433 cluster.getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons) < 0.1) {
434 continue;
435 }
436
437 const double dt99 = cluster.getDeltaTime99();
438
439 // Store cluster information
440 SelectedECLCluster selectedCluster;
441 selectedCluster.p4Lab = cUtil.Get4MomentumFromCluster(&cluster, clustervertex,
443 selectedCluster.p4CMS = boostrotate.rotateLabToCms() * selectedCluster.p4Lab; // was clustp4COM
444 selectedCluster.cluster = &cluster;
445 selectedCluster.clusterTime = time / dt99; // was clustT
446 selectedCluster.energyLab = selectedCluster.p4Lab.E();
447 selectedCluster.energyCMS = selectedCluster.p4CMS.E(); // was first of ECOMPair
448 selectedCluster.isTrack = cluster.isTrack(); // was tempCharge
449
450 selectedClusters.push_back(selectedCluster);
451
452 if (selectedCluster.energyCMS > m_E2min) {
453 calculationResult["nElow"] += 1;
454 }
455 if (selectedCluster.energyCMS > m_E0min) {
456 calculationResult["nEmedium"] += 1;
457 }
458 if (selectedCluster.energyCMS > m_Ehigh) {
459 calculationResult["nEhigh"] += 1;
460 }
461
462 if (selectedCluster.energyCMS > m_E0min and std::abs(selectedCluster.clusterTime) < 10) {
463 calculationResult["nVetoClust"] += 1;
464 }
465
466
467 //..Single cluster trigger objects use charge, Zernike moment, and thetaLab
468 const double thetaLab = selectedCluster.p4Lab.Theta() * TMath::RadToDeg();
469 const double zmva = cluster.getZernikeMVA();
470 const bool photon = zmva > 0.5 and not selectedCluster.isTrack;
471 const bool electron = zmva > 0.5 and selectedCluster.isTrack;
472
473 //..For 1 track radiative Bhabha control sample
474 if (selectedCluster.energyCMS > 2. and selectedCluster.isTrack) {
475 calculationResult["chrgClust2GeV"] += 1;
476 }
477 if (selectedCluster.energyCMS > 0.45 and not selectedCluster.isTrack) {
478 const bool isInAcceptance = 17. < thetaLab and thetaLab < 150.;
479 if (isInAcceptance) {calculationResult["neutClust045GeVAcc"] += 1;}
480 const bool isInBarrel = 30. < thetaLab and thetaLab < 130.;
481 if (isInBarrel) {calculationResult["neutClust045GeVBarrel"] += 1;}
482 }
483
484
485 // improved 3 cluster (4 cluster) trigger
486 const bool notInHighBackgroundEndcapRegion = 18.5 < thetaLab and thetaLab < 139.3;
487 if (selectedCluster.energyLab > m_EminLab and notInHighBackgroundEndcapRegion) {
488 calculationResult["nE180Lab"] += 1;
489 }
490
491 if (selectedCluster.energyLab > m_EminLab4Cluster and notInHighBackgroundEndcapRegion) {
492 calculationResult["nE300Lab"] += 1;
493 }
494
495 if (selectedCluster.energyLab > m_EminLab3Cluster and notInHighBackgroundEndcapRegion) {
496 calculationResult["nE500Lab"] += 1;
497 }
498
499 if (selectedCluster.energyCMS > m_Ehigh and notInHighBackgroundEndcapRegion) {
500 calculationResult["nE2000CMS"] += 1;
501 }
502
503 //..For two-photon fusion ALP trigger
504 if (selectedCluster.energyLab > 0.25) {
505 calculationResult["nE250Lab"] += 1;
506 }
507 if (selectedCluster.energyCMS > 4.) {
508 calculationResult["nE4000CMS"] += 1;
509 }
510
511 //..Single cluster triggers
512 if (selectedCluster.energyCMS > m_EsinglePhoton) {
513 calculationResult["nEsingleClust"] += 1;
514
515 const bool barrelRegion = thetaLab > 45 and thetaLab < 115;
516 const bool extendedBarrelRegion = thetaLab > 30 and thetaLab < 130;
517 const bool endcapRegion = (thetaLab > 22 and thetaLab < 45) or (thetaLab > 115 and thetaLab < 145);
518
519 if (photon and barrelRegion) {
520 calculationResult["nEsinglePhotonBarrel"] += 1;
521 }
522
523 if (photon and extendedBarrelRegion) {
524 calculationResult["nEsinglePhotonExtendedBarrel"] += 1;
525 }
526
527 if (electron and barrelRegion) {
528 calculationResult["nEsingleElectronBarrel"] += 1;
529 }
530
531 if (electron and extendedBarrelRegion) {
532 calculationResult["nEsingleElectronExtendedBarrel"] += 1;
533 }
534
535 if (photon and endcapRegion) {
536 calculationResult["nEsinglePhotonEndcap"] += 1;
537 }
538 }
539
540 if (selectedCluster.energyCMS > m_reducedEsinglePhoton) {
541 const bool reducedBarrelRegion = thetaLab > 44 and thetaLab < 98;
542
543 if (photon and reducedBarrelRegion) {
544 calculationResult["nReducedEsinglePhotonReducedBarrel"] += 1;
545 }
546 }
547
548 if (selectedCluster.energyCMS > m_Ehigh) {
549 // TODO: different definition!
550 const bool barrelRegion = thetaLab > 32 and thetaLab < 130;
551 const bool endcapRegion = (thetaLab > 22 and thetaLab < 32) or (thetaLab > 130 and thetaLab < 145);
552 const bool lowAngleRegion = thetaLab < 22 or thetaLab > 145;
553
554 if (not selectedCluster.isTrack and barrelRegion) {
555 calculationResult["n2GeVNeutBarrel"] += 1;
556 }
557 if (not selectedCluster.isTrack and endcapRegion) {
558 calculationResult["n2GeVNeutEndcap"] += 1;
559 }
560 if (selectedCluster.isTrack and not lowAngleRegion) {
561 calculationResult["n2GeVChrg"] += 1;
562 }
563 if (lowAngleRegion) {
564 calculationResult["nEhighLowAng"] += 1;
565 }
566 if (photon and barrelRegion) {
567 calculationResult["n2GeVPhotonBarrel"] += 1;
568 }
569 if (photon and endcapRegion) {
570 calculationResult["n2GeVPhotonEndcap"] += 1;
571 }
572 }
573 } // end of loop over clusters
574
575 //..Order clusters by CMS energy
576 std::sort(selectedClusters.begin(), selectedClusters.end(), [](const auto & lhs, const auto & rhs) {
577 return lhs.energyCMS > rhs.energyCMS;
578 });
579
580 // -- Bhabha and two-photon lepton preparation -- //
581 for (short charge : { -1, 1}) {
582 auto& maximumPtTrack = maximumPtTracks.at(charge);
583 if (not maximumPtTrack) {
584 continue;
585 }
586
587 // Single leg Bhabha selections
588 if (maximumPtTrack->clusterEnergySumCMS > 4.5) {
589 calculationResult["ee1leg"] = 1;
590 }
591
592 // single muon
593 if (maximumPtTrack->pCMS > 3 and maximumPtTrack->clusterEnergySumLab > 0. and maximumPtTrack->clusterEnergySumLab < 1.) {
594 calculationResult["singleMuon"] = 1;
595 }
596 }
597
598 // Bhabha selections using two tracks
599 if (maximumPtTracks.at(-1) and maximumPtTracks.at(1)) {
600 const MaximumPtTrack& negativeTrack = *maximumPtTracks.at(-1);
601 const MaximumPtTrack& positiveTrack = *maximumPtTracks.at(1);
602
603 double dphi = std::abs(negativeTrack.p4CMS.Phi() - positiveTrack.p4CMS.Phi()) * TMath::RadToDeg();
604 if (dphi > 180) {
605 dphi = 360 - dphi;
606 }
607
608 const double negativeClusterSum = negativeTrack.clusterEnergySumCMS;
609 const double negativeClusterSumLab = negativeTrack.clusterEnergySumLab;
610 const double negativeP = negativeTrack.pCMS;
611 const double positiveClusterSum = positiveTrack.clusterEnergySumCMS;
612 const double positiveClusterSumLab = positiveTrack.clusterEnergySumLab;
613 const double positiveP = positiveTrack.pCMS;
614
615 // two leg Bhabha
616 const double thetaSum = (negativeTrack.p4CMS.Theta() + positiveTrack.p4CMS.Theta()) * TMath::RadToDeg();
617 const double dthetaSum = std::abs(thetaSum - 180);
618 const auto back2back = dphi > 175 and dthetaSum < 15;
619 if (back2back and negativeClusterSum > 3 and positiveClusterSum > 3 and
620 (negativeClusterSum > 4.5 or positiveClusterSum > 4.5)) {
621 calculationResult["ee2leg"] = 1;
622 }
623
624 // one leg one track Bhabha
625 if (back2back and ((negativeClusterSum > 4.5 and positiveP > 3) or (positiveClusterSum > 4.5 and negativeP > 3))) {
626 calculationResult["ee1leg1trk"] = 1;
627 }
628
629
630 // one leg plus one electron
631 if ((negativeClusterSum > 4.5 and positiveClusterSum > 0.8 * positiveP) or
632 (positiveClusterSum > 4.5 and negativeClusterSum > 0.8 * negativeP)) {
633 calculationResult["ee1leg1e"] = 1;
634 }
635
636 // eeee, eemm, mu mu, and radiative Bhabha selection
637 const ROOT::Math::PxPyPzEVector p4Miss = p4ofCOM - negativeTrack.p4CMS - positiveTrack.p4CMS;
638 const double pmissTheta = p4Miss.Theta() * TMath::RadToDeg();
639 const double pmissp = p4Miss.P();
640 const double relMissAngle0 = B2Vector3D(negativeTrack.p4CMS.Vect()).Angle(B2Vector3D(p4Miss.Vect())) * TMath::RadToDeg();
641 const double relMissAngle1 = B2Vector3D(positiveTrack.p4CMS.Vect()).Angle(B2Vector3D(p4Miss.Vect())) * TMath::RadToDeg();
642
643 const bool electronEP = positiveClusterSum > 0.8 * positiveP or negativeClusterSum > 0.8 * negativeP;
644 const bool notMuonPair = negativeClusterSumLab > 1 or positiveClusterSumLab > 1;
645 const double highp = std::max(negativeP, positiveP);
646 const double lowp = std::min(negativeP, positiveP);
647 const bool lowEdep = negativeClusterSumLab < 0.5 and positiveClusterSumLab < 0.5;
648
649 //..Select two photon fusion lepton pairs
650 if (calculationResult["maximumPCMS"] < 2 and dphi > 160 and (pmissTheta < 25. or pmissTheta > 155.)) {
651 calculationResult["eexxSelect"] = 1;
652 if (electronEP) {
653 calculationResult["eeee"] = 1;
654 } else {
655 calculationResult["eemm"] = 1;
656 }
657 }
658
659 //..Veto two photon fusion lepton pairs
660 if ((pmissTheta < 20. or pmissTheta > 160.) and
661 ((calculationResult["maximumPCMS"] < 1.2 and dphi > 150.) or
662 (calculationResult["maximumPCMS"] < 2. and 175. < dphi))) {
663 calculationResult["eexx"] = 1;
664 }
665
666 //..Radiative Bhabhas with missing momentum in acceptance
667 if (negativeP > 1. and pmissTheta > 10. and pmissTheta < 170. and positiveP > 1. and dphi < 170. and pmissp > 1. and electronEP) {
668 if (calculationResult["nTrkLoose"] == 2 and calculationResult["nTrkTight"] >= 1) {
669 calculationResult["radBhabha"] = 1;
670 }
671 if (calculationResult["nTrkLooseB"] == 2 and calculationResult["nTrkTightB"] >= 1) {
672 calculationResult["radBhabhaB"] = 1;
673 }
674 }
675
676 //..ISR radiative Bhabha with missing momentum down beam pipe
677 if (negativeP > 2. and positiveP > 2. and 2 == calculationResult["nTrkLoose"] and
678 calculationResult["nTrkTight"] >= 1 and dphi > 175. and
679 (pmissTheta < 5. or pmissTheta > 175.) and electronEP) {
680 calculationResult["isrRadBhabha"] = 1;
681 }
682
683 //..Veto Bhabha with hard brem
684 if (highp > 4.5 and notMuonPair and pmissp > 1. and (relMissAngle0 < 5. or relMissAngle1 < 5.)) {
685 if (calculationResult["nTrkLoose"] == 2) { calculationResult["eeBrem"] = 1;}
686 if (calculationResult["nTrkLooseB"] >= 1) { calculationResult["eeBremB"] = 1;}
687 }
688
689 //..Veto e+e- -> mu+mu-
690 if (highp > 4.5 and lowEdep and thetaSum > 175. and thetaSum < 185. and dphi > 175.) {
691 if (calculationResult["nTrkLoose"] == 2) {calculationResult["muonPairV"] = 1;}
692 if (calculationResult["nTrkLooseB"] == 2) {calculationResult["muonPairVB"] = 1;}
693 }
694
695 //..Select e+e- -> mu+mu-
696 if (highp > 3. and lowp > 2.5 and dphi > 165. and
697 ((negativeClusterSumLab > 0. and negativeClusterSumLab < 1.) or
698 (positiveClusterSumLab > 0. and positiveClusterSumLab < 1.))) {
699 calculationResult["selectmumu"] = 1;
700 }
701 }
702
703 //..Filters and vetoes using two maximum-energy clusters
704 if (selectedClusters.size() >= 2) {
705 const SelectedECLCluster& firstCluster = selectedClusters[0];
706 const SelectedECLCluster& secondCluster = selectedClusters[1];
707
708 double dphi = std::abs(firstCluster.p4CMS.Phi() - secondCluster.p4CMS.Phi()) * TMath::RadToDeg();
709 if (dphi > 180) {
710 dphi = 360 - dphi;
711 }
712 double thetaSum = (firstCluster.p4CMS.Theta() + secondCluster.p4CMS.Theta()) * TMath::RadToDeg();
713 double dthetaSum = std::abs(thetaSum - 180);
714
715 //..Quantities for two-photon fusion triggers
716 calculationResult["dphiCmsClust"] = dphi;
717 for (int ic = 0; ic < 2; ic++) {
718 const double thetaLab = selectedClusters[ic].p4Lab.Theta() * TMath::RadToDeg();
719 const bool isInAcceptance = 17. < thetaLab and thetaLab < 150.;
720 const ECLCluster* cluster = selectedClusters[ic].cluster;
721 const double zmva = cluster->getZernikeMVA();
722 const bool photon = zmva > 0.5 and not selectedClusters[ic].isTrack;
723 if (isInAcceptance and photon) {calculationResult["nMaxEPhotonAcc"] += 1;}
724 }
725
726 const double firstEnergy = firstCluster.p4CMS.E();
727 const double secondEnergy = secondCluster.p4CMS.E();
728
729 const bool highEnergetic = firstEnergy > 3 and secondEnergy > 3 and (firstEnergy > 4.5 or secondEnergy > 4.5);
730
731 if (dphi > 160 and dphi < 178 and dthetaSum < 15 and highEnergetic) {
732 calculationResult["ee2clst"] = 1;
733 }
734
735 if (dphi > 178 and dthetaSum < 15 and highEnergetic) {
736 calculationResult["gg2clst"] = 1;
737 }
738
739 if ((calculationResult["ee2clst"] == 1 or calculationResult["gg2clst"] == 1) and
740 calculationResult["ee1leg"] == 1) {
741 calculationResult["ee1leg1clst"] = 1;
742 }
743
744 const double Elab0 = firstCluster.p4Lab.E();
745 const double Elab1 = secondCluster.p4Lab.E();
746
747 // gg and mumu accept triggers using ECL
748 if (firstEnergy > 2 and secondEnergy > 2) {
749 const double thetaLab0 = firstCluster.p4Lab.Theta() * TMath::RadToDeg();
750 const double thetaLab1 = secondCluster.p4Lab.Theta() * TMath::RadToDeg();
751
752 const bool barrel0 = 32. < thetaLab0 and thetaLab0 < 130.;
753 const bool barrel1 = 32. < thetaLab1 and thetaLab1 < 130.;
754 const bool oneClustersAbove4 = firstEnergy > 4 or secondEnergy > 4;
755 const bool oneIsNeutral = not firstCluster.isTrack or not secondCluster.isTrack;
756 const bool bothAreNeutral = not firstCluster.isTrack and not secondCluster.isTrack;
757 const bool oneIsBarrel = barrel0 or barrel1;
758 const bool dphiCutExtraLoose = dphi > 175;
759 const bool dphiCutLoose = dphi > 177;
760 const bool dphiCutTight = dphi > 177.5;
761
762 if (dphiCutExtraLoose and oneIsNeutral and oneIsBarrel) {
763 calculationResult["ggBarrelVL"] = 1;
764 }
765 if (oneClustersAbove4 and dphiCutLoose and oneIsNeutral and oneIsBarrel) {
766 calculationResult["ggBarrelLoose"] = 1;
767 }
768 if (oneClustersAbove4 and dphiCutTight and bothAreNeutral and oneIsBarrel) {
769 calculationResult["ggBarrelTight"] = 1;
770 }
771 if (dphiCutExtraLoose and oneIsNeutral and not oneIsBarrel) {
772 calculationResult["ggEndcapVL"] = 1;
773 }
774 if (oneClustersAbove4 and dphiCutLoose and oneIsNeutral and not oneIsBarrel) {
775 calculationResult["ggEndcapLoose"] = 1;
776 }
777 if (oneClustersAbove4 and dphiCutTight and bothAreNeutral and not oneIsBarrel) {
778 calculationResult["ggEndcapTight"] = 1;
779 }
780 }
781
782 const double minEnergy = std::min(Elab0, Elab1);
783 const double maxEnergy = std::max(Elab0, Elab1);
784 if (dphi > 155 and thetaSum > 165 and thetaSum < 195 and minEnergy > 0.15 and minEnergy < 0.5 and
785 maxEnergy > 0.15 and maxEnergy < 0.5) {
786 calculationResult["muonPairECL"] = 1;
787 }
788
789 //..diPhoton line
790 const double thetaLab0 = firstCluster.p4Lab.Theta() * TMath::RadToDeg();
791 const double thetaLab1 = secondCluster.p4Lab.Theta() * TMath::RadToDeg();
792 const bool inHieRegion0 = thetaLab0 > 26. and thetaLab0 < 130.;
793 const bool inHieRegion1 = thetaLab1 > 26. and thetaLab1 < 130.;
794 const bool firstIsNeutral = not firstCluster.isTrack;
795 const bool secondIsNeutral = not secondCluster.isTrack;
796
797 if (secondEnergy > 0.3 and inHieRegion0 and inHieRegion1 and firstIsNeutral and secondIsNeutral) {
798 const ROOT::Math::PxPyPzEVector ggP4CMS = firstCluster.p4CMS + secondCluster.p4CMS;
799 if (ggP4CMS.pt() > 1.) {calculationResult["ggHighPt"] = 1;}
800 }
801
802 } // end of two-clusters lines
803
804 // Bhabha accept triggers.
805 // Use theta_lab of negative track if available; otherwise cluster.
806 double thetaFlatten = 0;
807
808 // One leg, one cluster.
809 for (short charge : { -1, 1}) {
810 const auto& maximumPtTrack = maximumPtTracks.at(charge);
811 if (not maximumPtTrack) {
812 continue;
813 }
814
815 if (maximumPtTrack->clusterEnergySumCMS > 1.5) {
816 double invMass = 0.;
817 double tempFlatten = 0.;
818 for (const SelectedECLCluster& cluster : selectedClusters) {
819 double tempInvMass = (maximumPtTrack->p4Lab + cluster.p4Lab).M();
820 if (tempInvMass > invMass) {
821 invMass = tempInvMass;
822 if (charge == 1) {
823 tempFlatten = cluster.p4Lab.Theta() * TMath::RadToDeg();
824 }
825 }
826 }
827 if (charge == -1) {
828 tempFlatten = maximumPtTrack->p4Lab.Theta() * TMath::RadToDeg();
829 }
830 if (invMass > 5.29) {
831 calculationResult["selectee1leg1clst"] = 1;
832 thetaFlatten = tempFlatten;
833 }
834 }
835 }
836
837 // Two legs
838 if (maximumPtTracks.at(-1) and maximumPtTracks.at(1)) {
839 const MaximumPtTrack& negativeTrack = *maximumPtTracks.at(-1);
840 const MaximumPtTrack& positiveTrack = *maximumPtTracks.at(1);
841 const double invMass = (negativeTrack.p4Lab + positiveTrack.p4Lab).M();
842 if (invMass > 5.29 and (negativeTrack.clusterEnergySumCMS > 1.5 or positiveTrack.clusterEnergySumCMS > 1.5)) {
843 calculationResult["selectee1leg1trk"] = 1;
844 }
845
846 thetaFlatten = negativeTrack.p4Lab.Theta() * TMath::RadToDeg();
847
848 // Two tracks but (exactly) 1 high energy cluster
849 const bool missNegClust = positiveTrack.clusterEnergySumCMS > 4.5 and negativeTrack.clusterEnergySumCMS < 1.;
850 const bool missPosClust = negativeTrack.clusterEnergySumCMS > 4.5 and positiveTrack.clusterEnergySumCMS < 1.;
851 if ((invMass > 9.) and (missNegClust or missPosClust)) {
852 calculationResult["eeOneClust"] = 1;
853 }
854 }
855
856 if (calculationResult["selectee1leg1trk"] == 1 or calculationResult["selectee1leg1clst"] == 1) {
857 for (int iflat = 0; iflat < 9; iflat++) {
858 const std::string& eeFlatName = "eeFlat" + std::to_string(iflat);
859 calculationResult[eeFlatName] =
860 thetaFlatten >= flatBoundaries[iflat] and thetaFlatten < flatBoundaries[iflat + 1];
861 if (calculationResult[eeFlatName]) {
862 calculationResult["selectee"] = 1;
863 }
864 }
865 }
866
867 // Single tag pi0 / eta dedicated lines
868 if (calculationResult["nTrkLoose"] == 1 and calculationResult["maximumPCMS"] > 0.8 and selectedClusters.size() >= 2) {
869
870 decltype(selectedClusters) selectedSingleTagClusters(selectedClusters.size());
871 auto lastItem = std::copy_if(selectedClusters.begin(), selectedClusters.end(), selectedSingleTagClusters.begin(),
872 [](auto & cluster) {
873 const bool isNeutralCluster = not cluster.isTrack;
874 const bool hasEnoughEnergy = cluster.energyLab > 0.1;
875 const double clusterThetaLab = cluster.p4Lab.Theta() * TMath::RadToDeg();
876 const bool isInAcceptance = 17 < clusterThetaLab and clusterThetaLab < 150.;
877 return isNeutralCluster and hasEnoughEnergy and isInAcceptance;
878 });
879 selectedSingleTagClusters.resize(std::distance(selectedSingleTagClusters.begin(), lastItem));
880
881 if (selectedSingleTagClusters.size() >= 2) { // One track and at least two clusters are found
882
883 const auto& track = maximumPtTracks.at(-1) ? *maximumPtTracks.at(-1) : *maximumPtTracks.at(1);
884 // in real signal, the pi0 decay daughters are always the two most-energetic neutral clusters.
885 const auto& firstCluster = selectedSingleTagClusters[0];
886 const auto& secondCluster = selectedSingleTagClusters[1];
887
888 const ROOT::Math::PxPyPzEVector trackP4CMS = track.p4CMS;
889 const ROOT::Math::PxPyPzEVector pi0P4CMS = firstCluster.p4CMS + secondCluster.p4CMS;
890
891 const bool passPi0ECMS = pi0P4CMS.E() > 1. and pi0P4CMS.E() < 0.525 * p4ofCOM.M();
892 const double thetaSumCMS = (pi0P4CMS.Theta() + trackP4CMS.Theta()) * TMath::RadToDeg();
893 const bool passThetaSum = thetaSumCMS < 170. or thetaSumCMS > 190.;
894
895 double dphiCMS = std::abs(trackP4CMS.Phi() - pi0P4CMS.Phi()) * TMath::RadToDeg();
896 if (dphiCMS > 180) {
897 dphiCMS = 360 - dphiCMS;
898 }
899 const bool passdPhi = dphiCMS > 160.;
900
901 if (passPi0ECMS and passThetaSum and passdPhi and pi0P4CMS.M() < 0.7) {
902 calculationResult["singleTagLowMass"] = 1;
903 } else if (passPi0ECMS and passThetaSum and passdPhi and pi0P4CMS.M() > 0.7) {
904 calculationResult["singleTagHighMass"] = 1;
905 }
906 }
907 }
908
909 //..Updated to use new track definitions
910 if (calculationResult["nTrkLooseB"] == 1 and calculationResult["maximumPCMSB"] > 0.8 and selectedClusters.size() >= 2) {
911
912 decltype(selectedClusters) selectedSingleTagClusters(selectedClusters.size());
913 auto lastItem = std::copy_if(selectedClusters.begin(), selectedClusters.end(), selectedSingleTagClusters.begin(),
914 [](auto & cluster) {
915 const bool isNeutralCluster = not cluster.isTrack;
916 const bool hasEnoughEnergy = cluster.energyLab > 0.1;
917 const double clusterThetaLab = cluster.p4Lab.Theta() * TMath::RadToDeg();
918 const bool isInAcceptance = 17 < clusterThetaLab and clusterThetaLab < 150.;
919 return isNeutralCluster and hasEnoughEnergy and isInAcceptance;
920 });
921 selectedSingleTagClusters.resize(std::distance(selectedSingleTagClusters.begin(), lastItem));
922
923 if (selectedSingleTagClusters.size() >= 2) { // One track and at least two clusters are found
924
925 const auto& track = maximumPtTracks.at(-1) ? *maximumPtTracks.at(-1) : *maximumPtTracks.at(1);
926 // in real signal, the pi0 decay daughters are always the two most-energetic neutral clusters.
927 const auto& firstCluster = selectedSingleTagClusters[0];
928 const auto& secondCluster = selectedSingleTagClusters[1];
929
930 const ROOT::Math::PxPyPzEVector trackP4CMS = track.p4CMS;
931 const ROOT::Math::PxPyPzEVector pi0P4CMS = firstCluster.p4CMS + secondCluster.p4CMS;
932
933 const bool passPi0ECMS = pi0P4CMS.E() > 1. and pi0P4CMS.E() < 0.525 * p4ofCOM.M();
934 const double thetaSumCMS = (pi0P4CMS.Theta() + trackP4CMS.Theta()) * TMath::RadToDeg();
935 const bool passThetaSum = thetaSumCMS < 170. or thetaSumCMS > 190.;
936
937 double dphiCMS = std::abs(trackP4CMS.Phi() - pi0P4CMS.Phi()) * TMath::RadToDeg();
938 if (dphiCMS > 180) {
939 dphiCMS = 360 - dphiCMS;
940 }
941 const bool passdPhi = dphiCMS > 160.;
942
943 if (passPi0ECMS and passThetaSum and passdPhi and pi0P4CMS.M() < 0.7) {
944 calculationResult["singleTagLowMassB"] = 1;
945 } else if (passPi0ECMS and passThetaSum and passdPhi and pi0P4CMS.M() > 0.7) {
946 calculationResult["singleTagHighMassB"] = 1;
947 }
948 }
949 }
950
951 //..Cosmic selection. Need exactly two tracks.
952 if (m_tracks.getEntries() == 2) {
953
954 const auto negTrack = maximumPtTracksWithoutZCut.at(-1);
955 const auto posTrack = maximumPtTracksWithoutZCut.at(1);
956
957 if (negTrack and posTrack) {
958
959 const double maxNegpT = negTrack->pT;
960 const double maxPospT = posTrack->pT;
961 const double maxClusterENeg = negTrack->clusterEnergySumLab;
962 const double maxClusterEPos = posTrack->clusterEnergySumLab;
963
964 const ROOT::Math::PxPyPzEVector& momentumLabNeg(negTrack->p4Lab);
965 const ROOT::Math::PxPyPzEVector& momentumLabPos(posTrack->p4Lab);
966
967 const double& z0Neg = negTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getZ0();
968 const double& d0Neg = negTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getD0();
969 const double& z0Pos = posTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getZ0();
970 const double& d0Pos = posTrack->track->getTrackFitResultWithClosestMass(Const::pion)->getD0();
971
972 // Select cosmic using these tracks
973 const bool goodMagneticRegion = (z0Neg<m_goodMagneticRegionZ0 or abs(d0Neg)>m_goodMagneticRegionD0)
974 and (z0Pos<m_goodMagneticRegionZ0 or abs(d0Pos)>m_goodMagneticRegionD0);
975 if (maxNegpT > m_cosmicMinPt and maxPospT > m_cosmicMinPt and maxClusterENeg < m_cosmicMaxClusterEnergy
976 and maxClusterEPos < m_cosmicMaxClusterEnergy and goodMagneticRegion) {
977 double dphiLab = std::abs(momentumLabNeg.Phi() - momentumLabPos.Phi()) * TMath::RadToDeg();
978 if (dphiLab > 180) {
979 dphiLab = 360 - dphiLab;
980 }
981
982 const double thetaSumLab = (momentumLabNeg.Theta() + momentumLabPos.Theta()) * TMath::RadToDeg();
983
984 constexpr double phiBackToBackTolerance = 2.;
985 constexpr double thetaBackToBackTolerance = 2.;
986 if ((180 - dphiLab) < phiBackToBackTolerance and std::abs(180 - thetaSumLab) < thetaBackToBackTolerance) {
987 calculationResult["cosmic"] = 1;
988 }
989 }
990 }
991 }
992
993 //..Displaced vertex
994 // See https://indico.belle2.org/event/8973/ and references therein
995 if (maximumPtTracksWithoutZCut.at(-1) and maximumPtTracksWithoutZCut.at(1)) {
996
997 //..Make particles of the two highest pt tracks (without IP cut)
998 const auto nTrack = maximumPtTracksWithoutZCut.at(-1)->track;
999 Particle* nParticle = new Particle(nTrack, Const::pion);
1000
1001 const auto pTrack = maximumPtTracksWithoutZCut.at(1)->track;
1002 Particle* pParticle = new Particle(pTrack, Const::pion);
1003
1004 //..Make a vertex of these two
1006 vertexFit->addParticle(nParticle);
1007 vertexFit->addParticle(pParticle);
1008 vertexFit->doFit();
1009
1010 //..Vertex properties
1011 const double chisq = vertexFit->getCHIsq();
1012 const int ndf = vertexFit->getNDF();
1013 const double vertexProb = TMath::Prob(chisq, ndf);
1014 const auto vertexLocation = vertexFit->getVertex();
1015 const double vertexXY = vertexLocation.perp();
1016 const double vertexTheta = vertexLocation.theta() * TMath::RadToDeg();
1017
1018 //..Angular differance of two tracks to reject cosmics
1019 // Tolerance could be reduced from 10 deg to 2 deg if needed for physics reasons,
1020 // for a 5% increase in the rate of selected displaced vertex triggers.
1021 // See https://gitlab.desy.de/belle2/software/basf2/-/merge_requests/1867
1022 const ROOT::Math::PxPyPzEVector& momentumLabNeg(maximumPtTracksWithoutZCut.at(-1)->p4Lab);
1023 const ROOT::Math::PxPyPzEVector& momentumLabPos(maximumPtTracksWithoutZCut.at(1)->p4Lab);
1024 const double thetaSumLab = (momentumLabNeg.Theta() + momentumLabPos.Theta()) * TMath::RadToDeg();
1025 double dPhiLab = std::abs(momentumLabNeg.Phi() - momentumLabPos.Phi()) * TMath::RadToDeg();
1026 if (dPhiLab > 180) {
1027 dPhiLab = 360 - dPhiLab;
1028 }
1029 const double backToBackTolerance = 10.; // degrees
1030 const bool backToBackLab = std::abs(thetaSumLab - 180.) < backToBackTolerance and std::abs(dPhiLab - 180.) < backToBackTolerance;
1031
1032 //..Select a displaced vertex
1033 const double minProbChiVertex = 0.005; // minimum probability chi square. Many backgrounds have bad chi sq.
1034 const double minXYVertex = 3.; // minimum xy of vertex (cm). Should pass track filters below this.
1035 const double maxXYVertex = 60.; // maximum xy. Insufficient CDC for good reconstruction beyond this.
1036 const double minThetaVertex = 30.; // Large Bhabha background at low angles.
1037 const double maxThetaVertex = 120.; // Large Bhabha background at low angles.
1038 if (vertexProb > minProbChiVertex and vertexXY > minXYVertex and vertexXY < maxXYVertex and vertexTheta > minThetaVertex
1039 and vertexTheta < maxThetaVertex
1040 and not backToBackLab) {calculationResult["displacedVertex"] = 1;}
1041
1042 //..Clean up
1043 delete nParticle;
1044 delete pParticle;
1045 delete vertexFit;
1046
1047 }
1048
1049}
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
Definition: B2Vector3.h:302
Class to provide momentum-related information from ECLClusters.
Definition: ClusterUtils.h:36
const ROOT::Math::PxPyPzEVector Get4MomentumFromCluster(const ECLCluster *cluster, ECLCluster::EHypothesisBit hypo)
Returns four momentum vector.
Definition: ClusterUtils.cc:25
const ROOT::Math::XYZVector GetIPPosition()
Returns default IP position from beam parameters.
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
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 ROOT::Math::LorentzRotation rotateLabToCms() const
Returns Lorentz transformation from Lab to CMS.
Class to store reconstructed particles.
Definition: Particle.h:75
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
@ TTYP_DPHY
delayed physics events for background
Definition: TRGSummary.h:65
@ TTYP_POIS
poisson random trigger
Definition: TRGSummary.h:73
@ TTYP_RAND
random trigger events
Definition: TRGSummary.h:67
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
ROOT::Math::PxPyPzEVector 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
virtual int getNDF(void) const
Get an NDF of the fit.
Definition: KFitBase.cc:114
enum KFitError::ECode addParticle(const Particle *particle)
Add a particle to the fitter.
Definition: KFitBase.cc:59
VertexFitKFit is a derived class from KFitBase to perform vertex-constraint kinematical fit.
Definition: VertexFitKFit.h:34
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode doFit(void)
Perform a vertex-constraint fit.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
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
ROOT::Math::PxPyPzEVector p4CMS
the 4 momentum in CMS system
double clusterEnergySumLab
the sum of related cluster energies in lab system
ROOT::Math::PxPyPzEVector p4Lab
the 4 momentum in lab system
double clusterEnergySumCMS
the sum of related cluster energies in CMS system
double pCMS
the momentum magnitude in CMS 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
ROOT::Math::PxPyPzEVector p4CMS
the 4 momentum in CMS system
ROOT::Math::PxPyPzEVector 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