CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PackedCandidate.cc
Go to the documentation of this file.
6 
8 using namespace logintpack;
9 
10 void pat::PackedCandidate::pack(bool unpackAfterwards) {
11  packedPt_ = MiniFloatConverter::float32to16(p4_.Pt());
12  packedEta_ = int16_t(std::round(p4_.Eta()/6.0f*std::numeric_limits<int16_t>::max()));
13  packedPhi_ = int16_t(std::round(p4_.Phi()/3.2f*std::numeric_limits<int16_t>::max()));
14  packedM_ = MiniFloatConverter::float32to16(p4_.M());
15  if (unpackAfterwards) unpack(); // force the values to match with the packed ones
16 }
17 
18 void pat::PackedCandidate::packVtx(bool unpackAfterwards) {
19  reco::VertexRef pvRef = vertexRef();
20  Point pv = pvRef.isNonnull() ? pvRef->position() : Point();
21  float dxPV = vertex_.X() - pv.X(), dyPV = vertex_.Y() - pv.Y(); //, rPV = std::hypot(dxPV, dyPV);
22  float s = std::sin(float(p4_.Phi())+dphi_), c = std::cos(float(p4_.Phi()+dphi_)); // not the fastest option, but we're in reduced precision already, so let's avoid more roundoffs
23  dxy_ = - dxPV * s + dyPV * c;
24  // if we want to go back to the full x,y,z we need to store also
25  // float dl = dxPV * c + dyPV * s;
26  // float xRec = - dxy_ * s + dl * c, yRec = dxy_ * c + dl * s;
27  float pzpt = p4_.Pz()/p4_.Pt();
28  dz_ = vertex_.Z() - pv.Z() - (dxPV*c + dyPV*s) * pzpt;
29  packedDxy_ = MiniFloatConverter::float32to16(dxy_*100);
30  packedDz_ = pvRef.isNonnull() ? MiniFloatConverter::float32to16(dz_*100) : int16_t(std::round(dz_/40.f*std::numeric_limits<int16_t>::max()));
31  packedDPhi_ = int16_t(std::round(dphi_/3.2f*std::numeric_limits<int16_t>::max()));
32  packedCovarianceDxyDxy_ = MiniFloatConverter::float32to16(dxydxy_*10000.);
33  packedCovarianceDxyDz_ = MiniFloatConverter::float32to16(dxydz_*10000.);
34  packedCovarianceDzDz_ = MiniFloatConverter::float32to16(dzdz_*10000.);
35 // packedCovarianceDxyDxy_ = pack8log(dxydxy_,-15,-1); // MiniFloatConverter::float32to16(dxydxy_*10000.);
36 // packedCovarianceDxyDz_ = pack8log(dxydz_,-20,-1); //MiniFloatConverter::float32to16(dxydz_*10000.);
37 // packedCovarianceDzDz_ = pack8log(dzdz_,-13,-1); //MiniFloatConverter::float32to16(dzdz_*10000.);
38 
39  packedCovarianceDptDpt_ = pack8logCeil(dptdpt_,-15,0);
40  packedCovarianceDetaDeta_ = pack8logCeil(detadeta_,-20,-5);
41  packedCovarianceDphiDphi_ = pack8logCeil(dphidphi_,-15,0);
42  packedCovarianceDphiDxy_ = pack8log(dphidxy_,-17,-4); // MiniFloatConverter::float32to16(dphidxy_*10000.);
43  packedCovarianceDlambdaDz_ = pack8log(dlambdadz_,-17,-4); // MiniFloatConverter::float32to16(dlambdadz_*10000.);
44 
45  /*packedCovarianceDptDpt_ = pack8logCeil(dptdpt_,-15,5,32);
46  packedCovarianceDetaDeta_ = pack8logCeil(detadeta_,-20,0,32);
47  packedCovarianceDphiDphi_ = pack8logCeil(dphidphi_,-15,5,32);
48  packedCovarianceDphiDxy_ = pack8log(dphidxy_,-17,-4); // MiniFloatConverter::float32to16(dphidxy_*10000.);
49  packedCovarianceDlambdaDz_ = pack8log(dlambdadz_,-17,-4); // MiniFloatConverter::float32to16(dlambdadz_*10000.);
50 
51 */
52 /* packedCovarianceDphiDxy_ = MiniFloatConverter::float32to16(dphidxy_*10000.);
53  packedCovarianceDlambdaDz_ = MiniFloatConverter::float32to16(dlambdadz_*10000.);
54  packedCovarianceDetaDeta_ = MiniFloatConverter::float32to16(detadeta_*10000.);
55  packedCovarianceDphiDphi_ = MiniFloatConverter::float32to16(dphidphi_*10000.);
56  packedCovarianceDptDpt_ = MiniFloatConverter::float32to16(dptdpt_*10000.);
57 */
58 // packedCovarianceDphiDxy_ = pack8log(dphidxy_,-17,-4); // MiniFloatConverter::float32to16(dphidxy_*10000.);
59 // packedCovarianceDlambdaDz_ = pack8log(dlambdadz_,-17,-4); // MiniFloatConverter::float32to16(dlambdadz_*10000.);
60 
61  if (unpackAfterwards) unpackVtx();
62 }
63 
65  float pt = MiniFloatConverter::float16to32(packedPt_);
66  double shift = (pt<1. ? 0.1*pt : 0.1/pt); // shift particle phi to break degeneracies in angular separations
67  double sign = ( ( int(pt*10) % 2 == 0 ) ? 1 : -1 ); // introduce a pseudo-random sign of the shift
68  double phi = int16_t(packedPhi_)*3.2f/std::numeric_limits<int16_t>::max() + sign*shift*3.2/std::numeric_limits<int16_t>::max();
69  p4_ = PolarLorentzVector(pt,
70  int16_t(packedEta_)*6.0f/std::numeric_limits<int16_t>::max(),
71  phi,
73  p4c_ = p4_;
74  unpacked_ = true;
75 }
77  reco::VertexRef pvRef = vertexRef();
78  dphi_ = int16_t(packedDPhi_)*3.2f/std::numeric_limits<int16_t>::max(),
79  dxy_ = MiniFloatConverter::float16to32(packedDxy_)/100.;
80  dz_ = pvRef.isNonnull() ? MiniFloatConverter::float16to32(packedDz_)/100. : int16_t(packedDz_)*40.f/std::numeric_limits<int16_t>::max();
81  Point pv = pvRef.isNonnull() ? pvRef->position() : Point();
82  float phi = p4_.Phi()+dphi_, s = std::sin(phi), c = std::cos(phi);
83  vertex_ = Point(pv.X() - dxy_ * s,
84  pv.Y() + dxy_ * c,
85  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
86 // dxydxy_ = unpack8log(packedCovarianceDxyDxy_,-15,-1);
87 // dxydz_ = unpack8log(packedCovarianceDxyDz_,-20,-1);
88 // dzdz_ = unpack8log(packedCovarianceDzDz_,-13,-1);
89  dphidxy_ = unpack8log(packedCovarianceDphiDxy_,-17,-4);
90  dlambdadz_ = unpack8log(packedCovarianceDlambdaDz_,-17,-4);
91  dptdpt_ = unpack8log(packedCovarianceDptDpt_,-15,0);
92  detadeta_ = unpack8log(packedCovarianceDetaDeta_,-20,-5);
93  dphidphi_ = unpack8log(packedCovarianceDphiDphi_,-15,0);
94 /*
95  dphidxy_ = unpack8log(packedCovarianceDphiDxy_,-17,-4);
96  dlambdadz_ = unpack8log(packedCovarianceDlambdaDz_,-17,-4);
97  dptdpt_ = unpack8log(packedCovarianceDptDpt_,-15,5,32);
98  detadeta_ = unpack8log(packedCovarianceDetaDeta_,-20,0,32);
99  dphidphi_ = unpack8log(packedCovarianceDphiDphi_,-15,5,32);
100 */
101 
102 /* dphidxy_ = MiniFloatConverter::float16to32(packedCovarianceDphiDxy_)/10000.;
103  dlambdadz_ = MiniFloatConverter::float16to32(packedCovarianceDlambdaDz_)/10000.;
104  dptdpt_ = MiniFloatConverter::float16to32(packedCovarianceDptDpt_)/10000.;
105  detadeta_ = MiniFloatConverter::float16to32(packedCovarianceDetaDeta_)/10000.;
106  dphidphi_ = MiniFloatConverter::float16to32(packedCovarianceDphiDphi_)/10000.;
107 */
108  dxydxy_ = MiniFloatConverter::float16to32(packedCovarianceDxyDxy_)/10000.;
109  dxydz_ =MiniFloatConverter::float16to32(packedCovarianceDxyDz_)/10000.;
110  dzdz_ =MiniFloatConverter::float16to32(packedCovarianceDzDz_)/10000.;
111 /* dphidxy_ = MiniFloatConverter::float16to32(packedCovarianceDphiDxy_)/10000.;
112  dlambdadz_ =MiniFloatConverter::float16to32(packedCovarianceDlambdaDz_)/10000.;
113 */
114  unpackedVtx_ = true;
115 }
116 
118 
119 
120 float pat::PackedCandidate::dxy(const Point &p) const {
121  maybeUnpackBoth();
122  return -(vertex_.X()-p.X()) * std::sin(float(p4_.Phi())) + (vertex_.Y()-p.Y()) * std::cos(float(p4_.Phi()));
123 }
124 float pat::PackedCandidate::dz(const Point &p) const {
125  maybeUnpackBoth();
126  return (vertex_.Z()-p.Z()) - ((vertex_.X()-p.X()) * std::cos(float(p4_.Phi())) + (vertex_.Y()-p.Y()) * std::sin(float(p4_.Phi()))) * p4_.Pz()/p4_.Pt();
127 }
128 
130  maybeUnpackBoth();
132 // m(0,0)=0.5e-4/pt()/pt(); //TODO: tune
133 // m(1,1)=6e-6; //TODO: tune
134 // m(2,2)=1.5e-5/pt()/pt(); //TODO: tune
135  m(0,0)=dptdpt_/pt()/pt(); //TODO: tune
136  m(1,1)=detadeta_; //TODO: tune
137  m(2,2)=dphidphi_/pt()/pt(); //TODO: tune
138  m(2,3)=dphidxy_;
139  m(3,2)=dphidxy_;
140  m(4,1)=dlambdadz_;
141  m(1,4)=dlambdadz_;
142  m(3,3)=dxydxy_;
143  m(3,4)=dxydz_;
144  m(4,3)=dxydz_;
145  m(4,4)=dzdz_;
146  math::RhoEtaPhiVector p3(p4_.pt(),p4_.eta(),phiAtVtx());
147  int numberOfPixelHits_ = packedHits_ & 0x7 ;
148  int numberOfHits_ = (packedHits_>>3) + numberOfPixelHits_;
149 
150  int ndof = numberOfHits_+numberOfPixelHits_-5;
151  reco::HitPattern hp, hpExpIn;
152  int i=0;
153  LostInnerHits innerLost = lostInnerHits();
154 
155  track_ = reco::Track(normalizedChi2_*ndof,ndof,vertex_,math::XYZVector(p3.x(),p3.y(),p3.z()),charge(),m,reco::TrackBase::undefAlgorithm,reco::TrackBase::loose);
156 
157  if(innerLost == validHitInFirstPixelBarrelLayer){
159  i++;
160  }
161  for(;i<numberOfPixelHits_; i++) {
162  track_.appendHitPattern(PXBDetId(i > 1 ? 3 : 2, 0, 0), TrackingRecHit::valid);
163  }
164 
165  for(;i<numberOfHits_;i++) {
166  track_.appendHitPattern(TIBDetId(1, 0, 0, 1, 1, 0), TrackingRecHit::valid);
167  }
168 
169  switch (innerLost) {
170  case validHitInFirstPixelBarrelLayer:
171  break;
172  case noLostInnerHits:
173  break;
174  case oneLostInnerHit:
175  track_.appendHitPattern(PXBDetId(1, 0, 0), TrackingRecHit::missing_inner);
176  break;
177  case moreLostInnerHits:
178  track_.appendHitPattern(PXBDetId(1, 0, 0), TrackingRecHit::missing_inner);
179  track_.appendHitPattern(PXBDetId(2, 0, 0), TrackingRecHit::missing_inner);
180  break;
181  };
182 
183  if (trackHighPurity()) track_.setQuality(reco::TrackBase::highPurity);
184 
185  unpackedTrk_ = true;
186 }
187 
189 
191  throw cms::Exception("Invalid Reference")
192  << "this Candidate has no master clone reference."
193  << "Can't call masterClone() method.\n";
194 }
195 
197  return false;
198 }
199 
201  return false;
202 }
203 
204 
206  throw cms::Exception("Invalid Reference")
207  << "this Candidate has no master clone ptr."
208  << "Can't call masterClonePtr() method.\n";
209 }
210 
212  return 0;
213 }
214 
216  return 0;
217 }
218 
220  return p4() == o.p4() && vertex() == o.vertex() && charge() == o.charge();
221 // return p4() == o.p4() && charge() == o.charge();
222 }
223 
225  return 0;
226 }
227 
229  return 0;
230 }
231 
234  << "This Candidate type does not implement daughter(std::string). "
235  << "Please use CompositeCandidate or NamedCompositeCandidate.\n";
236 }
237 
240  << "This Candidate type does not implement daughter(std::string). "
241  << "Please use CompositeCandidate or NamedCompositeCandidate.\n";
242 }
243 
244 
245 
247  return 0;
248 }
249 
251  return 0;
252 }
253 
255  return 0;
256 }
257 
259  return 0;
260 }
261 
264  << "reco::ConcreteCandidate does not implement vertex covariant matrix.\n";
265 }
266 
269  << "reco::ConcreteCandidate does not implement vertex covariant matrix.\n";
270 }
271 
272 
273 bool pat::PackedCandidate::longLived() const {return false;}
274 
275 bool pat::PackedCandidate::massConstraint() const {return false;}
276 
277 // puppiweight
278 void pat::PackedCandidate::setPuppiWeight(float p, float p_nolep) {
279  // Set both weights at once to avoid misconfigured weights if called in the wrong order
280  packedPuppiweight_ = pack8logClosed((p-0.5)*2,-2,0,64);
281  packedPuppiweightNoLepDiff_ = pack8logClosed((p_nolep-0.5)*2,-2,0,64) - packedPuppiweight_;
282 }
283 
284 float pat::PackedCandidate::puppiWeight() const { return unpack8logClosed(packedPuppiweight_,-2,0,64)/2. + 0.5;}
285 
286 float pat::PackedCandidate::puppiWeightNoLep() const { return unpack8logClosed(packedPuppiweightNoLepDiff_+packedPuppiweight_,-2,0,64)/2. + 0.5;}
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]
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
virtual size_t numberOfMothers() const
number of mothers
virtual bool hasMasterClonePtr() const
float puppiWeightNoLep() const
Weight from full PUPPI.
void setPuppiWeight(float p, float p_nolep=0.0)
size_t size_type
Definition: Candidate.h:30
double unpack8log(int8_t i, double lmin, double lmax, uint8_t base=128)
Definition: liblogintpack.h:47
int8_t pack8log(double x, double lmin, double lmax, uint8_t base=128)
Definition: liblogintpack.h:20
virtual void fillVertexCovariance(CovarianceMatrix &v) const
fill SMatrix
double sign(double x)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual double vertexNdof() const
virtual bool overlap(const reco::Candidate &) const
check overlap with another Candidate
std::pair< double, double > Point
Definition: CaloEllipse.h:18
void unpackVtx() const
RhoEtaPhiVectorD RhoEtaPhiVector
spatial vector with cylindrical internal representation using pseudorapidity
Definition: Vector3D.h:32
static float float16to32(uint16_t h)
Definition: libminifloat.h:10
int8_t pack8logCeil(double x, double lmin, double lmax, uint8_t base=128)
Definition: liblogintpack.h:8
virtual const reco::Candidate * daughter(size_type) const
return daughter at a given position (throws an exception)
int8_t pack8logClosed(double x, double lmin, double lmax, uint8_t base=128)
Definition: liblogintpack.h:34
susybsm::HSCParticleRefProd hp
Definition: classes.h:27
double p4[4]
Definition: TauolaWrapper.h:92
static uint16_t float32to16(float x)
Definition: libminifloat.h:15
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual const reco::CandidateBaseRef & masterClone() const
int j
Definition: DBlmapReader.cc:9
bool appendHitPattern(const TrackingRecHit &hit)
append a single hit to the HitPattern
Definition: TrackBase.h:431
double f[11][100]
virtual const Point & vertex() const =0
vertex position
virtual double vertexNormalizedChi2() const
chi-squared divided by n.d.o.f.
virtual int charge() const =0
electric charge
virtual const reco::CandidatePtr & masterClonePtr() const
LostInnerHits
Enumerator specifying the.
virtual bool massConstraint() const
do mass constraint?
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
Geom::Phi< T > phi() const
virtual size_t numberOfDaughters() const
number of daughters
CovarianceMatrix vertexCovariance() const
return SMatrix
double unpack8logClosed(int8_t i, double lmin, double lmax, uint8_t base=128)
reverse of pack8logClosed
Definition: liblogintpack.h:57
virtual bool longLived() const
is long lived?
virtual const reco::Candidate * mother(size_type) const
return mother at a given position (throws an exception)
void packVtx(bool unpackAfterwards=true)
static unsigned int const shift
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
void pack(bool unpackAfterwards=true)
virtual float dxy() const
dxy with respect to the PV ref
virtual double vertexChi2() const
chi-squares
virtual ~PackedCandidate()
destructor
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77
virtual bool hasMasterClone() const
double p3[4]
Definition: TauolaWrapper.h:91
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:39