Belle II Software development
FilterCalculator Class Reference

Implementation of a calculator used in the SoftwareTriggerModule to fill a SoftwareTriggerObject for doing HLT cuts. More...

#include <FilterCalculator.h>

Inheritance diagram for FilterCalculator:
SoftwareTriggerCalculation

Public Member Functions

 FilterCalculator ()
 Set the default names for the store object particle lists.
 
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.
 
void writeDebugOutput (const std::unique_ptr< TTree > &debugOutputTTree)
 Function to write out debug output into the given TTree.
 
void addDebugOutput (const StoreObjPtr< SoftwareTriggerVariables > &storeObject, const std::string &prefix)
 Function to write out debug output into the given StoreObject.
 
const SoftwareTriggerObject & fillInCalculations ()
 Main function of this class: calculate the needed variables using the overwritten doCalculation function and write out the values into the results object (with their names).
 

Private Attributes

StoreArray< Trackm_tracks
 Store Array of the tracks to be used.
 
StoreArray< ECLClusterm_eclClusters
 Store Array of the ecl clusters to be used.
 
StoreObjPtr< TRGSummarym_l1Trigger
 Store Object with the trigger result.
 
StoreArray< CDCTriggerUnpacker::NNBitStreamm_bitsNN
 Store Object with the trigger NN bits.
 
double m_looseTrkZ0 = 10 * Unit::cm
 which Z0 defines a loose track
 
double m_tightTrkZ0 = 2 * Unit::cm
 which Z0 defines a tight track
 
double m_E2min = 0.2
 which CMS energy defines nElow
 
double m_E0min = 0.3
 which CMS energy defines nEmedium
 
double m_Ehigh = 2
 which CMS energy defines nEhigh
 
double m_EminLab = 0.18
 which lab energy defines nE180Lab
 
double m_EminLab4Cluster = 0.3
 which lab energy defines nE300Lab
 
double m_EminLab3Cluster = 0.5
 which lab energy defines nE500Lab
 
double m_EsinglePhoton = 1
 which CMS energy defines nEsingleClust
 
double m_reducedEsinglePhoton = 0.5
 which CMS energy defines nReducedEsingle clusters
 
double m_cosmicMinPt = 0.5 * Unit::GeV
 which LAB pt defines a cosmic
 
double m_cosmicMaxClusterEnergy = 1.0 * Unit::GeV
 which LAB cluster energy vetoes a cosmic candidate
 
double m_goodMagneticRegionZ0 = 57.
 maximum z0 for well understood magnetic field (cm)
 
double m_goodMagneticRegionD0 = 26.5
 minimum d0 for well understood magnetic field, if z0 is large (cm)
 
DBObjPtr< HLTPrefilterParametersm_hltPrefilterParameters
 Objects relevant to HLTPrefilter monitoring HLTprefilterParameters Database OjbPtr.
 
HLTPrefilter::TimingCutState m_timingPrefilter
 Helper instance for timing based prefilter.
 
HLTPrefilter::CDCECLCutState m_cdceclPrefilter
 Helper instance for CDC-ECL occupancy based prefilter.
 
SoftwareTriggerObject m_calculationResult
 Internal storage of the result of the calculation.
 
bool m_debugPrepared = false
 Flag to not add the branches twice to the TTree.
 

Detailed Description

Implementation of a calculator used in the SoftwareTriggerModule to fill a SoftwareTriggerObject for doing HLT cuts.

This calculator exports variables needed for the trigger HLT part of the path ( = filtering out events)

This class implements the two main functions requireStoreArrays and doCalculation of the SoftwareTriggerCalculation class.

Definition at line 34 of file FilterCalculator.h.

Constructor & Destructor Documentation

◆ FilterCalculator()

Set the default names for the store object particle lists.

Definition at line 73 of file FilterCalculator.cc.

73 : m_bitsNN("CDCTriggerNNBits")
74{
75
76}
StoreArray< CDCTriggerUnpacker::NNBitStream > m_bitsNN
Store Object with the trigger NN bits.

Member Function Documentation

◆ addDebugOutput()

void addDebugOutput ( const StoreObjPtr< SoftwareTriggerVariables > & storeObject,
const std::string & prefix )
inherited

Function to write out debug output into the given StoreObject.

Needs an already prefilled calculationResult for this (probably using the fillInCalculations function). All added variables are prefixed with the given prefix string.

Definition at line 34 of file SoftwareTriggerCalculation.cc.

35 {
36 for (auto& identifierWithValue : m_calculationResult) {
37 const std::string& identifier = identifierWithValue.first;
38 const double value = identifierWithValue.second;
39
40 storeObject->append(prefix + "_" + identifier, value);
41 }
42 }

