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  * $Date: 2011/01/24 18:27:40 $
6  * $Revision: 1.4 $
7  */
8 
11 #include "TLorentzVector.h"
12 
13 #include "CLHEP/Units/defs.h"
14 #include "CLHEP/Units/PhysicalConstants.h"
15 
17 
18 using namespace edm;
19 
21  _wmanager(iPSet),
22  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
23  _flavor(iPSet.getParameter<int>("decaysTo")),
24  _name(iPSet.getParameter<std::string>("name"))
25 {
26  dbe = 0;
28 }
29 
31 
33 {
34  if(dbe){
36  std::string folderName = "Generator/DrellYan";
37  folderName+=_name;
38  dbe->setCurrentFolder(folderName.c_str());
39 
40  // Number of analyzed events
41  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
42 
43  //Kinematics
44  Zmass = dbe->book1D("Zmass","inv. Mass Z", 70 ,0,140);
45  ZmassPeak = dbe->book1D("ZmassPeak","inv. Mass Z", 80 ,80 ,100);
46  Zpt = dbe->book1D("Zpt","Z pt",100,0,200);
47  ZptLog = dbe->book1D("ZptLog","log(Z pt)",100,0.,5.);
48  Zrap = dbe->book1D("Zrap", "Z y", 100, -5, 5);
49  Zdaughters = dbe->book1D("Zdaughters", "Z daughters", 60, -30, 30);
50 
51  dilep_mass = dbe->book1D("dilep_mass","inv. Mass dilepton", 70 ,0,140);
52  dilep_massPeak = dbe->book1D("dilep_massPeak","inv. Mass dilepton", 80 ,80 ,100);
53  dilep_pt = dbe->book1D("dilep_pt","dilepton pt",100,0,200);
54  dilep_ptLog = dbe->book1D("dilep_ptLog","log(dilepton pt)",100,0.,5.);
55  dilep_rap = dbe->book1D("dilep_rap", "dilepton y", 100, -5, 5);
56 
57  gamma_energy = dbe->book1D("gamma_energy", "photon energy in Z rest frame", 200, 0., 100.);
58  cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton", "cos_theta_gamma_lepton in Z rest frame", 200, -1, 1);
59 
60  leadpt = dbe->book1D("leadpt","leading lepton pt", 200, 0., 200.);
61  secpt = dbe->book1D("secpt","second lepton pt", 200, 0., 200.);
62  leadeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);
63  seceta = dbe->book1D("seceta","second lepton eta", 100, -5., 5.);
64 
65  }
66  return;
67 }
68 
70 void DrellYanValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
71 {
73  iSetup.getData( fPDGTable );
74  return;
75 }
76 void DrellYanValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
78 {
79 
80  // we *DO NOT* rely on a Z entry in the particle listings!
81 
84  iEvent.getByLabel(hepmcCollection_, evt);
85 
86  //Get EVENT
87  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
88 
89  double weight = _wmanager.weight(iEvent);
90 
91  //std::cout << "weight: " << weight << std::endl;
92 
93  nEvt->Fill(0.5,weight);
94 
95  std::vector<const HepMC::GenParticle*> allproducts;
96 
97  //requires status 1 for leptons and neutrinos (except tau)
98  int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? 1 : 3;
99 
100  bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;
101 
102  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
103  if (vetotau) {
104  if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
105  }
106  if((*iter)->status()==requiredstatus) {
107  if(abs((*iter)->pdg_id())==_flavor)
108  allproducts.push_back(*iter);
109  }
110  }
111 
112  //nothing to do if we don't have 2 particles
113  if (allproducts.size() < 2) return;
114 
115  //sort them in pt
116  std::sort(allproducts.begin(), allproducts.end(), HepMCValidationHelper::sortByPt);
117 
118  //get the first element and the first following element with opposite charge
119  std::vector<const HepMC::GenParticle*> products;
120  products.push_back(allproducts.front());
121  const HepPDT::ParticleData* PData1 = fPDGTable->particle(HepPDT::ParticleID(allproducts.front()->pdg_id()));
122  double charge1 = PData1->charge();
123  for (unsigned int i = 1; i < allproducts.size(); ++i ){
124  const HepPDT::ParticleData* PData2 = fPDGTable->particle(HepPDT::ParticleID(allproducts[i]->pdg_id()));
125  double charge2 = PData2->charge();
126  if (charge1*charge2 < 0) products.push_back(allproducts[i]);
127  }
128 
129  //if we did not find any opposite charge pair there is nothing to do
130  if (products.size() < 2) return;
131 
132  assert(products[0]->momentum().perp() >= products[1]->momentum().perp());
133 
134  //leading lepton with pt > 20.
135  if (products[0]->momentum().perp() < 20.) return;
136 
137  //assemble FourMomenta
138  TLorentzVector lep1(products[0]->momentum().x(), products[0]->momentum().y(), products[0]->momentum().z(), products[0]->momentum().t());
139  TLorentzVector lep2(products[1]->momentum().x(), products[1]->momentum().y(), products[1]->momentum().z(), products[1]->momentum().t());
140  TLorentzVector dilepton_mom = lep1 + lep2;
141  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
142 
143  //mass > 60.
144  if (dilepton_mom.M() < 60.) return;
145 
146  //find possible qed fsr photons
147  std::vector<const HepMC::GenParticle*> fsrphotons;
148  HepMCValidationHelper::findFSRPhotons(products, myGenEvent, 0.1, fsrphotons);
149 
150  Zdaughters->Fill(products[0]->pdg_id(),weight);
151  Zdaughters->Fill(products[1]->pdg_id(),weight);
152 
153  std::vector<TLorentzVector> gammasMomenta;
154  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
155  TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t());
156  dilepton_andphoton_mom += phomom;
157  Zdaughters->Fill(fsrphotons[ipho]->pdg_id(),weight);
158  gammasMomenta.push_back(phomom);
159  }
160  //Fill Z histograms
161  Zmass->Fill(dilepton_andphoton_mom.M(),weight);
162  ZmassPeak->Fill(dilepton_andphoton_mom.M(),weight);
163  Zpt->Fill(dilepton_andphoton_mom.Pt(),weight);
164  ZptLog->Fill(log10(dilepton_andphoton_mom.Pt()),weight);
165  Zrap->Fill(dilepton_andphoton_mom.Rapidity(),weight);
166 
167  //Fill dilepton histograms
168  dilep_mass->Fill(dilepton_mom.M(),weight);
169  dilep_massPeak->Fill(dilepton_mom.M(),weight);
170  dilep_pt->Fill(dilepton_mom.Pt(),weight);
171  dilep_ptLog->Fill(log10(dilepton_mom.Pt()),weight);
172  dilep_rap->Fill(dilepton_mom.Rapidity(),weight);
173 
174  //Fill lepton histograms
175  leadpt->Fill(lep1.Pt(),weight);
176  secpt->Fill(lep2.Pt(),weight);
177  leadeta->Fill(lep1.Eta(),weight);
178  seceta->Fill(lep2.Eta(),weight);
179 
180  //boost everything in the Z frame
181  TVector3 boost = dilepton_andphoton_mom.BoostVector();
182  boost*=-1.;
183  lep1.Boost(boost);
184  lep2.Boost(boost);
185  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
186  gammasMomenta[ipho].Boost(boost);
187  }
188  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
189 
190  //fill gamma histograms
191  if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
192  gamma_energy->Fill(gammasMomenta.front().E(),weight);
193  double dphi = lep1.DeltaR(gammasMomenta.front()) < lep2.DeltaR(gammasMomenta.front()) ?
194  lep1.DeltaPhi(gammasMomenta.front()) : lep2.DeltaPhi(gammasMomenta.front());
195  cos_theta_gamma_lepton->Fill(cos(dphi),weight);
196  }
197 
198 }//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:717
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 &)
#define abs(x)
Definition: mlp_lapack.h:159
MonitorElement * Zrap
MonitorElement * leadpt
double double double z
void getData(T &iHolder) const
Definition: EventSetup.h:67
WeightManager _wmanager
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;.
MonitorElement * Zdaughters
MonitorElement * dilep_pt
HepPDT::ParticleData ParticleData
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
MonitorElement * ZptLog
virtual void beginJob()
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
std::string _name
decay flavor name
MonitorElement * seceta
T perp() const
Magnitude of transverse component.
MonitorElement * dilep_mass
MonitorElement * ZmassPeak
MonitorElement * dilep_ptLog
x
Definition: VDTMath.h:216
double weight(const edm::Event &)
DrellYanValidation(const edm::ParameterSet &)
MonitorElement * Zpt
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
Definition: Run.h:33