CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DrellYanValidation.cc
Go to the documentation of this file.
1 /*class DrellYanValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  *
5  */
6 
10 #include "TLorentzVector.h"
11 
12 #include "CLHEP/Units/defs.h"
13 #include "CLHEP/Units/PhysicalConstants.h"
14 
16 
17 using namespace edm;
18 
20  wmanager_(iPSet,consumesCollector()),
21  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
22  _flavor(iPSet.getParameter<int>("decaysTo")),
23  _name(iPSet.getParameter<std::string>("name"))
24 {
25  dbe = 0;
27 
28  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
29 }
30 
32 
34 {
35  if(dbe){
37  std::string folderName = "Generator/DrellYan";
38  folderName+=_name;
39  dbe->setCurrentFolder(folderName.c_str());
40 
41  // Number of analyzed events
42  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
43 
44  //Kinematics
45  Zmass = dbe->book1D("Zmass","inv. Mass Z", 70 ,0,140);
46  ZmassPeak = dbe->book1D("ZmassPeak","inv. Mass Z", 80 ,80 ,100);
47  Zpt = dbe->book1D("Zpt","Z pt",100,0,200);
48  ZptLog = dbe->book1D("ZptLog","log(Z pt)",100,0.,5.);
49  Zrap = dbe->book1D("Zrap", "Z y", 100, -5, 5);
50  Zdaughters = dbe->book1D("Zdaughters", "Z daughters", 60, -30, 30);
51 
52  dilep_mass = dbe->book1D("dilep_mass","inv. Mass dilepton", 70 ,0,140);
53  dilep_massPeak = dbe->book1D("dilep_massPeak","inv. Mass dilepton", 80 ,80 ,100);
54  dilep_pt = dbe->book1D("dilep_pt","dilepton pt",100,0,200);
55  dilep_ptLog = dbe->book1D("dilep_ptLog","log(dilepton pt)",100,0.,5.);
56  dilep_rap = dbe->book1D("dilep_rap", "dilepton y", 100, -5, 5);
57 
58  gamma_energy = dbe->book1D("gamma_energy", "photon energy in Z rest frame", 200, 0., 100.);
59  cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton", "cos_theta_gamma_lepton in Z rest frame", 200, -1, 1);
60 
61  leadpt = dbe->book1D("leadpt","leading lepton pt", 200, 0., 200.);
62  secpt = dbe->book1D("secpt","second lepton pt", 200, 0., 200.);
63  leadeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);
64  seceta = dbe->book1D("seceta","second lepton eta", 100, -5., 5.);
65 
66  }
67  return;
68 }
69 
71 void DrellYanValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
72 {
74  iSetup.getData( fPDGTable );
75  return;
76 }
77 void DrellYanValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
79 {
80 
81  // we *DO NOT* rely on a Z entry in the particle listings!
82 
85  iEvent.getByToken(hepmcCollectionToken_, evt);
86 
87  //Get EVENT
88  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
89 
90  double weight = wmanager_.weight(iEvent);
91 
92  //std::cout << "weight: " << weight << std::endl;
93 
94  nEvt->Fill(0.5,weight);
95 
96  std::vector<const HepMC::GenParticle*> allproducts;
97 
98  //requires status 1 for leptons and neutrinos (except tau)
99  int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? 1 : 3;
100 
101  bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;
102 
103  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
104  if (vetotau) {
105  if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
106  }
107  if((*iter)->status()==requiredstatus) {
108  if(abs((*iter)->pdg_id())==_flavor)
109  allproducts.push_back(*iter);
110  }
111  }
112 
113  //nothing to do if we don't have 2 particles
114  if (allproducts.size() < 2) return;
115 
116  //sort them in pt
117  std::sort(allproducts.begin(), allproducts.end(), HepMCValidationHelper::sortByPt);
118 
119  //get the first element and the first following element with opposite charge
120  std::vector<const HepMC::GenParticle*> products;
121  products.push_back(allproducts.front());
122  const HepPDT::ParticleData* PData1 = fPDGTable->particle(HepPDT::ParticleID(allproducts.front()->pdg_id()));
123  double charge1 = PData1->charge();
124  for (unsigned int i = 1; i < allproducts.size(); ++i ){
125  const HepPDT::ParticleData* PData2 = fPDGTable->particle(HepPDT::ParticleID(allproducts[i]->pdg_id()));
126  double charge2 = PData2->charge();
127  if (charge1*charge2 < 0) products.push_back(allproducts[i]);
128  }
129 
130  //if we did not find any opposite charge pair there is nothing to do
131  if (products.size() < 2) return;
132 
133  assert(products[0]->momentum().perp() >= products[1]->momentum().perp());
134 
135  //leading lepton with pt > 20.
136  if (products[0]->momentum().perp() < 20.) return;
137 
138  //assemble FourMomenta
139  TLorentzVector lep1(products[0]->momentum().x(), products[0]->momentum().y(), products[0]->momentum().z(), products[0]->momentum().t());
140  TLorentzVector lep2(products[1]->momentum().x(), products[1]->momentum().y(), products[1]->momentum().z(), products[1]->momentum().t());
141  TLorentzVector dilepton_mom = lep1 + lep2;
142  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
143 
144  //mass > 60.
145  if (dilepton_mom.M() < 60.) return;
146 
147  //find possible qed fsr photons
148  std::vector<const HepMC::GenParticle*> fsrphotons;
149  HepMCValidationHelper::findFSRPhotons(products, myGenEvent, 0.1, fsrphotons);
150 
151  Zdaughters->Fill(products[0]->pdg_id(),weight);
152  Zdaughters->Fill(products[1]->pdg_id(),weight);
153 
154  std::vector<TLorentzVector> gammasMomenta;
155  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
156  TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t());
157  dilepton_andphoton_mom += phomom;
158  Zdaughters->Fill(fsrphotons[ipho]->pdg_id(),weight);
159  gammasMomenta.push_back(phomom);
160  }
161  //Fill Z histograms
162  Zmass->Fill(dilepton_andphoton_mom.M(),weight);
163  ZmassPeak->Fill(dilepton_andphoton_mom.M(),weight);
164  Zpt->Fill(dilepton_andphoton_mom.Pt(),weight);
165  ZptLog->Fill(log10(dilepton_andphoton_mom.Pt()),weight);
166  Zrap->Fill(dilepton_andphoton_mom.Rapidity(),weight);
167 
168  //Fill dilepton histograms
169  dilep_mass->Fill(dilepton_mom.M(),weight);
170  dilep_massPeak->Fill(dilepton_mom.M(),weight);
171  dilep_pt->Fill(dilepton_mom.Pt(),weight);
172  dilep_ptLog->Fill(log10(dilepton_mom.Pt()),weight);
173  dilep_rap->Fill(dilepton_mom.Rapidity(),weight);
174 
175  //Fill lepton histograms
176  leadpt->Fill(lep1.Pt(),weight);
177  secpt->Fill(lep2.Pt(),weight);
178  leadeta->Fill(lep1.Eta(),weight);
179  seceta->Fill(lep2.Eta(),weight);
180 
181  //boost everything in the Z frame
182  TVector3 boost = dilepton_andphoton_mom.BoostVector();
183  boost*=-1.;
184  lep1.Boost(boost);
185  lep2.Boost(boost);
186  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
187  gammasMomenta[ipho].Boost(boost);
188  }
189  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
190 
191  //fill gamma histograms
192  if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
193  gamma_energy->Fill(gammasMomenta.front().E(),weight);
194  double dphi = lep1.DeltaR(gammasMomenta.front()) < lep2.DeltaR(gammasMomenta.front()) ?
195  lep1.DeltaPhi(gammasMomenta.front()) : lep2.DeltaPhi(gammasMomenta.front());
196  cos_theta_gamma_lepton->Fill(cos(dphi),weight);
197  }
198 
199 }//analyze
MonitorElement * secpt
MonitorElement * gamma_energy
int i
Definition: DBlmapReader.cc:9
MonitorElement * dilep_rap
virtual void analyze(const edm::Event &, const edm::EventSetup &)
void findFSRPhotons(const std::vector< const HepMC::GenParticle * > &leptons, const std::vector< const HepMC::GenParticle * > &all, double deltaR, std::vector< const HepMC::GenParticle * > &photons)
MonitorElement * Zmass
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:873
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
edm::InputTag hepmcCollection_
MonitorElement * cos_theta_gamma_lepton
MonitorElement * nEvt
MonitorElement * leadeta
bool sortByPt(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
virtual void endRun(const edm::Run &, const edm::EventSetup &)
MonitorElement * Zrap
MonitorElement * leadpt
float float float z
void getData(T &iHolder) const
Definition: EventSetup.h:67
void Fill(long long x)
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
int _flavor
decay flavor
int iEvent
Definition: GenABIO.cc:243
MonitorElement * dilep_massPeak
ESProducts< T, S > products(const T &i1, const S &i2)
Definition: ESProducts.h:189
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
DQMStore * dbe
ME&#39;s &quot;container&quot;.
WeightManager wmanager_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * Zdaughters
MonitorElement * dilep_pt
HepPDT::ParticleData ParticleData
MonitorElement * ZptLog
virtual void beginJob()
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
std::string _name
decay flavor name
MonitorElement * seceta
T perp() const
Magnitude of transverse component.
MonitorElement * dilep_mass
MonitorElement * ZmassPeak
MonitorElement * dilep_ptLog
Definition: DDAxes.h:10
int weight
Definition: histoStyle.py:50
double weight(const edm::Event &)
DrellYanValidation(const edm::ParameterSet &)
MonitorElement * Zpt
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:585
Definition: Run.h:41