◆ doCalculation()

void doCalculation ( SoftwareTriggerObject & calculationResult)
overridevirtual

Actually write out the variables into the map.

< number of loose tracks

< number of tight tracks

< Bhabha, 2 tracks with accompanying ECL info

< number of clusters with E*>m_Emedium (higher threshold)

< number of clusters with E*>m_Elow (lower threshold)

< number of clusters with E*>m_Ehigh (2 GeV)

< number of clusters with Elab > m_EminLab outside of high background endcap region

< number of clusters with Elab > m_EminLab4Cluster outside of high background endcap region

< number of clusters with Elab > m_EminLab3Cluster outside of high background endcap region

< number of clusters with Ecms > m_Ehigh outside of high background endcap region

< number of clusters with E*>4 GeV

< neutral clusters with Elab>250 MeV (anywhere)

< Neutral, zmva>0.5, in [17,150], max 2 energy clusters

< dphi* between 2 max E clusters

< net charge of loose tracks

< maximum p* of loose tracks (GeV/c)

< maximum pLab of loose tracks (GeV/c)

< Bhabha, 2 tracks, at least 1 of which has ECL info

< number of clusters with E*>m_Ehigh, low angles

< clusters with E*> m_EsinglePhoton (1 GeV)

< neutral clusters with E*> 1 GeV in [45,115]

< neutral clusters with E*> 1 GeV in [32,130]

< neutral clusters with E*> 1 GeV, not barrel or low

< charged clusters with E*> 1 GeV in [45,115]

< charged clusters with E*> 1 GeV in [32,130]

< charged clusters with E*> 0.5 GeV in [44,98]

< clusters with E>m_Emedium and |t|/dt99 < 10

< charged clusters with E*>2 GeV

< neutral clusters with E*>0.45 GeV in [17,150]

< neutral clusters with E*>0.45 GeV in [30,130]

< track + ECL info consistent with Bhabha

< Bhabha, 2 ECL clusters, one of which is charged

< Bhabha, 1 track has high energy cluster, other is electron

< Bhabha, 2 ECL clusters

< gg veto, 2 ECL clusters

< number of looseB tracks

< number of tightB tracks

< maximum pcms of looseB tracks (GeV/c)

< net charge of looseB tracks

< muon pair veto counting looseB tracks

< Bhabha veto hard Brem one leg

< single tag pi0 or eta

< single tag higher mass

< radiative Bhabha, recoil p in detector

< number of "C" tracks

< maximum pcms of "C" tracks GeV/c

< pCms of maximum p negatively charged track GeV/c

< clusterE lab of max p negatively charged track GeV

< pCms of maximum p positively charged track GeV/c

< clusterE lab of max p positively charged track GeV

< delta phi cms of max p positive and negative tracks deg

< Events in the injection strip

< Events with high CDC-ECL occupancy

Implements SoftwareTriggerCalculation.

Definition at line 86 of file FilterCalculator.cc.

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

◆ fillInCalculations()

const SoftwareTriggerObject & fillInCalculations ( )
inherited

Main function of this class: calculate the needed variables using the overwritten doCalculation function and write out the values into the results object (with their names).

Please make sure to override (or clear) the variables! Otherwise it can happen that their old values are still in the object.

What variables exactly are added to the result depends on the implementation details of the class.

Definition at line 44 of file SoftwareTriggerCalculation.cc.

45 {
46 const unsigned int sizeBeforeCheck = m_calculationResult.size();
47 doCalculation(m_calculationResult);
48
49 if (m_calculationResult.size() != sizeBeforeCheck and sizeBeforeCheck > 0) {
50 B2WARNING("Calculator added more variables (" << m_calculationResult.size() <<
51 ") than there were before (" << sizeBeforeCheck << "). Probably something strange is going on!");
52 }
53
54 return m_calculationResult;
55 }

◆ requireStoreArrays()

void requireStoreArrays ( )
overridevirtual

Require the particle list. We do not need more here.

Implements SoftwareTriggerCalculation.

Definition at line 78 of file FilterCalculator.cc.

79{
80 m_tracks.isRequired();
81 m_eclClusters.isRequired();
82 m_l1Trigger.isOptional();
83 m_bitsNN.isOptional();
84}

◆ writeDebugOutput()

void writeDebugOutput ( const std::unique_ptr< TTree > & debugOutputTTree)
inherited

Function to write out debug output into the given TTree.

Needs an already prefilled calculationResult for this (probably using the fillInCalculations function).

Definition at line 19 of file SoftwareTriggerCalculation.cc.

