CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DQMExample_Step1.cc
Go to the documentation of this file.
3 
4 
5 // Geometry
9 #include "TLorentzVector.h"
10 
11 
12 #include <iostream>
13 #include <iomanip>
14 #include <stdio.h>
15 #include <string>
16 #include <sstream>
17 #include <math.h>
18 
19 //
20 // -------------------------------------- Constructor --------------------------------------------
21 //
23 {
24  edm::LogInfo("DQMExample_Step1") << "Constructor DQMExample_Step1::DQMExample_Step1 " << std::endl;
25 
26  // Get parameters from configuration file
27  theElectronCollection_ = consumes<reco::GsfElectronCollection>(ps.getParameter<edm::InputTag>("electronCollection"));
28  theCaloJetCollection_ = consumes<reco::CaloJetCollection>(ps.getParameter<edm::InputTag>("caloJetCollection"));
29  thePfMETCollection_ = consumes<reco::PFMETCollection>(ps.getParameter<edm::InputTag>("pfMETCollection"));
30  theConversionCollection_ = consumes<reco::ConversionCollection>(ps.getParameter<edm::InputTag>("conversionsCollection"));
31  thePVCollection_ = consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("PVCollection"));
32  theBSCollection_ = consumes<reco::BeamSpot>(ps.getParameter<edm::InputTag>("beamSpotCollection"));
33  triggerEvent_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("TriggerEvent"));
34  triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
35  triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
36  triggerPath_ = ps.getParameter<std::string>("TriggerPath");
37 
38 
39  // cuts:
40  ptThrL1_ = ps.getUntrackedParameter<double>("PtThrL1");
41  ptThrL2_ = ps.getUntrackedParameter<double>("PtThrL2");
42  ptThrJet_ = ps.getUntrackedParameter<double>("PtThrJet");
43  ptThrMet_ = ps.getUntrackedParameter<double>("PtThrMet");
44 
45 }
46 
47 //
48 // -- Destructor
49 //
51 {
52  edm::LogInfo("DQMExample_Step1") << "Destructor DQMExample_Step1::~DQMExample_Step1 " << std::endl;
53 }
54 
55 //
56 // -------------------------------------- beginJob --------------------------------------------
57 //
59 {
60  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::beginJob " << std::endl;
61 }
62 //
63 // -------------------------------------- beginRun --------------------------------------------
64 //
66 {
67  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::beginRun" << std::endl;
68 }
69 //
70 // -------------------------------------- bookHistos --------------------------------------------
71 //
73 {
74  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::bookHistograms" << std::endl;
75 
76  //book at beginRun
77  bookHistos(ibooker_);
78 }
79 //
80 // -------------------------------------- beginLuminosityBlock --------------------------------------------
81 //
83  edm::EventSetup const& context)
84 {
85  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::beginLuminosityBlock" << std::endl;
86 }
87 
88 
89 //
90 // -------------------------------------- Analyze --------------------------------------------
91 //
93 {
94  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::analyze" << std::endl;
95 
96 
97  //-------------------------------
98  //--- Vertex Info
99  //-------------------------------
101  e.getByToken(thePVCollection_, vertexHandle);
102  if ( !vertexHandle.isValid() )
103  {
104  edm::LogError ("DQMClientExample") << "invalid collection: vertex" << "\n";
105  return;
106  }
107 
108  int vertex_number = vertexHandle->size();
109  reco::VertexCollection::const_iterator v = vertexHandle->begin();
110 
111  math::XYZPoint PVPoint(-999, -999, -999);
112  if(vertex_number != 0)
113  PVPoint = math::XYZPoint(v->position().x(), v->position().y(), v->position().z());
114 
115  PVPoint_=PVPoint;
116 
117  //-------------------------------
118  //--- MET
119  //-------------------------------
120  edm::Handle<reco::PFMETCollection> pfMETCollection;
121  e.getByToken(thePfMETCollection_, pfMETCollection);
122  if ( !pfMETCollection.isValid() )
123  {
124  edm::LogError ("DQMClientExample") << "invalid collection: MET" << "\n";
125  return;
126  }
127  //-------------------------------
128  //--- Electrons
129  //-------------------------------
130  edm::Handle<reco::GsfElectronCollection> electronCollection;
131  e.getByToken(theElectronCollection_, electronCollection);
132  if ( !electronCollection.isValid() )
133  {
134  edm::LogError ("DQMClientExample") << "invalid collection: electrons" << "\n";
135  return;
136  }
137 
138  float nEle=0;
139  int posEle=0, negEle=0;
140  const reco::GsfElectron* ele1 = NULL;
141  const reco::GsfElectron* ele2 = NULL;
142  for (reco::GsfElectronCollection::const_iterator recoElectron=electronCollection->begin(); recoElectron!=electronCollection->end(); ++recoElectron)
143  {
144  //decreasing pT
145  if( MediumEle(e,eSetup,*recoElectron) )
146  {
147  if(!ele1 && recoElectron->pt() > ptThrL1_)
148  ele1 = &(*recoElectron);
149 
150  else if(!ele2 && recoElectron->pt() > ptThrL2_)
151  ele2 = &(*recoElectron);
152 
153  }
154 
155  if(recoElectron->charge()==1)
156  posEle++;
157  else if(recoElectron->charge()==-1)
158  negEle++;
159 
160  } // end of loop over electrons
161 
162  nEle = posEle+negEle;
163 
164  //-------------------------------
165  //--- Jets
166  //-------------------------------
167  edm::Handle<reco::CaloJetCollection> caloJetCollection;
168  e.getByToken (theCaloJetCollection_,caloJetCollection);
169  if ( !caloJetCollection.isValid() )
170  {
171  edm::LogError ("DQMClientExample") << "invalid collection: jets" << "\n";
172  return;
173  }
174 
175  int nJet = 0;
176  const reco::CaloJet* jet1 = NULL;
177  const reco::CaloJet* jet2 = NULL;
178 
179  for (reco::CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin(); i_calojet != caloJetCollection->end(); ++i_calojet)
180  {
181  //remove jet-ele matching
182  if(ele1)
183  if (Distance(*i_calojet,*ele1) < 0.3) continue;
184 
185  if(ele2)
186  if (Distance(*i_calojet,*ele2) < 0.3) continue;
187 
188  if (i_calojet->pt() < ptThrJet_) continue;
189 
190  nJet++;
191 
192  if (!jet1)
193  jet1 = &(*i_calojet);
194 
195  else if (!jet2)
196  jet2 = &(*i_calojet);
197  }
198 
199  // ---------------------------
200  // ---- Analyze Trigger Event
201  // ---------------------------
202 
203  //check what is in the menu
205  e.getByToken(triggerResults_,hltresults);
206 
207  if(!hltresults.isValid())
208  {
209  edm::LogError ("DQMClientExample") << "invalid collection: TriggerResults" << "\n";
210  return;
211  }
212 
213  bool hasFired = false;
214  const edm::TriggerNames& trigNames = e.triggerNames(*hltresults);
215  unsigned int numTriggers = trigNames.size();
216 
217  for( unsigned int hltIndex=0; hltIndex<numTriggers; ++hltIndex )
218  {
219  if (trigNames.triggerName(hltIndex)==triggerPath_ && hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex))
220  hasFired = true;
221  }
222 
223 
224 
225  //access the trigger event
227  e.getByToken(triggerEvent_, triggerEvent);
228  if( triggerEvent.failedToGet() )
229  {
230  edm::LogError ("DQMClientExample") << "invalid collection: TriggerEvent" << "\n";
231  return;
232  }
233 
234 
235  reco::Particle* ele1_HLT = NULL;
236  int nEle_HLT = 0;
237 
238  size_t filterIndex = triggerEvent->filterIndex( triggerFilter_ );
239  trigger::TriggerObjectCollection triggerObjects = triggerEvent->getObjects();
240  if( !(filterIndex >= triggerEvent->sizeFilters()) )
241  {
242  const trigger::Keys& keys = triggerEvent->filterKeys( filterIndex );
243  std::vector<reco::Particle> triggeredEle;
244 
245  for( size_t j = 0; j < keys.size(); ++j )
246  {
247  trigger::TriggerObject foundObject = triggerObjects[keys[j]];
248  if( abs( foundObject.particle().pdgId() ) != 11 ) continue; //make sure that it is an electron
249 
250  triggeredEle.push_back( foundObject.particle() );
251  ++nEle_HLT;
252  }
253 
254  if( triggeredEle.size() >= 1 )
255  ele1_HLT = &(triggeredEle.at(0));
256  }
257 
258  //-------------------------------
259  //--- Fill the histos
260  //-------------------------------
261 
262  //vertex
263  h_vertex_number -> Fill( vertex_number );
264 
265  //met
266  h_pfMet -> Fill( pfMETCollection->begin()->et() );
267 
268  //multiplicities
269  h_eMultiplicity->Fill(nEle);
270  h_jMultiplicity->Fill(nJet);
271  h_eMultiplicity_HLT->Fill(nEle_HLT);
272 
273  //leading not matched
274  if(ele1)
275  {
276  h_ePt_leading->Fill(ele1->pt());
277  h_eEta_leading->Fill(ele1->eta());
278  h_ePhi_leading->Fill(ele1->phi());
279  }
280  if(ele1_HLT)
281  {
282  h_ePt_leading_HLT->Fill(ele1_HLT->pt());
283  h_eEta_leading_HLT->Fill(ele1_HLT->eta());
284  h_ePhi_leading_HLT->Fill(ele1_HLT->phi());
285  }
286  //leading Jet
287  if(jet1)
288  {
289  h_jPt_leading->Fill(jet1->pt());
290  h_jEta_leading->Fill(jet1->eta());
291  h_jPhi_leading->Fill(jet1->phi());
292  }
293 
294 
295  //fill only when the trigger candidate mathes with the reco one
296  if( ele1 && ele1_HLT && deltaR(*ele1_HLT,*ele1) < 0.3 && hasFired==true )
297  {
298  h_ePt_leading_matched->Fill(ele1->pt());
299  h_eEta_leading_matched->Fill(ele1->eta());
300  h_ePhi_leading_matched->Fill(ele1->phi());
301 
302  h_ePt_leading_HLT_matched->Fill(ele1_HLT->pt());
303  h_eEta_leading_HLT_matched->Fill(ele1_HLT->eta());
304  h_ePhi_leading_HLT_matched->Fill(ele1_HLT->phi());
305 
306  h_ePt_diff->Fill(ele1->pt()-ele1_HLT->pt());
307  }
308 }
309 //
310 // -------------------------------------- endLuminosityBlock --------------------------------------------
311 //
313 {
314  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::endLuminosityBlock" << std::endl;
315 }
316 
317 
318 //
319 // -------------------------------------- endRun --------------------------------------------
320 //
322 {
323  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::endRun" << std::endl;
324 }
325 
326 //
327 // -------------------------------------- endJob --------------------------------------------
328 //
330 {
331  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::endJob" << std::endl;
332 }
333 
334 
335 //
336 // -------------------------------------- book histograms --------------------------------------------
337 //
339 {
340  ibooker_.cd();
341  ibooker_.setCurrentFolder("Physics/TopTest");
342 
343  h_vertex_number = ibooker_.book1D("Vertex_number", "Number of event vertices in collection", 40,-0.5, 39.5 );
344 
345  h_pfMet = ibooker_.book1D("pfMet", "Pf Missing E_{T}; GeV" , 20, 0.0 , 100);
346 
347  h_eMultiplicity = ibooker_.book1D("NElectrons","# of electrons per event",10,0.,10.);
348  h_ePt_leading_matched = ibooker_.book1D("ElePt_leading_matched","Pt of leading electron",50,0.,100.);
349  h_eEta_leading_matched = ibooker_.book1D("EleEta_leading_matched","Eta of leading electron",50,-5.,5.);
350  h_ePhi_leading_matched = ibooker_.book1D("ElePhi_leading_matched","Phi of leading electron",50,-3.5,3.5);
351 
352  h_ePt_leading = ibooker_.book1D("ElePt_leading","Pt of leading electron",50,0.,100.);
353  h_eEta_leading = ibooker_.book1D("EleEta_leading","Eta of leading electron",50,-5.,5.);
354  h_ePhi_leading = ibooker_.book1D("ElePhi_leading","Phi of leading electron",50,-3.5,3.5);
355 
356  h_jMultiplicity = ibooker_.book1D("NJets","# of electrons per event",10,0.,10.);
357  h_jPt_leading = ibooker_.book1D("JetPt_leading","Pt of leading Jet",150,0.,300.);
358  h_jEta_leading = ibooker_.book1D("JetEta_leading","Eta of leading Jet",50,-5.,5.);
359  h_jPhi_leading = ibooker_.book1D("JetPhi_leading","Phi of leading Jet",50,-3.5,3.5);
360 
361  h_eMultiplicity_HLT = ibooker_.book1D("NElectrons_HLT","# of electrons per event @HLT",10,0.,10.);
362  h_ePt_leading_HLT = ibooker_.book1D("ElePt_leading_HLT","Pt of leading electron @HLT",50,0.,100.);
363  h_eEta_leading_HLT = ibooker_.book1D("EleEta_leading_HLT","Eta of leading electron @HLT",50,-5.,5.);
364  h_ePhi_leading_HLT = ibooker_.book1D("ElePhi_leading_HLT","Phi of leading electron @HLT",50,-3.5,3.5);
365 
366  h_ePt_leading_HLT_matched = ibooker_.book1D("ElePt_leading_HLT_matched","Pt of leading electron @HLT",50,0.,100.);
367  h_eEta_leading_HLT_matched = ibooker_.book1D("EleEta_leading_HLT_matched","Eta of leading electron @HLT",50,-5.,5.);
368  h_ePhi_leading_HLT_matched = ibooker_.book1D("ElePhi_leading_HLT_matched","Phi of leading electron @HLT",50,-3.5,3.5);
369 
370  h_ePt_diff = ibooker_.book1D("ElePt_diff_matched","pT(RECO) - pT(HLT) for mathed candidates",100,-10,10.);
371 
372  ibooker_.cd();
373 
374 }
375 
376 
377 //
378 // -------------------------------------- functions --------------------------------------------
379 //
381  return deltaR(c1,c2);
382 }
383 
385  return deltaPhi(c1.p4().phi(),c2.p4().phi());
386 }
387 
388 // This always returns only a positive deltaPhi
389 double DQMExample_Step1::calcDeltaPhi(double phi1, double phi2) {
390  double deltaPhi = phi1 - phi2;
391  if (deltaPhi < 0) deltaPhi = -deltaPhi;
392  if (deltaPhi > 3.1415926) {
393  deltaPhi = 2 * 3.1415926 - deltaPhi;
394  }
395  return deltaPhi;
396 }
397 
398 //
399 // -------------------------------------- electronID --------------------------------------------
400 //
402 {
403 
404  //********* CONVERSION TOOLS
406  iEvent.getByToken(theConversionCollection_, conversions_h);
407 
408  bool isMediumEle = false;
409 
410  float pt = electron.pt();
411  float eta = electron.eta();
412 
413  int isEB = electron.isEB();
414  float sigmaIetaIeta = electron.sigmaIetaIeta();
415  float DetaIn = electron.deltaEtaSuperClusterTrackAtVtx();
416  float DphiIn = electron.deltaPhiSuperClusterTrackAtVtx();
417  float HOverE = electron.hadronicOverEm();
418  float ooemoop = (1.0/electron.ecalEnergy() - electron.eSuperClusterOverP()/electron.ecalEnergy());
419 
420  int mishits = electron.gsfTrack()->trackerExpectedHitsInner().numberOfHits();
421  int nAmbiguousGsfTracks = electron.ambiguousGsfTracksSize();
422 
423  reco::GsfTrackRef eleTrack = electron.gsfTrack() ;
424  float dxy = eleTrack->dxy(PVPoint_);
425  float dz = eleTrack->dz (PVPoint_);
426 
428  iEvent.getByToken(theBSCollection_, BSHandle);
429  const reco::BeamSpot BS = *BSHandle;
430 
431  bool isConverted = ConversionTools::hasMatchedConversion(electron, conversions_h, BS.position());
432 
433  // default
434  if( (pt > 12.) && (fabs(eta) < 2.5) &&
435  ( ( (isEB == 1) && (fabs(DetaIn) < 0.004) ) || ( (isEB == 0) && (fabs(DetaIn) < 0.007) ) ) &&
436  ( ( (isEB == 1) && (fabs(DphiIn) < 0.060) ) || ( (isEB == 0) && (fabs(DphiIn) < 0.030) ) ) &&
437  ( ( (isEB == 1) && (sigmaIetaIeta < 0.010) ) || ( (isEB == 0) && (sigmaIetaIeta < 0.030) ) ) &&
438  ( ( (isEB == 1) && (HOverE < 0.120) ) || ( (isEB == 0) && (HOverE < 0.100) ) ) &&
439  ( ( (isEB == 1) && (fabs(ooemoop) < 0.050) ) || ( (isEB == 0) && (fabs(ooemoop) < 0.050) ) ) &&
440  ( ( (isEB == 1) && (fabs(dxy) < 0.020) ) || ( (isEB == 0) && (fabs(dxy) < 0.020) ) ) &&
441  ( ( (isEB == 1) && (fabs(dz) < 0.100) ) || ( (isEB == 0) && (fabs(dz) < 0.100) ) ) &&
442  ( ( (isEB == 1) && (!isConverted) ) || ( (isEB == 0) && (!isConverted) ) ) &&
443  ( mishits == 0 ) &&
444  ( nAmbiguousGsfTracks == 0 )
445  )
446  isMediumEle=true;
447 
448  return isMediumEle;
449 }
double calcDeltaPhi(double phi1, double phi2)
MonitorElement * h_ePhi_leading_matched
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:204
MonitorElement * h_jMultiplicity
math::XYZPoint PVPoint_
bool MediumEle(const edm::Event &iEvent, const edm::EventSetup &iESetup, const reco::GsfElectron &electron)
edm::EDGetTokenT< reco::GsfElectronCollection > theElectronCollection_
Jets made from CaloTowers.
Definition: CaloJet.h:29
float eSuperClusterOverP() const
Definition: GsfElectron.h:243
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void cd(void)
Definition: DQMStore.cc:266
void endLuminosityBlock(edm::LuminosityBlock const &lumi, edm::EventSetup const &eSetup)
reco::Particle particle(reco::Particle::Charge q=0, const reco::Particle::Point &vertex=reco::Particle::Point(0, 0, 0), int status=0, bool integerCharge=true) const
Definition: TriggerObject.h:69
MonitorElement * h_ePhi_leading_HLT_matched
#define NULL
Definition: scimark2.h:8
MonitorElement * h_ePhi_leading
MonitorElement * h_eEta_leading_HLT
MonitorElement * h_jPhi_leading
edm::EDGetTokenT< reco::ConversionCollection > theConversionCollection_
Strings::size_type size() const
Definition: TriggerNames.cc:39
edm::EDGetTokenT< trigger::TriggerEvent > triggerEvent_
T eta() const
double Distance(const reco::Candidate &c1, const reco::Candidate &c2)
bool isEB() const
Definition: GsfElectron.h:347
double eta() const
momentum pseudorapidity
Definition: Particle.cc:168
DQMExample_Step1(const edm::ParameterSet &ps)
void Fill(long long x)
MonitorElement * h_eEta_leading
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< reco::BeamSpot > theBSCollection_
MonitorElement * h_ePt_leading_HLT
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:22
double phi() const
momentum azimuthal angle
Definition: Particle.cc:159
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:247
int iEvent
Definition: GenABIO.cc:230
MonitorElement * h_ePt_diff
MonitorElement * h_ePt_leading_HLT_matched
float sigmaIetaIeta() const
Definition: GsfElectron.h:399
float hadronicOverEm() const
Definition: GsfElectron.h:440
MonitorElement * h_jEta_leading
void bookHistos(DQMStore::IBooker &)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:250
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:112
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
MonitorElement * h_eEta_leading_HLT_matched
MonitorElement * h_vertex_number
MonitorElement * h_ePt_leading
bool isValid() const
Definition: HandleBase.h:76
double pt() const
transverse momentum
Definition: Particle.cc:155
static bool hasMatchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
MonitorElement * h_eMultiplicity
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void analyze(edm::Event const &e, edm::EventSetup const &eSetup)
MonitorElement * h_pfMet
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:81
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
bool failedToGet() const
Definition: HandleBase.h:80
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
void endRun(edm::Run const &run, edm::EventSetup const &eSetup)
MonitorElement * h_jPt_leading
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< size_type > Keys
MonitorElement * h_ePt_leading_matched
virtual ~DQMExample_Step1()
GsfTrackRefVector::size_type ambiguousGsfTracksSize() const
Definition: GsfElectron.h:633
int pdgId() const
PDG identifier.
Definition: Particle.cc:246
edm::EDGetTokenT< reco::PFMETCollection > thePfMETCollection_
edm::InputTag triggerFilter_
float ecalEnergy() const
Definition: GsfElectron.h:772
MonitorElement * h_ePhi_leading_HLT
double DistancePhi(const reco::Candidate &c1, const reco::Candidate &c2)
std::string triggerPath_
const Point & position() const
position
Definition: BeamSpot.h:62
void beginLuminosityBlock(edm::LuminosityBlock const &lumi, edm::EventSetup const &eSetup)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * h_eMultiplicity_HLT
virtual float pt() const GCC11_FINAL
transverse momentum
MonitorElement * h_eEta_leading_matched
edm::EDGetTokenT< reco::VertexCollection > thePVCollection_
edm::EDGetTokenT< reco::CaloJetCollection > theCaloJetCollection_
Definition: Run.h:41
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:183
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector