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