CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Zto2lFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Zto2lFilter
4 // Class: Zto2lFilter
5 //
13 //
14 // Original Author: Aruna Nayak
15 // Created: Thu Aug 23 11:37:45 CEST 2007
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
23 
24 #include <vector>
26 #include "TLorentzVector.h"
27 
28 //
29 // constants, enums and typedefs
30 //
31 
32 using namespace std;
33 using namespace edm;
34 //
35 // static data member definitions
36 //
37 
38 //
39 // constructors and destructor
40 //
42 {
43  //now do what ever initialization is needed
44  fLabel_ = iConfig.getUntrackedParameter("moduleLabel",std::string("generator"));
45  maxEtaLepton_ = iConfig.getUntrackedParameter<double>("MaxEtaLepton");
46  minInvariantMass_ = iConfig.getUntrackedParameter<double>("MindiLeptonInvariantMass");
47 
48 }
49 
50 
52 {
53 
54  // do anything here that needs to be done at desctruction time
55  // (e.g. close files, deallocate resources etc.)
56 
57 }
58 
59 
60 //
61 // member functions
62 //
63 
64 // ------------ method called on each new Event ------------
65 bool
67 {
68  using namespace edm;
69 
70  bool accept = false;
71 
72  Handle<HepMCProduct> EvtHandle ;
73  iEvent.getByLabel( fLabel_, "unsmeared", EvtHandle ) ;
74  const HepMC::GenEvent* evt = EvtHandle->GetEvent();
75 
76  vector<TLorentzVector> Lepton; Lepton.clear();
77  for (HepMC::GenEvent::particle_const_iterator p = evt->particles_begin();
78  p != evt->particles_end(); ++p) {
79  if((*p)->status()==3){
80  if ( abs((*p)->pdg_id()) == 11 || abs((*p)->pdg_id()) == 13 || abs((*p)->pdg_id()) == 15 ){
81  if(fabs((*p)->momentum().eta()) < maxEtaLepton_){
82  TLorentzVector LeptP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
83  Lepton.push_back(LeptP);
84  }
85  }
86  }
87  }
88  if(Lepton.size() == 2){
89  if((Lepton[0]+Lepton[1]).M() > minInvariantMass_ )accept = true;
90  }
91  //delete evt;
92  return accept;
93 }
94 
95 // ------------ method called once each job just before starting event loop ------------
96 
97 // ------------ method called once each job just after ending the event loop ------------
98 void
100 }
101 
102 //define this as a plug-in
103 //DEFINE_FWK_MODULE(Zto2lFilter);
Zto2lFilter(const edm::ParameterSet &)
Definition: Zto2lFilter.cc:41
T getUntrackedParameter(std::string const &, T const &) const
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
virtual void endJob()
Definition: Zto2lFilter.cc:99
int iEvent
Definition: GenABIO.cc:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
virtual bool filter(edm::Event &, const edm::EventSetup &)
Definition: Zto2lFilter.cc:66