Belle II Software release-09-00-14
SkimSampleCalculator Class Reference

Implementation of a calculator used in the SoftwareTriggerModule to fill a SoftwareTriggerObject for selecting particles for skimming and data quality monitoring. More...

#include <SkimSampleCalculator.h>

Inheritance diagram for SkimSampleCalculator:
Collaboration diagram for SkimSampleCalculator:

Public Member Functions

 SkimSampleCalculator ()
 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

StoreObjPtr< ParticleListm_pionParticles
 Internal storage of the tracks as particles.
 
StoreObjPtr< ParticleListm_gammaParticles
 Internal storage of the ECL clusters as particles.
 
StoreObjPtr< ParticleListm_pionHadParticles
 Internal storage of the tracks as particles (definition for hadronb).
 
StoreObjPtr< ParticleListm_pionTauParticles
 Internal storage of the tracks as particles (definition for tau skims).
 
StoreObjPtr< ParticleListm_KsParticles
 Internal storage of the K_S0's.
 
StoreObjPtr< ParticleListm_LambdaParticles
 Internal storage of the Lambda0's.
 
StoreObjPtr< ParticleListm_DstParticles
 Internal storage of the D*'s.
 
StoreObjPtr< ParticleListm_offIpParticles
 Internal storage of the tracks for alignment calibration.
 
std::string m_filterL1TrgNN = ""
 HLT filter line for the TRG skim.
 
StoreObjPtr< ParticleListm_BpParticles
 Internal storage of the B+'s.
 
StoreObjPtr< ParticleListm_BzParticles
 Internal storage of the B0's.
 
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 selecting particles for skimming and data quality monitoring.

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

Definition at line 28 of file SkimSampleCalculator.h.

Constructor & Destructor Documentation

◆ SkimSampleCalculator()

Set the default names for the store object particle lists.

Definition at line 42 of file SkimSampleCalculator.cc.

42 :
43 m_pionParticles("pi+:skim"), m_gammaParticles("gamma:skim"), m_pionHadParticles("pi+:hadb"), m_pionTauParticles("pi+:tau"),
44 m_KsParticles("K_S0:merged"), m_LambdaParticles("Lambda0:merged"), m_DstParticles("D*+:d0pi"), m_offIpParticles("pi+:offip"),
45 m_filterL1TrgNN("software_trigger_cut&filter&L1_trigger_nn_info"),
46 m_BpParticles("B+:BtoCharmForHLT"), m_BzParticles("B0:BtoCharmForHLT")
47{
48
49}
StoreObjPtr< ParticleList > m_pionParticles
Internal storage of the tracks as particles.
StoreObjPtr< ParticleList > m_gammaParticles
Internal storage of the ECL clusters as particles.
StoreObjPtr< ParticleList > m_BzParticles
Internal storage of the B0's.
StoreObjPtr< ParticleList > m_LambdaParticles
Internal storage of the Lambda0's.
StoreObjPtr< ParticleList > m_pionHadParticles
Internal storage of the tracks as particles (definition for hadronb).
StoreObjPtr< ParticleList > m_KsParticles
Internal storage of the K_S0's.
std::string m_filterL1TrgNN
HLT filter line for the TRG skim.
StoreObjPtr< ParticleList > m_DstParticles
Internal storage of the D*'s.
StoreObjPtr< ParticleList > m_pionTauParticles
Internal storage of the tracks as particles (definition for tau skims).
StoreObjPtr< ParticleList > m_BpParticles
Internal storage of the B+'s.
StoreObjPtr< ParticleList > m_offIpParticles
Internal storage of the tracks for alignment calibration.

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 }
SoftwareTriggerObject m_calculationResult
Internal storage of the result of the calculation.

◆ doCalculation()

void doCalculation ( SoftwareTriggerObject &  calculationResult)
overridevirtual

Actually write out the variables into the map.

Implements SoftwareTriggerCalculation.

Definition at line 66 of file SkimSampleCalculator.cc.

