CMS 3D CMS Logo

BPHMassSymSelect.h
Go to the documentation of this file.
1 #ifndef HeavyFlavorAnalysis_SpecificDecay_BPHMassSymSelect_h
2 #define HeavyFlavorAnalysis_SpecificDecay_BPHMassSymSelect_h
3 
13 //----------------------
14 // Base Class Headers --
15 //----------------------
17 
18 //------------------------------------
19 // Collaborating Class Declarations --
20 //------------------------------------
23 
24 //---------------
25 // C++ Headers --
26 //---------------
27 #include <string>
28 
29 // ---------------------
30 // -- Class Interface --
31 // ---------------------
32 
34 public:
38  : nPos(np), nNeg(nn), mSel(ms) {}
39 
40  // deleted copy constructor and assignment operator
41  BPHMassSymSelect(const BPHMassSymSelect& x) = delete;
42  BPHMassSymSelect& operator=(const BPHMassSymSelect& x) = delete;
43 
46  ~BPHMassSymSelect() override {}
47 
50  bool accept(const BPHDecayMomentum& cand) const override {
52  if (mSel->accept(cand))
53  return true;
54 
55  const reco::Candidate* pp = cand.getDaug(nPos);
56  const reco::Candidate* np = cand.getDaug(nNeg);
57 
58  reco::Candidate* pc = cand.originalReco(pp)->clone();
59  reco::Candidate* nc = cand.originalReco(np)->clone();
60 
61  pc->setMass(np->p4().mass());
62  nc->setMass(pp->p4().mass());
63  const reco::Candidate::LorentzVector s4 = pc->p4() + nc->p4();
64  double mass = s4.mass();
65 
66  delete pc;
67  delete nc;
68  return ((mass >= mSel->getMassMin()) && (mass <= mSel->getMassMax()));
69  }
70 
71 private:
75 };
76 
77 #endif
reco::Candidate::setMass
virtual void setMass(double m)=0
set particle mass
np
int np
Definition: AMPTWrapper.h:43
BPHMassSymSelect::~BPHMassSymSelect
~BPHMassSymSelect() override
Definition: BPHMassSymSelect.h:46
BPHMassCuts::getMassMin
double getMassMin() const
get current mass cuts
Definition: BPHMassCuts.h:58
BPHMassSymSelect::nNeg
std::string nNeg
Definition: BPHMassSymSelect.h:73
BPHMassSymSelect::nPos
std::string nPos
Definition: BPHMassSymSelect.h:72
BPHMomentumSelect
Definition: BPHMomentumSelect.h:29
BPHDecayMomentum.h
DDAxes::x
BPHMassSymSelect
Definition: BPHMassSymSelect.h:33
BPHMassSymSelect::mSel
const BPHMassSelect * mSel
Definition: BPHMassSymSelect.h:74
BPHMassSymSelect::operator=
BPHMassSymSelect & operator=(const BPHMassSymSelect &x)=delete
BPHMassSelect
Definition: BPHMassSelect.h:31
BPHDecayMomentum
Definition: BPHDecayMomentum.h:35
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
BPHMomentumSelect.h
BPHMassSelect.h
cand
Definition: decayParser.h:32
BPHMassSymSelect::accept
bool accept(const BPHDecayMomentum &cand) const override
select particle
Definition: BPHMassSymSelect.h:51
reco::Candidate
Definition: Candidate.h:27
groupFilesInBlocks.nn
nn
Definition: groupFilesInBlocks.py:150
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
reco::Candidate::p4
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
BPHMassSelect::accept
bool accept(const BPHDecayMomentum &cand) const override
select particle
Definition: BPHMassSelect.h:48
BPHMassSymSelect::BPHMassSymSelect
BPHMassSymSelect(const std::string &np, const std::string &nn, const BPHMassSelect *ms)
Definition: BPHMassSymSelect.h:37
createTree.pp
pp
Definition: createTree.py:17
reco::Candidate::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36