CMS 3D CMS Logo

PackedCandidate.cc
Go to the documentation of this file.
6 
8 using namespace logintpack;
9 
12 
13 void pat::PackedCandidate::pack(bool unpackAfterwards) {
14  packedPt_ = MiniFloatConverter::float32to16(p4_.load()->Pt());
15  packedEta_ = int16_t(std::round(p4_.load()->Eta()/6.0f*std::numeric_limits<int16_t>::max()));
16  packedPhi_ = int16_t(std::round(p4_.load()->Phi()/3.2f*std::numeric_limits<int16_t>::max()));
17  packedM_ = MiniFloatConverter::float32to16(p4_.load()->M());
18  if (unpackAfterwards) {
19  delete p4_.exchange(nullptr);
20  delete p4c_.exchange(nullptr);
21  unpack(); // force the values to match with the packed ones
22  }
23 }
24 
25 void pat::PackedCandidate::packVtx(bool unpackAfterwards) {
26  reco::VertexRef pvRef = vertexRef();
27  Point pv = pvRef.isNonnull() ? pvRef->position() : Point();
28  float dxPV = vertex_.load()->X() - pv.X(), dyPV = vertex_.load()->Y() - pv.Y(); //, rPV = std::hypot(dxPV, dyPV);
29  float s = std::sin(float(p4_.load()->Phi())+dphi_), c = std::cos(float(p4_.load()->Phi()+dphi_)); // not the fastest option, but we're in reduced precision already, so let's avoid more roundoffs
30  dxy_ = - dxPV * s + dyPV * c;
31  // if we want to go back to the full x,y,z we need to store also
32  // float dl = dxPV * c + dyPV * s;
33  // float xRec = - dxy_ * s + dl * c, yRec = dxy_ * c + dl * s;
34  float pzpt = p4_.load()->Pz()/p4_.load()->Pt();
35  dz_ = vertex_.load()->Z() - pv.Z() - (dxPV*c + dyPV*s) * pzpt;
36  packedDxy_ = MiniFloatConverter::float32to16(dxy_*100);
37  packedDz_ = pvRef.isNonnull() ? MiniFloatConverter::float32to16(dz_*100) : int16_t(std::round(dz_/40.f*std::numeric_limits<int16_t>::max()));
38  packedDPhi_ = int16_t(std::round(dphi_/3.2f*std::numeric_limits<int16_t>::max()));
39  packedDEta_ = MiniFloatConverter::float32to16(deta_);
40  packedDTrkPt_ = MiniFloatConverter::float32to16(dtrkpt_);
41 
42  if (unpackAfterwards) {
43  delete vertex_.exchange(nullptr);
44  unpackVtx();
45  }
46 }
47 
49  float pt = MiniFloatConverter::float16to32(packedPt_);
50  double shift = (pt<1. ? 0.1*pt : 0.1/pt); // shift particle phi to break degeneracies in angular separations
51  double sign = ( ( int(pt*10) % 2 == 0 ) ? 1 : -1 ); // introduce a pseudo-random sign of the shift
52  double phi = int16_t(packedPhi_)*3.2f/std::numeric_limits<int16_t>::max() + sign*shift*3.2/std::numeric_limits<int16_t>::max();
53  auto p4 = std::make_unique<PolarLorentzVector>(pt,
54  int16_t(packedEta_)*6.0f/std::numeric_limits<int16_t>::max(),
55  phi,
57  auto p4c = std::make_unique<LorentzVector>( *p4 );
58  PolarLorentzVector* expectp4= nullptr;
59  if( p4_.compare_exchange_strong(expectp4,p4.get()) ) {
60  p4.release();
61  }
62 
63  //p4c_ works as the guard for unpacking so it
64  // must be set last
65  LorentzVector* expectp4c = nullptr;
66  if(p4c_.compare_exchange_strong(expectp4c, p4c.get()) ) {
67  p4c.release();
68  }
69 }
70 
72  packedCovariance_.dptdpt = packCovarianceElement(m,0,0);
73  packedCovariance_.detadeta = packCovarianceElement(m,1,1);
74  packedCovariance_.dphidphi = packCovarianceElement(m,2,2);
75  packedCovariance_.dxydxy =packCovarianceElement(m,3,3);
76  packedCovariance_.dzdz = packCovarianceElement(m,4,4);
77  packedCovariance_.dxydz = packCovarianceElement(m,3,4);
78  packedCovariance_.dlambdadz = packCovarianceElement(m,1,4);
79  packedCovariance_.dphidxy = packCovarianceElement(m,2,3);
80  //unpack afterwards
81  if(unpackAfterwards) unpackCovariance();
82 }
83 
85  const CovarianceParameterization & p=covarianceParameterization();
86  if(p.isValid())
87  {
88  auto m = std::make_unique<reco::TrackBase::CovarianceMatrix>() ;
89  for(int i=0;i<5;i++)
90  for(int j=0;j<5;j++){
91  (*m)(i,j)=0;
92  }
93  unpackCovarianceElement(*m,packedCovariance_.dptdpt,0,0);
94  unpackCovarianceElement(*m,packedCovariance_.detadeta,1,1);
95  unpackCovarianceElement(*m,packedCovariance_.dphidphi,2,2);
96  unpackCovarianceElement(*m,packedCovariance_.dxydxy,3,3);
97  unpackCovarianceElement(*m,packedCovariance_.dzdz,4,4);
98  unpackCovarianceElement(*m,packedCovariance_.dxydz,3,4);
99  unpackCovarianceElement(*m,packedCovariance_.dlambdadz,1,4);
100  unpackCovarianceElement(*m,packedCovariance_.dphidxy,2,3);
101  reco::TrackBase::CovarianceMatrix* expected = nullptr;
102  if( m_.compare_exchange_strong(expected,m.get()) ) {
103  m.release();
104  }
105 
106  } else {
108  << "You do not have a valid track parameters file loaded. "
109  << "Please check that the release version is compatible with your input data"
110  <<"or avoid accessing track parameter uncertainties. ";
111  }
112 }
113 
115  reco::VertexRef pvRef = vertexRef();
116  dphi_ = int16_t(packedDPhi_)*3.2f/std::numeric_limits<int16_t>::max(),
117  deta_ = MiniFloatConverter::float16to32(packedDEta_);
118  dtrkpt_ = MiniFloatConverter::float16to32(packedDTrkPt_);
119  dxy_ = MiniFloatConverter::float16to32(packedDxy_)/100.;
120  dz_ = pvRef.isNonnull() ? MiniFloatConverter::float16to32(packedDz_)/100. : int16_t(packedDz_)*40.f/std::numeric_limits<int16_t>::max();
121  Point pv = pvRef.isNonnull() ? pvRef->position() : Point();
122  float phi = p4_.load()->Phi()+dphi_, s = std::sin(phi), c = std::cos(phi);
123  auto vertex = std::make_unique<Point>(pv.X() - dxy_ * s,
124  pv.Y() + dxy_ * c,
125  pv.Z() + dz_ ); // for our choice of using the PCA to the PV, by definition the remaining term -(dx*cos(phi) + dy*sin(phi))*(pz/pt) is zero
126 
127 
128 
129 
130  Point* expected = nullptr;
131  if( vertex_.compare_exchange_strong(expected,vertex.get()) ) {
132  vertex.release();
133  }
134 }
135 
137  delete p4_.load();
138  delete p4c_.load();
139  delete vertex_.load();
140  delete track_.load();
141  delete m_.load();
142 }
143 
144 
145 float pat::PackedCandidate::dxy(const Point &p) const {
146  maybeUnpackBoth();
147  const float phi = float(p4_.load()->Phi())+dphi_;
148  return -(vertex_.load()->X()-p.X()) * std::sin(phi) + (vertex_.load()->Y()-p.Y()) * std::cos(phi);
149 }
150 float pat::PackedCandidate::dz(const Point &p) const {
151  maybeUnpackBoth();
152  const float phi = float(p4_.load()->Phi())+dphi_;
153  const float pzpt = deta_ ? std::sinh(etaAtVtx()) : p4_.load()->Pz()/p4_.load()->Pt();
154  return (vertex_.load()->Z()-p.Z()) - ((vertex_.load()->X()-p.X()) * std::cos(phi) + (vertex_.load()->Y()-p.Y()) * std::sin(phi)) * pzpt;
155 }
156 
158  maybeUnpackBoth();
159  math::RhoEtaPhiVector p3(ptTrk(),etaAtVtx(),phiAtVtx());
160  maybeUnpackCovariance();
161  int numberOfStripLayers = stripLayersWithMeasurement(), numberOfPixelLayers = pixelLayersWithMeasurement();
162  int numberOfPixelHits = this->numberOfPixelHits();
163  int numberOfHits = this->numberOfHits();
164 
165  int ndof = numberOfHits+numberOfPixelHits-5;
166  LostInnerHits innerLost = lostInnerHits();
167 
168  auto track = std::make_unique<reco::Track>(normalizedChi2_*ndof,ndof,*vertex_,math::XYZVector(p3.x(),p3.y(),p3.z()),charge(),*(m_.load()),reco::TrackBase::undefAlgorithm,reco::TrackBase::loose);
169  int i=0;
170  if ( firstHit_ == 0) { //Backward compatible
171  if(innerLost == validHitInFirstPixelBarrelLayer){
172  track->appendTrackerHitPattern(PixelSubdetector::PixelBarrel, 1, 0, TrackingRecHit::valid);
173  i=1;
174  }
175  } else {
176  track->appendHitPattern(firstHit_,TrackingRecHit::valid);
177  }
178 
179  if(firstHit_!=0 && reco::HitPattern::pixelHitFilter(firstHit_)) i=1;
180 
181  // add hits to match the number of laters and validHitInFirstPixelBarrelLayer
182  if(innerLost == validHitInFirstPixelBarrelLayer){
183  // then to encode the number of layers, we add more hits on distinct layers (B2, B3, B4, F1, ...)
184  for(; i<numberOfPixelLayers; i++) {
185  if (i <= 3) {
186  track->appendTrackerHitPattern(PixelSubdetector::PixelBarrel, i+1, 0, TrackingRecHit::valid);
187  } else {
188  track->appendTrackerHitPattern(PixelSubdetector::PixelEndcap, i-3, 0, TrackingRecHit::valid);
189  }
190  }
191  } else {
192  // to encode the information on the layers, we add one valid hits per layer but skipping PXB1
193  int iOffset=0;
194  if(firstHit_!=0 && reco::HitPattern::pixelHitFilter(firstHit_)) {
195  iOffset=reco::HitPattern::getLayer(firstHit_);
197  } else {iOffset=1; }
198  for(;i<numberOfPixelLayers; i++) {
199  if (i+iOffset <= 2 ) { track->appendTrackerHitPattern(PixelSubdetector::PixelBarrel, i+iOffset+1, 0, TrackingRecHit::valid); }
200  else {track->appendTrackerHitPattern(PixelSubdetector::PixelEndcap, i+iOffset-3+1, 0, TrackingRecHit::valid); }
201 
202  }
203  }
204  // add extra hits (overlaps, etc), all on the first layer with a hit - to avoid increasing the layer count
205  for(;i<numberOfPixelHits; i++) {
206  if(firstHit_ !=0 && reco::HitPattern::pixelHitFilter(firstHit_)) {
207  track->appendTrackerHitPattern(reco::HitPattern::getSubStructure(firstHit_), reco::HitPattern::getLayer(firstHit_), 0, TrackingRecHit::valid);
208  } else {
209  track->appendTrackerHitPattern(PixelSubdetector::PixelBarrel, (innerLost == validHitInFirstPixelBarrelLayer ? 1 : 2), 0, TrackingRecHit::valid);
210  }
211  }
212  // now start adding strip layers, putting one hit on each layer so that the hitPattern.stripLayersWithMeasurement works.
213  // we don't know what the layers where, so we just start with TIB (4 layers), then TOB (6 layers), then TEC (9)
214  // and then TID(3), so that we can get a number of valid strip layers up to 4+6+9+3
215  if(firstHit_!=0 && reco::HitPattern::stripHitFilter(firstHit_)) i+=1;
216  int slOffset=0;
217  if(firstHit_!=0 && reco::HitPattern::stripHitFilter(firstHit_)) {
218  slOffset=reco::HitPattern::getLayer(firstHit_)-1;
219  if(reco::HitPattern::getSubStructure(firstHit_)==StripSubdetector::TID) slOffset+=4;
220  if(reco::HitPattern::getSubStructure(firstHit_)==StripSubdetector::TOB) slOffset+=7;
221  if(reco::HitPattern::getSubStructure(firstHit_)==StripSubdetector::TEC) slOffset+=13;
222  }
223  for(int sl=slOffset; sl < numberOfStripLayers+slOffset; ++sl, ++i) {
224  if (sl < 4) track->appendTrackerHitPattern(StripSubdetector::TIB, sl +1, 1, TrackingRecHit::valid);
225  else if (sl < 4+3) track->appendTrackerHitPattern(StripSubdetector::TID, (sl- 4)+1, 1, TrackingRecHit::valid);
226  else if (sl < 7+6) track->appendTrackerHitPattern(StripSubdetector::TOB, (sl-7)+1, 1, TrackingRecHit::valid);
227  else if (sl < 13+9) track->appendTrackerHitPattern(StripSubdetector::TEC, (sl-13)+1, 1, TrackingRecHit::valid);
228  else break; // wtf?
229  }
230  // finally we account for extra strip hits beyond the one-per-layer added above. we put them all on TIB1,
231  // to avoid incrementing the number of layersWithMeasurement.
232  for(;i<numberOfHits;i++) {
233  if(reco::HitPattern::stripHitFilter(firstHit_)) {
234  track->appendTrackerHitPattern(reco::HitPattern::getSubStructure(firstHit_), reco::HitPattern::getLayer(firstHit_), 1, TrackingRecHit::valid);
235  } else {
236  track->appendTrackerHitPattern(StripSubdetector::TIB, 1, 1, TrackingRecHit::valid);
237  }
238  }
239 
240 
241  switch (innerLost) {
242  case validHitInFirstPixelBarrelLayer:
243  break;
244  case noLostInnerHits:
245  break;
246  case oneLostInnerHit:
247  track->appendTrackerHitPattern(PixelSubdetector::PixelBarrel, 1, 0, TrackingRecHit::missing_inner);
248  break;
249  case moreLostInnerHits:
250  track->appendTrackerHitPattern(PixelSubdetector::PixelBarrel, 1, 0, TrackingRecHit::missing_inner);
251  track->appendTrackerHitPattern(PixelSubdetector::PixelBarrel, 2, 0, TrackingRecHit::missing_inner);
252  break;
253  };
254 
255  if (trackHighPurity()) track->setQuality(reco::TrackBase::highPurity);
256 
257  reco::Track* expected = nullptr;
258  if( track_.compare_exchange_strong(expected,track.get()) ) {
259  track.release();
260  }
261 
262 
263 }
264 
266 
268  throw cms::Exception("Invalid Reference")
269  << "this Candidate has no master clone reference."
270  << "Can't call masterClone() method.\n";
271 }
272 
274  return false;
275 }
276 
278  return false;
279 }
280 
281 
283  throw cms::Exception("Invalid Reference")
284  << "this Candidate has no master clone ptr."
285  << "Can't call masterClonePtr() method.\n";
286 }
287 
289  return 0;
290 }
291 
293  return 0;
294 }
295 
297  return p4() == o.p4() && vertex() == o.vertex() && charge() == o.charge();
298 // return p4() == o.p4() && charge() == o.charge();
299 }
300 
302  return nullptr;
303 }
304 
306  return nullptr;
307 }
308 
311  << "This Candidate type does not implement daughter(std::string). "
312  << "Please use CompositeCandidate or NamedCompositeCandidate.\n";
313 }
314 
317  << "This Candidate type does not implement daughter(std::string). "
318  << "Please use CompositeCandidate or NamedCompositeCandidate.\n";
319 }
320 
321 
322 
324  return nullptr;
325 }
326 
328  return 0;
329 }
330 
332  return 0;
333 }
334 
336  return 0;
337 }
338 
339 double pat::PackedCandidate::vertexCovariance(int i, int j) const {
341  << "reco::ConcreteCandidate does not implement vertex covariant matrix.\n";
342 }
343 
346  << "reco::ConcreteCandidate does not implement vertex covariant matrix.\n";
347 }
348 
349 
350 bool pat::PackedCandidate::longLived() const {return false;}
351 
352 bool pat::PackedCandidate::massConstraint() const {return false;}
353 
354 // puppiweight
355 void pat::PackedCandidate::setPuppiWeight(float p, float p_nolep) {
356  // Set both weights at once to avoid misconfigured weights if called in the wrong order
357  packedPuppiweight_ = pack8logClosed((p-0.5)*2,-2,0,64);
358  packedPuppiweightNoLepDiff_ = pack8logClosed((p_nolep-0.5)*2,-2,0,64) - packedPuppiweight_;
359 }
360 
361 float pat::PackedCandidate::puppiWeight() const { return unpack8logClosed(packedPuppiweight_,-2,0,64)/2. + 0.5;}
362 
363 float pat::PackedCandidate::puppiWeightNoLep() const { return unpack8logClosed(packedPuppiweightNoLepDiff_+packedPuppiweight_,-2,0,64)/2. + 0.5;}
364 
367  rawCaloFraction_ = std::numeric_limits<uint8_t>::max(); // Set to overflow value
368  else
369  rawCaloFraction_ = 100*p;
370 }
371 
373  hcalFraction_ = 100*p;
374 }
375 
377  isIsolatedChargedHadron_ = p;
378 }
379 
380 void pat::PackedCandidate::setDTimeAssociatedPV(float aTime, float aTimeError) {
381  if (aTime == 0 && aTimeError == 0) {
382  packedTime_ = 0; packedTimeError_ = 0;
383  } else if (aTimeError == 0) {
384  packedTimeError_ = 0;
385  packedTime_ = packTimeNoError(aTime);
386  } else {
387  packedTimeError_ = packTimeError(aTimeError);
388  aTimeError = unpackTimeError(packedTimeError_); // for reproducibility
389  packedTime_ = packTimeWithError(aTime, aTimeError);
390  }
391 }
392 
394 uint8_t pat::PackedCandidate::packTimeError(float timeError) {
395  if (timeError <= 0) return 0;
396  // log-scale packing.
397  // for MIN_TIMEERROR = 0.002, EXPO_TIMEERROR = 5:
398  // minimum value 0.002 = 2ps (packed as 1)
399  // maximum value 0.5 ns (packed as 255)
400  // constant *relative* precision of about 2%
401  return std::max<uint8_t>( std::min(std::round(std::ldexp(std::log2(timeError/MIN_TIMEERROR), +EXPO_TIMEERROR)), 255.f), 1);
402 }
403 float pat::PackedCandidate::unpackTimeError(uint8_t timeError) {
404  return timeError > 0 ? MIN_TIMEERROR * std::exp2(std::ldexp(float(timeError),-EXPO_TIMEERROR)) : -1.0f;
405 }
407  if (time == 0) return 0.f;
408  return (time > 0 ? MIN_TIME_NOERROR : -MIN_TIME_NOERROR) * std::exp2(std::ldexp(float(std::abs(time)),-EXPO_TIME_NOERROR));
409 }
411  // encoding in log scale to store times in a large range with few bits.
412  // for MIN_TIME_NOERROR = 0.0002 and EXPO_TIME_NOERROR = 6:
413  // smallest non-zero time = 0.2 ps (encoded as +/-1)
414  // one BX, +/- 12.5 ns, is fully covered with 11 bits (+/- 1023)
415  // 12 bits cover by far any plausible value (+/-2047 corresponds to about +/- 0.8 ms!)
416  // constant *relative* ~1% precision
417  if (std::abs(time) < MIN_TIME_NOERROR) return 0; // prevent underflows
418  float fpacked = std::ldexp(std::log2(std::abs(time/MIN_TIME_NOERROR)),+EXPO_TIME_NOERROR);
419  return (time > 0 ? +1 : -1)*std::min(std::round(fpacked), 2047.f);
420 }
421 float pat::PackedCandidate::unpackTimeWithError(int16_t time, uint8_t timeError) {
422  if (time % 2 == 0) {
423  // no overflow: drop rightmost bit and unpack in units of timeError
424  return std::ldexp(unpackTimeError(timeError), EXPO_TIME_WITHERROR) * float(time/2);
425  } else {
426  // overflow: drop rightmost bit, unpack using the noError encoding
428  }
429 }
430 int16_t pat::PackedCandidate::packTimeWithError(float time, float timeError) {
431  // Encode in units of timeError * 2^EXPO_TIME_WITHERROR (~1.6% if EXPO_TIME_WITHERROR = -6)
432  // the largest value that can be stored in 14 bits + sign bit + overflow bit is about 260 sigmas
433  // values larger than that will be stored using the no-timeError packing (with less precision).
434  // overflows of these kinds should happen only for particles that are late arriving, out-of-time,
435  // or mis-reconstructed, as timeError is O(20ps) and the beam spot witdth is O(200ps)
436  float fpacked = std::round(time/std::ldexp(timeError, EXPO_TIME_WITHERROR));
437  if (std::abs(fpacked) < 16383.f) { // 16383 = (2^14 - 1) = largest absolute value for a signed 15 bit integer
438  return int16_t(fpacked) * 2; // make it even, and fit in a signed 16 bit int
439  } else {
440  int16_t packed = packTimeNoError(time); // encode
441  return packed * 2 + (time > 0 ? +1 : -1); // make it odd, to signal that there was an overlow
442  }
443 }
444 
float puppiWeight() const
Set both weights at once (with option for only full PUPPI)
virtual float dz(size_t ipv=0) const
dz with respect to the PV[ipv]
static uint32_t getLayer(uint16_t pattern)
Definition: HitPattern.h:700
const reco::Candidate * daughter(size_type) const override
return daughter at a given position (throws an exception)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
static bool pixelHitFilter(uint16_t pattern)
Definition: HitPattern.h:564
void fillVertexCovariance(CovarianceMatrix &v) const override
fill SMatrix
float puppiWeightNoLep() const
Weight from full PUPPI.
void setPuppiWeight(float p, float p_nolep=0.0)
size_t size_type
Definition: Candidate.h:30
CovarianceMatrix vertexCovariance() const override
return SMatrix
bool hasMasterClonePtr() const override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double vertexChi2() const override
chi-squares
bool longLived() const override
is long lived?
static float unpackTimeNoError(int16_t time)
double vertexNdof() const override
void setRawCaloFraction(float p)
Weight from PUPPI removing leptons.
const reco::CandidatePtr & masterClonePtr() const override
static int16_t packTimeNoError(float time)
static int16_t packTimeWithError(float time, float timeError)
std::pair< double, double > Point
Definition: CaloEllipse.h:18
static uint8_t packTimeError(float timeError)
static to allow unit testing
RhoEtaPhiVectorD RhoEtaPhiVector
spatial vector with cylindrical internal representation using pseudorapidity
Definition: Vector3D.h:32
static float float16to32(uint16_t h)
Definition: libminifloat.h:12
bool massConstraint() const override
do mass constraint?
size_t numberOfMothers() const override
number of mothers
int8_t pack8logClosed(double x, double lmin, double lmax, uint8_t base=128)
void setHcalFraction(float p)
Raw ECAL+HCAL energy over candidate energy for isolated charged hadrons.
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
double p4[4]
Definition: TauolaWrapper.h:92
static uint16_t float32to16(float x)
Definition: libminifloat.h:17
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static bool stripHitFilter(uint16_t pattern)
Definition: HitPattern.h:595
static CovarianceParameterization covarianceParameterization_
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double vertexNormalizedChi2() const override
chi-squared divided by n.d.o.f.
double f[11][100]
~PackedCandidate() override
destructor
T min(T a, T b)
Definition: MathUtil.h:58
static uint32_t getSubStructure(uint16_t pattern)
Definition: HitPattern.h:691
LostInnerHits
Enumerator specifying the.
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const reco::Candidate * mother(size_type) const override
return mother at a given position (throws an exception)
static std::once_flag covariance_load_flag
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void unpackCovariance() const
void setDTimeAssociatedPV(float aTime, float aTimeError=0)
set time measurement
virtual int charge() const =0
electric charge
double unpack8logClosed(int8_t i, double lmin, double lmax, uint8_t base=128)
reverse of pack8logClosed
size_t numberOfDaughters() const override
number of daughters
void packVtx(bool unpackAfterwards=true)
static unsigned int const shift
void setIsIsolatedChargedHadron(bool p)
Fraction of Ecal and Hcal for HF and neutral hadrons and isolated charged hadrons.
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
virtual const Point & vertex() const =0
vertex position
void pack(bool unpackAfterwards=true)
void packCovariance(const reco::TrackBase::CovarianceMatrix &m, bool unpackAfterwards=true)
const reco::CandidateBaseRef & masterClone() const override
virtual float dxy() const
dxy with respect to the PV ref
bool overlap(const reco::Candidate &) const override
check overlap with another Candidate
static float unpackTimeWithError(int16_t time, uint8_t timeError)
bool hasMasterClone() const override
static float unpackTimeError(uint8_t timeError)
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77
double p3[4]
Definition: TauolaWrapper.h:91
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:39