CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MCZll.cc
Go to the documentation of this file.
1 
3 
5 #include <iostream>
6 
8 
9 using namespace edm;
10 using namespace std;
11 
12 
14  token_(consumes<edm::HepMCProduct>(edm::InputTag(iConfig.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared"))),
15  nEvents_(0),
16  nAccepted_(0)
17 {
18  leptonFlavour_ = iConfig.getUntrackedParameter<int>("leptonFlavour",11);
19  leptonPtMin_ = iConfig.getUntrackedParameter<double>("leptonPtMin",5.);
20  leptonPtMax_ = iConfig.getUntrackedParameter<double>("leptonPtMax",99999.);
21  leptonEtaMin_ = iConfig.getUntrackedParameter<double>("leptonEtaMin",0.);
22  leptonEtaMax_ = iConfig.getUntrackedParameter<double>("leptonEtaMax",2.7);
23  zMassRange_.first = iConfig.getUntrackedParameter<double>("zMassMin",60.);
24  zMassRange_.second = iConfig.getUntrackedParameter<double>("zMassMax",120.);
25  filter_ = iConfig.getUntrackedParameter<bool>("filter",true);
26  std::ostringstream str;
27  str << "=========================================================\n"
28  << "Filter MCZll being constructed with parameters: "
29  << "\nleptonFlavour " << leptonFlavour_
30  << "\nleptonPtMin " << leptonPtMin_
31  << "\nleptonPtMax " << leptonPtMax_
32  << "\nleptonEtaMin " << leptonEtaMin_
33  << "\nleptonEtaMax " << leptonEtaMax_
34  << "\nzMassMin " << zMassRange_.first
35  << "\nzMassMax " << zMassRange_.second
36  << "\n=========================================================" ;
37  edm::LogVerbatim("MCZllInfo") << str.str() ;
38  if (filter_)
39  produces<HepMCProduct>();
40 }
41 
42 
44 {
45 
46  // do anything here that needs to be done at desctruction time
47  // (e.g. close files, deallocate resources etc.)
48 
49 }
50 
52 {
53  edm::LogVerbatim("MCZllInfo") << "================MCZll report========================================\n"
54  << "Events read " << nEvents_ << " Events accepted " << nAccepted_ << "\nEfficiency " << ((double)nAccepted_)/((double)nEvents_)
55  << "\n====================================================================" << std::endl;
56 }
57 
58 // ------------ method called to skim the data ------------
60 {
61  std::auto_ptr<HepMCProduct> bare_product(new HepMCProduct());
62 
63  nEvents_++;
64  using namespace edm;
65  bool accepted = false;
67  iEvent.getByToken(token_, evt);
68  HepMC::GenEvent * myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
69  HepMC::GenEvent * zEvent = new HepMC::GenEvent();
70 
71  if (myGenEvent->signal_process_id() != 1)
72  {
73  delete myGenEvent;
74  delete zEvent;
75  return false;
76  }
77 
78 
79  //found a prompt Z
80 
81  for ( HepMC::GenEvent::particle_iterator p = myGenEvent->particles_begin();
82  p != myGenEvent->particles_end(); ++p )
83  {
84  if ( !accepted && ( (*p)->pdg_id() == 23 ) && (*p)->status() == 3 )
85  {
86  accepted=true;
87  HepMC::GenVertex* zVertex = new HepMC::GenVertex();
89  zVertex->add_particle_in(myZ);
90  // std::cout << (*p)->momentum().invariantMass() << std::endl;
91  if ((*p)->momentum().m() < zMassRange_.first || (*p)->momentum().m() > zMassRange_.second)
92  accepted = false;
93  std::vector<HepMC::GenParticle*> children;
94  HepMC::GenVertex* outVertex=(*p)->end_vertex();
95  for(HepMC::GenVertex::particles_out_const_iterator iter = outVertex->particles_out_const_begin();
96  iter != outVertex->particles_out_const_end(); iter++)
97  children.push_back(*iter);
98  std::vector<HepMC::GenParticle*>::const_iterator aDaughter;
99  for (aDaughter = children.begin();aDaughter != children.end();aDaughter++)
100  {
101  HepMC::GenParticle* myDa= new HepMC::GenParticle(*(*aDaughter));
102  zVertex->add_particle_out(myDa);
103  if ((*aDaughter)->status() == 2)
104  continue;
105  // (*aDaughter)->print();
106 
107  if (! (abs((*aDaughter)->pdg_id()) == abs(leptonFlavour_)) )
108  accepted = false;
109  // std::cout << (*aDaughter)->momentum().perp() << " " << (*aDaughter)->momentum().eta() << std::endl;
110  if ((*aDaughter)->momentum().perp() < leptonPtMin_)
111  accepted = false;
112  if ((*aDaughter)->momentum().perp() > leptonPtMax_)
113  accepted = false;
114  if (fabs((*aDaughter)->momentum().eta()) > leptonEtaMax_)
115  accepted = false;
116  if (fabs((*aDaughter)->momentum().eta()) < leptonEtaMin_)
117  accepted = false;
118  }
119  zEvent->add_vertex( zVertex );
120  if (accepted)
121  break;
122 
123  }
124 
125  }
126 
127 
128  if (accepted)
129  {
130  if(zEvent)
131  bare_product->addHepMCData(zEvent);
132  if (filter_)
133  iEvent.put(bare_product);
134  nAccepted_++;
135  // std::cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<< std::endl;
136  LogDebug("MCZll") << "Event " << iEvent.id().event() << " accepted" << std::endl;
137  // std::cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<< std::endl;
138  // myGenEvent->print();
139  delete myGenEvent;
140  return true;
141  }
142 
143  delete myGenEvent;
144  delete zEvent;
145  return false;
146 
147 
148 }
149 
#define LogDebug(id)
edm::EDGetTokenT< edm::HepMCProduct > token_
Definition: MCZll.h:50
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
bool filter_
Definition: MCZll.h:59
double leptonPtMin_
Definition: MCZll.h:52
double leptonEtaMin_
Definition: MCZll.h:54
int leptonFlavour_
Definition: MCZll.h:51
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
double leptonEtaMax_
Definition: MCZll.h:55
virtual void endJob()
Definition: MCZll.cc:51
int iEvent
Definition: GenABIO.cc:230
~MCZll()
Definition: MCZll.cc:43
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double leptonPtMax_
Definition: MCZll.h:53
unsigned int nEvents_
Definition: MCZll.h:57
edm::EventID id() const
Definition: EventBase.h:60
MCZll(const edm::ParameterSet &)
Definition: MCZll.cc:13
unsigned int nAccepted_
Definition: MCZll.h:58
std::pair< double, double > zMassRange_
Definition: MCZll.h:56
virtual bool filter(edm::Event &, const edm::EventSetup &)
Definition: MCZll.cc:59