20 {
21 if (not m_debugPrepared) {
22 for (auto& identifierWithValue : m_calculationResult) {
23 const std::string& identifier = identifierWithValue.first;
24 double& value = identifierWithValue.second;
25
26 debugOutputTTree->Branch(identifier.c_str(), &value);
27 }
28 m_debugPrepared = true;
29 }
30
31 debugOutputTTree->Fill();
32 }

Member Data Documentation

◆ m_bitsNN

Store Object with the trigger NN bits.

Definition at line 53 of file FilterCalculator.h.

◆ m_calculationResult

SoftwareTriggerObject m_calculationResult
privateinherited

Internal storage of the result of the calculation.

Definition at line 74 of file SoftwareTriggerCalculation.h.

◆ m_cdceclPrefilter

HLTPrefilter::CDCECLCutState m_cdceclPrefilter
private

Helper instance for CDC-ECL occupancy based prefilter.

Definition at line 90 of file FilterCalculator.h.

◆ m_cosmicMaxClusterEnergy

double m_cosmicMaxClusterEnergy = 1.0 * Unit::GeV
private

which LAB cluster energy vetoes a cosmic candidate

Definition at line 78 of file FilterCalculator.h.

◆ m_cosmicMinPt

double m_cosmicMinPt = 0.5 * Unit::GeV
private

which LAB pt defines a cosmic

Definition at line 76 of file FilterCalculator.h.

◆ m_debugPrepared

bool m_debugPrepared = false
privateinherited

Flag to not add the branches twice to the TTree.

Definition at line 76 of file SoftwareTriggerCalculation.h.

◆ m_E0min

double m_E0min = 0.3
private

which CMS energy defines nEmedium

Definition at line 62 of file FilterCalculator.h.

◆ m_E2min

double m_E2min = 0.2
private

which CMS energy defines nElow

Definition at line 60 of file FilterCalculator.h.

◆ m_eclClusters

StoreArray<ECLCluster> m_eclClusters
private

Store Array of the ecl clusters to be used.

Definition at line 49 of file FilterCalculator.h.

◆ m_Ehigh

double m_Ehigh = 2
private

which CMS energy defines nEhigh

Definition at line 64 of file FilterCalculator.h.

◆ m_EminLab

double m_EminLab = 0.18
private

which lab energy defines nE180Lab

Definition at line 66 of file FilterCalculator.h.

◆ m_EminLab3Cluster

double m_EminLab3Cluster = 0.5
private

which lab energy defines nE500Lab

Definition at line 70 of file FilterCalculator.h.

◆ m_EminLab4Cluster

double m_EminLab4Cluster = 0.3
private

which lab energy defines nE300Lab

Definition at line 68 of file FilterCalculator.h.

◆ m_EsinglePhoton

double m_EsinglePhoton = 1
private

which CMS energy defines nEsingleClust

Definition at line 72 of file FilterCalculator.h.

◆ m_goodMagneticRegionD0

double m_goodMagneticRegionD0 = 26.5
private

minimum d0 for well understood magnetic field, if z0 is large (cm)

Definition at line 82 of file FilterCalculator.h.

◆ m_goodMagneticRegionZ0

double m_goodMagneticRegionZ0 = 57.
private

maximum z0 for well understood magnetic field (cm)

Definition at line 80 of file FilterCalculator.h.

◆ m_hltPrefilterParameters

DBObjPtr<HLTPrefilterParameters> m_hltPrefilterParameters
private

Objects relevant to HLTPrefilter monitoring HLTprefilterParameters Database OjbPtr.

HLT prefilter parameters

Definition at line 86 of file FilterCalculator.h.

◆ m_l1Trigger

StoreObjPtr<TRGSummary> m_l1Trigger
private

Store Object with the trigger result.

Definition at line 51 of file FilterCalculator.h.

◆ m_looseTrkZ0

double m_looseTrkZ0 = 10 * Unit::cm
private

which Z0 defines a loose track

Definition at line 56 of file FilterCalculator.h.

◆ m_reducedEsinglePhoton

double m_reducedEsinglePhoton = 0.5
private

which CMS energy defines nReducedEsingle clusters

Definition at line 74 of file FilterCalculator.h.

◆ m_tightTrkZ0

double m_tightTrkZ0 = 2 * Unit::cm
private

which Z0 defines a tight track

Definition at line 58 of file FilterCalculator.h.

◆ m_timingPrefilter

HLTPrefilter::TimingCutState m_timingPrefilter
private

Helper instance for timing based prefilter.

Definition at line 88 of file FilterCalculator.h.

◆ m_tracks

StoreArray<Track> m_tracks
private

Store Array of the tracks to be used.

Definition at line 47 of file FilterCalculator.h.


The documentation for this class was generated from the following files: