CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QcdPhotonsDQM.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author Michael B. Anderson, University of Wisconsin Madison
5  */
6 
8 
12 
17 
19 
21 
22 // Physics Objects
25 
26 // Vertex
28 
29 // For removing ECAL Spikes
34 
37 
38 // geometry
39 //#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
40 //#include "Geometry/Records/interface/IdealGeometryRecord.h"
41 //#include "Geometry/Records/interface/CaloGeometryRecord.h"
42 //#include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
43 //#include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
44 //#include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
45 //#include "Geometry/EcalAlgo/interface/EcalPreshowerGeometry.h"
46 
47 // Math stuff
50 
51 #include <vector>
52 
53 #include <string>
54 #include <cmath>
55 using namespace std;
56 using namespace edm;
57 using namespace reco;
58 
60  // Get parameters from configuration file
61  theTriggerPathToPass_ = parameters.getParameter<string>("triggerPathToPass");
62  thePlotTheseTriggersToo_ =
63  parameters.getParameter<vector<string> >("plotTheseTriggersToo");
64  theJetCollectionLabel_ = parameters.getParameter<InputTag>("jetCollection");
65  trigTagToken_ = consumes<edm::TriggerResults>(
66  parameters.getUntrackedParameter<edm::InputTag>("trigTag"));
67  thePhotonCollectionToken_ = consumes<reco::PhotonCollection>(
68  parameters.getParameter<InputTag>("photonCollection"));
69  theJetCollectionToken_ = consumes<edm::View<reco::Jet> >(
70  parameters.getParameter<InputTag>("jetCollection"));
71  theVertexCollectionToken_ = consumes<reco::VertexCollection>(
72  parameters.getParameter<InputTag>("vertexCollection"));
73  theMinJetPt_ = parameters.getParameter<double>("minJetPt");
74  theMinPhotonEt_ = parameters.getParameter<double>("minPhotonEt");
75  theRequirePhotonFound_ = parameters.getParameter<bool>("requirePhotonFound");
76  thePlotPhotonMaxEt_ = parameters.getParameter<double>("plotPhotonMaxEt");
77  thePlotPhotonMaxEta_ = parameters.getParameter<double>("plotPhotonMaxEta");
78  thePlotJetMaxEta_ = parameters.getParameter<double>("plotJetMaxEta");
79  theBarrelRecHitTag_ = parameters.getParameter<InputTag>("barrelRecHitTag");
80  theEndcapRecHitTag_ = parameters.getParameter<InputTag>("endcapRecHitTag");
81  theBarrelRecHitToken_ = consumes<EcalRecHitCollection>(
82  parameters.getParameter<InputTag>("barrelRecHitTag"));
83  theEndcapRecHitToken_ = consumes<EcalRecHitCollection>(
84  parameters.getParameter<InputTag>("endcapRecHitTag"));
85 
86  // coverity says...
87  h_deltaEt_photon_jet = 0;
88  h_deltaPhi_jet_jet2 = 0;
89  h_deltaPhi_photon_jet = 0;
90  h_deltaPhi_photon_jet2 = 0;
91  h_deltaR_jet_jet2 = 0;
92  h_deltaR_photon_jet2 = 0;
93  h_jet2_eta = 0;
94  h_jet2_pt = 0;
95  h_jet2_ptOverPhotonEt = 0;
96  h_jet_count = 0;
97  h_jet_eta = 0;
98  h_jet_pt = 0;
99  h_photon_count_bar = 0;
100  h_photon_count_end = 0;
101  h_photon_et = 0;
102  h_photon_et_beforeCuts = 0;
103  h_photon_et_jetco = 0;
104  h_photon_et_jetcs = 0;
105  h_photon_et_jetfo = 0;
106  h_photon_et_jetfs = 0;
107  h_photon_eta = 0;
108  h_triggers_passed = 0;
109 
110 }
111 
113 
114 
116  edm::Run const &, edm::EventSetup const & ){
117 
118  logTraceName = "QcdPhotonAnalyzer";
119 
120  LogTrace(logTraceName) << "Parameters initialization";
121 
122  ibooker.setCurrentFolder("Physics/QcdPhotons"); // Use folder with name of PAG
123 
124  std::stringstream aStringStream;
125  std::string aString;
126  aStringStream << theMinJetPt_;
127  aString = aStringStream.str();
128 
129  // Monitor of triggers passed
130  int numOfTriggersToMonitor = thePlotTheseTriggersToo_.size();
131  h_triggers_passed = ibooker.book1D("triggers_passed",
132  "Events passing these trigger paths", numOfTriggersToMonitor,
133  0, numOfTriggersToMonitor);
134  for (int i = 0; i < numOfTriggersToMonitor; i++) {
135  h_triggers_passed->setBinLabel(i + 1, thePlotTheseTriggersToo_[i]);
136  }
137 
138  // Keep the number of plots and number of bins to a minimum!
139  h_photon_et_beforeCuts = ibooker.book1D("photon_et_beforeCuts",
140  "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
141  h_photon_et = ibooker.book1D("photon_et",
142  "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
143  h_photon_eta = ibooker.book1D("photon_eta",
144  "#gamma with highest E_{T};#eta(#gamma)", 40,
145  -thePlotPhotonMaxEta_, thePlotPhotonMaxEta_);
146  h_photon_count_bar = ibooker.book1D("photon_count_bar",
147  "Number of #gamma's passing selection (Barrel);Number of #gamma's", 8,
148  -0.5, 7.5);
149  h_photon_count_end = ibooker.book1D("photon_count_end",
150  "Number of #gamma's passing selection (Endcap);Number of #gamma's", 8,
151  -0.5, 7.5);
152  h_jet_pt = ibooker.book1D("jet_pt",
153  "Jet with highest p_{T} (from " + theJetCollectionLabel_.label() +
154  ");p_{T}(1^{st} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
155  h_jet_eta = ibooker.book1D("jet_eta",
156  "Jet with highest p_{T} (from " + theJetCollectionLabel_.label() +
157  ");#eta(1^{st} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
158  h_deltaPhi_photon_jet = ibooker.book1D("deltaPhi_photon_jet",
159  "#Delta#phi between Highest E_{T} #gamma and jet;#Delta#phi(#gamma,1^{st} jet)",
160  20, 0, 3.1415926);
161  h_deltaPhi_jet_jet2 = ibooker.book1D("deltaPhi_jet_jet2",
162  "#Delta#phi between Highest E_{T} jet and 2^{nd} "
163  "jet;#Delta#phi(1^{st} jet,2^{nd} jet)",
164  20, 0, 3.1415926);
165  h_deltaEt_photon_jet = ibooker.book1D("deltaEt_photon_jet",
166  "(E_{T}(#gamma)-p_{T}(jet))/E_{T}(#gamma) when #Delta#phi(#gamma,1^{st} "
167  "jet) > 2.8;#DeltaE_{T}(#gamma,1^{st} jet)/E_{T}(#gamma)",
168  20, -1.0, 1.0);
169  h_jet_count = ibooker.book1D("jet_count",
170  "Number of " + theJetCollectionLabel_.label() + " (p_{T} > " + aString +
171  " GeV);Number of Jets", 8, -0.5, 7.5);
172  h_jet2_pt = ibooker.book1D("jet2_pt",
173  "Jet with 2^{nd} highest p_{T} (from " + theJetCollectionLabel_.label() +
174  ");p_{T}(2^{nd} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
175  h_jet2_eta = ibooker.book1D("jet2_eta",
176  "Jet with 2^{nd} highest p_{T} (from " + theJetCollectionLabel_.label() +
177  ");#eta(2^{nd} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
178  h_jet2_ptOverPhotonEt = ibooker.book1D("jet2_ptOverPhotonEt",
179  "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} Jet)/E_{T}(#gamma)",
180  20, 0.0, 4.0);
181  h_deltaPhi_photon_jet2 = ibooker.book1D("deltaPhi_photon_jet2",
182  "#Delta#phi between Highest E_{T} #gamma and 2^{nd} "
183  "highest jet;#Delta#phi(#gamma,2^{nd} jet)", 20, 0, 3.1415926);
184  h_deltaR_jet_jet2 = ibooker.book1D("deltaR_jet_jet2",
185  "#DeltaR between Highest Jet and 2^{nd} Highest;#DeltaR(1^{st} jet,2^{nd} jet)",
186  30, 0, 6.0);
187  h_deltaR_photon_jet2 = ibooker.book1D("deltaR_photon_jet2",
188  "#DeltaR between Highest E_{T} #gamma and 2^{nd} "
189  "jet;#DeltaR(#gamma, 2^{nd} jet)", 30, 0, 6.0);
190 
191  // Photon Et for different jet configurations
192  Float_t bins_et[] = {15, 20, 30, 50, 80};
193  int num_bins_et = 4;
194  h_photon_et_jetcs = ibooker.book1D("photon_et_jetcs",
195  "#gamma with highest E_{T} (#eta(jet)<1.45, "
196  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
197  h_photon_et_jetco = ibooker.book1D("photon_et_jetco",
198  "#gamma with highest E_{T} (#eta(jet)<1.45, "
199  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
200  h_photon_et_jetfs = ibooker.book1D("photon_et_jetfs",
201  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
202  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
203  h_photon_et_jetfo = ibooker.book1D("photon_et_jetfo",
204  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
205  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
206 
207  auto setSumw2 = [](MonitorElement* me) {
208  if (me->getTH1F()->GetSumw2N() == 0) {
209  me->getTH1F()->Sumw2();
210  }
211  };
212 
213  setSumw2(h_photon_et_jetcs);
214  setSumw2(h_photon_et_jetco);
215  setSumw2(h_photon_et_jetfs);
216  setSumw2(h_photon_et_jetfo);
217 }
218 
219 void QcdPhotonsDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
220  LogTrace(logTraceName) << "Analysis of event # ";
221 
223  // Did event pass HLT paths?
224  Handle<TriggerResults> HLTresults;
225  iEvent.getByToken(trigTagToken_, HLTresults);
226  if (!HLTresults.isValid()) {
227  // LogWarning("") << ">>> TRIGGER collection does not exist !!!";
228  return;
229  }
230  const edm::TriggerNames& trigNames = iEvent.triggerNames(*HLTresults);
231 
232  bool passed_HLT = false;
233 
234  // See if event passed trigger paths
235  // increment that bin in the trigger plot
236  for (unsigned int i = 0; i < thePlotTheseTriggersToo_.size(); i++) {
237  passed_HLT = false;
238  for (unsigned int ti = 0; (ti < trigNames.size()) && !passed_HLT; ++ti) {
239  size_t pos = trigNames.triggerName(ti).find(thePlotTheseTriggersToo_[i]);
240  if (pos == 0) passed_HLT = HLTresults->accept(ti);
241  }
242  if (passed_HLT) h_triggers_passed->Fill(i);
243  }
244 
245  // grab photons
246  Handle<PhotonCollection> photonCollection;
247  iEvent.getByToken(thePhotonCollectionToken_, photonCollection);
248 
249  // If photon collection is empty, exit
250  if (!photonCollection.isValid()) return;
251 
252  // Quit if the event did not pass the HLT path we care about
253  passed_HLT = false;
254  {
255  // bool found=false;
256  for (unsigned int ti = 0; ti < trigNames.size(); ++ti) {
257  size_t pos = trigNames.triggerName(ti).find(theTriggerPathToPass_);
258  if (pos == 0) {
259  passed_HLT = HLTresults->accept(ti);
260  // found=true;
261  break;
262  }
263  }
264 
265  // Assumption: reco photons are ordered by Et
266  for (PhotonCollection::const_iterator recoPhoton =
267  photonCollection->begin();
268  recoPhoton != photonCollection->end(); recoPhoton++) {
269  // stop looping over photons once we get to too low Et
270  if (recoPhoton->et() < theMinPhotonEt_) break;
271 
272  h_photon_et_beforeCuts->Fill(recoPhoton->et());
273  break; // leading photon only
274  }
275 
276  if (!passed_HLT) {
277  return;
278  }
279  }
280 
282 
283  // std::cout << "\tpassed main trigger (" << theTriggerPathToPass_ << ")" <<
284  // std::endl;
285 
287  // Does event have valid vertex?
288  // Get the primary event vertex
289  Handle<VertexCollection> vertexHandle;
290  iEvent.getByToken(theVertexCollectionToken_, vertexHandle);
291  VertexCollection vertexCollection = *(vertexHandle.product());
292  // double vtx_ndof = -1.0;
293  // double vtx_z = 0.0;
294  // bool vtx_isFake = true;
295  // if (vertexCollection.size()>0) {
296  // vtx_ndof = vertexCollection.begin()->ndof();
297  // vtx_z = vertexCollection.begin()->z();
298  // vtx_isFake = false;
299  //}
300  // if (vtx_isFake || fabs(vtx_z)>15 || vtx_ndof<4) return;
301 
302  int nvvertex = 0;
303  for (unsigned int i = 0; i < vertexCollection.size(); ++i) {
304  if (vertexCollection[i].isValid()) nvvertex++;
305  }
306  if (nvvertex == 0) return;
307 
309 
310  // std::cout << "\tpassed vertex selection" << std::endl;
311 
313  // Did the event pass certain L1 Technical Trigger bits?
314  // It's probably beam halo
315  // TODO: ADD code
317 
318  // For finding spikes
319  Handle<EcalRecHitCollection> EBReducedRecHits;
320  iEvent.getByToken(theBarrelRecHitToken_, EBReducedRecHits);
321  Handle<EcalRecHitCollection> EEReducedRecHits;
322  iEvent.getByToken(theEndcapRecHitToken_, EEReducedRecHits);
323  EcalClusterLazyTools lazyTool(iEvent, iSetup, theBarrelRecHitToken_,
324  theEndcapRecHitToken_);
325 
326  // Find the highest et "decent" photon
327  float photon_et = -9.0;
328  float photon_eta = -9.0;
329  float photon_phi = -9.0;
330  bool photon_passPhotonID = false;
331  bool found_lead_pho = false;
332  int photon_count_bar = 0;
333  int photon_count_end = 0;
334  // False Assumption: reco photons are ordered by Et
335  // find the photon with highest et
336  auto pho_maxet = std::max_element(
337  photonCollection->begin(), photonCollection->end(),
339  const PhotonCollection::value_type& b) { return a.et() < b.et(); });
340  if (pho_maxet != photonCollection->end() &&
341  pho_maxet->et() >= theMinPhotonEt_) {
342  /*
343  // Ignore ECAL Spikes
344  const reco::CaloClusterPtr seed = pho_maxet->superCluster()->seed();
345  DetId id = lazyTool.getMaximum(*seed).first; // Cluster shape variables
346  // float time = -999., outOfTimeChi2 = -999., chi2 = -999.; // UNUSED
347  int flags=-1, severity = -1;
348  const EcalRecHitCollection & rechits = ( pho_maxet->isEB() ?
349  *EBReducedRecHits : *EEReducedRecHits);
350  EcalRecHitCollection::const_iterator it = rechits.find( id );
351  if( it != rechits.end() ) {
352  // time = it->time(); // UNUSED
353  // outOfTimeChi2 = it->outOfTimeChi2(); // UNUSED
354  // chi2 = it->chi2(); // UNUSED
355  flags = it->recoFlag();
356 
357  edm::ESHandle<EcalSeverityLevelAlgo> sevlv;
358  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
359  severity = sevlv->severityLevel( id, rechits);
360  }
361  bool isNotSpike = ((pho_maxet->isEB() && (severity!=3 && severity!=4 ) &&
362  (flags != 2) ) || pho_maxet->isEE());
363  if (!isNotSpike) continue; // move on to next photon
364  // END of determining ECAL Spikes
365  */
366 
367  bool pho_current_passPhotonID = false;
368  bool pho_current_isEB = pho_maxet->isEB();
369  bool pho_current_isEE = pho_maxet->isEE();
370 
371  if (pho_current_isEB && (pho_maxet->sigmaIetaIeta() < 0.01 ||
372  pho_maxet->hadronicOverEm() < 0.05)) {
373  // Photon object in barrel passes photon ID
374  pho_current_passPhotonID = true;
375  photon_count_bar++;
376  } else if (pho_current_isEE && (pho_maxet->hadronicOverEm() < 0.05)) {
377  // Photon object in endcap passes photon ID
378  pho_current_passPhotonID = true;
379  photon_count_end++;
380  }
381 
382  if (!found_lead_pho) {
383  found_lead_pho = true;
384  photon_passPhotonID = pho_current_passPhotonID;
385  photon_et = pho_maxet->et();
386  photon_eta = pho_maxet->eta();
387  photon_phi = pho_maxet->phi();
388  }
389  }
390 
391  // If user requires a photon to be found, but none is, return.
392  // theRequirePhotonFound should pretty much always be set to 'True'
393  // except when running on qcd monte carlo just to see the jets.
394  if (theRequirePhotonFound_ &&
395  (!photon_passPhotonID || photon_et < theMinPhotonEt_))
396  return;
397 
399  // Find the highest et jet
400  Handle<View<Jet> > jetCollection;
401  iEvent.getByToken(theJetCollectionToken_, jetCollection);
402  if (!jetCollection.isValid()) return;
403 
404  float jet_pt = -8.0;
405  float jet_eta = -8.0;
406  float jet_phi = -8.0;
407  int jet_count = 0;
408  float jet2_pt = -9.0;
409  float jet2_eta = -9.0;
410  float jet2_phi = -9.0;
411  // Assumption: jets are ordered by Et
412  for (unsigned int i_jet = 0; i_jet < jetCollection->size(); i_jet++) {
413  const Jet* jet = &jetCollection->at(i_jet);
414 
415  float jet_current_pt = jet->pt();
416 
417  // don't care about jets that overlap with the lead photon
418  if (deltaR(jet->eta(), jet->phi(), photon_eta, photon_phi) < 0.5) continue;
419  // stop looping over jets once we get to too low Et
420  if (jet_current_pt < theMinJetPt_) break;
421 
422  jet_count++;
423  if (jet_current_pt > jet_pt) {
424  jet2_pt = jet_pt; // 2nd highest jet get's et from current highest
425  jet2_eta = jet_eta;
426  jet2_phi = jet_phi;
427  jet_pt =
428  jet_current_pt; // current highest jet gets et from the new highest
429  jet_eta = jet->eta();
430  jet_phi = jet->phi();
431  } else if (jet_current_pt > jet2_pt) {
432  jet2_pt = jet_current_pt;
433  jet2_eta = jet->eta();
434  jet2_phi = jet->phi();
435  }
436  }
438 
440  // Fill histograms if a jet found
441  // NOTE: if a photon was required to be found, but wasn't
442  // we wouldn't have made it to this point in the code
443  if (jet_pt > 0.0) {
444 
445  // Photon Plots
446  h_photon_et->Fill(photon_et);
447  h_photon_eta->Fill(photon_eta);
448  h_photon_count_bar->Fill(photon_count_bar);
449  h_photon_count_end->Fill(photon_count_end);
450 
451  // Photon Et hists for different orientations to the jet
452  if (fabs(photon_eta) < 1.45 &&
453  photon_passPhotonID) { // Lead photon is in barrel
454  if (fabs(jet_eta) < 1.45) { // jet is in barrel
455  if (photon_eta * jet_eta > 0) {
456  h_photon_et_jetcs->Fill(photon_et);
457  } else {
458  h_photon_et_jetco->Fill(photon_et);
459  }
460  } else if (jet_eta > 1.55 && jet_eta < 2.5) { // jet is in endcap
461  if (photon_eta * jet_eta > 0) {
462  h_photon_et_jetfs->Fill(photon_et);
463  } else {
464  h_photon_et_jetfo->Fill(photon_et);
465  }
466  }
467  } // END of Lead Photon is in Barrel
468 
469  // Jet Plots
470  h_jet_pt->Fill(jet_pt);
471  h_jet_eta->Fill(jet_eta);
472  h_jet_count->Fill(jet_count);
473  h_deltaPhi_photon_jet->Fill(abs(deltaPhi(photon_phi, jet_phi)));
474  if (abs(deltaPhi(photon_phi, jet_phi)) > 2.8)
475  h_deltaEt_photon_jet->Fill((photon_et - jet_pt) / photon_et);
476 
477  // 2nd Highest Jet Plots
478  if (jet2_pt > 0.0) {
479  h_jet2_pt->Fill(jet2_pt);
480  h_jet2_eta->Fill(jet2_eta);
481  h_jet2_ptOverPhotonEt->Fill(jet2_pt / photon_et);
482  h_deltaPhi_photon_jet2->Fill(abs(deltaPhi(photon_phi, jet2_phi)));
483  h_deltaPhi_jet_jet2->Fill(abs(deltaPhi(jet_phi, jet2_phi)));
484  h_deltaR_jet_jet2->Fill(deltaR(jet_eta, jet_phi, jet2_eta, jet2_phi));
485  h_deltaR_photon_jet2->Fill(
486  deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi));
487  }
488  }
489  // End of Filling histograms
491 }
492 
493 // Local Variables:
494 // show-trailing-whitespace: t
495 // truncate-lines: t
496 // End:
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:208
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
Base class for all types of Jets.
Definition: Jet.h:20
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Strings::size_type size() const
Definition: TriggerNames.cc:39
tuple vertexCollection
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
QcdPhotonsDQM(const edm::ParameterSet &)
Constructor.
int iEvent
Definition: GenABIO.cc:230
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:75
Container::value_type value_type
#define LogTrace(id)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
T const * product() const
Definition: Handle.h:81
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
TH1F * getTH1F(void) const
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
virtual ~QcdPhotonsDQM()
Destructor.
void analyze(const edm::Event &, const edm::EventSetup &)
Get the analysis.
virtual double phi() const
momentum azimuthal angle
Definition: Run.h:41