CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PackedGenParticle.cc
Go to the documentation of this file.
4 
5 void pat::PackedGenParticle::pack(bool unpackAfterwards) {
7  packedY_ = int16_t(p4_.load()->Rapidity() / 6.0f * std::numeric_limits<int16_t>::max());
8  packedPhi_ = int16_t(p4_.load()->Phi() / 3.2f * std::numeric_limits<int16_t>::max());
10  if (unpackAfterwards) {
11  delete p4_.exchange(nullptr);
12  delete p4c_.exchange(nullptr);
13  unpack(); // force the values to match with the packed ones
14  }
15 }
16 
18  float y = int16_t(packedY_) * 6.0f / std::numeric_limits<int16_t>::max();
19  float pt = MiniFloatConverter::float16to32(packedPt_);
20  float m = MiniFloatConverter::float16to32(packedM_);
21  float pz = std::tanh(y) * std::sqrt((m * m + pt * pt) / (1. - std::tanh(y) * std::tanh(y)));
22  float eta = 0;
23  if (pt != 0.) {
24  eta = std::asinh(pz / pt);
25  }
26  double shift = (pt < 1. ? 0.1 * pt : 0.1 / pt); // shift particle phi to break degeneracies in angular separations
27  double sign = ((int(pt * 10) % 2 == 0) ? 1 : -1); // introduce a pseudo-random sign of the shift
28  double phi = int16_t(packedPhi_) * 3.2f / std::numeric_limits<int16_t>::max() +
29  sign * shift * 3.2 / std::numeric_limits<int16_t>::max();
30  auto p4 = std::make_unique<PolarLorentzVector>(pt, eta, phi, m);
31  PolarLorentzVector* expectp4 = nullptr;
32  if (p4_.compare_exchange_strong(expectp4, p4.get())) {
33  p4.release();
34  }
35  auto p4c = std::make_unique<LorentzVector>(*p4_);
36  LorentzVector* expectp4c = nullptr;
37  if (p4c_.compare_exchange_strong(expectp4c, p4c.get())) {
38  p4c.release();
39  }
40 }
41 
43  delete p4_.load();
44  delete p4c_.load();
45 }
46 
47 float pat::PackedGenParticle::dxy(const Point& p) const {
48  unpack();
49  return -(vertex_.X() - p.X()) * std::sin(float(p4_.load()->Phi())) +
50  (vertex_.Y() - p.Y()) * std::cos(float(p4_.load()->Phi()));
51 }
52 float pat::PackedGenParticle::dz(const Point& p) const {
53  unpack();
54  return (vertex_.Z() - p.X()) - ((vertex_.X() - p.X()) * std::cos(float(p4_.load()->Phi())) +
55  (vertex_.Y() - p.Y()) * std::sin(float(p4_.load()->Phi()))) *
56  p4_.load()->Pz() / p4_.load()->Pt();
57 }
58 
60 
62  throw cms::Exception("Invalid Reference") << "this Candidate has no master clone reference."
63  << "Can't call masterClone() method.\n";
64 }
65 
66 bool pat::PackedGenParticle::hasMasterClone() const { return false; }
67 
68 bool pat::PackedGenParticle::hasMasterClonePtr() const { return false; }
69 
71  throw cms::Exception("Invalid Reference") << "this Candidate has no master clone ptr."
72  << "Can't call masterClonePtr() method.\n";
73 }
74 
75 size_t pat::PackedGenParticle::numberOfDaughters() const { return 0; }
76 
78  if (motherRef().isNonnull())
79  return 1;
80  return 0;
81 }
82 
84  return p4() == o.p4() && vertex() == o.vertex() && charge() == o.charge();
85  // return p4() == o.p4() && charge() == o.charge();
86 }
87 
89 
90 const reco::Candidate* pat::PackedGenParticle::mother(size_type) const { return motherRef().get(); }
91 
94  << "This Candidate type does not implement daughter(std::string). "
95  << "Please use CompositeCandidate or NamedCompositeCandidate.\n";
96 }
97 
100  << "This Candidate type does not implement daughter(std::string). "
101  << "Please use CompositeCandidate or NamedCompositeCandidate.\n";
102 }
103 
105 
106 double pat::PackedGenParticle::vertexChi2() const { return 0; }
107 
108 double pat::PackedGenParticle::vertexNdof() const { return 0; }
109 
111 
114  << "reco::ConcreteCandidate does not implement vertex covariant matrix.\n";
115 }
116 
119  << "reco::ConcreteCandidate does not implement vertex covariant matrix.\n";
120 }
121 
122 bool pat::PackedGenParticle::isElectron() const { return false; }
123 
124 bool pat::PackedGenParticle::isMuon() const { return false; }
125 
126 bool pat::PackedGenParticle::isGlobalMuon() const { return false; }
127 
128 bool pat::PackedGenParticle::isStandAloneMuon() const { return false; }
129 
130 bool pat::PackedGenParticle::isTrackerMuon() const { return false; }
131 
132 bool pat::PackedGenParticle::isCaloMuon() const { return false; }
133 
134 bool pat::PackedGenParticle::isPhoton() const { return false; }
135 
136 bool pat::PackedGenParticle::isConvertedPhoton() const { return false; }
137 
138 bool pat::PackedGenParticle::isJet() const { return false; }
139 
140 bool pat::PackedGenParticle::longLived() const { return false; }
141 
142 bool pat::PackedGenParticle::massConstraint() const { return false; }
bool hasMasterClone() const override
bool isStandAloneMuon() const override
std::pair< unsigned int, unsigned int > unpack(cond::Time_t since)
size_t size_type
Definition: Candidate.h:29
bool isCaloMuon() const override
bool longLived() const override
is long lived?
double vertexNdof() const override
double sign(double x)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
size_t numberOfDaughters() const override
number of daughters
bool isJet() const override
static float float16to32(uint16_t h)
Definition: libminifloat.h:13
std::atomic< LorentzVector * > p4c_
void pack(bool unpackAfterwards=true)
double vertexNormalizedChi2() const override
chi-squared divided by n.d.o.f.
bool hasMasterClonePtr() const override
void fillVertexCovariance(CovarianceMatrix &v) const override
fill SMatrix
size_t numberOfMothers() const override
number of mothers
double vertexChi2() const override
chi-squares
T sqrt(T t)
Definition: SSEVec.h:19
static uint16_t float32to16(float x)
Definition: libminifloat.h:17
bool overlap(const reco::Candidate &) const override
check overlap with another Candidate
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
~PackedGenParticle() override
destructor
virtual const Point & vertex() const =0
vertex position
virtual int charge() const =0
electric charge
bool isElectron() const override
get a component
const reco::CandidatePtr & masterClonePtr() const override
bool massConstraint() const override
do mass constraint?
const reco::CandidateBaseRef & masterClone() const override
std::atomic< PolarLorentzVector * > p4_
the four vector
bool isConvertedPhoton() const override
const reco::Candidate * daughter(size_type) const override
return daughter at a given position (throws an exception)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
const reco::Candidate * mother(size_type) const override
return mother at a given position (throws an exception)
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
bool isGlobalMuon() const override
T get() const
get a component
Definition: Candidate.h:221
static unsigned int const shift
bool isTrackerMuon() const override
bool isPhoton() const override
virtual float dz() const
dz with respect to the PV ref
CovarianceMatrix vertexCovariance() const override
return SMatrix
bool isMuon() const override
virtual float dxy() const
dxy with respect to the PV ref
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38