CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
WValidation.cc
Go to the documentation of this file.
1 /*class WValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  *
5  */
6 
9 #include "TLorentzVector.h"
10 
11 #include "CLHEP/Units/defs.h"
12 #include "CLHEP/Units/PhysicalConstants.h"
13 
16 using namespace edm;
17 
19  wmanager_(iPSet,consumesCollector()),
20  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
21  _flavor(iPSet.getParameter<int>("decaysTo")),
22  _name(iPSet.getParameter<std::string>("name"))
23 {
24 
25  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
26 }
27 
29 
32  std::string folderName = "Generator/W";
33  folderName+=_name;
34  DQMHelper dqm(&i); i.setCurrentFolder(folderName.c_str());
35 
36  // Number of analyzed events
37  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1.,"bin","Number of Events");
38 
39  //Kinematics
40  Wmass = dqm.book1dHisto("Wmass","inv. Mass W", 70 ,0,140,"M_{T}^{W} (GeV)","Number of Events");
41  WmassPeak = dqm.book1dHisto("WmassPeak","inv. Mass W", 80 ,80 ,100,"M_{T}^{W} (GeV)","Number of Events");
42  Wpt = dqm.book1dHisto("Wpt","W pt",100,0,200,"P_{T}^{W} (GeV)","Number of Events");
43  WptLog = dqm.book1dHisto("WptLog","log(W pt)",100,0.,5.,"Log_{10}(P_{T}^{W}) (GeV)","Number of Events");
44  Wrap = dqm.book1dHisto("Wrap", "W y", 100, -5, 5,"Y^{W}","Number of Events");
45  Wdaughters = dqm.book1dHisto("Wdaughters", "W daughters", 60, -30, 30,"W daughters (PDG ID)","Number of Events");
46 
47  lepmet_mT = dqm.book1dHisto("lepmet_mT","lepton-met transverse mass", 70 ,0,140,"M_{T}^{Lepton_{T}+E_{T}^{Miss}} (GeV)","Number of Events");
48  lepmet_mTPeak = dqm.book1dHisto("lepmet_mTPeak","lepton-met transverse mass", 80 ,80 ,100,"M_{T}^{Lepton_{T}+E_{T}^{Miss}} (GeV)","Number of Events");
49  lepmet_pt = dqm.book1dHisto("lepmet_pt","lepton-met",100,0,200,"P_{T}^{Lepton_{T}+E_{T}^{Miss}} (GeV)","Number of Events");
50  lepmet_ptLog = dqm.book1dHisto("lepmet_ptLog","log(lepton-met pt)",100,0.,5.,"log_{10}(P_{T}^{Lepton_{T}+E_{T}^{Miss}}) (Log_{10}(GeV))","Number of Events");
51 
52  gamma_energy = dqm.book1dHisto("gamma_energy", "photon energy in W rest frame", 200, 0., 100.,"E_{#gamma}^{W rest-frame}","Number of Events");
53  cos_theta_gamma_lepton = dqm.book1dHisto("cos_theta_gamma_lepton", "cos_theta_gamma_lepton in W rest frame", 200, -1, 1,"cos(#theta_{#gamma-lepton}^{W rest-frame})","Number of Events");
54 
55  leppt = dqm.book1dHisto("leadpt","lepton pt", 200, 0., 200.,"P_{t}^{Lead-Lepton} (GeV)","Number of Events");
56  met = dqm.book1dHisto("met","met", 200, 0., 200.,"E_{T}^{Miss} (GeV)","Number of Events");
57  lepeta = dqm.book1dHisto("leadeta","leading lepton eta", 100, -5., 5.,"#eta^{Lead-Lepton}","Number of Events");
58 
59  return;
60 }
61 
63  c.getData( fPDGTable );
64 }
65 
67 {
68 
69  // we *DO NOT* rely on a Z entry in the particle listings!
70 
73  iEvent.getByToken(hepmcCollectionToken_, evt);
74 
75  //Get EVENT
76  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
77 
78  double weight = wmanager_.weight(iEvent);
79 
80  nEvt->Fill(0.5,weight);
81 
82  std::vector<const HepMC::GenParticle*> allleptons;
83  std::vector<const HepMC::GenParticle*> allneutrinos;
84 
85  //requires status 1 for leptons and neutrinos (except tau)
86  int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) ==13 ) ? 1 : 3;
87 
88  bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;
89 
90  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
91  if (vetotau) {
92  if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
93  }
94  if((*iter)->status()==requiredstatus) {
95  //@todo: improve this selection
96  if((*iter)->pdg_id()==_flavor)
97  allleptons.push_back(*iter);
98  else if (abs((*iter)->pdg_id()) == abs(_flavor)+1)
99  allneutrinos.push_back(*iter);
100  }
101  }
102 
103  //nothing to do if we don't have 2 particles
104  if (allleptons.empty() || allneutrinos.empty()) return;
105 
106  //sort them in pt
107  std::sort(allleptons.begin(), allleptons.end(), HepMCValidationHelper::sortByPt);
108  std::sort(allneutrinos.begin(), allneutrinos.end(), HepMCValidationHelper::sortByPt);
109 
110  //get the first lepton and the first neutrino, and check that one is particle one is antiparticle (product of pdgids < 0)
111  std::vector<const HepMC::GenParticle*> products;
112  if (allleptons.front()->pdg_id() * allneutrinos.front()->pdg_id() > 0) return;
113 
114  //require at least 20 GeV on the lepton
115  if (allleptons.front()->momentum().perp() < 20. || allneutrinos.front()->momentum().perp() < 20. ) return;
116 
117  //find possible qed fsr photons
118  std::vector<const HepMC::GenParticle*> selectedLepton;
119  selectedLepton.push_back(allleptons.front());
120  std::vector<const HepMC::GenParticle*> fsrphotons;
121  HepMCValidationHelper::findFSRPhotons(selectedLepton, myGenEvent, 0.1, fsrphotons);
122 
123  Wdaughters->Fill(allleptons.front()->pdg_id(),weight);
124  Wdaughters->Fill(allneutrinos.front()->pdg_id(),weight);
125 
126  //assemble FourMomenta
127  TLorentzVector lep1(allleptons[0]->momentum().x(), allleptons[0]->momentum().y(), allleptons[0]->momentum().z(), allleptons[0]->momentum().t());
128  TLorentzVector lep2(allneutrinos[0]->momentum().x(), allneutrinos[0]->momentum().y(), allneutrinos[0]->momentum().z(), allneutrinos[0]->momentum().t());
129  TLorentzVector dilepton_mom = lep1 + lep2;
130  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
131  std::vector<TLorentzVector> gammasMomenta;
132  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
133  TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t());
134  dilepton_andphoton_mom += phomom;
135  Wdaughters->Fill(fsrphotons[ipho]->pdg_id(),weight);
136  gammasMomenta.push_back(phomom);
137  }
138  //Fill "true" W histograms
139  Wmass->Fill(dilepton_andphoton_mom.M(),weight);
140  WmassPeak->Fill(dilepton_andphoton_mom.M(),weight);
141  Wpt->Fill(dilepton_andphoton_mom.Pt(),weight);
142  WptLog->Fill(log10(dilepton_andphoton_mom.Pt()),weight);
143  Wrap->Fill(dilepton_andphoton_mom.Rapidity(),weight);
144 
145  TLorentzVector met_mom = HepMCValidationHelper::genMet(myGenEvent, -3., 3.);
146  TLorentzVector lep1T(lep1.Px(), lep1.Py(), 0., lep1.Et());
147  TLorentzVector lepmet_mom = lep1T + met_mom;
148  //Fill lepmet histograms
149  lepmet_mT->Fill(lepmet_mom.M(),weight);
150  lepmet_mTPeak->Fill(lepmet_mom.M(),weight);
151  lepmet_pt->Fill(lepmet_mom.Pt(),weight);
152  lepmet_ptLog->Fill(log10(lepmet_mom.Pt()),weight);
153 
154  //Fill lepton histograms
155  leppt->Fill(lep1.Pt(),weight);
156  lepeta->Fill(lep1.Eta(),weight);
157  met->Fill(met_mom.Pt(),weight);
158 
159  //boost everything in the W frame
160  TVector3 boost = dilepton_andphoton_mom.BoostVector();
161  boost*=-1.;
162  lep1.Boost(boost);
163  lep2.Boost(boost);
164  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
165  gammasMomenta[ipho].Boost(boost);
166  }
167  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
168 
169  //fill gamma histograms
170  if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
171  gamma_energy->Fill(gammasMomenta.front().E(),weight);
172  double dphi = lep1.DeltaR(gammasMomenta.front());
173  cos_theta_gamma_lepton->Fill(cos(dphi),weight);
174  }
175 
176 
177 }//analyze
WValidation(const edm::ParameterSet &)
Definition: WValidation.cc:18
edm::InputTag hepmcCollection_
Definition: WValidation.h:46
int i
Definition: DBlmapReader.cc:9
MonitorElement * leppt
Definition: WValidation.h:54
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 * Wrap
Definition: WValidation.h:52
MonitorElement * nEvt
Definition: WValidation.h:51
WeightManager wmanager_
Definition: WValidation.h:45
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
std::string _name
decay flavor name
Definition: WValidation.h:60
bool sortByPt(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
MonitorElement * lepeta
Definition: WValidation.h:54
virtual ~WValidation()
Definition: WValidation.cc:28
MonitorElement * book1dHisto(std::string name, std::string title, int n, double xmin, double xmax, std::string xaxis, std::string yaxis)
Definition: DQMHelper.cc:12
void getData(T &iHolder) const
Definition: EventSetup.h:79
MonitorElement * lepmet_mT
Definition: WValidation.h:53
void Fill(long long x)
MonitorElement * lepmet_pt
Definition: WValidation.h:53
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
Definition: WValidation.h:62
ESProducts< T, S > products(const T &i1, const S &i2)
Definition: ESProducts.h:189
tuple lep1
print &#39;MRbb(1b)&#39;,event.mr_bb
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * gamma_energy
Definition: WValidation.h:55
virtual void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
Definition: WValidation.cc:30
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
Definition: WValidation.h:49
virtual void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
Definition: WValidation.cc:62
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
MonitorElement * WptLog
Definition: WValidation.h:52
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
MonitorElement * lepmet_mTPeak
Definition: WValidation.h:53
MonitorElement * WmassPeak
Definition: WValidation.h:52
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: WValidation.cc:66
MonitorElement * met
Definition: WValidation.h:54
MonitorElement * Wpt
Definition: WValidation.h:52
MonitorElement * lepmet_ptLog
Definition: WValidation.h:53
MonitorElement * Wdaughters
Definition: WValidation.h:52
double weight(const edm::Event &)
MonitorElement * Wmass
Definition: WValidation.h:52
MonitorElement * cos_theta_gamma_lepton
Definition: WValidation.h:55
Definition: Run.h:43
int _flavor
decay flavor
Definition: WValidation.h:58