CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BPHHistoSpecificDecay.cc
Go to the documentation of this file.
1 
3 
5 
7 
9 
12 
16 
18 
20 #include "TMath.h"
21 #include "Math/VectorUtil.h"
22 #include "TVector3.h"
23 
24 #include <TH1.h>
25 #include <TTree.h>
26 #include <TFile.h>
27 
28 #include <set>
29 #include <string>
30 #include <iostream>
31 #include <fstream>
32 #include <cmath>
33 
34 using namespace std;
35 
36 #define SET_LABEL(NAME, PSET) (NAME = PSET.getParameter<string>(#NAME))
37 // SET_LABEL(xyz,ps);
38 // is equivalent to
39 // xyz = ps.getParameter<string>( "xyx" );
40 
41 #define CAT3(A, B, C) A##B##C
42 #define STRING_NX(A) #A
43 #define STRING(A) STRING_NX(A)
44 #define CHK_TRIG(RESULTS, NAMES, INDEX, PATH) \
45  if (NAMES[INDEX].find(STRING(CAT3(HLT_, PATH, _v))) == 0) { \
46  flag_##PATH = RESULTS->accept(INDEX); \
47  continue; \
48  }
49 // CHK_TRIG( trigResults, names, i, xyz );
50 // is equivalent to
51 // if ( names[i].find( "HLT_xyz_v" ) == 0 ) { flag_xyz = trigResults->accept( i ); break; }
52 
54 public:
55  static double cAlpha(const reco::Vertex* pvtx, const reco::Vertex* svtx, float px, float py) {
56  TVector3 disp(svtx->x() - pvtx->x(), svtx->y() - pvtx->y(), 0);
57  TVector3 cmom(px, py, 0);
58  return disp.Dot(cmom) / (disp.Perp() * cmom.Perp());
59  }
60  static double cAlpha(const reco::Vertex* pvtx, const reco::Vertex* svtx, const TVector3& cmom) {
61  TVector3 disp(svtx->x() - pvtx->x(), svtx->y() - pvtx->y(), 0);
62  return disp.Dot(cmom) / (disp.Perp() * cmom.Perp());
63  }
64  static void dist2D(const reco::Vertex* pvtx,
65  const reco::Vertex* svtx,
66  float px,
67  float py,
68  float mass,
69  double& ctauPV,
70  double& ctauErrPV) {
71  dist2D(pvtx, svtx, px, py, cAlpha(pvtx, svtx, px, py), mass, ctauPV, ctauErrPV);
72  return;
73  }
74  static void dist2D(const reco::Vertex* pvtx,
75  const reco::Vertex* svtx,
76  float px,
77  float py,
78  double cosAlpha,
79  float mass,
80  double& ctauPV,
81  double& ctauErrPV) {
82  TVector3 cmom(px, py, 0);
83  AlgebraicVector3 vmom(px, py, 0);
84  VertexDistanceXY vdistXY;
85  Measurement1D distXY = vdistXY.distance(*svtx, *pvtx);
86  ctauPV = distXY.value() * cosAlpha * mass / cmom.Perp();
87  GlobalError sve = svtx->error();
88  GlobalError pve = pvtx->error();
89  AlgebraicSymMatrix33 vXYe = sve.matrix() + pve.matrix();
90  ctauErrPV = sqrt(ROOT::Math::Similarity(vmom, vXYe)) * mass / cmom.Perp2();
91  return;
92  }
93 };
94 
95 class BPHUserData {
96 public:
97  template <class T>
98  static const T* get(const pat::CompositeCandidate& cand, const string& name) {
99  if (cand.hasUserData(name))
100  return cand.userData<T>(name);
101  return nullptr;
102  }
103  static float get(const pat::CompositeCandidate& cand, const string& name, float d = 0.0) {
104  if (cand.hasUserFloat(name))
105  return cand.userFloat(name);
106  return d;
107  }
108  template <class T>
109  static const T* getByRef(const pat::CompositeCandidate& cand, const string& name) {
110  if (cand.hasUserData(name)) {
111  typedef edm::Ref<std::vector<T>> objRef;
112  const objRef* ref = cand.userData<objRef>(name);
113  if (ref == nullptr)
114  return nullptr;
115  if (ref->isNull())
116  return nullptr;
117  return ref->get();
118  }
119  return nullptr;
120  }
121 };
122 
124 public:
125  static vector<const reco::Candidate*> get(const pat::CompositeCandidate& cand, float massMin, float massMax) {
126  int i;
127  int n = cand.numberOfDaughters();
128  vector<const reco::Candidate*> v;
129  v.reserve(n);
130  for (i = 0; i < n; ++i) {
131  const reco::Candidate* dptr = cand.daughter(i);
132  float mass = dptr->mass();
133  if ((mass > massMin) && (mass < massMax))
134  v.push_back(dptr);
135  }
136  return v;
137  }
138 };
139 
141 public:
142  BPHSoftMuonSelect(int cutTrackerLayers = 5,
143  int cutPixelLayers = 0,
144  float maxDxy = 0.3,
145  float maxDz = 20.0,
146  bool goodMuon = true,
147  bool highPurity = true)
148  : cutTL(cutTrackerLayers), cutPL(cutPixelLayers), maxXY(maxDxy), maxZ(maxDz), gM(goodMuon), hP(highPurity) {}
150  bool accept(const reco::Candidate& cand, const reco::Vertex* pv) const {
151  const pat::Muon* p = dynamic_cast<const pat::Muon*>(&cand);
152  if (p == nullptr)
153  return false;
155  return false;
156  if (p->innerTrack()->hitPattern().trackerLayersWithMeasurement() <= cutTL)
157  return false;
158  if (p->innerTrack()->hitPattern().pixelLayersWithMeasurement() <= cutPL)
159  return false;
160  if (hP && !p->innerTrack()->quality(reco::TrackBase::highPurity))
161  return false;
162  if (pv == nullptr)
163  return true;
164  const reco::Vertex::Point& pos = pv->position();
165  if (fabs(p->innerTrack()->dxy(pos)) >= maxXY)
166  return false;
167  if (fabs(p->innerTrack()->dz(pos)) >= maxZ)
168  return false;
169  return true;
170  }
171 
172 private:
173  const reco::Vertex* pv;
174  int cutTL;
175  int cutPL;
176  float maxXY;
177  float maxZ;
178  bool gM;
179  bool hP;
180 };
181 
183 public:
184  BPHDaughterSelect(float ptMinLoose,
185  float ptMinTight,
186  float etaMaxLoose,
187  float etaMaxTight,
188  const BPHSoftMuonSelect* softMuonselector = nullptr)
189  : pLMin(ptMinLoose), pTMin(ptMinTight), eLMax(etaMaxLoose), eTMax(etaMaxTight), sms(softMuonselector) {}
190  ~BPHDaughterSelect() override {}
191  bool accept(const pat::CompositeCandidate& cand, const reco::Vertex* pv = nullptr) const override {
192  return accept(cand, pLMin, pTMin, eLMax, eTMax, pv, sms);
193  }
194  static bool accept(const pat::CompositeCandidate& cand,
195  float ptMinLoose,
196  float ptMinTight,
197  float etaMaxLoose,
198  float etaMaxTight,
199  const reco::Vertex* pv = nullptr,
200  const BPHSoftMuonSelect* softMuonselector = nullptr) {
201  const reco::Candidate* dptr0 = cand.daughter(0);
202  const reco::Candidate* dptr1 = cand.daughter(1);
203  if (dptr0 == nullptr)
204  return false;
205  if (dptr1 == nullptr)
206  return false;
207  float pt0 = dptr0->pt();
208  float pt1 = dptr1->pt();
209  if ((pt0 < ptMinLoose) || (pt1 < ptMinLoose))
210  return false;
211  if ((pt0 < ptMinTight) && (pt1 < ptMinTight))
212  return false;
213  float eta0 = fabs(dptr0->eta());
214  float eta1 = fabs(dptr1->eta());
215  if ((etaMaxLoose > 0) && ((eta0 > etaMaxLoose) || (eta1 > etaMaxLoose)))
216  return false;
217  if ((etaMaxTight > 0) && ((eta0 > etaMaxTight) && (eta1 > etaMaxTight)))
218  return false;
219  if (softMuonselector != nullptr) {
220  const reco::Vertex* pvtx = BPHUserData::getByRef<reco::Vertex>(cand, "primaryVertex");
221  if (pvtx == nullptr)
222  return false;
223  if (!softMuonselector->accept(*dptr0, pvtx))
224  return false;
225  if (!softMuonselector->accept(*dptr1, pvtx))
226  return false;
227  }
228  return true;
229  }
230 
231 private:
232  float pLMin;
233  float pTMin;
234  float eLMax;
235  float eTMax;
237 };
238 
240 public:
242  float massMin, float massMax, float ptMin = -1.0, float etaMax = -1.0, float rapidityMax = -1.0)
243  : mMin(massMin), mMax(massMax), pMin(ptMin), eMax(etaMax), yMax(rapidityMax) {}
245  bool accept(const pat::CompositeCandidate& cand, const reco::Vertex* pv = nullptr) const override {
246  if (((mMin > 0) && (mMax < 0)) || ((mMin < 0) && (mMax > 0)) || ((mMin > 0) && (mMax > 0) && (mMin < mMax))) {
247  float mass = cand.mass();
248  if (mass < mMin)
249  return false;
250  if ((mMax > 0) && (mass > mMax))
251  return false;
252  }
253  if (cand.pt() < pMin)
254  return false;
255  if ((eMax > 0) && (fabs(cand.eta()) > eMax))
256  return false;
257  if ((yMax > 0) && (fabs(cand.rapidity()) > yMax))
258  return false;
259  return true;
260  }
261 
262 private:
263  float mMin;
264  float mMax;
265  float pMin;
266  float eMax;
267  float yMax;
268 };
269 
271 public:
272  BPHFittedBasicSelect(float massMin, float massMax, float ptMin = -1.0, float etaMax = -1.0, float rapidityMax = -1.0)
273  : mMin(massMin), mMax(massMax), pMin(ptMin), eMax(etaMax), yMax(rapidityMax) {}
274  ~BPHFittedBasicSelect() override {}
275  bool accept(const pat::CompositeCandidate& cand, const reco::Vertex* pv = nullptr) const override {
276  if (!cand.hasUserFloat("fitMass"))
277  return false;
278  float mass = cand.userFloat("fitMass");
279  if (((mMin > 0) && (mMax < 0)) || ((mMin < 0) && (mMax > 0)) || ((mMin > 0) && (mMax > 0) && (mMin < mMax))) {
280  if (mass < mMin)
281  return false;
282  if ((mMax > 0) && (mass > mMax))
283  return false;
284  }
285  const Vector3DBase<float, GlobalTag>* fmom = BPHUserData::get<Vector3DBase<float, GlobalTag>>(cand, "fitMomentum");
286  if (fmom == nullptr)
287  return false;
288  if (pMin > 0) {
289  if (fmom->transverse() < pMin)
290  return false;
291  }
292  if (eMax > 0) {
293  if (fabs(fmom->eta()) > eMax)
294  return false;
295  }
296  if (yMax > 0) {
297  float x = fmom->x();
298  float y = fmom->y();
299  float z = fmom->z();
300  float e = sqrt((x * x) + (y * y) + (z * z) + (mass * mass));
301  float r = log((e + z) / (e - z)) / 2;
302  if (fabs(r) > yMax)
303  return false;
304  }
305  return true;
306  }
307 
308 private:
309  float mMin;
310  float mMax;
311  float pMin;
312  float eMax;
313  float yMax;
314 };
315 
317 public:
318  BPHGenericVertexSelect(char vType, float probMin, float cosMin = -2.0, float sigMin = -1.0, char dMode = 'r')
319  : type(vType), pMin(probMin), cMin(cosMin), sMin(sigMin), mode(dMode) {}
321  bool accept(const pat::CompositeCandidate& cand, const reco::Vertex* pvtx) const override {
322  if (pvtx == nullptr)
323  return false;
324  const reco::Vertex* svtx = nullptr;
325  float px;
326  float py;
327  float mass;
328  switch (type) {
329  case 'c':
330  svtx = BPHUserData::get<reco::Vertex>(cand, "vertex");
331  px = cand.px();
332  py = cand.py();
333  mass = cand.mass();
334  break;
335  case 'f':
336  svtx = BPHUserData::get<reco::Vertex>(cand, "fitVertex");
337  {
338  const Vector3DBase<float, GlobalTag>* fmom =
339  BPHUserData::get<Vector3DBase<float, GlobalTag>>(cand, "fitMomentum");
340  if (fmom == nullptr)
341  return false;
342  px = fmom->x();
343  py = fmom->y();
344  }
345  if (!cand.hasUserFloat("fitMass"))
346  return false;
347  mass = cand.userFloat("fitMass");
348  break;
349  default:
350  return false;
351  }
352  if (svtx == nullptr)
353  return false;
354  if (pMin > 0) {
355  if (ChiSquaredProbability(svtx->chi2(), svtx->ndof()) < pMin)
356  return false;
357  }
358  if ((cMin > -1.0) || (sMin > 0)) {
359  float cosAlpha = VertexAnalysis::cAlpha(pvtx, svtx, px, py);
360  if (cosAlpha < cMin)
361  return false;
362  if (sMin < 0)
363  return true;
364  double ctauPV;
365  double ctauErrPV;
366  VertexAnalysis::dist2D(pvtx, svtx, px, py, cosAlpha, mass, ctauPV, ctauErrPV);
367  float dTest;
368  switch (mode) {
369  case 'a':
370  case 'd':
371  dTest = ctauPV;
372  break;
373  case 'r':
374  case 's':
375  default:
376  dTest = ctauPV / ctauErrPV;
377  break;
378  }
379  if (dTest < sMin)
380  return false;
381  }
382  return true;
383  }
384 
385 private:
386  char type;
387  float pMin;
388  float cMin;
389  float sMin;
390  char mode;
391 };
392 
394  useTrig = (!SET_LABEL(trigResultsLabel, ps).empty());
395  useOnia = (!SET_LABEL(oniaCandsLabel, ps).empty());
396  useSd = (!SET_LABEL(sdCandsLabel, ps).empty());
397  useSs = (!SET_LABEL(ssCandsLabel, ps).empty());
398  useBu = (!SET_LABEL(buCandsLabel, ps).empty());
399  useBd = (!SET_LABEL(bdCandsLabel, ps).empty());
400  useBs = (!SET_LABEL(bsCandsLabel, ps).empty());
401  useK0 = (!SET_LABEL(k0CandsLabel, ps).empty());
402  useL0 = (!SET_LABEL(l0CandsLabel, ps).empty());
403  useB0 = (!SET_LABEL(b0CandsLabel, ps).empty());
404  useLb = (!SET_LABEL(lbCandsLabel, ps).empty());
405  useBc = (!SET_LABEL(bcCandsLabel, ps).empty());
406  useX3872 = (!SET_LABEL(x3872CandsLabel, ps).empty());
407  if (useTrig)
408  consume<edm::TriggerResults>(trigResultsToken, trigResultsLabel);
409  if (useOnia)
410  consume<vector<pat::CompositeCandidate>>(oniaCandsToken, oniaCandsLabel);
411  if (useSd)
412  consume<vector<pat::CompositeCandidate>>(sdCandsToken, sdCandsLabel);
413  if (useSs)
414  consume<vector<pat::CompositeCandidate>>(ssCandsToken, ssCandsLabel);
415  if (useBu)
416  consume<vector<pat::CompositeCandidate>>(buCandsToken, buCandsLabel);
417  if (useBd)
418  consume<vector<pat::CompositeCandidate>>(bdCandsToken, bdCandsLabel);
419  if (useBs)
420  consume<vector<pat::CompositeCandidate>>(bsCandsToken, bsCandsLabel);
421  if (useK0)
422  consume<vector<pat::CompositeCandidate>>(k0CandsToken, k0CandsLabel);
423  if (useL0)
424  consume<vector<pat::CompositeCandidate>>(l0CandsToken, l0CandsLabel);
425  if (useB0)
426  consume<vector<pat::CompositeCandidate>>(b0CandsToken, b0CandsLabel);
427  if (useLb)
428  consume<vector<pat::CompositeCandidate>>(lbCandsToken, lbCandsLabel);
429  if (useBc)
430  consume<vector<pat::CompositeCandidate>>(bcCandsToken, bcCandsLabel);
431  if (useX3872)
432  consume<vector<pat::CompositeCandidate>>(x3872CandsToken, x3872CandsLabel);
433 
434  static const BPHSoftMuonSelect* sms = new BPHSoftMuonSelect;
435 
437 
438  double phiIMassMin = 0.85;
439  double phiIMassMax = 1.20;
440  double phiIPtMin = 18.0;
441  double phiIEtaMax = -1.0;
442  double phiIYMax = -1.0;
443  double phiBMassMin = 0.85;
444  double phiBMassMax = 1.20;
445  double phiBPtMin = 14.0;
446  double phiBEtaMax = -1.0;
447  double phiBYMax = 1.25;
448  double jPsiIMassMin = 2.80;
449  double jPsiIMassMax = 3.40;
450  double jPsiIPtMin = 25.0;
451  double jPsiIEtaMax = -1.0;
452  double jPsiIYMax = -1.0;
453  double jPsiBMassMin = 2.80;
454  double jPsiBMassMax = 3.40;
455  double jPsiBPtMin = 20.0;
456  double jPsiBEtaMax = -1.0;
457  double jPsiBYMax = 1.25;
458  double psi2IMassMin = 3.40;
459  double psi2IMassMax = 4.00;
460  double psi2IPtMin = 18.0;
461  double psi2IEtaMax = -1.0;
462  double psi2IYMax = -1.0;
463  double psi2BMassMin = 3.40;
464  double psi2BMassMax = 4.00;
465  double psi2BPtMin = 10.0;
466  double psi2BEtaMax = -1.0;
467  double psi2BYMax = 1.25;
468  double upsIMassMin = 8.50;
469  double upsIMassMax = 11.0;
470  double upsIPtMin = 15.0;
471  double upsIEtaMax = -1.0;
472  double upsIYMax = -1.0;
473  double upsBMassMin = 8.50;
474  double upsBMassMax = 11.0;
475  double upsBPtMin = 12.0;
476  // double upsBEtaMax = 1.5; // 2017
477  // double upsBYMax = -1.0; // 2017
478  double upsBEtaMax = -1.0; // 2018
479  double upsBYMax = 1.4; // 2018
480 
481  double oniaProbMin = 0.005;
482  double oniaCosMin = -2.0;
483  double oniaSigMin = -1.0;
484 
485  double oniaMuPtMinLoose = 2.0;
486  double oniaMuPtMinTight = -1.0;
487  double oniaMuEtaMaxLoose = -1.0;
488  double oniaMuEtaMaxTight = -1.0;
489 
490  phiIBasicSelect = new BPHCompositeBasicSelect(phiIMassMin, phiIMassMax, phiIPtMin, phiIEtaMax, phiIYMax);
491  phiBBasicSelect = new BPHCompositeBasicSelect(phiBMassMin, phiBMassMax, phiBPtMin, phiBEtaMax, phiBYMax);
492  jPsiIBasicSelect = new BPHCompositeBasicSelect(jPsiIMassMin, jPsiIMassMax, jPsiIPtMin, jPsiIEtaMax, jPsiIYMax);
493  jPsiBBasicSelect = new BPHCompositeBasicSelect(jPsiBMassMin, jPsiBMassMax, jPsiBPtMin, jPsiBEtaMax, jPsiBYMax);
494  psi2IBasicSelect = new BPHCompositeBasicSelect(psi2IMassMin, psi2IMassMax, psi2IPtMin, psi2IEtaMax, psi2IYMax);
495  psi2BBasicSelect = new BPHCompositeBasicSelect(psi2BMassMin, psi2BMassMax, psi2BPtMin, psi2BEtaMax, psi2BYMax);
496  upsIBasicSelect = new BPHCompositeBasicSelect(upsIMassMin, upsIMassMax, upsIPtMin, upsIEtaMax, upsIYMax);
497  upsBBasicSelect = new BPHCompositeBasicSelect(upsBMassMin, upsBMassMax, upsBPtMin, upsBEtaMax, upsBYMax);
498  // oniaVertexSelect = new BPHCompositeVertexSelect(
499  oniaVertexSelect = new BPHGenericVertexSelect('c', oniaProbMin, oniaCosMin, oniaSigMin);
500  oniaDaughterSelect =
501  new BPHDaughterSelect(oniaMuPtMinLoose, oniaMuPtMinTight, oniaMuEtaMaxLoose, oniaMuEtaMaxTight, sms);
502 
503  double npJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
504  double npJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
505  double npJPsiPtMin = 8.0;
506  double npJPsiEtaMax = -1.0;
507  double npJPsiYMax = -1.0;
508  double npMuPtMinLoose = 4.0;
509  double npMuPtMinTight = -1.0;
510  double npMuEtaMaxLoose = 2.2;
511  double npMuEtaMaxTight = -1.0;
512 
513  npJPsiBasicSelect = new BPHCompositeBasicSelect(npJPsiMassMin, npJPsiMassMax, npJPsiPtMin, npJPsiEtaMax, npJPsiYMax);
514  npJPsiDaughterSelect = new BPHDaughterSelect(npMuPtMinLoose, npMuPtMinTight, npMuEtaMaxLoose, npMuEtaMaxTight, sms);
515 
517 
518  double buIMassMin = 0.0;
519  double buIMassMax = 999999.0;
520  double buIPtMin = 27.0;
521  double buIEtaMax = -1.0;
522  double buIYMax = -1.0;
523  double buIJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
524  double buIJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
525  double buIJPsiPtMin = 25.0;
526  double buIJPsiEtaMax = -1.0;
527  double buIJPsiYMax = -1.0;
528  double buIProbMin = 0.15;
529  double buICosMin = -2.0;
530  double buISigMin = -1.0;
531  // double buIMuPtMinLoose = -1.0;
532  // double buIMuPtMinTight = -1.0;
533  // double buIMuEtaMaxLoose = -1.0;
534  // double buIMuEtaMaxTight = -1.0;
535 
536  buIKPtMin = 2.0;
537 
538  buIBasicSelect = new BPHFittedBasicSelect(buIMassMin, buIMassMax, buIPtMin, buIEtaMax, buIYMax);
539  buIJPsiBasicSelect =
540  new BPHCompositeBasicSelect(buIJPsiMassMin, buIJPsiMassMax, buIJPsiPtMin, buIJPsiEtaMax, buIJPsiYMax);
541  buIVertexSelect = new BPHGenericVertexSelect('f', buIProbMin, buICosMin, buISigMin);
542  buIJPsiDaughterSelect = nullptr;
543  // buIJPsiDaughterSelect = new BPHDaughterSelect(
544  // buIMuPtMinLoose , buIMuPtMinTight ,
545  // buIMuEtaMaxLoose, buMuEtaMaxTight, sms );
546 
547  double buDMassMin = 0.0;
548  double buDMassMax = 999999.0;
549  double buDPtMin = 10.0;
550  double buDEtaMax = -1.0;
551  double buDYMax = -1.0;
552  double buDJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
553  double buDJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
554  double buDJPsiPtMin = 8.0;
555  double buDJPsiEtaMax = -1.0;
556  double buDJPsiYMax = -1.0;
557  double buDProbMin = 0.10;
558  double buDCosMin = 0.99;
559  double buDSigMin = 3.0;
560  // double buDMuPtMinLoose = -1.0;
561  // double buDMuPtMinTight = -1.0;
562  // double buDMuEtaMaxLoose = -1.0;
563  // double buDMuEtaMaxTight = -1.0;
564 
565  buDKPtMin = 1.6;
566 
567  buDBasicSelect = new BPHFittedBasicSelect(buDMassMin, buDMassMax, buDPtMin, buDEtaMax, buDYMax);
568  buDJPsiBasicSelect =
569  new BPHCompositeBasicSelect(buDJPsiMassMin, buDJPsiMassMax, buDJPsiPtMin, buDJPsiEtaMax, buDJPsiYMax);
570  buDVertexSelect = new BPHGenericVertexSelect('f', buDProbMin, buDCosMin, buDSigMin);
571  buDJPsiDaughterSelect = nullptr;
572  // buDJPsiDaughterSelect = new BPHDaughterSelect(
573  // buDMuPtMinLoose , buDMuPtMinTight ,
574  // buDMuEtaMaxLoose, buDMuEtaMaxTight, sms );
575 
577 
578  double bdIMassMin = 0.0;
579  double bdIMassMax = 999999.0;
580  double bdIPtMin = 27.0;
581  double bdIEtaMax = -1.0;
582  double bdIYMax = -1.0;
583  double bdIJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
584  double bdIJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
585  double bdIJPsiPtMin = 25.0;
586  double bdIJPsiEtaMax = -1.0;
587  double bdIJPsiYMax = -1.0;
588  double bdIKx0MassMin = BPHParticleMasses::kx0Mass - 0.100;
589  double bdIKx0MassMax = BPHParticleMasses::kx0Mass + 0.100;
590  double bdIKx0PtMin = -1.0;
591  double bdIKx0EtaMax = -1.0;
592  double bdIKx0YMax = -1.0;
593  double bdIProbMin = 0.15;
594  double bdICosMin = -2.0;
595  double bdISigMin = -1.0;
596  // double bdIMuPtMinLoose = -1.0;
597  // double bdIMuPtMinTight = -1.0;
598  // double bdIMuEtaMaxLoose = -1.0;
599  // double bdIMuEtaMaxTight = -1.0;
600 
601  bdIBasicSelect = new BPHFittedBasicSelect(bdIMassMin, bdIMassMax, bdIPtMin, bdIEtaMax, bdIYMax);
602  bdIJPsiBasicSelect =
603  new BPHCompositeBasicSelect(bdIJPsiMassMin, bdIJPsiMassMax, bdIJPsiPtMin, bdIJPsiEtaMax, bdIJPsiYMax);
604  bdIKx0BasicSelect = new BPHCompositeBasicSelect(bdIKx0MassMin, bdIKx0MassMax, bdIKx0PtMin, bdIKx0EtaMax, bdIKx0YMax);
605  bdIVertexSelect = new BPHGenericVertexSelect('f', bdIProbMin, bdICosMin, bdISigMin);
606  bdIJPsiDaughterSelect = nullptr;
607  // bdIJPsiDaughterSelect = new BPHDaughterSelect(
608  // bdIMuPtMinLoose , bdIMuPtMinTight ,
609  // bdIMuEtaMaxLoose, bdIMuEtaMaxTight, sms );
610 
611  double bdDMassMin = 0.0;
612  double bdDMassMax = 999999.0;
613  double bdDPtMin = 10.0;
614  double bdDEtaMax = -1.0;
615  double bdDYMax = -1.0;
616  double bdDJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
617  double bdDJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
618  double bdDJPsiPtMin = 8.0;
619  double bdDJPsiEtaMax = -1.0;
620  double bdDJPsiYMax = -1.0;
621  double bdDKx0MassMin = BPHParticleMasses::kx0Mass - 0.100;
622  double bdDKx0MassMax = BPHParticleMasses::kx0Mass + 0.100;
623  double bdDKx0PtMin = -1.0;
624  double bdDKx0EtaMax = -1.0;
625  double bdDKx0YMax = -1.0;
626  double bdDProbMin = 0.10;
627  double bdDCosMin = 0.99;
628  double bdDSigMin = 3.0;
629  // double bdDMuPtMinLoose = -1.0;
630  // double bdDMuPtMinTight = -1.0;
631  // double bdDMuEtaMaxLoose = -1.0;
632  // double bdDMuEtaMaxTight = -1.0;
633 
634  bdDBasicSelect = new BPHFittedBasicSelect(bdDMassMin, bdDMassMax, bdDPtMin, bdDEtaMax, bdDYMax);
635  bdDJPsiBasicSelect =
636  new BPHCompositeBasicSelect(bdDJPsiMassMin, bdDJPsiMassMax, bdDJPsiPtMin, bdDJPsiEtaMax, bdDJPsiYMax);
637  bdDKx0BasicSelect = new BPHCompositeBasicSelect(bdDKx0MassMin, bdDKx0MassMax, bdDKx0PtMin, bdDKx0EtaMax, bdDKx0YMax);
638  bdDVertexSelect = new BPHGenericVertexSelect('f', bdDProbMin, bdDCosMin, bdDSigMin);
639  bdDJPsiDaughterSelect = nullptr;
640  // bdDJPsiDaughterSelect = new BPHDaughterSelect(
641  // bdDMuPtMinLoose , bdDMuPtMinTight ,
642  // bdDMuEtaMaxLoose, bdDMuEtaMaxTight, sms );
643 
645 
646  double bsIMassMin = 0.0;
647  double bsIMassMax = 999999.0;
648  double bsIPtMin = 27.0;
649  double bsIEtaMax = -1.0;
650  double bsIYMax = -1.0;
651  double bsIJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
652  double bsIJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
653  double bsIJPsiPtMin = 25.0;
654  double bsIJPsiEtaMax = -1.0;
655  double bsIJPsiYMax = -1.0;
656  double bsIPhiMassMin = BPHParticleMasses::phiMass - 0.010;
657  double bsIPhiMassMax = BPHParticleMasses::phiMass + 0.010;
658  double bsIPhiPtMin = -1.0;
659  double bsIPhiEtaMax = -1.0;
660  double bsIPhiYMax = -1.0;
661  double bsIProbMin = 0.15;
662  double bsICosMin = -2.0;
663  double bsISigMin = -1.0;
664  // double bsIMuPtMinLoose = -1.0;
665  // double bsIMuPtMinTight = -1.0;
666  // double bsIMuEtaMaxLoose = -1.0;
667  // double bsIMuEtaMaxTight = -1.0;
668 
669  bsIBasicSelect = new BPHFittedBasicSelect(bsIMassMin, bsIMassMax, bsIPtMin, bsIEtaMax, bsIYMax);
670  bsIJPsiBasicSelect =
671  new BPHCompositeBasicSelect(bsIJPsiMassMin, bsIJPsiMassMax, bsIJPsiPtMin, bsIJPsiEtaMax, bsIJPsiYMax);
672  bsIPhiBasicSelect = new BPHCompositeBasicSelect(bsIPhiMassMin, bsIPhiMassMax, bsIPhiPtMin, bsIPhiEtaMax, bsIPhiYMax);
673  bsIVertexSelect = new BPHGenericVertexSelect('f', bsIProbMin, bsICosMin, bsISigMin);
674  bsIJPsiDaughterSelect = nullptr;
675  // bsIJPsiDaughterSelect = new BPHDaughterSelect(
676  // bsIMuPtMinLoose , bsIMuPtMinTight ,
677  // bsIMuEtaMaxLoose, bsIMuEtaMaxTight, sms );
678 
679  double bsDMassMin = 0.0;
680  double bsDMassMax = 999999.0;
681  double bsDPtMin = 10.0;
682  double bsDEtaMax = -1.0;
683  double bsDYMax = -1.0;
684  double bsDJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
685  double bsDJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
686  double bsDJPsiPtMin = 8.0;
687  double bsDJPsiEtaMax = -1.0;
688  double bsDJPsiYMax = -1.0;
689  double bsDPhiMassMin = BPHParticleMasses::phiMass - 0.010;
690  double bsDPhiMassMax = BPHParticleMasses::phiMass + 0.010;
691  double bsDPhiPtMin = -1.0;
692  double bsDPhiEtaMax = -1.0;
693  double bsDPhiYMax = -1.0;
694  double bsDProbMin = 0.10;
695  double bsDCosMin = 0.99;
696  double bsDSigMin = 3.0;
697  // double bsDMuPtMinLoose = -1.0;
698  // double bsDMuPtMinTight = -1.0;
699  // double bsDMuEtaMaxLoose = -1.0;
700  // double bsDMuEtaMaxTight = -1.0;
701 
702  bsDBasicSelect = new BPHFittedBasicSelect(bsDMassMin, bsDMassMax, bsDPtMin, bsDEtaMax, bsDYMax);
703  bsDJPsiBasicSelect =
704  new BPHCompositeBasicSelect(bsDJPsiMassMin, bsDJPsiMassMax, bsDJPsiPtMin, bsDJPsiEtaMax, bsDJPsiYMax);
705  bsDPhiBasicSelect = new BPHCompositeBasicSelect(bsDPhiMassMin, bsDPhiMassMax, bsDPhiPtMin, bsDPhiEtaMax, bsDPhiYMax);
706  bsDVertexSelect = new BPHGenericVertexSelect('f', bsDProbMin, bsDCosMin, bsDSigMin);
707  bsDJPsiDaughterSelect = nullptr;
708  // bsDJPsiDaughterSelect = new BPHDaughterSelect(
709  // bsDMuPtMinLoose , bsDMuPtMinTight ,
710  // bsDMuEtaMaxLoose, bsDMuEtaMaxTight, sms );
711 
713 
714  double b0IMassMin = 0.0;
715  double b0IMassMax = 999999.0;
716  double b0IPtMin = 27.0;
717  double b0IEtaMax = -1.0;
718  double b0IYMax = -1.0;
719  double b0IJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
720  double b0IJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
721  double b0IJPsiPtMin = 25.0;
722  double b0IJPsiEtaMax = -1.0;
723  double b0IJPsiYMax = -1.0;
724  double b0IK0sMassMin = BPHParticleMasses::k0sMass - 0.010;
725  double b0IK0sMassMax = BPHParticleMasses::k0sMass + 0.010;
726  double b0IK0sPtMin = -1.0;
727  double b0IK0sEtaMax = -1.0;
728  double b0IK0sYMax = -1.0;
729  double b0IProbMin = 0.15;
730  double b0ICosMin = -2.0;
731  double b0ISigMin = -1.0;
732  // double b0IMuPtMinLoose = -1.0;
733  // double b0IMuPtMinTight = -1.0;
734  // double b0IMuEtaMaxLoose = -1.0;
735  // double b0IMuEtaMaxTight = -1.0;
736 
737  b0IBasicSelect = new BPHFittedBasicSelect(b0IMassMin, b0IMassMax, b0IPtMin, b0IEtaMax, b0IYMax);
738  b0IJPsiBasicSelect =
739  new BPHCompositeBasicSelect(b0IJPsiMassMin, b0IJPsiMassMax, b0IJPsiPtMin, b0IJPsiEtaMax, b0IJPsiYMax);
740  b0IK0sBasicSelect = new BPHFittedBasicSelect(b0IK0sMassMin, b0IK0sMassMax, b0IK0sPtMin, b0IK0sEtaMax, b0IK0sYMax);
741  b0IVertexSelect = new BPHGenericVertexSelect('f', b0IProbMin, b0ICosMin, b0ISigMin);
742  b0IJPsiDaughterSelect = nullptr;
743  // b0IJPsiDaughterSelect = new BPHDaughterSelect(
744  // b0IMuPtMinLoose , b0IMuPtMinTight ,
745  // b0IMuEtaMaxLoose, b0IMuEtaMaxTight, sms );
746 
747  double b0DMassMin = 0.0;
748  double b0DMassMax = 999999.0;
749  double b0DPtMin = 10.0;
750  double b0DEtaMax = -1.0;
751  double b0DYMax = -1.0;
752  double b0DJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
753  double b0DJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
754  double b0DJPsiPtMin = 8.0;
755  double b0DJPsiEtaMax = -1.0;
756  double b0DJPsiYMax = -1.0;
757  double b0DK0sMassMin = BPHParticleMasses::k0sMass - 0.010;
758  double b0DK0sMassMax = BPHParticleMasses::k0sMass + 0.010;
759  double b0DK0sPtMin = -1.0;
760  double b0DK0sEtaMax = -1.0;
761  double b0DK0sYMax = -1.0;
762  double b0DProbMin = 0.10;
763  double b0DCosMin = 0.99;
764  double b0DSigMin = 3.0;
765  // double b0DMuPtMinLoose = -1.0;
766  // double b0DMuPtMinTight = -1.0;
767  // double b0DMuEtaMaxLoose = -1.0;
768  // double b0DMuEtaMaxTight = -1.0;
769 
770  b0DBasicSelect = new BPHFittedBasicSelect(b0DMassMin, b0DMassMax, b0DPtMin, b0DEtaMax, b0DYMax);
771  b0DJPsiBasicSelect =
772  new BPHCompositeBasicSelect(b0DJPsiMassMin, b0DJPsiMassMax, b0DJPsiPtMin, b0DJPsiEtaMax, b0DJPsiYMax);
773  b0DK0sBasicSelect = new BPHFittedBasicSelect(b0DK0sMassMin, b0DK0sMassMax, b0DK0sPtMin, b0DK0sEtaMax, b0DK0sYMax);
774  b0DVertexSelect = new BPHGenericVertexSelect('f', b0DProbMin, b0DCosMin, b0DSigMin);
775  b0DJPsiDaughterSelect = nullptr;
776  // b0DJPsiDaughterSelect = new BPHDaughterSelect(
777  // b0DMuPtMinLoose , b0DMuPtMinTight ,
778  // b0DMuEtaMaxLoose, b0DMuEtaMaxTight, sms );
779 
781 
782  double lbIMassMin = 0.0;
783  double lbIMassMax = 999999.0;
784  double lbIPtMin = 27.0;
785  double lbIEtaMax = -1.0;
786  double lbIYMax = -1.0;
787  double lbIJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
788  double lbIJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
789  double lbIJPsiPtMin = 25.0;
790  double lbIJPsiEtaMax = -1.0;
791  double lbIJPsiYMax = -1.0;
792  double lbILambda0MassMin = BPHParticleMasses::lambda0Mass - 0.006;
793  double lbILambda0MassMax = BPHParticleMasses::lambda0Mass + 0.006;
794  double lbILambda0PtMin = -1.0;
795  double lbILambda0EtaMax = -1.0;
796  double lbILambda0YMax = -1.0;
797  double lbIProbMin = 0.10;
798  double lbICosMin = -2.0;
799  double lbISigMin = -1.0;
800  // double lbIMuPtMinLoose = -1.0;
801  // double lbIMuPtMinTight = -1.0;
802  // double lbIMuEtaMaxLoose = -1.0;
803  // double lbIMuEtaMaxTight = -1.0;
804 
805  lbIBasicSelect = new BPHFittedBasicSelect(lbIMassMin, lbIMassMax, lbIPtMin, lbIEtaMax, lbIYMax);
806  lbIJPsiBasicSelect =
807  new BPHCompositeBasicSelect(lbIJPsiMassMin, lbIJPsiMassMax, lbIJPsiPtMin, lbIJPsiEtaMax, lbIJPsiYMax);
808  lbILambda0BasicSelect =
809  new BPHFittedBasicSelect(lbILambda0MassMin, lbILambda0MassMax, lbILambda0PtMin, lbILambda0EtaMax, lbILambda0YMax);
810  lbIVertexSelect = new BPHGenericVertexSelect('f', lbIProbMin, lbICosMin, lbISigMin);
811  lbIJPsiDaughterSelect = nullptr;
812  // lbIJPsiDaughterSelect = new BPHDaughterSelect(
813  // lbIMuPtMinLoose , lbIMuPtMinTight ,
814  // lbIMuEtaMaxLoose, lbIMuEtaMaxTight, sms );
815 
816  double lbDMassMin = 0.0;
817  double lbDMassMax = 999999.0;
818  double lbDPtMin = 10.0;
819  double lbDEtaMax = -1.0;
820  double lbDYMax = -1.0;
821  double lbDJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
822  double lbDJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
823  double lbDJPsiPtMin = 8.0;
824  double lbDJPsiEtaMax = -1.0;
825  double lbDJPsiYMax = -1.0;
826  double lbDLambda0MassMin = BPHParticleMasses::lambda0Mass - 0.006;
827  double lbDLambda0MassMax = BPHParticleMasses::lambda0Mass + 0.006;
828  double lbDLambda0PtMin = -1.0;
829  double lbDLambda0EtaMax = -1.0;
830  double lbDLambda0YMax = -1.0;
831  double lbDProbMin = 0.10;
832  double lbDCosMin = 0.99;
833  double lbDSigMin = 3.0;
834  // double lbDMuPtMinLoose = -1.0;
835  // double lbDMuPtMinTight = -1.0;
836  // double lbDMuEtaMaxLoose = -1.0;
837  // double lbDMuEtaMaxTight = -1.0;
838 
839  lbDBasicSelect = new BPHFittedBasicSelect(lbDMassMin, lbDMassMax, lbDPtMin, lbDEtaMax, lbDYMax);
840  lbDJPsiBasicSelect =
841  new BPHCompositeBasicSelect(lbDJPsiMassMin, lbDJPsiMassMax, lbDJPsiPtMin, lbDJPsiEtaMax, lbDJPsiYMax);
842  lbDLambda0BasicSelect =
843  new BPHFittedBasicSelect(lbDLambda0MassMin, lbDLambda0MassMax, lbDLambda0PtMin, lbDLambda0EtaMax, lbDLambda0YMax);
844  lbDVertexSelect = new BPHGenericVertexSelect('f', lbDProbMin, lbDCosMin, lbDSigMin);
845  lbDJPsiDaughterSelect = nullptr;
846  // lbDJPsiDaughterSelect = new BPHDaughterSelect(
847  // lbDMuPtMinLoose , lbDMuPtMinTight ,
848  // lbDMuEtaMaxLoose, lbDMuEtaMaxTight, sms );
849 
851 
852  double bcIMassMin = 0.0;
853  double bcIMassMax = 999999.0;
854  double bcIPtMin = 27.0;
855  double bcIEtaMax = -1.0;
856  double bcIYMax = -1.0;
857  double bcIJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
858  double bcIJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
859  double bcIJPsiPtMin = 25.0;
860  double bcIJPsiEtaMax = -1.0;
861  double bcIJPsiYMax = -1.0;
862  double bcIJPsiProbMin = 0.005;
863  double bcIProbMin = 0.10;
864  double bcICosMin = -2.0;
865  double bcISigMin = -1.0;
866  double bcIDistMin = 0.01;
867  // double bcIMuPtMinLoose = -1.0;
868  // double bcIMuPtMinTight = -1.0;
869  // double bcIMuEtaMaxLoose = -1.0;
870  // double bcIMuEtaMaxTight = -1.0;
871 
872  bcIPiPtMin = 3.5;
873 
874  bcIBasicSelect = new BPHFittedBasicSelect(bcIMassMin, bcIMassMax, bcIPtMin, bcIEtaMax, bcIYMax);
875  bcIJPsiBasicSelect =
876  new BPHCompositeBasicSelect(bcIJPsiMassMin, bcIJPsiMassMax, bcIJPsiPtMin, bcIJPsiEtaMax, bcIJPsiYMax);
877  bcIJPsiVertexSelect = new BPHGenericVertexSelect('c', bcIJPsiProbMin);
878  bcIVertexSelect = new BPHGenericVertexSelect('f', bcIProbMin, bcICosMin, bcISigMin, bcIDistMin);
879  bcIJPsiDaughterSelect = nullptr;
880  // bcIJPsiDaughterSelect = new BPHDaughterSelect(
881  // bcIMuPtMinLoose , bcIMuPtMinTight ,
882  // bcIMuEtaMaxLoose, bcMuEtaMaxTight, sms );
883 
884  double bcDMassMin = 0.0;
885  double bcDMassMax = 999999.0;
886  double bcDPtMin = 8.0;
887  double bcDEtaMax = -1.0;
888  double bcDYMax = -1.0;
889  double bcDJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
890  double bcDJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
891  double bcDJPsiPtMin = 7.0;
892  double bcDJPsiEtaMax = -1.0;
893  double bcDJPsiYMax = -1.0;
894  double bcDJPsiProbMin = 0.005;
895  double bcDProbMin = 0.10;
896  double bcDCosMin = 0.99;
897  double bcDSigMin = 3.0;
898  // double bcDMuPtMinLoose = -1.0;
899  // double bcDMuPtMinTight = -1.0;
900  // double bcDMuEtaMaxLoose = -1.0;
901  // double bcDMuEtaMaxTight = -1.0;
902 
903  bcJPsiDcaMax = 0.5;
904  bcDPiPtMin = 3.5;
905 
906  bcDBasicSelect = new BPHFittedBasicSelect(bcDMassMin, bcDMassMax, bcDPtMin, bcDEtaMax, bcDYMax);
907  bcDJPsiBasicSelect =
908  new BPHCompositeBasicSelect(bcDJPsiMassMin, bcDJPsiMassMax, bcDJPsiPtMin, bcDJPsiEtaMax, bcDJPsiYMax);
909  bcDJPsiVertexSelect = new BPHGenericVertexSelect('c', bcDJPsiProbMin);
910  bcDVertexSelect = new BPHGenericVertexSelect('f', bcDProbMin, bcDCosMin, bcDSigMin);
911  bcDJPsiDaughterSelect = nullptr;
912  // bcDJPsiDaughterSelect = new BPHDaughterSelect(
913  // bcDMuPtMinLoose , bcDMuPtMinTight ,
914  // bcDMuEtaMaxLoose, bcDMuEtaMaxTight, sms );
915 
917 
918  double x3872IMassMin = 0.0;
919  double x3872IMassMax = 999999.0;
920  double x3872IPtMin = 27.0;
921  double x3872IEtaMax = -1.0;
922  double x3872IYMax = -1.0;
923  double x3872IJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
924  double x3872IJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
925  double x3872IJPsiPtMin = 25.0;
926  double x3872IJPsiEtaMax = -1.0;
927  double x3872IJPsiYMax = -1.0;
928  double x3872IJPsiProbMin = 0.10;
929  double x3872IProbMin = 0.10;
930  double x3872ICosMin = -2.0;
931  double x3872ISigMin = -1.0;
932  double x3872IDistMin = 0.01;
933  // double x3872IMuPtMinLoose = -1.0;
934  // double x3872IMuPtMinTight = -1.0;
935  // double x3872IMuEtaMaxLoose = -1.0;
936  // double x3872IMuEtaMaxTight = -1.0;
937 
938  x3872JPsiDcaMax = 0.5;
939  x3872IPiPtMin = 1.2;
940 
941  x3872IBasicSelect = new BPHFittedBasicSelect(x3872IMassMin, x3872IMassMax, x3872IPtMin, x3872IEtaMax, x3872IYMax);
942  x3872IJPsiBasicSelect = new BPHCompositeBasicSelect(
943  x3872IJPsiMassMin, x3872IJPsiMassMax, x3872IJPsiPtMin, x3872IJPsiEtaMax, x3872IJPsiYMax);
944  x3872IJPsiVertexSelect = new BPHGenericVertexSelect('c', x3872IJPsiProbMin);
945  x3872IVertexSelect = new BPHGenericVertexSelect('f', x3872IProbMin, x3872ICosMin, x3872ISigMin, x3872IDistMin);
946  x3872IJPsiDaughterSelect = nullptr;
947  // x3872IJPsiDaughterSelect = new BPHDaughterSelect(
948  // x3872IMuPtMinLoose , x3872IMuPtMinTight,
949  // x3872IMuEtaMaxLoose, x3872MuEtaMaxTight,
950  // sms );
951 
952  double x3872DMassMin = 0.0;
953  double x3872DMassMax = 999999.0;
954  double x3872DPtMin = 8.0;
955  double x3872DEtaMax = -1.0;
956  double x3872DYMax = -1.0;
957  double x3872DJPsiMassMin = BPHParticleMasses::jPsiMass - 0.150;
958  double x3872DJPsiMassMax = BPHParticleMasses::jPsiMass + 0.150;
959  double x3872DJPsiPtMin = 7.0;
960  double x3872DJPsiEtaMax = -1.0;
961  double x3872DJPsiYMax = -1.0;
962  double x3872DJPsiProbMin = 0.10;
963  double x3872DProbMin = 0.10;
964  double x3872DCosMin = 0.99;
965  double x3872DSigMin = 3.0;
966  // double x3872DMuPtMinLoose = -1.0;
967  // double x3872DMuPtMinTight = -1.0;
968  // double x3872DMuEtaMaxLoose = -1.0;
969  // double x3872DMuEtaMaxTight = -1.0;
970 
971  x3872DPiPtMin = 1.2;
972 
973  x3872DBasicSelect = new BPHFittedBasicSelect(x3872DMassMin, x3872DMassMax, x3872DPtMin, x3872DEtaMax, x3872DYMax);
974  x3872DJPsiBasicSelect = new BPHCompositeBasicSelect(
975  x3872DJPsiMassMin, x3872DJPsiMassMax, x3872DJPsiPtMin, x3872DJPsiEtaMax, x3872DJPsiYMax);
976  x3872DJPsiVertexSelect = new BPHGenericVertexSelect('c', x3872DJPsiProbMin);
977  x3872DVertexSelect = new BPHGenericVertexSelect('f', x3872DProbMin, x3872DCosMin, x3872DSigMin);
978  x3872DJPsiDaughterSelect = nullptr;
979  // x3872DJPsiDaughterSelect = new BPHDaughterSelect(
980  // x3872DMuPtMinLoose , x3872DMuPtMinTight ,
981  // x3872DMuEtaMaxLoose, x3872DMuEtaMaxTight,,
982  // sms );
983 }
984 
986  delete phiIBasicSelect;
987  delete phiBBasicSelect;
988  delete jPsiIBasicSelect;
989  delete jPsiBBasicSelect;
990  delete psi2IBasicSelect;
991  delete psi2BBasicSelect;
992  delete upsIBasicSelect;
993  delete upsBBasicSelect;
994  delete oniaVertexSelect;
995  delete oniaDaughterSelect;
996 
997  delete npJPsiBasicSelect;
998  delete npJPsiDaughterSelect;
999 
1000  delete buIBasicSelect;
1001  delete buIJPsiBasicSelect;
1002  delete buIVertexSelect;
1003  delete buIJPsiDaughterSelect;
1004  delete buDBasicSelect;
1005  delete buDJPsiBasicSelect;
1006  delete buDVertexSelect;
1007  delete buDJPsiDaughterSelect;
1008 
1009  delete bdIBasicSelect;
1010  delete bdIJPsiBasicSelect;
1011  delete bdIKx0BasicSelect;
1012  delete bdIVertexSelect;
1013  delete bdIJPsiDaughterSelect;
1014  delete bdDBasicSelect;
1015  delete bdDJPsiBasicSelect;
1016  delete bdDKx0BasicSelect;
1017  delete bdDVertexSelect;
1018  delete bdDJPsiDaughterSelect;
1019 
1020  delete bsIBasicSelect;
1021  delete bsIJPsiBasicSelect;
1022  delete bsIPhiBasicSelect;
1023  delete bsIVertexSelect;
1024  delete bsIJPsiDaughterSelect;
1025  delete bsDBasicSelect;
1026  delete bsDJPsiBasicSelect;
1027  delete bsDPhiBasicSelect;
1028  delete bsDVertexSelect;
1029  delete bsDJPsiDaughterSelect;
1030 
1031  delete b0IBasicSelect;
1032  delete b0IJPsiBasicSelect;
1033  delete b0IK0sBasicSelect;
1034  delete b0IVertexSelect;
1035  delete b0IJPsiDaughterSelect;
1036  delete b0DBasicSelect;
1037  delete b0DJPsiBasicSelect;
1038  delete b0DK0sBasicSelect;
1039  delete b0DVertexSelect;
1040  delete b0DJPsiDaughterSelect;
1041 
1042  delete lbIBasicSelect;
1043  delete lbIJPsiBasicSelect;
1044  delete lbILambda0BasicSelect;
1045  delete lbIVertexSelect;
1046  delete lbIJPsiDaughterSelect;
1047  delete lbDBasicSelect;
1048  delete lbDJPsiBasicSelect;
1049  delete lbDLambda0BasicSelect;
1050  delete lbDVertexSelect;
1051  delete lbDJPsiDaughterSelect;
1052 
1053  delete bcIBasicSelect;
1054  delete bcIJPsiBasicSelect;
1055  delete bcIJPsiVertexSelect;
1056  delete bcIVertexSelect;
1057  delete bcIJPsiDaughterSelect;
1058  delete bcDBasicSelect;
1059  delete bcDJPsiBasicSelect;
1060  delete bcDJPsiVertexSelect;
1061  delete bcDVertexSelect;
1062  delete bcDJPsiDaughterSelect;
1063 
1064  delete x3872IBasicSelect;
1065  delete x3872IJPsiBasicSelect;
1066  delete x3872IJPsiVertexSelect;
1067  delete x3872IVertexSelect;
1068  delete x3872IJPsiDaughterSelect;
1069  delete x3872DBasicSelect;
1070  delete x3872DJPsiBasicSelect;
1071  delete x3872DJPsiVertexSelect;
1072  delete x3872DVertexSelect;
1073  delete x3872DJPsiDaughterSelect;
1074 }
1075 
1078  desc.add<string>("trigResultsLabel", "");
1079  desc.add<string>("oniaCandsLabel", "");
1080  desc.add<string>("sdCandsLabel", "");
1081  desc.add<string>("ssCandsLabel", "");
1082  desc.add<string>("buCandsLabel", "");
1083  desc.add<string>("bdCandsLabel", "");
1084  desc.add<string>("bsCandsLabel", "");
1085  desc.add<string>("k0CandsLabel", "");
1086  desc.add<string>("l0CandsLabel", "");
1087  desc.add<string>("b0CandsLabel", "");
1088  desc.add<string>("lbCandsLabel", "");
1089  desc.add<string>("bcCandsLabel", "");
1090  desc.add<string>("x3872CandsLabel", "");
1091  descriptions.add("process.bphHistoSpecificDecay", desc);
1092  return;
1093 }
1094 
1096  createHisto("massDIPhi", 50, 0.90, 1.15); // Phi mass inclusive
1097  createHisto("massTIPhi", 50, 0.90, 1.15); // Phi mass inclusive
1098  createHisto("massDBPhi", 50, 0.90, 1.15); // Phi mass barrel
1099  createHisto("massTBPhi", 50, 0.90, 1.15); // Phi mass barrel
1100  createHisto("massDIJPsi", 35, 2.95, 3.30); // JPsi mass inclusive
1101  createHisto("massTIJPsi", 35, 2.95, 3.30); // JPsi mass inclusive
1102  createHisto("massDBJPsi", 35, 2.95, 3.30); // JPsi mass barrel
1103  createHisto("massTBJPsi", 35, 2.95, 3.30); // JPsi mass barrel
1104  createHisto("massDIPsi2", 60, 3.40, 4.00); // Psi2 mass inclusive
1105  createHisto("massTIPsi2", 60, 3.40, 4.00); // Psi2 mass inclusive
1106  createHisto("massDBPsi2", 60, 3.40, 4.00); // Psi2 mass barrel
1107  createHisto("massTBPsi2", 60, 3.40, 4.00); // Psi2 mass barrel
1108  createHisto("massDIUps123", 115, 8.70, 11.0); // Ups mass inclusive
1109  createHisto("massTIUps123", 115, 8.70, 11.0); // Ups mass inclusive
1110  createHisto("massDBUps123", 115, 8.70, 11.0); // Ups mass barrel
1111  createHisto("massTBUps123", 115, 8.70, 11.0); // Ups mass barrel
1112  createHisto("massDIBu", 100, 5.00, 6.00); // Bu mass inclusive
1113  createHisto("massTIBu", 100, 5.00, 6.00); // Bu mass inclusive
1114  createHisto("massDDBu", 100, 5.00, 6.00); // Bu mass displaced
1115  createHisto("massTDBu", 100, 5.00, 6.00); // Bu mass displaced
1116  createHisto("massDIBd", 100, 5.00, 6.00); // Bd mass inclusive
1117  createHisto("massTIBd", 100, 5.00, 6.00); // Bd mass inclusive
1118  createHisto("massDDBd", 100, 5.00, 6.00); // Bd mass displaced
1119  createHisto("massTDBd", 100, 5.00, 6.00); // Bd mass displaced
1120  createHisto("massDIBs", 100, 5.00, 6.00); // Bs mass inclusive
1121  createHisto("massTIBs", 100, 5.00, 6.00); // Bs mass inclusive
1122  createHisto("massDDBs", 100, 5.00, 6.00); // Bs mass displaced
1123  createHisto("massTDBs", 100, 5.00, 6.00); // Bs mass displaced
1124  createHisto("massDIBc", 100, 6.00, 7.00); // Bc mass inclusive
1125  createHisto("massTIBc", 100, 6.00, 7.00); // Bc mass inclusive
1126  createHisto("massDDBc", 100, 6.00, 7.00); // Bc mass displaced
1127  createHisto("massTDBc", 100, 6.00, 7.00); // Bc mass displaced
1128  createHisto("massDIX3872", 40, 3.60, 4.00); // X3872 mass inclusive
1129  createHisto("massTIX3872", 40, 3.60, 4.00); // X3872 mass inclusive
1130  createHisto("massDDX3872", 40, 3.60, 4.00); // X3872 mass displaced
1131  createHisto("massTDX3872", 40, 3.60, 4.00); // X3872 mass displaced
1132  createHisto("mfitDIBu", 100, 5.00, 6.00); // Bu mass, with constr.
1133  createHisto("mfitTIBu", 100, 5.00, 6.00); // Bu mass, with constr.
1134  createHisto("mfitDDBu", 100, 5.00, 6.00); // Bu mass, with constr.
1135  createHisto("mfitTDBu", 100, 5.00, 6.00); // Bu mass, with constr.
1136  createHisto("mfitDIBd", 100, 5.00, 6.00); // Bd mass, with constr.
1137  createHisto("mfitTIBd", 100, 5.00, 6.00); // Bd mass, with constr.
1138  createHisto("mfitDDBd", 100, 5.00, 6.00); // Bd mass, with constr.
1139  createHisto("mfitTDBd", 100, 5.00, 6.00); // Bd mass, with constr.
1140  createHisto("mfitDIBs", 100, 5.00, 6.00); // Bs mass, with constr.
1141  createHisto("mfitTIBs", 100, 5.00, 6.00); // Bs mass, with constr.
1142  createHisto("mfitDDBs", 100, 5.00, 6.00); // Bs mass, with constr.
1143  createHisto("mfitTDBs", 100, 5.00, 6.00); // Bs mass, with constr.
1144  createHisto("mfitDIBc", 100, 6.00, 7.00); // Bc mass, with constr.
1145  createHisto("mfitTIBc", 100, 6.00, 7.00); // Bc mass, with constr.
1146  createHisto("mfitDDBc", 100, 6.00, 7.00); // Bc mass, with constr.
1147  createHisto("mfitTDBc", 100, 6.00, 7.00); // Bc mass, with constr.
1148  createHisto("mfitDIX3872", 40, 3.60, 4.00); // X3872 mass, with constr.
1149  createHisto("mfitTIX3872", 40, 3.60, 4.00); // X3872 mass, with constr.
1150  createHisto("mfitDDX3872", 40, 3.60, 4.00); // X3872 mass, with constr.
1151  createHisto("mfitTDX3872", 40, 3.60, 4.00); // X3872 mass, with constr.
1152  createHisto("massDIBuJPsi", 35, 2.95, 3.30); // JPsi mass in Bu decay
1153  createHisto("massTIBuJPsi", 35, 2.95, 3.30); // JPsi mass in Bu decay
1154  createHisto("massDDBuJPsi", 35, 2.95, 3.30); // JPsi mass in Bu decay
1155  createHisto("massTDBuJPsi", 35, 2.95, 3.30); // JPsi mass in Bu decay
1156  createHisto("massDIBdJPsi", 35, 2.95, 3.30); // JPsi mass in Bd decay
1157  createHisto("massTIBdJPsi", 35, 2.95, 3.30); // JPsi mass in Bd decay
1158  createHisto("massDDBdJPsi", 35, 2.95, 3.30); // JPsi mass in Bd decay
1159  createHisto("massTDBdJPsi", 35, 2.95, 3.30); // JPsi mass in Bd decay
1160  createHisto("massDIBdKx0", 50, 0.80, 1.05); // Kx0 mass in Bd decay
1161  createHisto("massTIBdKx0", 50, 0.80, 1.05); // Kx0 mass in Bd decay
1162  createHisto("massDDBdKx0", 50, 0.80, 1.05); // Kx0 mass in Bd decay
1163  createHisto("massTDBdKx0", 50, 0.80, 1.05); // Kx0 mass in Bd decay
1164  createHisto("massDIBsJPsi", 35, 2.95, 3.30); // JPsi mass in Bs decay
1165  createHisto("massTIBsJPsi", 35, 2.95, 3.30); // JPsi mass in Bs decay
1166  createHisto("massDDBsJPsi", 35, 2.95, 3.30); // JPsi mass in Bs decay
1167  createHisto("massTDBsJPsi", 35, 2.95, 3.30); // JPsi mass in Bs decay
1168  createHisto("massDIBsPhi", 50, 1.01, 1.03); // Phi mass in Bs decay
1169  createHisto("massTIBsPhi", 50, 1.01, 1.03); // Phi mass in Bs decay
1170  createHisto("massDDBsPhi", 50, 1.01, 1.03); // Phi mass in Bs decay
1171  createHisto("massTDBsPhi", 50, 1.01, 1.03); // Phi mass in Bs decay
1172  createHisto("massDIBcJPsi", 35, 2.95, 3.30); // JPsi mass in Bc decay
1173  createHisto("massTIBcJPsi", 35, 2.95, 3.30); // JPsi mass in Bc decay
1174  createHisto("massDDBcJPsi", 35, 2.95, 3.30); // JPsi mass in Bc decay
1175  createHisto("massTDBcJPsi", 35, 2.95, 3.30); // JPsi mass in Bc decay
1176  createHisto("massDIX3JPsi", 35, 2.95, 3.30); // JPsi mass in X3872 decay
1177  createHisto("massTIX3JPsi", 35, 2.95, 3.30); // JPsi mass in X3872 decay
1178  createHisto("massDDX3JPsi", 35, 2.95, 3.30); // JPsi mass in X3872 decay
1179  createHisto("massTDX3JPsi", 35, 2.95, 3.30); // JPsi mass in X3872 decay
1180  createHisto("massDK0s", 50, 0.40, 0.60); // K0s mass
1181  createHisto("mfitDK0s", 50, 0.40, 0.60); // K0s mass
1182  createHisto("massDLambda0", 60, 1.00, 1.30); // Lambda0 mass
1183  createHisto("mfitDLambda0", 60, 1.00, 1.30); // Lambda0 mass
1184  createHisto("massDIB0", 50, 5.00, 6.00); // B0 mass inclusive
1185  createHisto("massTIB0", 50, 5.00, 6.00); // B0 mass inclusive
1186  createHisto("massDDB0", 50, 5.00, 6.00); // B0 mass displaced
1187  createHisto("massTDB0", 50, 5.00, 6.00); // B0 mass displaced
1188  createHisto("mfitDIB0", 50, 5.00, 6.00); // B0 mass, with constr.
1189  createHisto("mfitTIB0", 50, 5.00, 6.00); // B0 mass, with constr.
1190  createHisto("mfitDDB0", 50, 5.00, 6.00); // B0 mass, with constr.
1191  createHisto("mfitTDB0", 50, 5.00, 6.00); // B0 mass, with constr.
1192  createHisto("massDIB0JPsi", 35, 2.95, 3.30); // JPsi mass in B0 decay
1193  createHisto("massTIB0JPsi", 35, 2.95, 3.30); // JPsi mass in B0 decay
1194  createHisto("massDDB0JPsi", 35, 2.95, 3.30); // JPsi mass in B0 decay
1195  createHisto("massTDB0JPsi", 35, 2.95, 3.30); // JPsi mass in B0 decay
1196  createHisto("massDIB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1197  createHisto("massTIB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1198  createHisto("massDDB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1199  createHisto("massTDB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1200  createHisto("mfitDIB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1201  createHisto("mfitTIB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1202  createHisto("mfitDDB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1203  createHisto("mfitTDB0K0s", 50, 0.40, 0.60); // K0s mass in B0 decay
1204  createHisto("massDILambdab", 25, 5.00, 6.00); // Lambdab mass inclusive
1205  createHisto("massTILambdab", 25, 5.00, 6.00); // Lambdab mass inclusive
1206  createHisto("massDDLambdab", 25, 5.00, 6.00); // Lambdab mass displaced
1207  createHisto("massTDLambdab", 25, 5.00, 6.00); // Lambdab mass displaced
1208  createHisto("mfitDILambdab", 25, 5.00, 6.00); // Lambdab mass, with constr.
1209  createHisto("mfitTILambdab", 25, 5.00, 6.00); // Lambdab mass, with constr.
1210  createHisto("mfitDDLambdab", 25, 5.00, 6.00); // Lambdab mass, with constr.
1211  createHisto("mfitTDLambdab", 25, 5.00, 6.00); // Lambdab mass, with constr.
1212  createHisto("massDILbJPsi", 35, 2.95, 3.30); // JPsi mass in Lambdab decay
1213  createHisto("massTILbJPsi", 35, 2.95, 3.30); // JPsi mass in Lambdab decay
1214  createHisto("massDDLbJPsi", 35, 2.95, 3.30); // JPsi mass in Lambdab decay
1215  createHisto("massTDLbJPsi", 35, 2.95, 3.30); // JPsi mass in Lambdab decay
1216  createHisto("massDILbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1217  createHisto("massTILbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1218  createHisto("massDDLbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1219  createHisto("massTDLbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1220  createHisto("mfitDILbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1221  createHisto("mfitTILbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1222  createHisto("mfitDDLbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1223  createHisto("mfitTDLbL0", 60, 1.00, 1.30); // L0 mass in Lambdab decay
1224 
1225  createHisto("massFull", 200, 2.00, 12.0); // Full onia mass
1226 
1227  createHisto("ctauDIJPsi", 60, -0.05, 0.25); // JPsi ctau inclusive
1228  createHisto("ctauTIJPsi", 60, -0.05, 0.25); // JPsi ctau inclusive
1229  createHisto("ctauDBJPsi", 60, -0.05, 0.25); // JPsi ctau barrel
1230  createHisto("ctauTBJPsi", 60, -0.05, 0.25); // JPsi ctau barrel
1231  createHisto("ctauDIBu", 60, -0.05, 0.25); // Bu ctau inclusive
1232  createHisto("ctauTIBu", 60, -0.05, 0.25); // Bu ctau inclusive
1233  createHisto("ctauDDBu", 60, -0.05, 0.25); // Bu ctau displaced
1234  createHisto("ctauTDBu", 60, -0.05, 0.25); // Bu ctau displaced
1235  createHisto("ctauDIBd", 60, -0.05, 0.25); // Bd ctau inclusive
1236  createHisto("ctauTIBd", 60, -0.05, 0.25); // Bd ctau inclusive
1237  createHisto("ctauDDBd", 60, -0.05, 0.25); // Bd ctau displaced
1238  createHisto("ctauTDBd", 60, -0.05, 0.25); // Bd ctau displaced
1239  createHisto("ctauDIBs", 60, -0.05, 0.25); // Bs ctau inclusive
1240  createHisto("ctauTIBs", 60, -0.05, 0.25); // Bs ctau inclusive
1241  createHisto("ctauDDBs", 60, -0.05, 0.25); // Bs ctau displaced
1242  createHisto("ctauTDBs", 60, -0.05, 0.25); // Bs ctau displaced
1243  createHisto("ctauDIB0", 60, -0.05, 0.25); // B0 ctau inclusive
1244  createHisto("ctauTIB0", 60, -0.05, 0.25); // B0 ctau inclusive
1245  createHisto("ctauDDB0", 60, -0.05, 0.25); // B0 ctau displaced
1246  createHisto("ctauTDB0", 60, -0.05, 0.25); // B0 ctau displaced
1247  createHisto("ctauDILambdab", 60, -0.05, 0.25); // Lambdab ctau inclusive
1248  createHisto("ctauTILambdab", 60, -0.05, 0.25); // Lambdab ctau inclusive
1249  createHisto("ctauDDLambdab", 60, -0.05, 0.25); // Lambdab ctau displaced
1250  createHisto("ctauTDLambdab", 60, -0.05, 0.25); // Lambdab ctau displaced
1251 
1252  recoName = new string;
1253  tree = fs->make<TTree>("BPHReco", "BPHReco");
1254  b_runNumber = tree->Branch("runNumber", &runNumber, "runNumber/i");
1255  b_lumiSection = tree->Branch("lumiSection", &lumiSection, "lumiSection/i");
1256  b_eventNumber = tree->Branch("eventNumber", &eventNumber, "eventNumber/i");
1257  b_recoName = tree->Branch("recoName", &recoName, 8192, 99);
1258  b_recoMass = tree->Branch("recoMass", &recoMass, "recoMass/F");
1259  b_recoTime = tree->Branch("recoTime", &recoTime, "recoTime/F");
1260  b_recoErrT = tree->Branch("recoErrT", &recoErrT, "recoErrT/F");
1261 
1262  return;
1263 }
1264 
1266  // if ( ev.id().run () != 316239 ) return;
1267  // if ( ev.id().event() != 170736782 ) return;
1268  static map<string, ofstream*> ofMap;
1269  if (ofMap.empty()) {
1270  ofMap["BarPhi"] = nullptr;
1271  ofMap["IncJPsi"] = nullptr;
1272  ofMap["BarJPsi"] = nullptr;
1273  ofMap["IncPsi2"] = nullptr;
1274  ofMap["BarPsi2"] = nullptr;
1275  ofMap["BarUpsilon123"] = nullptr;
1276  ofMap["InclusiveBu"] = nullptr;
1277  ofMap["DisplacedBu"] = nullptr;
1278  ofMap["InclusiveBd"] = nullptr;
1279  ofMap["DisplacedBd"] = nullptr;
1280  ofMap["InclusiveBs"] = nullptr;
1281  ofMap["DisplacedBs"] = nullptr;
1282  ofMap["K0s"] = nullptr;
1283  ofMap["Lambda0"] = nullptr;
1284  ofMap["InclusiveB0"] = nullptr;
1285  ofMap["DisplacedB0"] = nullptr;
1286  ofMap["InclusiveLambdab"] = nullptr;
1287  ofMap["DisplacedLambdab"] = nullptr;
1288  ofMap["InclusiveBc"] = nullptr;
1289  ofMap["DisplacedBc"] = nullptr;
1290  ofMap["InclusiveX3872"] = nullptr;
1291  ofMap["DisplacedX3872"] = nullptr;
1292  map<string, ofstream*>::iterator iter = ofMap.begin();
1293  map<string, ofstream*>::iterator iend = ofMap.end();
1294  string name = "list";
1295  while (iter != iend) {
1296  iter->second = new ofstream(name + iter->first);
1297  ++iter;
1298  }
1299  }
1300 
1301  // event number
1302  runNumber = ev.id().run();
1303  lumiSection = ev.id().luminosityBlock();
1304  eventNumber = ev.id().event();
1305 
1306  // get magnetic field
1308  es.get<IdealMagneticFieldRecord>().get(magneticField);
1309 
1310  // get object collections
1311  // collections are got through "BPHTokenWrapper" interface to allow
1312  // uniform access in different CMSSW versions
1313 
1315 
1317  const edm::TriggerNames* trigNames = nullptr;
1318  if (useTrig) {
1319  trigResultsToken.get(ev, trigResults);
1320  if (trigResults.isValid())
1321  trigNames = &ev.triggerNames(*trigResults);
1322  }
1323 
1324  bool flag_Dimuon25_Jpsi = false;
1325  bool flag_Dimuon20_Jpsi_Barrel_Seagulls = false;
1326  bool flag_Dimuon14_Phi_Barrel_Seagulls = false;
1327  bool flag_Dimuon18_PsiPrime = false;
1328  bool flag_Dimuon10_PsiPrime_Barrel_Seagulls = false;
1329  bool flag_Dimuon12_Upsilon_eta1p5 = false;
1330  bool flag_Dimuon12_Upsilon_y1p4 = false;
1331  bool flag_DoubleMu4_JpsiTrk_Displaced = false;
1332  if (trigNames != nullptr) {
1333  const edm::TriggerNames::Strings& names = trigNames->triggerNames();
1334  int iObj;
1335  int nObj = names.size();
1336  for (iObj = 0; iObj < nObj; ++iObj) {
1337  // cout << names[iObj] << endl;
1338  CHK_TRIG(trigResults, names, iObj, Dimuon25_Jpsi)
1339  CHK_TRIG(trigResults, names, iObj, Dimuon20_Jpsi_Barrel_Seagulls)
1340  CHK_TRIG(trigResults, names, iObj, Dimuon14_Phi_Barrel_Seagulls)
1341  CHK_TRIG(trigResults, names, iObj, Dimuon18_PsiPrime)
1342  CHK_TRIG(trigResults, names, iObj, Dimuon10_PsiPrime_Barrel_Seagulls)
1343  CHK_TRIG(trigResults, names, iObj, Dimuon12_Upsilon_eta1p5)
1344  CHK_TRIG(trigResults, names, iObj, Dimuon12_Upsilon_y1p4)
1345  CHK_TRIG(trigResults, names, iObj, DoubleMu4_JpsiTrk_Displaced)
1346  }
1347  }
1348 
1349  // cout << "Dimuon25_Jpsi "
1350  // << ( flag_Dimuon25_Jpsi ? 'A' : 'R' ) << endl;
1351  // cout << "Dimuon20_Jpsi_Barrel_Seagulls "
1352  // << ( flag_Dimuon20_Jpsi_Barrel_Seagulls ? 'A' : 'R' ) << endl;
1353  // cout << "Dimuon14_Phi_Barrel_Seagulls "
1354  // << ( flag_Dimuon14_Phi_Barrel_Seagulls ? 'A' : 'R' ) << endl;
1355  // cout << "Dimuon18_PsiPrime "
1356  // << ( flag_Dimuon18_PsiPrime ? 'A' : 'R' ) << endl;
1357  // cout << "Dimuon10_PsiPrime_Barrel_Seagulls "
1358  // << ( flag_Dimuon10_PsiPrime_Barrel_Seagulls ? 'A' : 'R' ) << endl;
1359  // cout << "Dimuon12_Upsilon_eta1p5 "
1360  // << ( flag_Dimuon12_Upsilon_eta1p5 ? 'A' : 'R' ) << endl;
1361  // cout << "DoubleMu4_JpsiTrk_Displaced "
1362  // << ( flag_DoubleMu4_JpsiTrk_Displaced ? 'A' : 'R' ) << endl;
1363 
1365 
1367  int iqo;
1368  int nqo = 0;
1369  if (useOnia) {
1370  oniaCandsToken.get(ev, oniaCands);
1371  nqo = oniaCands->size();
1372  }
1373 
1374  for (iqo = 0; iqo < nqo; ++iqo) {
1375  LogTrace("DataDump") << "*********** quarkonium " << iqo << "/" << nqo << " ***********";
1376  const pat::CompositeCandidate& cand = oniaCands->at(iqo);
1377  if (!oniaVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(cand, "primaryVertex")))
1378  continue;
1379  if (!oniaDaughterSelect->accept(cand))
1380  continue;
1381  fillHisto("Full", cand, 'c');
1382  if (phiBBasicSelect->accept(cand)) {
1383  fillHisto("DBPhi", cand, 'c');
1384  if (flag_Dimuon14_Phi_Barrel_Seagulls)
1385  fillHisto("TBPhi", cand, 'c');
1386  if (flag_Dimuon14_Phi_Barrel_Seagulls)
1387  *ofMap["BarPhi"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1388  << cand.mass() << endl;
1389  }
1390  if (jPsiIBasicSelect->accept(cand)) {
1391  fillHisto("DIJPsi", cand, 'c');
1392  if (flag_Dimuon25_Jpsi)
1393  fillHisto("TIJPsi", cand, 'c');
1394  if (flag_Dimuon25_Jpsi)
1395  *ofMap["IncJPsi"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1396  << cand.mass() << endl;
1397  }
1398  if (jPsiBBasicSelect->accept(cand)) {
1399  fillHisto("DBJPsi", cand, 'c');
1400  if (flag_Dimuon20_Jpsi_Barrel_Seagulls)
1401  fillHisto("TBJPsi", cand, 'c');
1402  if (flag_Dimuon20_Jpsi_Barrel_Seagulls)
1403  *ofMap["BarJPsi"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1404  << cand.mass() << endl;
1405  }
1406  if (psi2IBasicSelect->accept(cand)) {
1407  fillHisto("DIPsi2", cand, 'c');
1408  if (flag_Dimuon18_PsiPrime)
1409  fillHisto("TIPsi2", cand, 'c');
1410  if (flag_Dimuon18_PsiPrime)
1411  *ofMap["IncPsi2"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1412  << cand.mass() << endl;
1413  }
1414  if (psi2BBasicSelect->accept(cand)) {
1415  fillHisto("DBPsi2", cand, 'c');
1416  if (flag_Dimuon10_PsiPrime_Barrel_Seagulls)
1417  fillHisto("TBPsi2", cand, 'c');
1418  if (flag_Dimuon10_PsiPrime_Barrel_Seagulls)
1419  *ofMap["BarPsi2"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1420  << cand.mass() << endl;
1421  }
1422  if (upsBBasicSelect->accept(cand)) {
1423  fillHisto("DBUps123", cand, 'c');
1424  if (flag_Dimuon12_Upsilon_eta1p5 || flag_Dimuon12_Upsilon_y1p4)
1425  fillHisto("TBUps123", cand, 'c');
1426  if (flag_Dimuon12_Upsilon_eta1p5 || flag_Dimuon12_Upsilon_y1p4)
1427  *ofMap["BarUpsilon123"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1428  << cand.mass() << endl;
1429  }
1430  }
1431 
1433 
1435  int ibu;
1436  int nbu = 0;
1437  if (useBu) {
1438  buCandsToken.get(ev, buCands);
1439  nbu = buCands->size();
1440  }
1441 
1442  for (ibu = 0; ibu < nbu; ++ibu) {
1443  LogTrace("DataDump") << "*********** Bu " << ibu << "/" << nbu << " ***********";
1444  const pat::CompositeCandidate& cand = buCands->at(ibu);
1445  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1446  LogTrace("DataDump") << "JPsi: " << jPsi;
1447  if (jPsi == nullptr)
1448  continue;
1449  if (!npJPsiBasicSelect->accept(*jPsi))
1450  continue;
1451  if (!npJPsiDaughterSelect->accept(*jPsi))
1452  continue;
1453  const reco::Candidate* kptr = BPHDaughters::get(cand, 0.49, 0.50).front();
1454  if (kptr == nullptr)
1455  continue;
1456  if (buIBasicSelect->accept(cand) && buIJPsiBasicSelect->accept(*jPsi) &&
1457  // buIJPsiDaughterSelect->accept( *jPsi ) &&
1458  buIVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1459  (kptr->pt() > buIKPtMin)) {
1460  fillHisto("DIBu", cand, 'f');
1461  fillHisto("DIBuJPsi", *jPsi, 'c');
1462  if (flag_Dimuon25_Jpsi) {
1463  fillHisto("TIBu", cand, 'f');
1464  fillHisto("TIBuJPsi", *jPsi, 'c');
1465  *ofMap["InclusiveBu"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1466  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1467  }
1468  }
1469  if (buDBasicSelect->accept(cand) && buDJPsiBasicSelect->accept(*jPsi) &&
1470  // buDJPsiDaughterSelect->accept( *jPsi ) &&
1471  buDVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1472  (kptr->pt() > buDKPtMin)) {
1473  fillHisto("DDBu", cand, 'f');
1474  fillHisto("DDBuJPsi", *jPsi, 'c');
1475  if (flag_DoubleMu4_JpsiTrk_Displaced) {
1476  fillHisto("TDBu", cand, 'f');
1477  fillHisto("TDBuJPsi", *jPsi, 'c');
1478  *ofMap["DisplacedBu"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1479  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1480  }
1481  }
1482  }
1483 
1485 
1487  int ibd;
1488  int nbd = 0;
1489  if (useBd) {
1490  bdCandsToken.get(ev, bdCands);
1491  nbd = bdCands->size();
1492  }
1493 
1494  for (ibd = 0; ibd < nbd; ++ibd) {
1495  LogTrace("DataDump") << "*********** Bd " << ibd << "/" << nbd << " ***********";
1496  const pat::CompositeCandidate& cand = bdCands->at(ibd);
1497  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1498  LogTrace("DataDump") << "JPsi: " << jPsi;
1499  if (jPsi == nullptr)
1500  continue;
1501  if (!npJPsiBasicSelect->accept(*jPsi))
1502  continue;
1503  if (!npJPsiDaughterSelect->accept(*jPsi))
1504  continue;
1505  const pat::CompositeCandidate* kx0 = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToKx0");
1506  LogTrace("DataDump") << "Kx0: " << kx0;
1507  if (kx0 == nullptr)
1508  continue;
1509  if (bdIBasicSelect->accept(cand) && bdIJPsiBasicSelect->accept(*jPsi) && bdIKx0BasicSelect->accept(*kx0) &&
1510  // bdIJPsiDaughterSelect->accept( *jPsi ) &&
1511  bdIVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1512  fillHisto("DIBd", cand, 'f');
1513  fillHisto("DIBdJPsi", *jPsi, 'c');
1514  fillHisto("DIBdKx0", *kx0, 'c');
1515  if (flag_Dimuon25_Jpsi) {
1516  fillHisto("TIBd", cand, 'f');
1517  fillHisto("TIBdJPsi", *jPsi, 'c');
1518  fillHisto("TIBdKx0", *kx0, 'c');
1519  *ofMap["InclusiveBd"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1520  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1521  }
1522  }
1523  if (bdDBasicSelect->accept(cand) && bdDJPsiBasicSelect->accept(*jPsi) && bdDKx0BasicSelect->accept(*kx0) &&
1524  // bdDJPsiDaughterSelect->accept( *jPsi ) &&
1525  bdDVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1526  fillHisto("DDBd", cand, 'f');
1527  fillHisto("DDBdJPsi", *jPsi, 'c');
1528  fillHisto("DDBdKx0", *kx0, 'c');
1529  if (flag_DoubleMu4_JpsiTrk_Displaced) {
1530  fillHisto("TDBd", cand, 'f');
1531  fillHisto("TDBdJPsi", *jPsi, 'c');
1532  fillHisto("TDBdKx0", *kx0, 'c');
1533  *ofMap["DisplacedBd"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1534  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1535  }
1536  }
1537  }
1538 
1540 
1542  int ibs;
1543  int nbs = 0;
1544  if (useBs) {
1545  bsCandsToken.get(ev, bsCands);
1546  nbs = bsCands->size();
1547  }
1548 
1549  for (ibs = 0; ibs < nbs; ++ibs) {
1550  LogTrace("DataDump") << "*********** Bs " << ibs << "/" << nbs << " ***********";
1551  const pat::CompositeCandidate& cand = bsCands->at(ibs);
1552  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1553  LogTrace("DataDump") << "JPsi: " << jPsi;
1554  if (jPsi == nullptr)
1555  continue;
1556  if (!npJPsiBasicSelect->accept(*jPsi))
1557  continue;
1558  if (!npJPsiDaughterSelect->accept(*jPsi))
1559  continue;
1560  const pat::CompositeCandidate* phi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToPhi");
1561  LogTrace("DataDump") << "Phi: " << phi;
1562  if (phi == nullptr)
1563  continue;
1564  if (bsIBasicSelect->accept(cand) && bsIJPsiBasicSelect->accept(*jPsi) && bsIPhiBasicSelect->accept(*phi) &&
1565  // bsIJPsiDaughterSelect->accept( *jPsi ) &&
1566  bsIVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1567  fillHisto("DIBs", cand, 'f');
1568  fillHisto("DIBsJPsi", *jPsi, 'c');
1569  fillHisto("DIBsPhi", *phi, 'c');
1570  if (flag_Dimuon25_Jpsi) {
1571  fillHisto("TIBs", cand, 'f');
1572  fillHisto("TIBsJPsi", *jPsi, 'c');
1573  fillHisto("TIBsPhi", *phi, 'c');
1574  *ofMap["InclusiveBs"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1575  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1576  }
1577  }
1578  if (bsDBasicSelect->accept(cand) && bsDJPsiBasicSelect->accept(*jPsi) && bsDPhiBasicSelect->accept(*phi) &&
1579  // bsDJPsiDaughterSelect->accept( *jPsi ) &&
1580  bsDVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1581  fillHisto("DDBs", cand, 'f');
1582  fillHisto("DDBsJPsi", *jPsi, 'c');
1583  fillHisto("DDBsPhi", *phi, 'c');
1584  if (flag_DoubleMu4_JpsiTrk_Displaced) {
1585  fillHisto("TDBs", cand, 'f');
1586  fillHisto("TDBsJPsi", *jPsi, 'c');
1587  fillHisto("TDBsPhi", *phi, 'c');
1588  *ofMap["DisplacedBs"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1589  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1590  }
1591  }
1592  }
1593 
1595 
1597  int ik0;
1598  int nk0 = 0;
1599  if (useK0) {
1600  k0CandsToken.get(ev, k0Cands);
1601  nk0 = k0Cands->size();
1602  }
1603 
1604  for (ik0 = 0; ik0 < nk0; ++ik0) {
1605  LogTrace("DataDump") << "*********** K0 " << ik0 << "/" << nk0 << " ***********";
1606  const pat::CompositeCandidate& cand = k0Cands->at(ik0);
1607  fillHisto("DK0s", cand, 'f');
1608  *ofMap["K0s"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1609  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1610  }
1611 
1613 
1615  int il0;
1616  int nl0 = 0;
1617  if (useL0) {
1618  l0CandsToken.get(ev, l0Cands);
1619  nl0 = l0Cands->size();
1620  }
1621 
1622  for (il0 = 0; il0 < nl0; ++il0) {
1623  LogTrace("DataDump") << "*********** Lambda0 " << il0 << "/" << nl0 << " ***********";
1624  const pat::CompositeCandidate& cand = l0Cands->at(il0);
1625  fillHisto("DLambda0", cand, 'f');
1626  *ofMap["Lambda0"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1627  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1628  }
1629 
1631 
1633  int ib0;
1634  int nb0 = 0;
1635  if (useB0) {
1636  b0CandsToken.get(ev, b0Cands);
1637  nb0 = b0Cands->size();
1638  }
1639 
1640  // cout << nb0 << ' ' << ev.id().run() << ' ' << ev.id().event();
1641  // if ( nb0 ) cout << " *************************";
1642  // cout << endl;
1643 
1644  for (ib0 = 0; ib0 < nb0; ++ib0) {
1645  LogTrace("DataDump") << "*********** B0 " << ib0 << "/" << nb0 << " ***********";
1646  const pat::CompositeCandidate& cand = b0Cands->at(ib0);
1647  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1648  LogTrace("DataDump") << "JPsi: " << jPsi;
1649  if (jPsi == nullptr)
1650  continue;
1651  if (!npJPsiBasicSelect->accept(*jPsi))
1652  continue;
1653  if (!npJPsiDaughterSelect->accept(*jPsi))
1654  continue;
1655  const pat::CompositeCandidate* k0s = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToK0s");
1656  LogTrace("DataDump") << "K0s: " << k0s;
1657  if (k0s == nullptr)
1658  continue;
1659  if (b0IBasicSelect->accept(cand) && b0IJPsiBasicSelect->accept(*jPsi) && b0IK0sBasicSelect->accept(*k0s) &&
1660  // b0IJPsiDaughterSelect->accept( *jPsi ) &&
1661  b0IVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1662  fillHisto("DIB0", cand, 'f');
1663  fillHisto("DIB0JPsi", *jPsi, 'c');
1664  fillHisto("DIB0K0s", *k0s, 'c');
1665  if (flag_Dimuon25_Jpsi) {
1666  fillHisto("TIB0", cand, 'f');
1667  fillHisto("TIB0JPsi", *jPsi, 'c');
1668  fillHisto("TIB0K0s", *k0s, 'c');
1669  *ofMap["InclusiveB0"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1670  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1671  }
1672  }
1673  if (b0DBasicSelect->accept(cand) && b0DJPsiBasicSelect->accept(*jPsi) && b0DK0sBasicSelect->accept(*k0s) &&
1674  // b0DJPsiDaughterSelect->accept( *jPsi ) &&
1675  b0DVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1676  fillHisto("DDB0", cand, 'f');
1677  fillHisto("DDB0JPsi", *jPsi, 'c');
1678  fillHisto("DDB0K0s", *k0s, 'c');
1679  if (flag_DoubleMu4_JpsiTrk_Displaced) {
1680  fillHisto("TDB0", cand, 'f');
1681  fillHisto("TDB0JPsi", *jPsi, 'c');
1682  fillHisto("TDB0K0s", *k0s, 'c');
1683  *ofMap["DisplacedB0"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1684  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1685  }
1686  }
1687  }
1688 
1690 
1692  int ilb;
1693  int nlb = 0;
1694  if (useLb) {
1695  lbCandsToken.get(ev, lbCands);
1696  nlb = lbCands->size();
1697  }
1698 
1699  for (ilb = 0; ilb < nlb; ++ilb) {
1700  LogTrace("DataDump") << "*********** Lambdab " << ilb << "/" << nlb << " ***********";
1701  const pat::CompositeCandidate& cand = lbCands->at(ilb);
1702  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1703  LogTrace("DataDump") << "JPsi: " << jPsi;
1704  if (jPsi == nullptr)
1705  continue;
1706  if (!npJPsiBasicSelect->accept(*jPsi))
1707  continue;
1708  if (!npJPsiDaughterSelect->accept(*jPsi))
1709  continue;
1710  const pat::CompositeCandidate* l0 = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToLambda0");
1711  LogTrace("DataDump") << "Lambda0: " << l0;
1712  if (l0 == nullptr)
1713  continue;
1714  if (lbIBasicSelect->accept(cand) && lbIJPsiBasicSelect->accept(*jPsi) && lbILambda0BasicSelect->accept(*l0) &&
1715  // lbIJPsiDaughterSelect->accept( *jPsi ) &&
1716  lbIVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1717  fillHisto("DILambdab", cand, 'f');
1718  fillHisto("DILbJPsi", *jPsi, 'c');
1719  fillHisto("DILbL0", *l0, 'c');
1720  if (flag_Dimuon25_Jpsi) {
1721  fillHisto("TILambdab", cand, 'f');
1722  fillHisto("TILbJPsi", *jPsi, 'c');
1723  fillHisto("TILbL0", *l0, 'c');
1724  *ofMap["InclusiveLambdab"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1725  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1726  }
1727  }
1728  if (lbDBasicSelect->accept(cand) && lbDJPsiBasicSelect->accept(*jPsi) && lbDLambda0BasicSelect->accept(*l0) &&
1729  // lbDJPsiDaughterSelect->accept( *jPsi ) &&
1730  lbDVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex"))) {
1731  fillHisto("DDLambdab", cand, 'f');
1732  fillHisto("DDLbJPsi", *jPsi, 'c');
1733  fillHisto("DDLbL0", *l0, 'c');
1734  if (flag_DoubleMu4_JpsiTrk_Displaced) {
1735  fillHisto("TDLambdab", cand, 'f');
1736  fillHisto("TDLbJPsi", *jPsi, 'c');
1737  fillHisto("TDLbL0", *l0, 'c');
1738  *ofMap["DisplacedLambdab"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1739  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1740  }
1741  }
1742  }
1743 
1745 
1747  int ibc;
1748  int nbc = 0;
1749  if (useBc) {
1750  bcCandsToken.get(ev, bcCands);
1751  nbc = bcCands->size();
1752  }
1753 
1754  for (ibc = 0; ibc < nbc; ++ibc) {
1755  LogTrace("DataDump") << "*********** Bc " << ibc << "/" << nbc << " ***********";
1756  const pat::CompositeCandidate& cand = bcCands->at(ibc);
1757  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1758  LogTrace("DataDump") << "JPsi: " << jPsi;
1759  if (jPsi == nullptr)
1760  continue;
1761  // if ( BPHUserData::get( *jPsi, "dca", -1.0 ) < bcJPsiDcaMax ) continue;
1762  if (!npJPsiBasicSelect->accept(*jPsi))
1763  continue;
1764  if (!npJPsiDaughterSelect->accept(*jPsi))
1765  continue;
1766  const reco::Candidate* pptr = BPHDaughters::get(cand, 0.13, 0.14).front();
1767  if (pptr == nullptr)
1768  continue;
1769 
1770  if (bcIBasicSelect->accept(cand) && bcIJPsiBasicSelect->accept(*jPsi) &&
1771  // bcIJPsiDaughterSelect->accept( *jPsi ) &&
1772  bcIJPsiVertexSelect->accept(*jPsi, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1773  bcIVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1774  (pptr->pt() > bcIPiPtMin)) {
1775  fillHisto("DIBc", cand, 'f');
1776  fillHisto("DIBcJPsi", *jPsi, 'c');
1777  if (flag_Dimuon25_Jpsi) {
1778  fillHisto("TIBc", cand, 'f');
1779  fillHisto("TIBcJPsi", *jPsi, 'c');
1780  *ofMap["InclusiveBc"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1781  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1782  }
1783  }
1784  if (bcDBasicSelect->accept(cand) && bcDJPsiBasicSelect->accept(*jPsi) &&
1785  // bcDJPsiDaughterSelect->accept( *jPsi ) &&
1786  bcDVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1787  bcDVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1788  (pptr->pt() > bcDPiPtMin)) {
1789  fillHisto("DDBc", cand, 'f');
1790  fillHisto("DDBcJPsi", *jPsi, 'c');
1791  if (flag_DoubleMu4_JpsiTrk_Displaced) {
1792  fillHisto("TDBc", cand, 'f');
1793  fillHisto("TDBcJPsi", *jPsi, 'c');
1794  *ofMap["DisplacedBc"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1795  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1796  }
1797  }
1798  }
1799 
1801 
1803  int ix3872;
1804  int nx3872 = 0;
1805  if (useX3872) {
1806  x3872CandsToken.get(ev, x3872Cands);
1807  nx3872 = x3872Cands->size();
1808  }
1809 
1810  for (ix3872 = 0; ix3872 < nx3872; ++ix3872) {
1811  LogTrace("DataDump") << "*********** X3872 " << ix3872 << "/" << nx3872 << " ***********";
1812  const pat::CompositeCandidate& cand = x3872Cands->at(ix3872);
1813  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1814  LogTrace("DataDump") << "JPsi: " << jPsi;
1815  if (jPsi == nullptr)
1816  continue;
1817  // if ( BPHUserData::get( *jPsi, "dca", -1.0 ) < x3872JPsiDcaMax ) continue;
1818  if (!npJPsiBasicSelect->accept(*jPsi))
1819  continue;
1820  if (!npJPsiDaughterSelect->accept(*jPsi))
1821  continue;
1822  const reco::Candidate* ppt1 = BPHDaughters::get(cand, 0.13, 0.14)[0];
1823  const reco::Candidate* ppt2 = BPHDaughters::get(cand, 0.13, 0.14)[1];
1824  if (ppt1 == nullptr)
1825  continue;
1826  if (ppt2 == nullptr)
1827  continue;
1828  if (x3872IBasicSelect->accept(cand) && x3872IJPsiBasicSelect->accept(*jPsi) &&
1829  // x3872IJPsiDaughterSelect->accept( *jPsi ) &&
1830  x3872IJPsiVertexSelect->accept(*jPsi, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1831  x3872IVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1832  (ppt1->pt() > x3872IPiPtMin) && (ppt2->pt() > x3872IPiPtMin)) {
1833  fillHisto("DIX3872", cand, 'f');
1834  fillHisto("DIX3872JPsi", *jPsi, 'c');
1835  if (flag_Dimuon25_Jpsi) {
1836  fillHisto("TIX3872", cand, 'f');
1837  fillHisto("TIX3872JPsi", *jPsi, 'c');
1838  *ofMap["InclusiveX3872"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1839  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1840  }
1841  }
1842  if (x3872DBasicSelect->accept(cand) && x3872DJPsiBasicSelect->accept(*jPsi) &&
1843  // x3872DJPsiDaughterSelect->accept( *jPsi ) &&
1844  x3872DVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1845  x3872DVertexSelect->accept(cand, BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex")) &&
1846  (ppt1->pt() > x3872DPiPtMin) && (ppt2->pt() > x3872DPiPtMin)) {
1847  fillHisto("DDX3872", cand, 'f');
1848  fillHisto("DDX3872JPsi", *jPsi, 'c');
1849  if (flag_DoubleMu4_JpsiTrk_Displaced) {
1850  fillHisto("TDX3872", cand, 'f');
1851  fillHisto("TDX3872JPsi", *jPsi, 'c');
1852  *ofMap["DisplacedX3872"] << ev.id().run() << ' ' << ev.id().luminosityBlock() << ' ' << ev.id().event() << ' '
1853  << (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1) << endl;
1854  }
1855  }
1856  }
1857 
1858  return;
1859 }
1860 
1862  // tree->Write();
1863  return;
1864 }
1865 
1866 void BPHHistoSpecificDecay::fillHisto(const string& name, const pat::CompositeCandidate& cand, char svType) {
1867  *recoName = name;
1868  float mass = (cand.hasUserFloat("fitMass") ? cand.userFloat("fitMass") : -1);
1869  fillHisto("mass" + name, cand.mass());
1870  fillHisto("mfit" + name, mass);
1871  recoMass = mass;
1872  recoTime = -999999.0;
1873  recoErrT = -999999.0;
1874 
1875  const reco::Vertex* pvtx = BPHUserData::getByRef<reco::Vertex>(cand, "primaryVertex");
1876  if (pvtx == nullptr) {
1877  const pat::CompositeCandidate* jPsi = BPHUserData::getByRef<pat::CompositeCandidate>(cand, "refToJPsi");
1878  if (jPsi == nullptr)
1879  return;
1880  pvtx = BPHUserData::getByRef<reco::Vertex>(*jPsi, "primaryVertex");
1881  }
1882  // if ( pvtx == nullptr ) return;
1883 
1884  if (pvtx != nullptr) {
1885  const reco::Vertex* svtx = nullptr;
1886  if (svType == 'f')
1887  svtx = BPHUserData::get<reco::Vertex>(cand, "fitVertex");
1888  if (svtx == nullptr)
1889  svtx = BPHUserData::get<reco::Vertex>(cand, "vertex");
1890  // if ( svtx == nullptr ) return;
1891  if (svtx != nullptr) {
1892  float px;
1893  float py;
1894  const Vector3DBase<float, GlobalTag>* fmom =
1895  BPHUserData::get<Vector3DBase<float, GlobalTag>>(cand, "fitMomentum");
1896  if (fmom != nullptr) {
1897  px = fmom->x();
1898  py = fmom->y();
1899  } else {
1900  px = cand.px();
1901  py = cand.py();
1902  }
1903  double ctauPV2;
1904  double ctauErrPV2;
1905  VertexAnalysis::dist2D(pvtx, svtx, px, py, mass, ctauPV2, ctauErrPV2);
1906  recoTime = ctauPV2;
1907  recoErrT = ctauErrPV2;
1908  fillHisto("ctau" + name, ctauPV2);
1909  }
1910  }
1911  tree->Fill();
1912  return;
1913 }
1914 
1915 void BPHHistoSpecificDecay::fillHisto(const string& name, float x) {
1916  map<string, TH1F*>::iterator iter = histoMap.find(name);
1917  map<string, TH1F*>::iterator iend = histoMap.end();
1918  if (iter == iend)
1919  return;
1920  iter->second->Fill(x);
1921  return;
1922 }
1923 
1924 void BPHHistoSpecificDecay::createHisto(const string& name, int nbin, float hmin, float hmax) {
1925  histoMap[name] = fs->make<TH1F>(name.c_str(), name.c_str(), nbin, hmin, hmax);
1926  return;
1927 }
1928 
1930 
RunNumber_t run() const
Definition: EventID.h:38
Analysis-level particle class.
EventNumber_t event() const
Definition: EventID.h:40
BPHDaughterSelect(float ptMinLoose, float ptMinTight, float etaMaxLoose, float etaMaxTight, const BPHSoftMuonSelect *softMuonselector=nullptr)
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
#define CHK_TRIG(RESULTS, NAMES, INDEX, PATH)
static bool accept(const pat::CompositeCandidate &cand, float ptMinLoose, float ptMinTight, float etaMaxLoose, float etaMaxTight, const reco::Vertex *pv=nullptr, const BPHSoftMuonSelect *softMuonselector=nullptr)
const reco::Vertex * pv
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
static std::vector< std::string > checklist log
double pt() const final
transverse momentum
static void dist2D(const reco::Vertex *pvtx, const reco::Vertex *svtx, float px, float py, float mass, double &ctauPV, double &ctauErrPV)
BPHHistoSpecificDecay(const edm::ParameterSet &ps)
BPHCompositeBasicSelect(float massMin, float massMax, float ptMin=-1.0, float etaMax=-1.0, float rapidityMax=-1.0)
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
tuple maxDxy
Definition: gather_cfg.py:53
static const T * getByRef(const pat::CompositeCandidate &cand, const string &name)
std::vector< std::string > Strings
Definition: TriggerNames.h:58
virtual double pt() const =0
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const AlgebraicSymMatrix33 matrix() const
double y() const
y coordinate
Definition: Vertex.h:131
virtual double mass() const =0
mass
constexpr float ptMin
static void dist2D(const reco::Vertex *pvtx, const reco::Vertex *svtx, float px, float py, double cosAlpha, float mass, double &ctauPV, double &ctauErrPV)
bool hasUserFloat(const std::string &key) const
Return true if there is a user-defined float with a given name.
Definition: PATObject.h:373
bool accept(const pat::CompositeCandidate &cand, const reco::Vertex *pv=nullptr) const override
tuple magneticField
BPHGenericVertexSelect(char vType, float probMin, float cosMin=-2.0, float sigMin=-1.0, char dMode= 'r')
bool ev
static const double jPsiMass
static const double k0sMass
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:39
const Point & position() const
position
Definition: Vertex.h:127
Strings const & triggerNames() const
Definition: TriggerNames.cc:48
const std::string names[nVars_]
float userFloat(const std::string &key) const
Definition: PATObject.h:871
#define LogTrace(id)
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:275
tuple d
Definition: ztail.py:151
void analyze(const edm::Event &ev, const edm::EventSetup &es) override
static const double kx0Mass
bool accept(const pat::CompositeCandidate &cand, const reco::Vertex *pv=nullptr) const override
T sqrt(T t)
Definition: SSEVec.h:19
double rapidity() const final
rapidity
void createHisto(const std::string &name, int nbin, float hmin, float hmax)
double chi2() const
chi-squares
Definition: Vertex.h:116
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
float ChiSquaredProbability(double chiSquared, double nrDOF)
BPHFittedBasicSelect(float massMin, float massMax, float ptMin=-1.0, float etaMax=-1.0, float rapidityMax=-1.0)
const T * userData(const std::string &key) const
Returns user-defined data. Returns NULL if the data is not present, or not of type T...
Definition: PATObject.h:322
void fillHisto(const std::string &name, const pat::CompositeCandidate &cand, char svType)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
bool hasUserData(const std::string &key) const
Check if user data with a specific type is present.
Definition: PATObject.h:327
double ndof() const
Definition: Vertex.h:123
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
bool accept(const reco::Candidate &cand, const reco::Vertex *pv) const
static const double lambda0Mass
static const double phiMass
static double cAlpha(const reco::Vertex *pvtx, const reco::Vertex *svtx, float px, float py)
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:57
bool accept(const pat::CompositeCandidate &cand, const reco::Vertex *pvtx) const override
double x() const
x coordinate
Definition: Vertex.h:129
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double value() const
Definition: Measurement1D.h:25
Error error() const
return SMatrix
Definition: Vertex.h:163
reco::TrackRef innerTrack() const override
reference to Track reconstructed in the tracker only (reimplemented from reco::Muon) ...
Definition: Muon.h:72
edm::EventID id() const
Definition: EventBase.h:59
double mass() const final
mass
T get() const
Definition: EventSetup.h:88
ROOT::Math::SVector< double, 3 > AlgebraicVector3
BPHSoftMuonSelect(int cutTrackerLayers=5, int cutPixelLayers=0, float maxDxy=0.3, float maxDz=20.0, bool goodMuon=true, bool highPurity=true)
static vector< const reco::Candidate * > get(const pat::CompositeCandidate &cand, float massMin, float massMax)
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
tuple etaMax
Definition: Puppi_cff.py:46
long double T
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const BPHSoftMuonSelect * sms
Analysis-level muon class.
Definition: Muon.h:51
static double cAlpha(const reco::Vertex *pvtx, const reco::Vertex *svtx, const TVector3 &cmom)
#define SET_LABEL(NAME, PSET)
bool accept(const pat::CompositeCandidate &cand, const reco::Vertex *pv=nullptr) const override
virtual double eta() const =0
momentum pseudorapidity
double eta() const final
momentum pseudorapidity