67{
68 // Prefetch some later needed objects/values
69 const Particle* gammaWithMaximumRho = getElementWithMaximumRho<Particle>(m_gammaParticles);
70 const Particle* gammaWithSecondMaximumRho = getElementWithMaximumRhoBelow<Particle>(m_gammaParticles,
71 getRho(gammaWithMaximumRho));
72 const Particle* trackWithMaximumRho = getElementWithMaximumRho<Particle>(m_pionParticles);
73 const Particle* trackWithSecondMaximumRho = getElementWithMaximumRhoBelow<Particle>(m_pionParticles,
74 getRho(trackWithMaximumRho));
75
76 const double& rhoOfECLClusterWithMaximumRho = getRhoOfECLClusterWithMaximumRho(m_pionParticles, m_gammaParticles);
77 const double& rhoOfECLClusterWithSecondMaximumRho = getRhoOfECLClusterWithMaximumRhoBelow(m_pionParticles,
79 rhoOfECLClusterWithMaximumRho);
80
81 const double& rhoOfTrackWithMaximumRho = getRho(trackWithMaximumRho);
82 const double& rhoOfTrackWithSecondMaximumRho = getRho(trackWithSecondMaximumRho);
83 const double& rhoOfGammaWithMaximumRho = getRho(gammaWithMaximumRho);
84 const double& rhoOfGammaWithSecondMaximumRho = getRho(gammaWithSecondMaximumRho);
85
86 // Simple to calculate variables
87 // EC1CMSLE
88 calculationResult["EC1CMSLE"] = rhoOfECLClusterWithMaximumRho;
89
90 // EC2CMSLE
91 calculationResult["EC2CMSLE"] = rhoOfECLClusterWithSecondMaximumRho;
92
93 // EC12CMSLE
94 calculationResult["EC12CMSLE"] = rhoOfECLClusterWithMaximumRho + rhoOfECLClusterWithSecondMaximumRho;
95
96 // nTracksLE
97 calculationResult["nTracksLE"] = m_pionParticles->getListSize();
98
99 // nTracksTAU
100 calculationResult["nTracksTAU"] = m_pionTauParticles->getListSize();
101
102 // nGammasLE
103 calculationResult["nGammasLE"] = m_gammaParticles->getListSize();
104
105 // P1CMSBhabhaLE
106 calculationResult["P1CMSBhabhaLE"] = rhoOfTrackWithMaximumRho;
107
108 // P1CMSBhabhaLE/E_beam
109 calculationResult["P1OEbeamCMSBhabhaLE"] = rhoOfTrackWithMaximumRho / BeamEnergyCMS();
110
111 // P2CMSBhabhaLE
112 calculationResult["P2CMSBhabhaLE"] = rhoOfTrackWithSecondMaximumRho;
113
114 // P2CMSBhabhaLE/E_beam
115 calculationResult["P2OEbeamCMSBhabhaLE"] = rhoOfTrackWithSecondMaximumRho / BeamEnergyCMS();
116
117 // P12CMSBhabhaLE
118 calculationResult["P12CMSBhabhaLE"] = rhoOfTrackWithMaximumRho + rhoOfTrackWithSecondMaximumRho;
119
120 //G1CMSLE, the largest energy of gamma in CMS
121 calculationResult["G1CMSBhabhaLE"] = rhoOfGammaWithMaximumRho;
122 //G1OEbeamCMSLE, the largest energy of gamma in CMS over beam energy
123 calculationResult["G1OEbeamCMSBhabhaLE"] = rhoOfGammaWithMaximumRho / BeamEnergyCMS();
124
125 //G2CMSLE, the secondary largest energy of gamma in CMS
126 calculationResult["G2CMSBhabhaLE"] = rhoOfGammaWithSecondMaximumRho;
127 //G2OEbeamCMSLE, the largest energy of gamma in CMS over beam energy
128 calculationResult["G2OEbeamCMSBhabhaLE"] = rhoOfGammaWithSecondMaximumRho / BeamEnergyCMS();
129
130 //G12CMSLE, the secondary largest energy of gamma in CMS
131 calculationResult["G12CMSBhabhaLE"] = rhoOfGammaWithMaximumRho + rhoOfGammaWithSecondMaximumRho;
132 //G12CMSLE, the secondary largest energy of gamma in CMS over beam energy
133 calculationResult["G12OEbeamCMSBhabhaLE"] =
134 (rhoOfGammaWithMaximumRho + rhoOfGammaWithSecondMaximumRho) / BeamEnergyCMS();
135
136
137 // Medium hard to calculate variables
138 // ENeutralLE
139 if (gammaWithMaximumRho) {
140 calculationResult["ENeutralLE"] = getRho(gammaWithMaximumRho);
141 } else {
142 calculationResult["ENeutralLE"] = -1;
143 }
144
145 // nECLMatchTracksLE
146 const unsigned int numberOfTracksWithECLMatch = std::count_if(m_pionParticles->begin(), m_pionParticles->end(),
147 [](const Particle & particle) {
148 return particle.getECLCluster() != nullptr;
149 });
150 calculationResult["nECLMatchTracksLE"] = numberOfTracksWithECLMatch;
151
152 //nECLClustersLE
153 double neclClusters = -1.;
154 double eneclClusters = 0.;
155 StoreArray<ECLCluster> eclClusters;
156 ClusterUtils Cl;
157 double PzGamma = 0.;
158 double EsumGamma = 0.;
159 if (eclClusters.isValid()) {
160 const unsigned int numberOfECLClusters = std::count_if(eclClusters.begin(), eclClusters.end(),
161 [](const ECLCluster & eclcluster) {
162 return (eclcluster.hasHypothesis(
163 ECLCluster::EHypothesisBit::c_nPhotons)
164 and eclcluster.getEnergy(
165 ECLCluster::EHypothesisBit::c_nPhotons) > 0.1);
166 });
167 neclClusters = numberOfECLClusters;
168
169 for (int ncl = 0; ncl < eclClusters.getEntries(); ncl++) {
170 if (eclClusters[ncl]->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)
171 && eclClusters[ncl]->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons) > 0.1) {
172 eneclClusters += eclClusters[ncl]->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons);
173 if (!eclClusters[ncl]->getRelatedFrom<Track>()) {
174 ROOT::Math::PxPyPzEVector V4Gamma_CMS = PCmsLabTransform::labToCms(Cl.Get4MomentumFromCluster(eclClusters[ncl],
176 EsumGamma += V4Gamma_CMS.E();
177 PzGamma += V4Gamma_CMS.Pz();
178 }
179 }
180 }
181 }
182 calculationResult["nECLClustersLE"] = neclClusters;
183
184 int nb2bcc_PhiHigh = 0;
185 int nb2bcc_PhiLow = 0;
186 int nb2bcc_3D = 0;
187 ClusterUtils C;
188 for (int i = 0; i < eclClusters.getEntries() - 1; i++) {
189 if (!eclClusters[i]->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons))
190 continue;
191 ROOT::Math::PxPyPzEVector V4g1 = C.Get4MomentumFromCluster(eclClusters[i], ECLCluster::EHypothesisBit::c_nPhotons);
192 double Eg1 = V4g1.E();
193 for (int j = i + 1; j < eclClusters.getEntries(); j++) {
194 if (!eclClusters[j]->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons))
195 continue;
196 ROOT::Math::PxPyPzEVector V4g2 = C.Get4MomentumFromCluster(eclClusters[j], ECLCluster::EHypothesisBit::c_nPhotons);
197 double Eg2 = V4g2.E();
198 const B2Vector3D V3g1 = (PCmsLabTransform::labToCms(V4g1)).Vect();
199 const B2Vector3D V3g2 = (PCmsLabTransform::labToCms(V4g2)).Vect();
200 double Thetag1 = (PCmsLabTransform::labToCms(V4g1)).Theta() * TMath::RadToDeg();
201 double Thetag2 = (PCmsLabTransform::labToCms(V4g2)).Theta() * TMath::RadToDeg();
202 double deltphi = fabs(V3g1.DeltaPhi(V3g2) * TMath::RadToDeg());
203 double Tsum = Thetag1 + Thetag2;
204 if (deltphi > 170. && (Eg1 > 0.25 && Eg2 > 0.25)) nb2bcc_PhiHigh++;
205 if (deltphi > 170. && (Eg1 < 0.25 || Eg2 < 0.25)) nb2bcc_PhiLow++;
206 if (deltphi > 160. && (Tsum > 160. && Tsum < 200.)) nb2bcc_3D++;
207 }
208 }
209
210 calculationResult["nB2BCCPhiHighLE"] = nb2bcc_PhiHigh;
211 calculationResult["nB2BCCPhiLowLE"] = nb2bcc_PhiLow;
212 calculationResult["nB2BCC3DLE"] = nb2bcc_3D;
213
214
215 // AngleGTLE
216 double angleGTLE = -10.;
217 if (gammaWithMaximumRho) {
218 const B2Vector3D& V3g1 = gammaWithMaximumRho->getMomentum();
219 if (trackWithMaximumRho) {
220 const B2Vector3D& V3p1 = trackWithMaximumRho->getMomentum();
221 const double theta1 = V3g1.Angle(V3p1);
222 if (angleGTLE < theta1) angleGTLE = theta1;
223 }
224 if (trackWithSecondMaximumRho) {
225 const B2Vector3D& V3p2 = trackWithSecondMaximumRho->getMomentum();
226 const double theta2 = V3g1.Angle(V3p2);
227 if (angleGTLE < theta2) angleGTLE = theta2;
228 }
229 }
230
231 calculationResult["AngleGTLE"] = angleGTLE;
232
233 // AngleG1G2LE
234 double angleG1G2CMSLE = -10.;
235 if (gammaWithMaximumRho) {
236 const ROOT::Math::PxPyPzEVector& V4p1 = gammaWithMaximumRho->get4Vector();
237 if (gammaWithSecondMaximumRho) {
238 const ROOT::Math::PxPyPzEVector& V4p2 = gammaWithSecondMaximumRho->get4Vector();
239 const B2Vector3D V3p1 = (PCmsLabTransform::labToCms(V4p1)).Vect();
240 const B2Vector3D V3p2 = (PCmsLabTransform::labToCms(V4p2)).Vect();
241 angleG1G2CMSLE = V3p1.Angle(V3p2);
242 }
243 }
244
245 calculationResult["AngleG1G2CMSLE"] = angleG1G2CMSLE;
246
247 // maxAngleTTLE
248 double maxAngleTTLE = -10.;
249 int nJpsi = 0;
250 double Jpsi = 0.;
251 const double jPsiMasswindow = 0.11;
252 if (m_pionParticles->getListSize() >= 2) {
253 for (unsigned int i = 0; i < m_pionParticles->getListSize() - 1; i++) {
254 Particle* par1 = m_pionParticles->getParticle(i);
255 for (unsigned int j = i + 1; j < m_pionParticles->getListSize(); j++) {
256 Particle* par2 = m_pionParticles->getParticle(j);
257 ROOT::Math::PxPyPzEVector V4p1 = par1->get4Vector();
258 ROOT::Math::PxPyPzEVector V4p2 = par2->get4Vector();
259 ROOT::Math::PxPyPzEVector V4pSum = V4p1 + V4p2;
260 const auto chSum = par1->getCharge() + par2->getCharge();
261 const double mSum = V4pSum.M();
262 const double JpsidM = mSum - TDatabasePDG::Instance()->GetParticle(443)->Mass();
263 if (abs(JpsidM) < jPsiMasswindow && chSum == 0) nJpsi++;
264 const B2Vector3D V3p1 = (PCmsLabTransform::labToCms(V4p1)).Vect();
265 const B2Vector3D V3p2 = (PCmsLabTransform::labToCms(V4p2)).Vect();
266 const double temp = V3p1.Angle(V3p2);
267 if (maxAngleTTLE < temp) maxAngleTTLE = temp;
268 }
269 }
270 }
271
272 if (nJpsi != 0) Jpsi = 1;
273
274 calculationResult["maxAngleTTLE"] = maxAngleTTLE;
275 calculationResult["Jpsi"] = Jpsi;
276
277 //maxAngleGGLE
278 double maxAngleGGLE = -10.;
279 if (m_gammaParticles->getListSize() >= 2) {
280 for (unsigned int i = 0; i < m_gammaParticles->getListSize() - 1; i++) {
281 Particle* par1 = m_gammaParticles->getParticle(i);
282 for (unsigned int j = i + 1; j < m_gammaParticles->getListSize(); j++) {
283 Particle* par2 = m_gammaParticles->getParticle(j);
284 ROOT::Math::PxPyPzEVector V4p1 = par1->get4Vector();
285 ROOT::Math::PxPyPzEVector V4p2 = par2->get4Vector();
286 const B2Vector3D V3p1 = (PCmsLabTransform::labToCms(V4p1)).Vect();
287 const B2Vector3D V3p2 = (PCmsLabTransform::labToCms(V4p2)).Vect();
288 const double temp = V3p1.Angle(V3p2);
289 if (maxAngleGGLE < temp) maxAngleGGLE = temp;
290 }
291 }
292 }
293
294 calculationResult["maxAngleGGLE"] = maxAngleGGLE;
295
296 // nEidLE
297 const unsigned int nEidLE = std::count_if(m_pionParticles->begin(), m_pionParticles->end(),
298 [](const Particle & p) {
299 const double& momentum = p.getMomentumMagnitude();
300 const double& r_rho = getRho(&p);
301 const ECLCluster* eclTrack = p.getECLCluster();
302 if (eclTrack) {
303 const double& energyOverMomentum = eclTrack->getEnergy(
304 ECLCluster::EHypothesisBit::c_nPhotons) / momentum;
305 double r_rhotoebeam = r_rho / BeamEnergyCMS();
306 return (r_rhotoebeam) > 0.35 && energyOverMomentum > 0.8;
307 }
308 return false;
309 });
310
311 calculationResult["nEidLE"] = nEidLE;
312
313
314 // VisibleEnergyLE
315 const double visibleEnergyTracks = std::accumulate(m_pionParticles->begin(), m_pionParticles->end(), 0.0,
316 [](const double & visibleEnergy, const Particle & p) {
317 return visibleEnergy + p.getMomentumMagnitude();
318 });
319
320 const double visibleEnergyGammas = std::accumulate(m_gammaParticles->begin(), m_gammaParticles->end(), 0.0,
321 [](const double & visibleEnergy, const Particle & p) {
322 return visibleEnergy + p.getMomentumMagnitude();
323 });
324
325 calculationResult["VisibleEnergyLE"] = visibleEnergyTracks + visibleEnergyGammas;
326
327 // EtotLE
328 const double eTotTracks = std::accumulate(m_pionParticles->begin(), m_pionParticles->end(), 0.0,
329 [](const double & eTot, const Particle & p) {
330 const ECLCluster* eclCluster = p.getECLCluster();
331 if (eclCluster) {
332 const double eclEnergy = eclCluster->getEnergy(
333 ECLCluster::EHypothesisBit::c_nPhotons);
334 if (eclEnergy > 0.1) {
335 return eTot + eclCluster->getEnergy(
336 ECLCluster::EHypothesisBit::c_nPhotons);
337 }
338 }
339 return eTot;
340 });
341
342 const double eTotGammas = std::accumulate(m_gammaParticles->begin(), m_gammaParticles->end(), 0.0,
343 [](const double & eTot, const Particle & p) {
344 return eTot + p.getEnergy();
345 });
346 double Etot = eTotTracks + eTotGammas;
347 calculationResult["EtotLE"] = Etot;
348
349 //KLM inforamtion
350 // The clusters with the largest pentrate layers in KLM.
351 double numMaxLayerKLM = -1.;
352 double numSecMaxLayerKLM = -1.;
353 StoreArray<KLMCluster> klmClusters;
354 if (klmClusters.isValid()) {
355 for (const auto& klmCluster : klmClusters) {
356 double klmClusterLayer = klmCluster.getLayers();
357 if (numMaxLayerKLM < klmClusterLayer) {
358 numSecMaxLayerKLM = numMaxLayerKLM;
359 numMaxLayerKLM = klmClusterLayer;
360 } else if (numSecMaxLayerKLM < klmClusterLayer)
361 numSecMaxLayerKLM = klmClusterLayer;
362 }
363 }
364 calculationResult["N1KLMLayer"] = numMaxLayerKLM;
365 calculationResult["N2KLMLayer"] = numSecMaxLayerKLM;
366
367 //define bhabha_2trk, bhabha_1trk, eclbhabha
368 int charget1 = -10;
369 if (trackWithMaximumRho) charget1 = trackWithMaximumRho->getCharge();
370 int charget2 = -10;
371 if (trackWithSecondMaximumRho) charget2 = trackWithSecondMaximumRho->getCharge();
372
373 double Bhabha2Trk = 0.;
374 int ntrk_bha = m_pionParticles->getListSize();
375 double rp1ob = rhoOfTrackWithMaximumRho / BeamEnergyCMS();
376 double rp2ob = rhoOfTrackWithSecondMaximumRho / BeamEnergyCMS();
377 bool bhabha2trk_tag =
378 ntrk_bha >= 2 && maxAngleTTLE > 2.88 && nEidLE >= 1 && rp1ob > 0.35 && rp2ob > 0.35 && (Etot) > 4.0
379 && (abs(charget1) == 1 && abs(charget2) == 1 && (charget1 + charget2) == 0);
380 if (bhabha2trk_tag) Bhabha2Trk = 1;
381 calculationResult["Bhabha2Trk"] = Bhabha2Trk;
382
383 double Bhabha1Trk = 0.;
384 double rc1ob = rhoOfGammaWithMaximumRho / BeamEnergyCMS();
385 double rc2ob = rhoOfGammaWithSecondMaximumRho / BeamEnergyCMS();
386 bool bhabha1trk_tag = ntrk_bha == 1 && rp1ob > 0.35 && rc1ob > 0.35 && angleGTLE > 2.618;
387 if (bhabha1trk_tag) Bhabha1Trk = 1;
388 calculationResult["Bhabha1Trk"] = Bhabha1Trk;
389
390 double ggSel = 0.;
391 bool gg_tag = ntrk_bha <= 1 && nEidLE == 0 && rc1ob > 0.35 && rc2ob > 0.2 && Etot > 4.0 && maxAngleGGLE > 2.618;
392 if (gg_tag) ggSel = 1;
393 calculationResult["GG"] = ggSel;
394
395 // Bhabha skim with ECL information only (bhabhaecl)
396 double BhabhaECL = 0.;
397 ClusterUtils Cls;
398 for (int i = 0; i < eclClusters.getEntries() - 1; i++) {
399 if (!eclClusters[i]->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons))
400 continue;
401
402 ROOT::Math::PxPyPzEVector V4g1 = PCmsLabTransform::labToCms(Cls.Get4MomentumFromCluster(eclClusters[i],
404 double Eg1ob = V4g1.E() / (2 * BeamEnergyCMS());
405 for (int j = i + 1; j < eclClusters.getEntries(); j++) {
406 if (!eclClusters[j]->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons))
407 continue;
408 ROOT::Math::PxPyPzEVector V4g2 = PCmsLabTransform::labToCms(Cls.Get4MomentumFromCluster(eclClusters[j],
410 double Eg2ob = V4g2.E() / (2 * BeamEnergyCMS());
411 const B2Vector3D V3g1 = V4g1.Vect();
412 const B2Vector3D V3g2 = V4g2.Vect();
413 double Thetag1 = V4g1.Theta() * TMath::RadToDeg();
414 double Thetag2 = V4g2.Theta() * TMath::RadToDeg();
415 double deltphi = fabs(V3g1.DeltaPhi(V3g2) * TMath::RadToDeg());
416 double Tsum = Thetag1 + Thetag2;
417 if ((deltphi > 165. && deltphi < 178.5) && (Eg1ob > 0.4 && Eg2ob > 0.4 && (Eg1ob > 0.45 || Eg2ob > 0.45)) && (Tsum > 178.
418 && Tsum < 182.)) BhabhaECL = 1;
419 }
420 }
421 calculationResult["BhabhaECL"] = BhabhaECL;
422
423 // Bhabha skim (BhabhaCDC) for CDC dE/dx calib studies
424 double BhabhaCDC = 0.;
425 // Radiative Bhabha skim (radee)
426 double radee = 0.;
427 const double lowdEdxEdge = 0.70, highdEdxEdge = 1.30;
428 const double lowEoPEdge = 0.70, highEoPEdge = 1.30;
429
430 if (m_pionParticles->getListSize() == 2) {
431
432 //------------First track variables----------------
433 for (unsigned int i = 0; i < m_pionParticles->getListSize() - 1; i++) {
434
435 Particle* part1 = m_pionParticles->getParticle(i);
436 if (!part1) continue;
437
438 const auto chargep1 = part1->getCharge();
439 if (abs(chargep1) != 1) continue;
440
441 const ECLCluster* eclTrack1 = part1->getECLCluster();
442 if (!eclTrack1) continue;
443 if (!eclTrack1->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
444
445 const double& momentum1 = part1->getMomentumMagnitude();
446 const double& energyOverMomentum1 = eclTrack1->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons) / momentum1;
447 if (energyOverMomentum1 <= lowEoPEdge || energyOverMomentum1 >= highEoPEdge) continue;
448
449 const Track* track1 = part1->getTrack();
450 if (!track1) continue;
451
452 const TrackFitResult* trackFit1 = track1->getTrackFitResultWithClosestMass(Const::pion);
453 if (!trackFit1) continue;
454 if (trackFit1->getHitPatternCDC().getNHits() <= 0) continue;
455
456 //------------Second track variables----------------
457 for (unsigned int j = i + 1; j < m_pionParticles->getListSize(); j++) {
458
459 Particle* part2 = m_pionParticles->getParticle(j);
460 if (!part2) continue;
461
462 const auto chargep2 = part2->getCharge();
463 if (abs(chargep2) != 1 || (chargep1 + chargep2 != 0)) continue;
464
465 const ECLCluster* eclTrack2 = part2->getECLCluster();
466 if (!eclTrack2) continue;
467 if (!eclTrack2->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
468
469 const double& momentum2 = part2->getMomentumMagnitude();
470 const double& energyOverMomentum2 = eclTrack2->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons) / momentum2;
471 if (energyOverMomentum2 <= lowEoPEdge || energyOverMomentum2 >= highEoPEdge) continue;
472
473 const Track* track2 = part2->getTrack();
474 if (!track2) continue;
475
476 const TrackFitResult* trackFit2 = track2->getTrackFitResultWithClosestMass(Const::pion);
477 if (!trackFit2) continue;
478 if (trackFit2->getHitPatternCDC().getNHits() <= 0) continue;
479
480 BhabhaCDC = 1;
481
482 const CDCDedxTrack* dedxTrack1 = track1->getRelatedTo<CDCDedxTrack>();
483 if (!dedxTrack1) continue;
484
485 const CDCDedxTrack* dedxTrack2 = track2->getRelatedTo<CDCDedxTrack>();
486 if (!dedxTrack2) continue;
487
488 double p1_dedxnosat = dedxTrack1->getDedxNoSat();
489 double p2_dedxnosat = dedxTrack2->getDedxNoSat();
490
491 if ((p1_dedxnosat > lowdEdxEdge && p1_dedxnosat < highdEdxEdge) || (p2_dedxnosat > lowdEdxEdge
492 && p2_dedxnosat < highdEdxEdge)) radee = 1;
493
494 }
495 }
496 }
497
498 calculationResult["BhabhaCDC"] = BhabhaCDC;
499 calculationResult["Radee"] = radee;
500
501 // Dimuon skim (mumutight) taken from the offline skim + Radiative dimuon (radmumu)
502 double mumutight = 0.;
503 double eMumuTotGammas = 0.;
504 int nTracks = 0;
505 double radmumu = 0.;
506 int nGammas = m_gammaParticles->getListSize();
507
508 for (int t = 0; t < nGammas; t++) {
509 const Particle* part = m_gammaParticles->getParticle(t);
510 const auto& frame = ReferenceFrame::GetCurrent();
511 eMumuTotGammas += frame.getMomentum(part).E();
512 }
513
514 StoreArray<Track> tracks;
515 nTracks = tracks.getEntries();
517 const ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
518 const auto& fr = ReferenceFrame::GetCurrent();
519
520 if (m_pionParticles->getListSize() == 2) {
521
522 //------------First track variables----------------
523 for (unsigned int k = 0; k < m_pionParticles->getListSize() - 1; k++) {
524
525 Particle* part1 = m_pionParticles->getParticle(k);
526 if (!part1) continue;
527
528 const auto chargep1 = part1->getCharge();
529 if (abs(chargep1) != 1) continue;
530
531 const ECLCluster* eclTrack1 = part1->getECLCluster();
532 if (!eclTrack1) continue;
533 if (!eclTrack1->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
534
535 const Track* track1 = part1->getTrack();
536 if (!track1) continue;
537
538 const TrackFitResult* trackFit1 = track1->getTrackFitResultWithClosestMass(Const::pion);
539 if (!trackFit1) continue;
540
541 const ROOT::Math::PxPyPzEVector V4p1 = trackFit1->get4Momentum();
542 const B2Vector3D V3p1 = (PCmsLabTransform::labToCms(V4p1)).Vect();
543
544 const double p1MomLab = V4p1.P();
545 double highestP = p1MomLab;
546 const double p1CDChits = trackFit1->getHitPatternCDC().getNHits();
547 const PIDLikelihood* p1Pid = part1->getPIDLikelihood();
548 bool p1hasKLMid = 0;
549 if (p1Pid) p1hasKLMid = p1Pid->isAvailable(Const::KLM);
550 const double p1isInCDC = Variable::inCDCAcceptance(part1);
551 const double p1clusPhi = Variable::eclClusterPhi(part1);
552
553 const double Pp1 = V3p1.Mag();
554 const double Thetap1 = (V3p1).Theta() * TMath::RadToDeg();
555 const double Phip1 = (V3p1).Phi() * TMath::RadToDeg();
556
557 const double enECLTrack1 = eclTrack1->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons);
558
559 const bool goodTrk1 = enECLTrack1 > 0 && enECLTrack1 < 0.5 && p1CDChits > 0
560 && ((p1hasKLMid == 0 && enECLTrack1 < 0.25 && p1MomLab < 2.0) || p1hasKLMid == 1) && p1isInCDC == 1;
561
562 //------------Second track variables----------------
563 for (unsigned int l = k + 1; l < m_pionParticles->getListSize(); l++) {
564
565 Particle* part2 = m_pionParticles->getParticle(l);
566 if (!part2) continue;
567
568 const auto chargep2 = part2->getCharge();
569 if (abs(chargep2) != 1 || (chargep1 + chargep2 != 0)) continue;
570
571 const ECLCluster* eclTrack2 = part2->getECLCluster();
572 if (!eclTrack2) continue;
573 if (!eclTrack2->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
574
575 const Track* track2 = part2->getTrack();
576 if (!track2) continue;
577
578 const TrackFitResult* trackFit2 = track2->getTrackFitResultWithClosestMass(Const::pion);
579 if (!trackFit2) continue;
580
581 const ROOT::Math::PxPyPzEVector V4p2 = trackFit2->get4Momentum();
582 const B2Vector3D V3p2 = (PCmsLabTransform::labToCms(V4p2)).Vect();
583
584 const double p2MomLab = V4p2.P();
585 double lowestP = p2MomLab;
586 const double p2CDChits = trackFit2->getHitPatternCDC().getNHits();
587 const PIDLikelihood* p2Pid = part2->getPIDLikelihood();
588 bool p2hasKLMid = 0;
589 if (p2Pid) p2hasKLMid = p2Pid->isAvailable(Const::KLM);
590 const double p2isInCDC = Variable::inCDCAcceptance(part2);
591 const double p2clusPhi = Variable::eclClusterPhi(part2);
592
593 const double Pp2 = V3p2.Mag();
594 const double Thetap2 = (V3p2).Theta() * TMath::RadToDeg();
595 const double Phip2 = (V3p2).Phi() * TMath::RadToDeg();
596
597 const double acopPhi = fabs(180 - fabs(Phip1 - Phip2));
598 const double acopTheta = fabs(fabs(Thetap1 + Thetap2) - 180);
599
600 const double enECLTrack2 = eclTrack2->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons);
601
602 const bool goodTrk2 = enECLTrack2 > 0 && enECLTrack2 < 0.5 && p2CDChits > 0
603 && ((p2hasKLMid == 0 && enECLTrack2 < 0.25 && p2MomLab < 2.0) || p2hasKLMid == 1) && p2isInCDC == 1;
604
605 double eTotMumuTracks = enECLTrack1 + enECLTrack2;
606 double EMumutot = eTotMumuTracks + eMumuTotGammas;
607
608 bool mumutight_tag = enECLTrack1 < 0.5 && enECLTrack2 < 0.5 && EMumutot < 2 && acopPhi < 10 && acopTheta < 10 && nTracks == 2
609 && Pp1 > 0.5 && Pp2 > 0.5;
610
611 if (mumutight_tag) mumutight = 1;
612
613 if (p1MomLab < p2MomLab) {
614 lowestP = highestP;
615 highestP = p2MomLab;
616 }
617
618 double diffPhi = p1clusPhi - p2clusPhi;
619 if (fabs(diffPhi) > M_PI) {
620 if (diffPhi > M_PI) {
621 diffPhi = diffPhi - 2 * M_PI;
622 } else {
623 diffPhi = 2 * M_PI + diffPhi;
624 }
625 }
626
627 const double recoilP = fr.getMomentum(pIN - V4p1 - V4p2).P();
628
629 const bool radmumu_tag = nTracks < 4 && goodTrk1 == 1 && goodTrk2 == 1 && highestP > 1 && lowestP < 3 && (p1hasKLMid == 1
630 || p2hasKLMid == 1) && abs(diffPhi) >= 0.5 * M_PI && recoilP > 0.1 && (enECLTrack1 <= 0.25 || enECLTrack2 <= 0.25);
631
632 if (radmumu_tag) radmumu = 1;
633
634 }
635 }
636 }
637
638 calculationResult["MumuTight"] = mumutight;
639 calculationResult["Radmumu"] = radmumu;
640
641 //Retrieve variables for HadronB skims
642 double EsumPiHad = 0;
643 double PzPiHad = 0;
644 int nHadTracks = m_pionHadParticles->getListSize();
645 double hadronb = 0;
646 double hadronb1 = 0;
647 double hadronb2 = 0;
648 std::vector<ROOT::Math::XYZVector> m_pionHadv3;
649
650 for (int nPiHad = 0; nPiHad < nHadTracks; nPiHad++) {
651 Particle* parPiHad = m_pionHadParticles->getParticle(nPiHad);
652 ROOT::Math::PxPyPzEVector V4PiHad = PCmsLabTransform::labToCms(parPiHad->get4Vector());
653 m_pionHadv3.push_back(parPiHad->getMomentum());
654 EsumPiHad += V4PiHad.E();
655 PzPiHad += V4PiHad.Pz();
656 }
657
658 double visibleEnergyCMSnorm = (EsumPiHad + EsumGamma) / (BeamEnergyCMS() * 2.0);
659 double EsumCMSnorm = eneclClusters / (BeamEnergyCMS() * 2.0);
660 double PzTotCMSnorm = (PzPiHad + PzGamma) / (BeamEnergyCMS() * 2.0);
661
662 bool hadronb_tag = nHadTracks >= 3 && visibleEnergyCMSnorm > 0.2 && abs(PzTotCMSnorm) < 0.5 && neclClusters > 1
663 && EsumCMSnorm > 0.1 && EsumCMSnorm < 0.8;
664
665 if (hadronb_tag) {
666 hadronb = 1;
667 FoxWolfram fw(m_pionHadv3);
668 fw.calculateBasicMoments();
669 double R2 = fw.getR(2);
670 if (R2 < 0.4) hadronb1 = 1;
671 if (hadronb1 && nHadTracks >= 5) hadronb2 = 1;
672 }
673
674 calculationResult["HadronB"] = hadronb;
675 calculationResult["HadronB1"] = hadronb1;
676 calculationResult["HadronB2"] = hadronb2;
677
678 // nKshort
679 int nKshort = 0;
680 double Kshort = 0.;
681 const double KsMassLow = 0.468;
682 const double KsMassHigh = 0.528;
683
684 if (m_KsParticles.isValid()) {
685 for (unsigned int i = 0; i < m_KsParticles->getListSize(); i++) {
686 const Particle* mergeKsCand = m_KsParticles->getParticle(i);
687 const double isKsCandGood = Variable::goodBelleKshort(mergeKsCand);
688 const double KsCandMass = mergeKsCand->getMass();
689 if (KsCandMass > KsMassLow && KsCandMass < KsMassHigh && isKsCandGood == 1.) nKshort++;
690 }
691 }
692
693 if (nKshort != 0) Kshort = 1;
694
695 calculationResult["Kshort"] = Kshort;
696
697 // 4 leptons skim
698 int nFourLep = 0;
699 double fourLep = 0.;
700
701 const double visibleEnergyCMS = visibleEnergyCMSnorm * BeamEnergyCMS() * 2.0;
702 const unsigned int n_particles = m_pionHadParticles->getListSize();
703
704 if (n_particles >= 2) {
705 for (unsigned int i = 0; i < n_particles - 1; i++) {
706 Particle* par1 = m_pionHadParticles->getParticle(i);
707 for (unsigned int j = i + 1; j < n_particles; j++) {
708 Particle* par2 = m_pionHadParticles->getParticle(j);
709 const auto chSum = par1->getCharge() + par2->getCharge();
710 const ROOT::Math::PxPyPzEVector V4p1 = par1->get4Vector();
711 const ROOT::Math::PxPyPzEVector V4p2 = par2->get4Vector();
712 const double opAng = V4p1.Theta() + V4p2.Theta();
713 const ROOT::Math::PxPyPzEVector V4pSum = V4p1 + V4p2;
714 const ROOT::Math::PxPyPzEVector V4pSumCMS = PCmsLabTransform::labToCms(V4pSum);
715 const double ptCMS = V4pSumCMS.Pt();
716 const double pzCMS = V4pSumCMS.Pz();
717 const double mSum = V4pSum.M();
718
719 const bool fourLepCand = chSum == 0 && (V4p1.P() > 0.4 && V4p2.P() > 0.4) && cos(opAng) > -0.997 && ptCMS < 0.15 && abs(pzCMS) < 2.5
720 && mSum < 6;
721
722 if (fourLepCand) nFourLep++;
723 }
724 }
725 }
726
727 if (nFourLep != 0 && visibleEnergyCMS < 6) fourLep = 1;
728
729 calculationResult["FourLep"] = fourLep;
730
731 // nLambda
732 unsigned int nLambda = 0;
733
734 if (m_LambdaParticles.isValid()) {
735 for (unsigned int i = 0; i < m_LambdaParticles->getListSize(); i++) {
736 const Particle* mergeLambdaCand = m_LambdaParticles->getParticle(i);
737 const double flightDist = Variable::flightDistance(mergeLambdaCand);
738 const double flightDistErr = Variable::flightDistanceErr(mergeLambdaCand);
739 const double flightSign = flightDist / flightDistErr;
740 const Particle* protCand = mergeLambdaCand->getDaughter(0);
741 const Particle* pionCand = mergeLambdaCand->getDaughter(1);
742 const double protMom = protCand->getP();
743 const double pionMom = pionCand->getP();
744 const double asymPDaughters = (protMom - pionMom) / (protMom + pionMom);
745 if (flightSign > 10 && asymPDaughters > 0.41) nLambda++;
746 }
747 }
748
749 if (nLambda > 0) {
750 calculationResult["Lambda"] = 1;
751 } else {
752 calculationResult["Lambda"] = 0;
753 }
754
755 // nDstp
756 unsigned int nDstp1 = 0;
757 unsigned int nDstp2 = 0;
758 unsigned int nDstp3 = 0;
759 unsigned int nDstp4 = 0;
760
761 if (m_DstParticles.isValid() && (ntrk_bha >= 3 && Bhabha2Trk == 0)) {
762 for (unsigned int i = 0; i < m_DstParticles->getListSize(); i++) {
763 const Particle* allDstCand = m_DstParticles->getParticle(i);
764 const double dstDecID = allDstCand->getExtraInfo("decayModeID");
765 if (dstDecID == 1.) nDstp1++;
766 if (dstDecID == 2.) nDstp2++;
767 if (dstDecID == 3.) nDstp3++;
768 if (dstDecID == 4.) nDstp4++;
769 }
770 }
771
772
773 if (nDstp1 > 0) {
774 calculationResult["Dstp1"] = 1;
775 } else {
776 calculationResult["Dstp1"] = 0;
777 }
778
779 if (nDstp2 > 0) {
780 calculationResult["Dstp2"] = 1;
781 } else {
782 calculationResult["Dstp2"] = 0;
783 }
784
785 if (nDstp3 > 0) {
786 calculationResult["Dstp3"] = 1;
787 } else {
788 calculationResult["Dstp3"] = 0;
789 }
790
791 if (nDstp4 > 0) {
792 calculationResult["Dstp4"] = 1;
793 } else {
794 calculationResult["Dstp4"] = 0;
795 }
796
797 // nTracksOffIP
798 calculationResult["nTracksOffIP"] = m_offIpParticles->getListSize();
799
800 // Flag for events with Trigger B2Link information
801 calculationResult["NeuroTRG"] = 0;
802 calculationResult["GammaGammaFilter"] = 0;
803
805 if (filter_result.isValid()) {
806 const std::map<std::string, int>& nonPrescaledResults = filter_result->getNonPrescaledResults();
807 if (nonPrescaledResults.find(m_filterL1TrgNN) != nonPrescaledResults.end()) {
808 const bool hasNN = (filter_result->getNonPrescaledResult(m_filterL1TrgNN) == SoftwareTriggerCutResult::c_accept);
809 if (hasNN) calculationResult["NeuroTRG"] = 1;
810 }
811 const bool ggEndcap = (filter_result->getNonPrescaledResult("software_trigger_cut&filter&ggEndcapLoose") ==
813 const bool ggBarrel = (filter_result->getNonPrescaledResult("software_trigger_cut&filter&ggBarrelLoose") ==
815 if (ggEndcap || ggBarrel) calculationResult["GammaGammaFilter"] = 1;
816 }
817
818 //Dimuon skim with invariant mass cut allowing at most one track not to be associated with ECL clusters
819
820 double mumuHighMass = 0.;
821
822 if (trackWithMaximumRho && trackWithSecondMaximumRho) {
823 int hasClus = 0;
824 double eclE1 = 0., eclE2 = 0.;
825
826 const auto charge1 = trackWithMaximumRho->getCharge();
827 const auto charge2 = trackWithSecondMaximumRho->getCharge();
828 const auto chSum = charge1 + charge2;
829
830 const ECLCluster* eclTrack1 = trackWithMaximumRho->getECLCluster();
831 if (eclTrack1) {
832 hasClus++;
834 }
835
836 const ECLCluster* eclTrack2 = trackWithSecondMaximumRho->getECLCluster();
837 if (eclTrack2) {
838 hasClus++;
840 }
841 const ROOT::Math::PxPyPzEVector V4p1 = PCmsLabTransform::labToCms(trackWithMaximumRho->get4Vector());
842 const ROOT::Math::PxPyPzEVector V4p2 = PCmsLabTransform::labToCms(trackWithSecondMaximumRho->get4Vector());
843
844 const ROOT::Math::PxPyPzEVector V4pSum = V4p1 + V4p2;
845 const double mSum = V4pSum.M();
846
847 const double thetaSumCMS = (V4p1.Theta() + V4p2.Theta()) * TMath::RadToDeg();
848 const double phi1CMS = V4p1.Phi() * TMath::RadToDeg();
849 const double phi2CMS = V4p2.Phi() * TMath::RadToDeg();
850
851 double diffPhi = phi1CMS - phi2CMS;
852 if (fabs(diffPhi) > 180) {
853 if (diffPhi > 180) {
854 diffPhi = diffPhi - 2 * 180;
855 } else {
856 diffPhi = 2 * 180 + diffPhi;
857 }
858 }
859 const double delThetaCMS = fabs(fabs(thetaSumCMS) - 180);
860 const double delPhiCMS = fabs(180 - fabs(diffPhi));
861
862 const bool mumuHighMassCand = chSum == 0 && (mSum > 8. && mSum < 12.) && hasClus > 0 && eclE1 <= 1
863 && eclE2 <= 1 && delThetaCMS < 10 && delPhiCMS < 10;
864
865 if (mumuHighMassCand) mumuHighMass = 1;
866
867 }
868
869 calculationResult["MumuHighM"] = mumuHighMass;
870
871 // BtoCharm skims
872 calculationResult["Bp"] = 0;
873 calculationResult["Bz"] = 0;
874
875 if (m_BpParticles.isValid() && (ntrk_bha >= 3 && Bhabha2Trk == 0)) {
876 calculationResult["Bp"] = m_BpParticles->getListSize() >= 1;
877 }
878
879 if (m_BzParticles.isValid() && (ntrk_bha >= 3 && Bhabha2Trk == 0)) {
880 calculationResult["Bz"] = m_BzParticles->getListSize() >= 1;
881 }
882
883}
DataType Theta() const
The polar angle.
Definition: B2Vector3.h:153
DataType DeltaPhi(const B2Vector3< DataType > &v) const
returns phi in the interval [-PI,PI)
Definition: B2Vector3.h:228
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:159
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
Definition: B2Vector3.h:302
Debug output for CDCDedxPID module.
Definition: CDCDedxTrack.h:25
double getDedxNoSat() const
Get dE/dx truncated mean without the saturation correction for this track.
Definition: CDCDedxTrack.h:107
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
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
ECL cluster data.
Definition: ECLCluster.h:27
bool hasHypothesis(EHypothesisBit bitmask) const
Return if specific hypothesis bit is set.
Definition: ECLCluster.h:351
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Definition: ECLCluster.cc:23
@ c_nPhotons
CR is split into n photons (N1)
Class to calculate the Fox-Wolfram moments up to order 8.
Definition: FoxWolfram.h:28
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.
ROOT::Math::PxPyPzEVector getBeamFourMomentum() const
Returns LAB four-momentum of e+e-, i.e.
static ROOT::Math::PxPyPzMVector labToCms(const ROOT::Math::PxPyPzMVector &vec)
Transforms Lorentz vector into CM System.
Class to collect log likelihoods from TOP, ARICH, dEdx, ECL and KLM aimed for output to mdst includes...
Definition: PIDLikelihood.h:29
bool isAvailable(Const::PIDDetectorSet set=Const::PIDDetectorSet::set()) const
Check whether PID information is available for at least one of the detectors in a given set.
Class to store reconstructed particles.
Definition: Particle.h:75
const Track * getTrack() const
Returns the pointer to the Track object that was used to create this Particle (ParticleType == c_Trac...
Definition: Particle.cc:845
const ECLCluster * getECLCluster() const
Returns the pointer to the ECLCluster object that was used to create this Particle (if ParticleType =...
Definition: Particle.cc:891
const PIDLikelihood * getPIDLikelihood() const
Returns the pointer to the PIDLikelihood object that is related to the Track, which was used to creat...
Definition: Particle.cc:871
double getCharge(void) const
Returns particle charge.
Definition: Particle.cc:622
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:547
ROOT::Math::XYZVector getMomentum() const
Returns momentum vector.
Definition: Particle.h:560
double getMomentumMagnitude() const
Returns momentum magnitude.
Definition: Particle.h:569
double getP() const
Returns momentum magnitude (same as getMomentumMagnitude but with shorter name)
Definition: Particle.h:578
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:631
double getExtraInfo(const std::string &name) const
Return given value if set.
Definition: Particle.cc:1289
double getMass() const
Returns invariant mass (= nominal for FS particles)
Definition: Particle.h:507
static const ReferenceFrame & GetCurrent()
Get current rest frame.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:288
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
iterator end()
Return iterator to last entry +1.
Definition: StoreArray.h:320
iterator begin()
Return iterator to first entry.
Definition: StoreArray.h:318
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:111
Values of the result of a track fit with a given particle hypothesis.
ROOT::Math::PxPyPzEVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
@ c_accept
Accept this event.

◆ 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();
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
55 }
virtual void doCalculation(SoftwareTriggerObject &m_calculationResult)=0
Override this function to implement your calculation.

◆ requireStoreArrays()

void requireStoreArrays ( )
overridevirtual

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

Implements SoftwareTriggerCalculation.

Definition at line 51 of file SkimSampleCalculator.cc.

52{
53 m_pionParticles.isRequired();
54 m_gammaParticles.isRequired();
55 m_pionHadParticles.isRequired();
56 m_pionTauParticles.isRequired();
57 m_KsParticles.isOptional();
58 m_LambdaParticles.isOptional();
59 m_DstParticles.isOptional();
60 m_offIpParticles.isRequired();
61 m_BpParticles.isOptional();
62 m_BzParticles.isOptional();
63
64};

◆ 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 }
bool m_debugPrepared
Flag to not add the branches twice to the TTree.

Member Data Documentation

◆ m_BpParticles

StoreObjPtr<ParticleList> m_BpParticles
private

Internal storage of the B+'s.

Definition at line 59 of file SkimSampleCalculator.h.

◆ m_BzParticles

StoreObjPtr<ParticleList> m_BzParticles
private

Internal storage of the B0's.

Definition at line 61 of file SkimSampleCalculator.h.

◆ m_calculationResult

SoftwareTriggerObject m_calculationResult
privateinherited

Internal storage of the result of the calculation.

Definition at line 74 of file SoftwareTriggerCalculation.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_DstParticles

StoreObjPtr<ParticleList> m_DstParticles
private

Internal storage of the D*'s.

Definition at line 53 of file SkimSampleCalculator.h.

◆ m_filterL1TrgNN

std::string m_filterL1TrgNN = ""
private

HLT filter line for the TRG skim.

Definition at line 57 of file SkimSampleCalculator.h.

◆ m_gammaParticles

StoreObjPtr<ParticleList> m_gammaParticles
private

Internal storage of the ECL clusters as particles.

Definition at line 43 of file SkimSampleCalculator.h.

◆ m_KsParticles

StoreObjPtr<ParticleList> m_KsParticles
private

Internal storage of the K_S0's.

Definition at line 49 of file SkimSampleCalculator.h.

◆ m_LambdaParticles

StoreObjPtr<ParticleList> m_LambdaParticles
private

Internal storage of the Lambda0's.

Definition at line 51 of file SkimSampleCalculator.h.

◆ m_offIpParticles

StoreObjPtr<ParticleList> m_offIpParticles
private

Internal storage of the tracks for alignment calibration.

Definition at line 55 of file SkimSampleCalculator.h.

◆ m_pionHadParticles

StoreObjPtr<ParticleList> m_pionHadParticles
private

Internal storage of the tracks as particles (definition for hadronb).

Definition at line 45 of file SkimSampleCalculator.h.

◆ m_pionParticles

StoreObjPtr<ParticleList> m_pionParticles
private

Internal storage of the tracks as particles.

Definition at line 41 of file SkimSampleCalculator.h.

◆ m_pionTauParticles

StoreObjPtr<ParticleList> m_pionTauParticles
private

Internal storage of the tracks as particles (definition for tau skims).

Definition at line 47 of file SkimSampleCalculator.h.


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