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