CMS 3D CMS Logo

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 
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  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
25  fPDGTableToken = esConsumes<edm::Transition::BeginRun>();
26 }
27 
29 
31  fPDGTable = c.getHandle(fPDGTableToken);
32 }
33 
36  std::string folderName = "Generator/DrellYan";
37  folderName += _name;
38  DQMHelper dqm(&i);
39  i.setCurrentFolder(folderName);
40 
41  // Number of analyzed events
42  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "bin", "Number of Events");
43 
44  //Kinematics
45  Zmass = dqm.book1dHisto("Zmass", "inv. Mass Z", 70, 0, 140, "M_{Z} (GeV)", "Number of Events");
46  ZmassPeak = dqm.book1dHisto("ZmassPeak", "inv. Mass Z", 80, 80, 100, "M_{Z} (GeV)", "Number of Events");
47  Zpt = dqm.book1dHisto("Zpt", "Z pt", 100, 0, 200, "P_{t}^{Z} (GeV)", "Number of Events");
48  ZptLog =
49  dqm.book1dHisto("ZptLog", "log(Z pt)", 100, 0., 5., "log_{10}(P_{t}^{Z}) (log_{10}(GeV))", "Number of Events");
50  Zrap = dqm.book1dHisto("Zrap", "Z y", 100, -5, 5, "Y_{Z}", "Number of Events");
51  Zdaughters = dqm.book1dHisto("Zdaughters", "Z daughters", 60, -30, 30, "Z daughters (PDG ID)", "Number of Events");
52 
53  dilep_mass = dqm.book1dHisto("dilep_mass", "inv. Mass dilepton", 70, 0, 140, "M_{ll} (GeV)", "Number of Events");
55  dqm.book1dHisto("dilep_massPeak", "inv. Mass dilepton", 80, 80, 100, "M_{ll} (GeV)", "Number of Events");
56  dilep_pt = dqm.book1dHisto("dilep_pt", "dilepton pt", 100, 0, 200, "P_{t}^{ll} (GeV)", "Number of Events");
57  dilep_ptLog = dqm.book1dHisto(
58  "dilep_ptLog", "log(dilepton pt)", 100, 0., 5., "log_{10}(P_{t}^{ll}) (log_{10}(GeV))", "Number of Events");
59  dilep_rap = dqm.book1dHisto("dilep_rap", "dilepton y", 100, -5, 5, "Y_{ll}", "Number of Events");
60 
61  gamma_energy = dqm.book1dHisto("gamma_energy",
62  "photon energy in Z rest frame",
63  200,
64  0.,
65  100.,
66  "E_{#gamma}^{Z rest-frame} (GeV)",
67  "Number of Events");
68  cos_theta_gamma_lepton = dqm.book1dHisto("cos_theta_gamma_lepton",
69  "cos_theta_gamma_lepton in Z rest frame",
70  200,
71  -1,
72  1,
73  "cos(#theta_{#gamma-lepton}^{Z rest-frame})",
74  "Number of Events");
75 
76  leadpt = dqm.book1dHisto("leadpt", "leading lepton pt", 200, 0., 200., "P_{t}^{1st-lepton}", "Number of Events");
77  secpt = dqm.book1dHisto("secpt", "second lepton pt", 200, 0., 200., "P_{t}^{2nd-lepton}", "Number of Events");
78  leadeta = dqm.book1dHisto("leadeta", "leading lepton eta", 100, -5., 5., "#eta^{1st-lepton}", "Number of Events");
79  seceta = dqm.book1dHisto("seceta", "second lepton eta", 100, -5., 5., "#eta^{2nd-lepton}", "Number of Events");
80 
81  return;
82 }
83 
85  // we *DO NOT* rely on a Z entry in the particle listings!
86 
89  iEvent.getByToken(hepmcCollectionToken_, evt);
90 
91  //Get EVENT
92  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
93 
94  double weight = wmanager_.weight(iEvent);
95 
96  //std::cout << "weight: " << weight << std::endl;
97 
98  nEvt->Fill(0.5, weight);
99 
100  std::vector<const HepMC::GenParticle*> allproducts;
101 
102  //requires status 1 for leptons and neutrinos (except tau)
103  int requiredstatus = (std::abs(_flavor) == 11 || std::abs(_flavor) == 12 || std::abs(_flavor) == 13 ||
104  std::abs(_flavor) == 14 || std::abs(_flavor) == 16)
105  ? 1
106  : 3;
107 
108  bool vetotau =
109  true; //(std::abs(_flavor) == 11 || std::abs(_flavor) == 12 || std::abs(_flavor) ==13 || std::abs(_flavor) ==14 || std::abs(_flavor) ==16) ? true : false;
110 
111  for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
112  iter != myGenEvent->particles_end();
113  ++iter) {
114  if (vetotau) {
115  if ((*iter)->status() == 3 && std::abs((*iter)->pdg_id()) == 15)
116  return;
117  }
118  if ((*iter)->status() == requiredstatus) {
119  if (std::abs((*iter)->pdg_id()) == _flavor)
120  allproducts.push_back(*iter);
121  }
122  }
123 
124  //nothing to do if we don't have 2 particles
125  if (allproducts.size() < 2)
126  return;
127 
128  //sort them in pt
129  std::sort(allproducts.begin(), allproducts.end(), HepMCValidationHelper::sortByPt);
130 
131  //get the first element and the first following element with opposite charge
132  std::vector<const HepMC::GenParticle*> products;
133  products.push_back(allproducts.front());
134  const HepPDT::ParticleData* PData1 = fPDGTable->particle(HepPDT::ParticleID(allproducts.front()->pdg_id()));
135  double charge1 = PData1->charge();
136  for (unsigned int i = 1; i < allproducts.size(); ++i) {
137  const HepPDT::ParticleData* PData2 = fPDGTable->particle(HepPDT::ParticleID(allproducts[i]->pdg_id()));
138  double charge2 = PData2->charge();
139  if (charge1 * charge2 < 0)
140  products.push_back(allproducts[i]);
141  }
142 
143  //if we did not find any opposite charge pair there is nothing to do
144  if (products.size() < 2)
145  return;
146 
147  assert(products[0]->momentum().perp() >= products[1]->momentum().perp());
148 
149  //leading lepton with pt > 20.
150  if (products[0]->momentum().perp() < 20.)
151  return;
152 
153  //assemble FourMomenta
154  TLorentzVector lep1(products[0]->momentum().x(),
155  products[0]->momentum().y(),
156  products[0]->momentum().z(),
157  products[0]->momentum().t());
158  TLorentzVector lep2(products[1]->momentum().x(),
159  products[1]->momentum().y(),
160  products[1]->momentum().z(),
161  products[1]->momentum().t());
162  TLorentzVector dilepton_mom = lep1 + lep2;
163  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
164 
165  //mass > 60.
166  if (dilepton_mom.M() < 60.)
167  return;
168 
169  //find possible qed fsr photons
170  std::vector<const HepMC::GenParticle*> fsrphotons;
171  HepMCValidationHelper::findFSRPhotons(products, myGenEvent, 0.1, fsrphotons);
172 
175 
176  std::vector<TLorentzVector> gammasMomenta;
177  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho) {
178  TLorentzVector phomom(fsrphotons[ipho]->momentum().x(),
179  fsrphotons[ipho]->momentum().y(),
180  fsrphotons[ipho]->momentum().z(),
181  fsrphotons[ipho]->momentum().t());
182  dilepton_andphoton_mom += phomom;
183  Zdaughters->Fill(fsrphotons[ipho]->pdg_id(), weight);
184  gammasMomenta.push_back(phomom);
185  }
186  //Fill Z histograms
187  Zmass->Fill(dilepton_andphoton_mom.M(), weight);
188  ZmassPeak->Fill(dilepton_andphoton_mom.M(), weight);
189  Zpt->Fill(dilepton_andphoton_mom.Pt(), weight);
190  ZptLog->Fill(log10(dilepton_andphoton_mom.Pt()), weight);
191  Zrap->Fill(dilepton_andphoton_mom.Rapidity(), weight);
192 
193  //Fill dilepton histograms
194  dilep_mass->Fill(dilepton_mom.M(), weight);
195  dilep_massPeak->Fill(dilepton_mom.M(), weight);
196  dilep_pt->Fill(dilepton_mom.Pt(), weight);
197  dilep_ptLog->Fill(log10(dilepton_mom.Pt()), weight);
198  dilep_rap->Fill(dilepton_mom.Rapidity(), weight);
199 
200  //Fill lepton histograms
201  leadpt->Fill(lep1.Pt(), weight);
202  secpt->Fill(lep2.Pt(), weight);
203  leadeta->Fill(lep1.Eta(), weight);
204  seceta->Fill(lep2.Eta(), weight);
205 
206  //boost everything in the Z frame
207  TVector3 boost = dilepton_andphoton_mom.BoostVector();
208  boost *= -1.;
209  lep1.Boost(boost);
210  lep2.Boost(boost);
211  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho) {
212  gammasMomenta[ipho].Boost(boost);
213  }
214  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
215 
216  //fill gamma histograms
217  if (!gammasMomenta.empty() && dilepton_andphoton_mom.M() > 50.) {
218  gamma_energy->Fill(gammasMomenta.front().E(), weight);
219  double dphi = lep1.DeltaR(gammasMomenta.front()) < lep2.DeltaR(gammasMomenta.front())
220  ? lep1.DeltaPhi(gammasMomenta.front())
221  : lep2.DeltaPhi(gammasMomenta.front());
223  }
224 
225 } //analyze
MonitorElement * secpt
MonitorElement * gamma_energy
MonitorElement * dilep_rap
Definition: CLHEP.h:16
MonitorElement * Zmass
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > fPDGTableToken
void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
ESProducts< std::remove_reference_t< TArgs >... > products(TArgs &&... args)
Definition: ESProducts.h:128
edm::InputTag hepmcCollection_
MonitorElement * cos_theta_gamma_lepton
MonitorElement * nEvt
MonitorElement * leadeta
bool sortByPt(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
~DrellYanValidation() override
MonitorElement * Zrap
Definition: weight.py:1
assert(be >=bs)
MonitorElement * leadpt
void findFSRPhotons(const std::vector< const HepMC::GenParticle *> &leptons, const std::vector< const HepMC::GenParticle *> &all, double deltaR, std::vector< const HepMC::GenParticle *> &photons)
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
void Fill(long long x)
int _flavor
decay flavor
int iEvent
Definition: GenABIO.cc:224
MonitorElement * dilep_massPeak
void analyze(edm::Event const &, edm::EventSetup const &) override
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
WeightManager wmanager_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * Zdaughters
MonitorElement * dilep_pt
HepPDT::ParticleData ParticleData
T perp() const
Magnitude of transverse component.
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
MonitorElement * ZptLog
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
std::string _name
decay flavor name
MonitorElement * seceta
MonitorElement * dilep_mass
HLT enums.
MonitorElement * ZmassPeak
MonitorElement * dilep_ptLog
lep1
print &#39;MRbb(1b)&#39;,event.mr_bb
double weight(const edm::Event &)
DrellYanValidation(const edm::ParameterSet &)
MonitorElement * Zpt
Definition: DQMStore.h:18
Definition: Run.h:45