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 
35  public:
36 
40  const BPHMassSelect* ms ): nPos( np ), nNeg( nn ),
41  mSel( ms ) {}
42 
45  ~BPHMassSymSelect() override {}
46 
49  bool accept( const BPHDecayMomentum& cand ) const override {
51 
52  if ( mSel->accept( cand ) ) return true;
53 
54  const
55  reco::Candidate* pp = cand.getDaug( nPos );
56  const
57  reco::Candidate* np = cand.getDaug( nNeg );
58 
59  reco::Candidate* pc = cand.originalReco( pp )->clone();
60  reco::Candidate* nc = cand.originalReco( np )->clone();
61 
62  pc->setMass( np->p4().mass() );
63  nc->setMass( pp->p4().mass() );
64  const reco::Candidate::LorentzVector s4 = pc->p4() + nc->p4();
65  double mass = s4.mass();
66 
67  delete pc;
68  delete nc;
69  return ( ( mass > mSel->getMassMin() ) &&
70  ( mass < mSel->getMassMax() ) );
71 
72  }
73 
74  private:
75 
76  // private copy and assigment constructors
77  BPHMassSymSelect ( const BPHMassSymSelect& x ) = delete;
78  BPHMassSymSelect& operator=( const BPHMassSymSelect& x ) = delete;
79 
83 
84 };
85 
86 
87 #endif
88 
bool accept(const BPHDecayMomentum &cand) const override
select particle
Definition: BPHMassSelect.h:48
BPHMassSymSelect(const std::string &np, const std::string &nn, const BPHMassSelect *ms)
int np
Definition: AMPTWrapper.h:33
virtual Candidate * clone() const =0
returns a clone of the Candidate object
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
~BPHMassSymSelect() override
bool accept(const BPHDecayMomentum &cand) const override
select particle
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
double getMassMin() const
get current mass cuts
Definition: BPHMassCuts.h:54
susybsm::MuonSegment ms
Definition: classes.h:31
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
virtual const reco::Candidate * getDaug(const std::string &name) const
virtual void setMass(double m)=0
set particle mass
BPHMassSymSelect & operator=(const BPHMassSymSelect &x)=delete
const BPHMassSelect * mSel