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