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