CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/ElectroWeakAnalysis/ZMuMu/plugins/ZToMuMuFilter.cc

Go to the documentation of this file.
00001 /* \class ZToMuMuFilter
00002  *
00003  * \author Juan Alcaraz, CIEMAT
00004  *
00005  */
00006 #include "FWCore/Framework/interface/EDFilter.h"
00007 #include "FWCore/Utilities/interface/InputTag.h"
00008 
00009 class ZToMuMuFilter : public edm::EDFilter {
00010 public:
00011   ZToMuMuFilter(const edm::ParameterSet &);
00012 private:
00013   virtual bool filter(edm::Event&, const edm::EventSetup&);
00014   edm::InputTag zCands_, muIso1_, muIso2_;
00015   double ptMin_, etaMin_, etaMax_, massMin_, massMax_, isoMax_;
00016 };
00017 
00018 #include "DataFormats/Common/interface/Handle.h"
00019 #include "FWCore/Framework/interface/Event.h"
00020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00021 #include "DataFormats/Candidate/interface/Candidate.h"
00022 #include "DataFormats/Candidate/interface/CandAssociation.h"
00023 using namespace edm;
00024 using namespace std;
00025 using namespace reco;
00026 
00027 ZToMuMuFilter::ZToMuMuFilter( const ParameterSet & cfg ) :
00028   zCands_(cfg.getParameter<InputTag>("zCands")),
00029   muIso1_(cfg.getParameter<InputTag>("muonIsolations1")),
00030   muIso2_(cfg.getParameter<InputTag>("muonIsolations2")),
00031   ptMin_(cfg.getParameter<double>("ptMin")),
00032   etaMin_(cfg.getParameter<double>("etaMin")),
00033   etaMax_(cfg.getParameter<double>("etaMax")),
00034   massMin_(cfg.getParameter<double>("massMin")),
00035   massMax_(cfg.getParameter<double>("massMax")),
00036   isoMax_(cfg.getParameter<double>("isoMax")) {
00037 }
00038 
00039 bool ZToMuMuFilter::filter (Event & ev, const EventSetup &) {
00040   Handle<CandidateCollection> zCands;
00041   ev.getByLabel(zCands_, zCands);
00042   Handle<CandDoubleAssociations> muIso1;
00043   ev.getByLabel(muIso1_, muIso1);
00044   Handle<CandDoubleAssociations> muIso2;
00045   ev.getByLabel(muIso2_, muIso2);
00046   unsigned int nZ = zCands->size();
00047   if (nZ == 0) return false;
00048   for(unsigned int i = 0; i < nZ; ++ i) {
00049     const Candidate & zCand = (*zCands)[i];
00050     double zMass = zCand.mass();
00051     if (zMass < massMin_ || zMass > massMax_) return false;
00052     if(zCand.numberOfDaughters()!=2) return false;
00053     const Candidate * dau0 = zCand.daughter(0);
00054     const Candidate * dau1 = zCand.daughter(1);
00055     double pt0 = dau0->pt(), pt1 = dau1->pt();
00056     if (pt0 < ptMin_ || pt1 < ptMin_) return false;
00057     double eta0 = dau0->eta(), eta1 = dau1->eta();
00058     if(eta0 < etaMin_ || eta0 > etaMax_) return false;
00059     if(eta1 < etaMin_ || eta1 > etaMax_) return false;
00060     CandidateRef mu0 = dau0->masterClone().castTo<CandidateRef>();
00061     CandidateRef mu1 = dau1->masterClone().castTo<CandidateRef>();
00062     double iso0 = (*muIso1)[mu0];
00063     double iso1 = (*muIso2)[mu1];
00064     if (iso0 > isoMax_) return false;
00065     if (iso1 > isoMax_) return false;
00066   }
00067   return true;
00068 }
00069 
00070 #include "FWCore/Framework/interface/MakerMacros.h"
00071 
00072 DEFINE_FWK_MODULE( ZToMuMuFilter );