CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZToMuMuFilter.cc
Go to the documentation of this file.
1 /* \class ZToMuMuFilter
2  *
3  * \author Juan Alcaraz, CIEMAT
4  *
5  */
10 
11 class ZToMuMuFilter : public edm::EDFilter {
12 public:
14 private:
15  virtual bool filter(edm::Event&, const edm::EventSetup&) override;
20 };
21 
25 using namespace edm;
26 using namespace std;
27 using namespace reco;
28 
30  zCandsToken_(consumes<CandidateCollection>(cfg.getParameter<InputTag>("zCands"))),
31  muIso1Token_(consumes<CandDoubleAssociations>(cfg.getParameter<InputTag>("muonIsolations1"))),
32  muIso2Token_(consumes<CandDoubleAssociations>(cfg.getParameter<InputTag>("muonIsolations2"))),
33  ptMin_(cfg.getParameter<double>("ptMin")),
34  etaMin_(cfg.getParameter<double>("etaMin")),
35  etaMax_(cfg.getParameter<double>("etaMax")),
36  massMin_(cfg.getParameter<double>("massMin")),
37  massMax_(cfg.getParameter<double>("massMax")),
38  isoMax_(cfg.getParameter<double>("isoMax")) {
39 }
40 
43  ev.getByToken(zCandsToken_, zCands);
45  ev.getByToken(muIso1Token_, muIso1);
47  ev.getByToken(muIso2Token_, muIso2);
48  unsigned int nZ = zCands->size();
49  if (nZ == 0) return false;
50  for(unsigned int i = 0; i < nZ; ++ i) {
51  const Candidate & zCand = (*zCands)[i];
52  double zMass = zCand.mass();
53  if (zMass < massMin_ || zMass > massMax_) return false;
54  if(zCand.numberOfDaughters()!=2) return false;
55  const Candidate * dau0 = zCand.daughter(0);
56  const Candidate * dau1 = zCand.daughter(1);
57  double pt0 = dau0->pt(), pt1 = dau1->pt();
58  if (pt0 < ptMin_ || pt1 < ptMin_) return false;
59  double eta0 = dau0->eta(), eta1 = dau1->eta();
60  if(eta0 < etaMin_ || eta0 > etaMax_) return false;
61  if(eta1 < etaMin_ || eta1 > etaMax_) return false;
64  double iso0 = (*muIso1)[mu0];
65  double iso1 = (*muIso2)[mu1];
66  if (iso0 > isoMax_) return false;
67  if (iso1 > isoMax_) return false;
68  }
69  return true;
70 }
71 
73 
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
tuple cfg
Definition: looper.py:293
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual double pt() const =0
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual double mass() const =0
mass
bool ev
edm::EDGetTokenT< reco::CandidateCollection > zCandsToken_
virtual size_type numberOfDaughters() const =0
number of daughters
ZToMuMuFilter(const edm::ParameterSet &)
edm::EDGetTokenT< reco::CandDoubleAssociations > muIso1Token_
REF castTo() const
Definition: RefToBase.h:278
virtual bool filter(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< reco::CandDoubleAssociations > muIso2Token_
virtual double eta() const =0
momentum pseudorapidity
virtual const CandidateBaseRef & masterClone() const =0