Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-08-12 23:11:27

0001 #ifndef __DataFormats_PatCandidates_PackedCandidate_h__
0002 #define __DataFormats_PatCandidates_PackedCandidate_h__
0003 
0004 #include "DataFormats/Candidate/interface/Candidate.h"
0005 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0006 #include "DataFormats/Common/interface/Association.h"
0007 #include "DataFormats/Common/interface/RefVector.h"
0008 #include "DataFormats/Math/interface/deltaPhi.h"
0009 #include "DataFormats/PatCandidates/interface/CovarianceParameterization.h"
0010 #include "DataFormats/VertexReco/interface/Vertex.h"
0011 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0012 #include "FWCore/Utilities/interface/thread_safety_macros.h"
0013 #include <atomic>
0014 #include <mutex>
0015 
0016 /* #include "DataFormats/Math/interface/PtEtaPhiMass.h" */
0017 
0018 // forward declare testing structure
0019 class testPackedCandidate;
0020 
0021 namespace pat {
0022   class PackedCandidate : public reco::Candidate {
0023   public:
0024     /// collection of daughter candidates
0025     typedef reco::CandidateCollection daughters;
0026     /// Lorentz vector
0027     typedef math::XYZTLorentzVector LorentzVector;
0028     /// Lorentz vector
0029     typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
0030     /// point in the space
0031     typedef math::XYZPoint Point;
0032     /// point in the space
0033     typedef math::XYZVector Vector;
0034 
0035     typedef unsigned int index;
0036     /// default constructor
0037     PackedCandidate()
0038         : packedPt_(0),
0039           packedEta_(0),
0040           packedPhi_(0),
0041           packedM_(0),
0042           packedDxy_(0),
0043           packedDz_(0),
0044           packedDPhi_(0),
0045           packedDEta_(0),
0046           packedDTrkPt_(0),
0047           packedCovariance_(),
0048           packedPuppiweight_(0),
0049           packedPuppiweightNoLepDiff_(0),
0050           rawCaloFraction_(0),
0051           rawHcalFraction_(0),
0052           caloFraction_(0),
0053           hcalFraction_(0),
0054           packedTime_(0),
0055           packedTimeError_(0),
0056           isIsolatedChargedHadron_(false),
0057           p4_(new PolarLorentzVector(0, 0, 0, 0)),
0058           p4c_(new LorentzVector(0, 0, 0, 0)),
0059           vertex_(new Point(0, 0, 0)),
0060           dphi_(0),
0061           deta_(0),
0062           dtrkpt_(0),
0063           track_(nullptr),
0064           pdgId_(0),
0065           qualityFlags_(0),
0066           pvRefKey_(reco::VertexRef::invalidKey()),
0067           m_(nullptr),
0068           packedHits_(0),
0069           packedLayers_(0),
0070           normalizedChi2_(0),
0071           covarianceVersion_(0),
0072           covarianceSchema_(0),
0073           firstHit_(0) {}
0074 
0075     explicit PackedCandidate(const reco::Candidate &c,
0076                              const reco::VertexRefProd &pvRefProd,
0077                              reco::VertexRef::key_type pvRefKey)
0078         : packedPuppiweight_(0),
0079           packedPuppiweightNoLepDiff_(0),
0080           rawCaloFraction_(0),
0081           rawHcalFraction_(0),
0082           caloFraction_(0),
0083           hcalFraction_(0),
0084           packedTime_(0),
0085           packedTimeError_(0),
0086           isIsolatedChargedHadron_(false),
0087           p4_(new PolarLorentzVector(c.pt(), c.eta(), c.phi(), c.mass())),
0088           p4c_(new LorentzVector(*p4_)),
0089           vertex_(new Point(c.vertex())),
0090           dphi_(0),
0091           deta_(0),
0092           dtrkpt_(0),
0093           track_(nullptr),
0094           pdgId_(c.pdgId()),
0095           qualityFlags_(0),
0096           pvRefProd_(pvRefProd),
0097           pvRefKey_(pvRefKey),
0098           m_(nullptr),
0099           packedHits_(0),
0100           packedLayers_(0),
0101           normalizedChi2_(0),
0102           covarianceVersion_(0),
0103           covarianceSchema_(0),
0104           firstHit_(0) {
0105       packBoth();
0106     }
0107 
0108     explicit PackedCandidate(const PolarLorentzVector &p4,
0109                              const Point &vtx,
0110                              float trkPt,
0111                              float etaAtVtx,
0112                              float phiAtVtx,
0113                              int pdgId,
0114                              const reco::VertexRefProd &pvRefProd,
0115                              reco::VertexRef::key_type pvRefKey)
0116         : packedPuppiweight_(0),
0117           packedPuppiweightNoLepDiff_(0),
0118           rawCaloFraction_(0),
0119           rawHcalFraction_(0),
0120           caloFraction_(0),
0121           hcalFraction_(0),
0122           packedTime_(0),
0123           packedTimeError_(0),
0124           isIsolatedChargedHadron_(false),
0125           p4_(new PolarLorentzVector(p4)),
0126           p4c_(new LorentzVector(*p4_)),
0127           vertex_(new Point(vtx)),
0128           dphi_(reco::deltaPhi(phiAtVtx, p4_.load()->phi())),
0129           deta_(std::abs(etaAtVtx - p4_.load()->eta()) >= kMinDEtaToStore_ ? etaAtVtx - p4_.load()->eta() : 0.),
0130           dtrkpt_(std::abs(trkPt - p4_.load()->pt()) >= kMinDTrkPtToStore_ ? trkPt - p4_.load()->pt() : 0.),
0131           track_(nullptr),
0132           pdgId_(pdgId),
0133           qualityFlags_(0),
0134           pvRefProd_(pvRefProd),
0135           pvRefKey_(pvRefKey),
0136           m_(nullptr),
0137           packedHits_(0),
0138           packedLayers_(0),
0139           normalizedChi2_(0),
0140           covarianceVersion_(0),
0141           covarianceSchema_(0),
0142           firstHit_(0) {
0143       packBoth();
0144     }
0145 
0146     explicit PackedCandidate(const LorentzVector &p4,
0147                              const Point &vtx,
0148                              float trkPt,
0149                              float etaAtVtx,
0150                              float phiAtVtx,
0151                              int pdgId,
0152                              const reco::VertexRefProd &pvRefProd,
0153                              reco::VertexRef::key_type pvRefKey)
0154         : packedPuppiweight_(0),
0155           packedPuppiweightNoLepDiff_(0),
0156           rawCaloFraction_(0),
0157           rawHcalFraction_(0),
0158           caloFraction_(0),
0159           hcalFraction_(0),
0160           packedTime_(0),
0161           packedTimeError_(0),
0162           isIsolatedChargedHadron_(false),
0163           p4_(new PolarLorentzVector(p4.Pt(), p4.Eta(), p4.Phi(), p4.M())),
0164           p4c_(new LorentzVector(p4)),
0165           vertex_(new Point(vtx)),
0166           dphi_(reco::deltaPhi(phiAtVtx, p4_.load()->phi())),
0167           deta_(std::abs(etaAtVtx - p4_.load()->eta()) >= kMinDEtaToStore_ ? etaAtVtx - p4_.load()->eta() : 0.),
0168           dtrkpt_(std::abs(trkPt - p4_.load()->pt()) >= kMinDTrkPtToStore_ ? trkPt - p4_.load()->pt() : 0.),
0169           track_(nullptr),
0170           pdgId_(pdgId),
0171           qualityFlags_(0),
0172           pvRefProd_(pvRefProd),
0173           pvRefKey_(pvRefKey),
0174           m_(nullptr),
0175           packedHits_(0),
0176           packedLayers_(0),
0177           normalizedChi2_(0),
0178           covarianceVersion_(0),
0179           covarianceSchema_(0),
0180           firstHit_(0) {
0181       packBoth();
0182     }
0183 
0184     PackedCandidate(const PackedCandidate &iOther)
0185         : packedPt_(iOther.packedPt_),
0186           packedEta_(iOther.packedEta_),
0187           packedPhi_(iOther.packedPhi_),
0188           packedM_(iOther.packedM_),
0189           packedDxy_(iOther.packedDxy_),
0190           packedDz_(iOther.packedDz_),
0191           packedDPhi_(iOther.packedDPhi_),
0192           packedDEta_(iOther.packedDEta_),
0193           packedDTrkPt_(iOther.packedDTrkPt_),
0194           packedCovariance_(iOther.packedCovariance_),
0195           packedPuppiweight_(iOther.packedPuppiweight_),
0196           packedPuppiweightNoLepDiff_(iOther.packedPuppiweightNoLepDiff_),
0197           rawCaloFraction_(iOther.rawCaloFraction_),
0198           rawHcalFraction_(iOther.rawHcalFraction_),
0199           caloFraction_(iOther.caloFraction_),
0200           hcalFraction_(iOther.hcalFraction_),
0201           packedTime_(iOther.packedTime_),
0202           packedTimeError_(iOther.packedTimeError_),
0203           isIsolatedChargedHadron_(iOther.isIsolatedChargedHadron_),
0204           // Need to trigger unpacking in iOther
0205           p4_(new PolarLorentzVector(iOther.polarP4())),
0206           p4c_(new LorentzVector(iOther.p4())),
0207           vertex_((iOther.vertex_ ? new Point(iOther.vertex()) : nullptr)),
0208           dxy_(vertex_ ? iOther.dxy_ : 0),
0209           dz_(vertex_ ? iOther.dz_ : 0),
0210           dphi_(vertex_ ? iOther.dphi_ : 0),
0211           deta_(vertex_ ? iOther.deta_ : 0),
0212           dtrkpt_(vertex_ ? iOther.dtrkpt_ : 0),
0213           track_(iOther.track_ ? new reco::Track(*iOther.track_) : nullptr),
0214           pdgId_(iOther.pdgId_),
0215           qualityFlags_(iOther.qualityFlags_),
0216           pvRefProd_(iOther.pvRefProd_),
0217           pvRefKey_(iOther.pvRefKey_),
0218           m_(iOther.m_ ? new reco::TrackBase::CovarianceMatrix(*iOther.m_) : nullptr),
0219           packedHits_(iOther.packedHits_),
0220           packedLayers_(iOther.packedLayers_),
0221           normalizedChi2_(iOther.normalizedChi2_),
0222           covarianceVersion_(iOther.covarianceVersion_),
0223           covarianceSchema_(iOther.covarianceSchema_),
0224           firstHit_(iOther.firstHit_) {}
0225 
0226     PackedCandidate(PackedCandidate &&iOther)
0227         : packedPt_(iOther.packedPt_),
0228           packedEta_(iOther.packedEta_),
0229           packedPhi_(iOther.packedPhi_),
0230           packedM_(iOther.packedM_),
0231           packedDxy_(iOther.packedDxy_),
0232           packedDz_(iOther.packedDz_),
0233           packedDPhi_(iOther.packedDPhi_),
0234           packedDEta_(iOther.packedDEta_),
0235           packedDTrkPt_(iOther.packedDTrkPt_),
0236           packedCovariance_(iOther.packedCovariance_),
0237           packedPuppiweight_(iOther.packedPuppiweight_),
0238           packedPuppiweightNoLepDiff_(iOther.packedPuppiweightNoLepDiff_),
0239           rawCaloFraction_(iOther.rawCaloFraction_),
0240           rawHcalFraction_(iOther.rawHcalFraction_),
0241           caloFraction_(iOther.caloFraction_),
0242           hcalFraction_(iOther.hcalFraction_),
0243           packedTime_(iOther.packedTime_),
0244           packedTimeError_(iOther.packedTimeError_),
0245           isIsolatedChargedHadron_(iOther.isIsolatedChargedHadron_),
0246           p4_(iOther.p4_.exchange(nullptr)),
0247           p4c_(iOther.p4c_.exchange(nullptr)),
0248           vertex_(iOther.vertex_.exchange(nullptr)),
0249           dxy_(iOther.dxy_),
0250           dz_(iOther.dz_),
0251           dphi_(iOther.dphi_),
0252           deta_(iOther.deta_),
0253           dtrkpt_(iOther.dtrkpt_),
0254           track_(iOther.track_.exchange(nullptr)),
0255           pdgId_(iOther.pdgId_),
0256           qualityFlags_(iOther.qualityFlags_),
0257           pvRefProd_(iOther.pvRefProd_),
0258           pvRefKey_(iOther.pvRefKey_),
0259           m_(iOther.m_.exchange(nullptr)),
0260           packedHits_(iOther.packedHits_),
0261           packedLayers_(iOther.packedLayers_),
0262           normalizedChi2_(iOther.normalizedChi2_),
0263           covarianceVersion_(iOther.covarianceVersion_),
0264           covarianceSchema_(iOther.covarianceSchema_),
0265           firstHit_(iOther.firstHit_) {}
0266 
0267     PackedCandidate &operator=(const PackedCandidate &iOther) {
0268       if (this == &iOther) {
0269         return *this;
0270       }
0271       packedPt_ = iOther.packedPt_;
0272       packedEta_ = iOther.packedEta_;
0273       packedPhi_ = iOther.packedPhi_;
0274       packedM_ = iOther.packedM_;
0275       packedDxy_ = iOther.packedDxy_;
0276       packedDz_ = iOther.packedDz_;
0277       packedDPhi_ = iOther.packedDPhi_;
0278       packedDEta_ = iOther.packedDEta_;
0279       packedDTrkPt_ = iOther.packedDTrkPt_;
0280       packedCovariance_ = iOther.packedCovariance_;
0281       packedPuppiweight_ = iOther.packedPuppiweight_;
0282       packedPuppiweightNoLepDiff_ = iOther.packedPuppiweightNoLepDiff_;
0283       rawCaloFraction_ = iOther.rawCaloFraction_;
0284       rawHcalFraction_ = iOther.rawHcalFraction_;
0285       caloFraction_ = iOther.caloFraction_;
0286       hcalFraction_ = iOther.hcalFraction_;
0287       packedTime_ = iOther.packedTime_;
0288       packedTimeError_ = iOther.packedTimeError_;
0289       isIsolatedChargedHadron_ = iOther.isIsolatedChargedHadron_;
0290       // Need to trigger unpacking in iOther
0291       if (p4_) {
0292         *p4_ = iOther.polarP4();
0293       } else {
0294         p4_.store(new PolarLorentzVector(iOther.polarP4()));
0295       }
0296       if (p4c_) {
0297         *p4c_ = iOther.p4();
0298       } else {
0299         p4c_.store(new LorentzVector(iOther.p4()));
0300       }
0301       if (vertex_) {
0302         *vertex_ = iOther.vertex();
0303       } else {
0304         vertex_.store(new Point(iOther.vertex()));
0305       }
0306       dxy_ = iOther.dxy_;
0307       dz_ = iOther.dz_;
0308       dphi_ = iOther.dphi_;
0309       deta_ = iOther.deta_;
0310       dtrkpt_ = iOther.dtrkpt_;
0311 
0312       if (!iOther.track_) {
0313         delete track_.exchange(nullptr);
0314       } else {
0315         if (!track_) {
0316           track_.store(new reco::Track(*iOther.track_));
0317         } else {
0318           *track_ = *(iOther.track_);
0319         }
0320       }
0321 
0322       pdgId_ = iOther.pdgId_;
0323       qualityFlags_ = iOther.qualityFlags_;
0324       pvRefProd_ = iOther.pvRefProd_;
0325       pvRefKey_ = iOther.pvRefKey_;
0326       if (!iOther.m_) {
0327         delete m_.exchange(nullptr);
0328       } else {
0329         if (!m_) {
0330           m_.store(new reco::Track::CovarianceMatrix(*iOther.m_));
0331         } else {
0332           *m_ = *(iOther.m_);
0333         }
0334       }
0335 
0336       packedHits_ = iOther.packedHits_;
0337       packedLayers_ = iOther.packedLayers_;
0338       normalizedChi2_ = iOther.normalizedChi2_;
0339       covarianceVersion_ = iOther.covarianceVersion_;
0340       covarianceSchema_ = iOther.covarianceSchema_;
0341       firstHit_ = iOther.firstHit_;
0342       return *this;
0343     }
0344 
0345     PackedCandidate &operator=(PackedCandidate &&iOther) {
0346       if (this == &iOther) {
0347         return *this;
0348       }
0349       packedPt_ = iOther.packedPt_;
0350       packedEta_ = iOther.packedEta_;
0351       packedPhi_ = iOther.packedPhi_;
0352       packedM_ = iOther.packedM_;
0353       packedDxy_ = iOther.packedDxy_;
0354       packedDz_ = iOther.packedDz_;
0355       packedDPhi_ = iOther.packedDPhi_;
0356       packedDEta_ = iOther.packedDEta_;
0357       packedDTrkPt_ = iOther.packedDTrkPt_;
0358       packedCovariance_ = iOther.packedCovariance_;
0359       packedPuppiweight_ = iOther.packedPuppiweight_;
0360       packedPuppiweightNoLepDiff_ = iOther.packedPuppiweightNoLepDiff_;
0361       rawCaloFraction_ = iOther.rawCaloFraction_;
0362       rawHcalFraction_ = iOther.rawHcalFraction_;
0363       caloFraction_ = iOther.caloFraction_;
0364       hcalFraction_ = iOther.hcalFraction_;
0365       packedTime_ = iOther.packedTime_;
0366       packedTimeError_ = iOther.packedTimeError_;
0367       isIsolatedChargedHadron_ = iOther.isIsolatedChargedHadron_;
0368       delete p4_.exchange(iOther.p4_.exchange(nullptr));
0369       delete p4c_.exchange(iOther.p4c_.exchange(nullptr));
0370       delete vertex_.exchange(iOther.vertex_.exchange(nullptr));
0371       dxy_ = iOther.dxy_;
0372       dz_ = iOther.dz_;
0373       dphi_ = iOther.dphi_;
0374       deta_ = iOther.deta_;
0375       dtrkpt_ = iOther.dtrkpt_;
0376       delete track_.exchange(iOther.track_.exchange(nullptr));
0377       pdgId_ = iOther.pdgId_;
0378       qualityFlags_ = iOther.qualityFlags_;
0379       pvRefProd_ = iOther.pvRefProd_;
0380       pvRefKey_ = iOther.pvRefKey_;
0381       delete m_.exchange(iOther.m_.exchange(nullptr));
0382       packedHits_ = iOther.packedHits_;
0383       packedLayers_ = iOther.packedLayers_;
0384       normalizedChi2_ = iOther.normalizedChi2_;
0385       covarianceVersion_ = iOther.covarianceVersion_;
0386       covarianceSchema_ = iOther.covarianceSchema_;
0387       firstHit_ = iOther.firstHit_;
0388       return *this;
0389     }
0390 
0391     /// destructor
0392     ~PackedCandidate() override;
0393     /// number of daughters
0394     size_t numberOfDaughters() const override;
0395     /// return daughter at a given position (throws an exception)
0396     const reco::Candidate *daughter(size_type) const override;
0397     /// number of mothers
0398     size_t numberOfMothers() const override;
0399     /// return mother at a given position (throws an exception)
0400     const reco::Candidate *mother(size_type) const override;
0401     /// return daughter at a given position (throws an exception)
0402     reco::Candidate *daughter(size_type) override;
0403     /// return daughter with a specified role name
0404     reco::Candidate *daughter(const std::string &s) override;
0405     /// return daughter with a specified role name
0406     const reco::Candidate *daughter(const std::string &s) const override;
0407     /// return the number of source Candidates
0408     /// ( the candidates used to construct this Candidate)
0409     size_t numberOfSourceCandidatePtrs() const override { return 0; }
0410     /// return a Ptr to one of the source Candidates
0411     /// ( the candidates used to construct this Candidate)
0412     reco::CandidatePtr sourceCandidatePtr(size_type i) const override { return reco::CandidatePtr(); }
0413 
0414     /// electric charge
0415     int charge() const override {
0416       switch (abs(pdgId_)) {
0417         case 211:
0418           return (pdgId_ > 0) - (pdgId_ < 0);
0419         case 11:
0420           return (-1) * (pdgId_ > 0) + (pdgId_ < 0);  // e
0421         case 13:
0422           return (-1) * (pdgId_ > 0) + (pdgId_ < 0);  // mu
0423         case 15:
0424           return (-1) * (pdgId_ > 0) + (pdgId_ < 0);  // tau
0425         case 24:
0426           return (pdgId_ > 0) - (pdgId_ < 0);  // W
0427         default:
0428           return 0;  // FIXME: charge is not defined
0429       }
0430     }
0431     /// set electric charge
0432     void setCharge(int charge) override {}
0433     /// electric charge
0434     int threeCharge() const override { return charge() * 3; }
0435     /// set electric charge
0436     void setThreeCharge(int threecharge) override {}
0437     /// four-momentum Lorentz vecto r
0438     const LorentzVector &p4() const override {
0439       if (!p4c_)
0440         unpack();
0441       return *p4c_;
0442     }
0443     /// four-momentum Lorentz vector
0444     const PolarLorentzVector &polarP4() const override {
0445       if (!p4c_)
0446         unpack();
0447       return *p4_;
0448     }
0449     /// spatial momentum vector
0450     Vector momentum() const override {
0451       if (!p4c_)
0452         unpack();
0453       return p4c_.load()->Vect();
0454     }
0455     /// boost vector to boost a Lorentz vector
0456     /// to the particle center of mass system
0457     Vector boostToCM() const override {
0458       if (!p4c_)
0459         unpack();
0460       return p4c_.load()->BoostToCM();
0461     }
0462     /// magnitude of momentum vector
0463     double p() const override {
0464       if (!p4c_)
0465         unpack();
0466       return p4c_.load()->P();
0467     }
0468     /// energy
0469     double energy() const override {
0470       if (!p4c_)
0471         unpack();
0472       return p4c_.load()->E();
0473     }
0474     /// transverse energy
0475     double et() const override { return (pt() <= 0) ? 0 : p4c_.load()->Et(); }
0476     /// transverse energy squared (use this for cuts)!
0477     double et2() const override { return (pt() <= 0) ? 0 : p4c_.load()->Et2(); }
0478     /// mass
0479     double mass() const override {
0480       if (!p4c_)
0481         unpack();
0482       return p4_.load()->M();
0483     }
0484     /// mass squared
0485     double massSqr() const override {
0486       if (!p4c_)
0487         unpack();
0488       auto m = p4_.load()->M();
0489       return m * m;
0490     }
0491 
0492     /// transverse mass
0493     double mt() const override {
0494       if (!p4c_)
0495         unpack();
0496       return p4_.load()->Mt();
0497     }
0498     /// transverse mass squared
0499     double mtSqr() const override {
0500       if (!p4c_)
0501         unpack();
0502       return p4_.load()->Mt2();
0503     }
0504     /// x coordinate of momentum vector
0505     double px() const override {
0506       if (!p4c_)
0507         unpack();
0508       return p4c_.load()->Px();
0509     }
0510     /// y coordinate of momentum vector
0511     double py() const override {
0512       if (!p4c_)
0513         unpack();
0514       return p4c_.load()->Py();
0515     }
0516     /// z coordinate of momentum vector
0517     double pz() const override {
0518       if (!p4c_)
0519         unpack();
0520       return p4c_.load()->Pz();
0521     }
0522     /// transverse momentum
0523     double pt() const override {
0524       if (!p4c_)
0525         unpack();
0526       return p4_.load()->Pt();
0527     }
0528     /// momentum azimuthal angle
0529     double phi() const override {
0530       if (!p4c_)
0531         unpack();
0532       return p4_.load()->Phi();
0533     }
0534 
0535     /// pt from the track (normally identical to pt())
0536     virtual double ptTrk() const {
0537       maybeUnpackBoth();
0538       return p4_.load()->pt() + dtrkpt_;
0539     }
0540     /// momentum azimuthal angle from the track (normally identical to phi())
0541     virtual float phiAtVtx() const {
0542       maybeUnpackBoth();
0543       float ret = p4_.load()->Phi() + dphi_;
0544       while (ret > float(M_PI))
0545         ret -= 2 * float(M_PI);
0546       while (ret < -float(M_PI))
0547         ret += 2 * float(M_PI);
0548       return ret;
0549     }
0550     /// eta from the track (normally identical to eta())
0551     virtual float etaAtVtx() const {
0552       maybeUnpackBoth();
0553       return p4_.load()->eta() + deta_;
0554     }
0555 
0556     /// momentum polar angle
0557     double theta() const override {
0558       if (!p4c_)
0559         unpack();
0560       return p4_.load()->Theta();
0561     }
0562     /// momentum pseudorapidity
0563     double eta() const override {
0564       if (!p4c_)
0565         unpack();
0566       return p4_.load()->Eta();
0567     }
0568     /// rapidity
0569     double rapidity() const override {
0570       if (!p4c_)
0571         unpack();
0572       return p4_.load()->Rapidity();
0573     }
0574     /// rapidity
0575     double y() const override {
0576       if (!p4c_)
0577         unpack();
0578       return p4_.load()->Rapidity();
0579     }
0580     /// set 4-momentum
0581     void setP4(const LorentzVector &p4) override {
0582       maybeUnpackBoth();  // changing px,py,pz changes also mapping between dxy,dz
0583                           // and x,y,z
0584       dphi_ += polarP4().Phi() - p4.Phi();
0585       deta_ += polarP4().Eta() - p4.Eta();
0586       dtrkpt_ += polarP4().Pt() - p4.Pt();
0587       *p4_ = PolarLorentzVector(p4.Pt(), p4.Eta(), p4.Phi(), p4.M());
0588       packBoth();
0589     }
0590     /// set 4-momentum
0591     void setP4(const PolarLorentzVector &p4) override {
0592       maybeUnpackBoth();  // changing px,py,pz changes also mapping between dxy,dz
0593                           // and x,y,z
0594       dphi_ += polarP4().Phi() - p4.Phi();
0595       deta_ += polarP4().Eta() - p4.Eta();
0596       dtrkpt_ += polarP4().Pt() - p4.Pt();
0597       *p4_ = p4;
0598       packBoth();
0599     }
0600     /// set particle mass
0601     void setMass(double m) override {
0602       if (!p4c_)
0603         unpack();
0604       *p4_ = PolarLorentzVector(p4_.load()->Pt(), p4_.load()->Eta(), p4_.load()->Phi(), m);
0605       pack();
0606     }
0607     void setPz(double pz) override {
0608       maybeUnpackBoth();  // changing px,py,pz changes also mapping between dxy,dz
0609                           // and x,y,z
0610       *p4c_ = LorentzVector(p4c_.load()->Px(), p4c_.load()->Py(), pz, p4c_.load()->E());
0611       dphi_ += polarP4().Phi() - p4c_.load()->Phi();
0612       deta_ += polarP4().Eta() - p4c_.load()->Eta();
0613       dtrkpt_ += polarP4().Pt() - p4c_.load()->Pt();
0614       *p4_ = PolarLorentzVector(p4c_.load()->Pt(), p4c_.load()->Eta(), p4c_.load()->Phi(), p4c_.load()->M());
0615       packBoth();
0616     }
0617     /// set impact parameters covariance
0618 
0619     // Note: mask is also the maximum value
0620     enum trackHitShiftsAndMasks { trackPixelHitsMask = 7, trackStripHitsMask = 31, trackStripHitsShift = 3 };
0621 
0622     // set number of tracker hits and layers
0623     virtual void setHits(const reco::Track &tk) {
0624       // first we count the number of layers with hits
0625       int numberOfPixelLayers_ = tk.hitPattern().pixelLayersWithMeasurement();
0626       if (numberOfPixelLayers_ > trackPixelHitsMask)
0627         numberOfPixelLayers_ = trackPixelHitsMask;
0628       int numberOfStripLayers_ = tk.hitPattern().stripLayersWithMeasurement();
0629       if (numberOfStripLayers_ > trackStripHitsMask)
0630         numberOfStripLayers_ = trackStripHitsMask;
0631       packedLayers_ = (numberOfPixelLayers_ & trackPixelHitsMask) | (numberOfStripLayers_ << trackStripHitsShift);
0632       // now we count number of additional hits, beyond the one-per-layer implied
0633       // by packedLayers_
0634       int numberOfPixelHits_ = tk.hitPattern().numberOfValidPixelHits() - numberOfPixelLayers_;
0635       if (numberOfPixelHits_ > trackPixelHitsMask)
0636         numberOfPixelHits_ = trackPixelHitsMask;
0637       int numberOfStripHits_ =
0638           tk.hitPattern().numberOfValidHits() - numberOfPixelHits_ - numberOfPixelLayers_ - numberOfStripLayers_;
0639       if (numberOfStripHits_ > trackStripHitsMask)
0640         numberOfStripHits_ = trackStripHitsMask;
0641 
0642       packedHits_ = (numberOfPixelHits_ & trackPixelHitsMask) | (numberOfStripHits_ << trackStripHitsShift);
0643     }
0644 
0645     virtual void setTrackProperties(const reco::Track &tk,
0646                                     const reco::Track::CovarianceMatrix &covariance,
0647                                     int quality,
0648                                     int covarianceVersion) {
0649       covarianceVersion_ = covarianceVersion;
0650       covarianceSchema_ = quality;
0651       normalizedChi2_ = tk.normalizedChi2();
0652       setHits(tk);
0653       maybeUnpackBoth();
0654       packBoth();
0655       packCovariance(covariance, false);
0656     }
0657 
0658     // set track properties using quality and covarianceVersion to define the
0659     // level of details in the cov. matrix
0660     virtual void setTrackProperties(const reco::Track &tk, int quality, int covarianceVersion) {
0661       setTrackProperties(tk, tk.covariance(), quality, covarianceVersion);
0662     }
0663 
0664     void setTrackPropertiesLite(unsigned int covSchema,
0665                                 unsigned int covarianceVersion,
0666                                 unsigned int nHits,
0667                                 unsigned int nPixelHits) {
0668       covarianceVersion_ = covarianceVersion;
0669       covarianceSchema_ = covSchema;
0670       packedHits_ =
0671           (nPixelHits & trackPixelHitsMask) | (((nHits - nPixelHits) & trackStripHitsMask) << trackStripHitsShift);
0672     }
0673 
0674     int numberOfPixelHits() const { return (packedHits_ & trackPixelHitsMask) + pixelLayersWithMeasurement(); }
0675     int numberOfHits() const {
0676       return (packedHits_ >> trackStripHitsShift) + stripLayersWithMeasurement() + numberOfPixelHits();
0677     }
0678     int pixelLayersWithMeasurement() const { return packedLayers_ & trackPixelHitsMask; }
0679     int stripLayersWithMeasurement() const { return (packedLayers_ >> trackStripHitsShift); }
0680     int trackerLayersWithMeasurement() const { return stripLayersWithMeasurement() + pixelLayersWithMeasurement(); }
0681     virtual void setCovarianceVersion(int v) { covarianceVersion_ = v; }
0682     int covarianceVersion() const { return covarianceVersion_; }
0683     int covarianceSchema() const { return covarianceSchema_; }
0684 
0685     /// vertex position
0686     const Point &vertex() const override {
0687       maybeUnpackBoth();
0688       return *vertex_;
0689     }  //{ if (fromPV_) return Point(0,0,0); else return Point(0,0,100); }
0690     /// x coordinate of vertex position
0691     double vx() const override {
0692       maybeUnpackBoth();
0693       return vertex_.load()->X();
0694     }  //{ return 0; }
0695     /// y coordinate of vertex position
0696     double vy() const override {
0697       maybeUnpackBoth();
0698       return vertex_.load()->Y();
0699     }  //{ return 0; }
0700     /// z coordinate of vertex position
0701     double vz() const override {
0702       maybeUnpackBoth();
0703       return vertex_.load()->Z();
0704     }  //{ if (fromPV_) return 0; else return 100; }
0705     /// set vertex
0706     void setVertex(const Point &vertex) override {
0707       maybeUnpackBoth();
0708       *vertex_ = vertex;
0709       packVtx();
0710     }
0711 
0712     /// This refers to the association to PV=ipv. >=PVLoose corresponds to JME
0713     /// definition, >=PVTight to isolation definition
0714     enum PVAssoc { NoPV = 0, PVLoose = 1, PVTight = 2, PVUsedInFit = 3 };
0715     const PVAssoc fromPV(size_t ipv = 0) const {
0716       reco::VertexRef pvRef = vertexRef();
0717       if (pvAssociationQuality() == UsedInFitTight and pvRef.key() == ipv)
0718         return PVUsedInFit;
0719       if (pvRef.key() == ipv or abs(pdgId()) == 13 or abs(pdgId()) == 11)
0720         return PVTight;
0721       if (pvAssociationQuality() == CompatibilityBTag and std::abs(dzAssociatedPV()) > std::abs(dz(ipv)))
0722         return PVTight;  // it is not closest, but at least prevents the B
0723                          // assignment stealing
0724       if (pvAssociationQuality() < UsedInFitLoose or pvRef->ndof() < 4.0)
0725         return PVLoose;
0726       return NoPV;
0727     }
0728 
0729     /// The following contains information about how the association to the PV,
0730     /// given in vertexRef, is obtained.
0731     ///
0732     enum PVAssociationQuality {
0733       NotReconstructedPrimary = 0,
0734       OtherDeltaZ = 1,
0735       CompatibilityBTag = 4,
0736       CompatibilityDz = 5,
0737       UsedInFitLoose = 6,
0738       UsedInFitTight = 7
0739     };
0740     const PVAssociationQuality pvAssociationQuality() const {
0741       return PVAssociationQuality((qualityFlags_ & assignmentQualityMask) >> assignmentQualityShift);
0742     }
0743     void setAssociationQuality(PVAssociationQuality q) {
0744       qualityFlags_ =
0745           (qualityFlags_ & ~assignmentQualityMask) | ((q << assignmentQualityShift) & assignmentQualityMask);
0746     }
0747 
0748     const reco::VertexRef vertexRef() const { return reco::VertexRef(pvRefProd_, pvRefKey_); }
0749 
0750     /// dxy with respect to the PV ref
0751     virtual float dxy() const {
0752       maybeUnpackBoth();
0753       return dxy_;
0754     }
0755 
0756     /// dz with respect to the PV[ipv]
0757     virtual float dz(size_t ipv = 0) const {
0758       maybeUnpackBoth();
0759       return dz_ + (*pvRefProd_)[pvRefKey_].position().z() - (*pvRefProd_)[ipv].position().z();
0760     }
0761     /// dz with respect to the PV ref
0762     virtual float dzAssociatedPV() const {
0763       maybeUnpackBoth();
0764       return dz_;
0765     }
0766 
0767     /// dxy with respect to another point
0768     virtual float dxy(const Point &p) const;
0769     /// dz  with respect to another point
0770     virtual float dz(const Point &p) const;
0771 
0772     /// uncertainty on dz
0773     float dzError() const override {
0774       maybeUnpackCovariance();
0775       return sqrt((*m_.load())(4, 4));
0776     }
0777     /// uncertainty on dxy
0778     float dxyError() const override {
0779       maybeUnpackCovariance();
0780       return sqrt((*m_.load())(3, 3));
0781     }
0782 
0783     /// Return reference to a pseudo track made with candidate kinematics,
0784     /// parameterized error for eta,phi,pt and full IP covariance
0785     virtual const reco::Track &pseudoTrack() const {
0786       if (!track_)
0787         unpackTrk();
0788       return *track_;
0789     }
0790 
0791     /// return a pointer to the track if present. otherwise, return a null pointer
0792     const reco::Track *bestTrack() const override {
0793       if (packedHits_ != 0 || packedLayers_ != 0) {
0794         maybeUnpackTrack();
0795         return track_.load();
0796       } else
0797         return nullptr;
0798     }
0799     /// Return true if a bestTrack can be extracted from this Candidate
0800     bool hasTrackDetails() const { return (packedHits_ != 0 || packedLayers_ != 0); }
0801     /// Return true if the original candidate had a track associated
0802     /// even if the PackedCandidate has no track
0803     bool fromTrackCandidate() const { return (packedDz_ != 0 || (packedDxy_ != 0 && packedDxy_ != 32768)); }
0804     /// true if the track had the highPurity quality bit
0805     bool trackHighPurity() const { return (qualityFlags_ & trackHighPurityMask) >> trackHighPurityShift; }
0806     /// set to true if the track had the highPurity quality bit
0807     void setTrackHighPurity(bool highPurity) {
0808       qualityFlags_ =
0809           (qualityFlags_ & ~trackHighPurityMask) | ((highPurity << trackHighPurityShift) & trackHighPurityMask);
0810     }
0811 
0812     /// Enumerator specifying the
0813     enum LostInnerHits {
0814       validHitInFirstPixelBarrelLayer = -1,
0815       noLostInnerHits = 0,  // it could still not have a hit in the first layer,
0816                             // e.g. if it crosses an inactive sensor
0817       oneLostInnerHit = 1,
0818       moreLostInnerHits = 2
0819     };
0820     LostInnerHits lostInnerHits() const {
0821       return LostInnerHits(int16_t((qualityFlags_ & lostInnerHitsMask) >> lostInnerHitsShift) - 1);
0822     }
0823     void setLostInnerHits(LostInnerHits hits) {
0824       int lost = hits;
0825       if (lost > 2)
0826         lost = 2;  // protection against misuse
0827       lost++;      // shift so it's 0 .. 3 instead of (-1) .. 2
0828       qualityFlags_ = (qualityFlags_ & ~lostInnerHitsMask) | ((lost << lostInnerHitsShift) & lostInnerHitsMask);
0829     }
0830 
0831     /// Set first hit from HitPattern
0832     void setFirstHit(uint16_t pattern) { firstHit_ = pattern; }
0833     /// Return first hit from HitPattern for tracks with high level details
0834     uint16_t firstHit() const { return firstHit_; }
0835 
0836     void setMuonID(bool isStandAlone, bool isGlobal) {
0837       int16_t muonFlags = isStandAlone | (2 * isGlobal);
0838       qualityFlags_ = (qualityFlags_ & ~muonFlagsMask) | ((muonFlags << muonFlagsShift) & muonFlagsMask);
0839     }
0840 
0841     void setGoodEgamma(bool isGoodEgamma = true) {
0842       int16_t egFlags = (isGoodEgamma << egammaFlagsShift) & egammaFlagsMask;
0843       qualityFlags_ = (qualityFlags_ & ~egammaFlagsMask) | egFlags;
0844     }
0845 
0846     /// PDG identifier
0847     int pdgId() const override { return pdgId_; }
0848     // set PDG identifier
0849     void setPdgId(int pdgId) override { pdgId_ = pdgId; }
0850     /// status word
0851     int status() const override { return qualityFlags_; } /*FIXME*/
0852     /// set status word
0853     void setStatus(int status) override {} /*FIXME*/
0854     /// long lived flag
0855     static const unsigned int longLivedTag = 0; /*FIXME*/
0856     /// set long lived flag
0857     void setLongLived() override {} /*FIXME*/
0858     /// is long lived?
0859     bool longLived() const override;
0860     /// do mass constraint flag
0861     static const unsigned int massConstraintTag = 0; /*FIXME*/
0862     /// set mass constraint flag
0863     void setMassConstraint() override {} /*FIXME*/
0864     /// do mass constraint?
0865     bool massConstraint() const override;
0866 
0867     /// returns a clone of the Candidate object
0868     PackedCandidate *clone() const override { return new PackedCandidate(*this); }
0869 
0870     /// chi-squares
0871     double vertexChi2() const override;
0872     /** Number of degrees of freedom
0873    *  Meant to be Double32_t for soft-assignment fitters:
0874    *  tracks may contribute to the vertex with fractional weights.
0875    *  The ndof is then = to the sum of the track weights.
0876    *  see e.g. CMS NOTE-2006/032, CMS NOTE-2004/002
0877    */
0878     double vertexNdof() const override;
0879     /// chi-squared divided by n.d.o.f.
0880     double vertexNormalizedChi2() const override;
0881     /// (i, j)-th element of error matrix, i, j = 0, ... 2
0882     double vertexCovariance(int i, int j) const override;
0883     /// return SMatrix
0884     CovarianceMatrix vertexCovariance() const override {
0885       CovarianceMatrix m;
0886       fillVertexCovariance(m);
0887       return m;
0888     }
0889     /// fill SMatrix
0890     void fillVertexCovariance(CovarianceMatrix &v) const override;
0891     /// returns true if this candidate has a reference to a master clone.
0892     /// This only happens if the concrete Candidate type is ShallowCloneCandidate
0893     bool hasMasterClone() const override;
0894     /// returns ptr to master clone, if existing.
0895     /// Throws an exception unless the concrete Candidate type is
0896     /// ShallowCloneCandidate
0897     const reco::CandidateBaseRef &masterClone() const override;
0898     /// returns true if this candidate has a ptr to a master clone.
0899     /// This only happens if the concrete Candidate type is
0900     /// ShallowClonePtrCandidate
0901     bool hasMasterClonePtr() const override;
0902     /// returns ptr to master clone, if existing.
0903     /// Throws an exception unless the concrete Candidate type is
0904     /// ShallowClonePtrCandidate
0905 
0906     const reco::CandidatePtr &masterClonePtr() const override;
0907 
0908     /// cast master clone reference to a concrete type
0909     template <typename Ref>
0910     Ref masterRef() const {
0911       return masterClone().template castTo<Ref>();
0912     }
0913 
0914     bool isElectron() const override { return false; }
0915     bool isMuon() const override { return false; }
0916     bool isStandAloneMuon() const override { return ((qualityFlags_ & muonFlagsMask) >> muonFlagsShift) & 1; }
0917     bool isGlobalMuon() const override { return ((qualityFlags_ & muonFlagsMask) >> muonFlagsShift) & 2; }
0918     bool isTrackerMuon() const override { return false; }
0919     bool isCaloMuon() const override { return false; }
0920     bool isPhoton() const override { return false; }
0921     bool isConvertedPhoton() const override { return false; }
0922     bool isJet() const override { return false; }
0923     bool isGoodEgamma() const { return (qualityFlags_ & egammaFlagsMask) != 0; }
0924 
0925     // puppiweights
0926     void setPuppiWeight(float p,
0927                         float p_nolep = 0.0);  /// Set both weights at once (with
0928                                                /// option for only full PUPPI)
0929     float puppiWeight() const;                 /// Weight from full PUPPI
0930     float puppiWeightNoLep() const;            /// Weight from PUPPI removing leptons
0931 
0932     // for the neutral fractions
0933     void setRawCaloFraction(float p);  /// Set the raw ECAL+HCAL energy over candidate
0934                                        /// energy for isolated charged hadrons
0935     float rawCaloFraction() const {
0936       return (rawCaloFraction_ / 100.);
0937     }                                  /// Raw ECAL+HCAL energy over candidate energy for isolated charged hadrons
0938     void setRawHcalFraction(float p);  /// Set the fraction of Hcal needed isolated charged hadrons
0939     float rawHcalFraction() const {
0940       return (rawHcalFraction_ / 100.);
0941     }                               /// Fraction of Hcal for isolated charged hadrons
0942     void setCaloFraction(float p);  /// Set the fraction of ECAL+HCAL energy over candidate energy
0943     float caloFraction() const {
0944       return (caloFraction_ / 100.);
0945     }                               /// Fraction of ECAL+HCAL energy over candidate energy
0946     void setHcalFraction(float p);  /// Set the fraction of Hcal needed for HF,
0947                                     /// neutral hadrons, and charged particles
0948     float hcalFraction() const {
0949       return (hcalFraction_ / 100.);
0950     }  /// Fraction of Hcal for HF, neutral hadrons, and charged particles
0951 
0952     // isolated charged hadrons
0953     void setIsIsolatedChargedHadron(bool p);  /// Set isolation (as in particle flow, i.e. at calorimeter
0954                                               /// surface rather than at PV) flat for charged hadrons
0955     bool isIsolatedChargedHadron() const {
0956       return isIsolatedChargedHadron_;
0957     }  /// Flag isolation (as in particle flow, i.e. at calorimeter surface rather
0958     /// than at PV) flag for charged hadrons
0959 
0960     struct PackedCovariance {
0961       PackedCovariance()
0962           : dxydxy(0), dxydz(0), dzdz(0), dlambdadz(0), dphidxy(0), dptdpt(0), detadeta(0), dphidphi(0) {}
0963       // 3D IP covariance
0964       uint16_t dxydxy;
0965       uint16_t dxydz;
0966       uint16_t dzdz;
0967       // other IP relevant elements
0968       uint16_t dlambdadz;
0969       uint16_t dphidxy;
0970       // other diag elements
0971       uint16_t dptdpt;
0972       uint16_t detadeta;
0973       uint16_t dphidphi;
0974     };
0975 
0976     /// time (wrt nominal zero of the collision)
0977     virtual float time() const { return vertexRef()->t() + dtimeAssociatedPV(); }
0978     /// dtime with respect to the PV[ipv]
0979     virtual float dtime(size_t ipv = 0) const {
0980       return dtimeAssociatedPV() + (*pvRefProd_)[pvRefKey_].t() - (*pvRefProd_)[ipv].t();
0981     }
0982     /// dtime with respect to the PV ref
0983     virtual float dtimeAssociatedPV() const {
0984       if (packedTime_ == 0)
0985         return 0.f;
0986       if (packedTimeError_ > 0)
0987         return unpackTimeWithError(packedTime_, packedTimeError_);
0988       else
0989         return unpackTimeNoError(packedTime_);
0990     }
0991     /// time measurement uncertainty (-1 if not available)
0992     virtual float timeError() const { return unpackTimeError(packedTimeError_); }
0993     /// set time measurement
0994     void setDTimeAssociatedPV(float aTime, float aTimeError = 0);
0995     /// set time measurement
0996     void setTime(float aTime, float aTimeError = 0) { setDTimeAssociatedPV(aTime - vertexRef()->t(), aTimeError); }
0997 
0998   private:
0999     void unpackCovarianceElement(reco::TrackBase::CovarianceMatrix &m, uint16_t packed, int i, int j) const {
1000       m(i, j) = covarianceParameterization().unpack(
1001           packed, covarianceSchema_, i, j, pt(), eta(), numberOfHits(), numberOfPixelHits());
1002     }
1003     uint16_t packCovarianceElement(const reco::TrackBase::CovarianceMatrix &m, int i, int j) const {
1004       return covarianceParameterization().pack(
1005           m(i, j), covarianceSchema_, i, j, pt(), eta(), numberOfHits(), numberOfPixelHits());
1006     }
1007 
1008   protected:
1009     friend class ::testPackedCandidate;
1010     static constexpr float kMinDEtaToStore_ = 0.001;
1011     static constexpr float kMinDTrkPtToStore_ = 0.001;
1012 
1013     uint16_t packedPt_, packedEta_, packedPhi_, packedM_;
1014     uint16_t packedDxy_, packedDz_, packedDPhi_, packedDEta_, packedDTrkPt_;
1015     PackedCovariance packedCovariance_;
1016 
1017     void pack(bool unpackAfterwards = true);
1018     void unpack() const;
1019     void packVtx(bool unpackAfterwards = true);
1020     void unpackVtx() const;
1021     void packCovariance(const reco::TrackBase::CovarianceMatrix &m, bool unpackAfterwards = true);
1022     void unpackCovariance() const;
1023     void maybeUnpackBoth() const {
1024       if (!p4c_)
1025         unpack();
1026       if (!vertex_)
1027         unpackVtx();
1028     }
1029     void maybeUnpackTrack() const {
1030       if (!track_)
1031         unpackTrk();
1032     }
1033     void maybeUnpackCovariance() const {
1034       if (!m_)
1035         unpackCovariance();
1036     }
1037     void packBoth() {
1038       pack(false);
1039       packVtx(false);
1040       delete p4_.exchange(nullptr);
1041       delete p4c_.exchange(nullptr);
1042       delete vertex_.exchange(nullptr);
1043       unpack();
1044       unpackVtx();
1045     }  // do it this way, so that we don't loose precision on the angles before
1046     // computing dxy,dz
1047     void unpackTrk() const;
1048 
1049     uint8_t packedPuppiweight_;
1050     int8_t packedPuppiweightNoLepDiff_;  // storing the DIFFERENCE of (all - "no
1051                                          // lep") for compression optimization
1052     uint8_t rawCaloFraction_;
1053     int8_t rawHcalFraction_;
1054     uint8_t caloFraction_;
1055     int8_t hcalFraction_;
1056     int16_t packedTime_;
1057     uint8_t packedTimeError_;
1058 
1059     bool isIsolatedChargedHadron_;
1060 
1061     /// the four vector
1062     mutable std::atomic<PolarLorentzVector *> p4_;
1063     mutable std::atomic<LorentzVector *> p4c_;
1064     /// vertex position
1065     mutable std::atomic<Point *> vertex_;
1066     CMS_THREAD_GUARD(vertex_) mutable float dxy_, dz_, dphi_, deta_, dtrkpt_;
1067     /// reco::Track
1068     mutable std::atomic<reco::Track *> track_;
1069     /// PDG identifier
1070     int pdgId_;
1071     uint16_t qualityFlags_;
1072     /// Use these to build a Ref to primary vertex
1073     reco::VertexRefProd pvRefProd_;
1074     reco::VertexRef::key_type pvRefKey_;
1075 
1076     /// IP covariance
1077     mutable std::atomic<reco::TrackBase::CovarianceMatrix *> m_;
1078     uint8_t packedHits_,
1079         packedLayers_;  // packedLayers_ -> layers with valid hits; packedHits_ ->
1080                         // extra hits beyond the one-per-layer implied by
1081                         // packedLayers_
1082 
1083     /// track quality information
1084     uint8_t normalizedChi2_;
1085     uint16_t covarianceVersion_;
1086     uint16_t covarianceSchema_;
1087     CMS_THREAD_SAFE static CovarianceParameterization covarianceParameterization_;
1088     // static std::atomic<CovarianceParameterization*>
1089     // covarianceParameterization_;
1090     static std::once_flag covariance_load_flag;
1091     const CovarianceParameterization &covarianceParameterization() const {
1092       if (!hasTrackDetails())
1093         throw edm::Exception(edm::errors::InvalidReference,
1094                              "Trying to access covariance matrix for a "
1095                              "PackedCandidate for which it's not available. "
1096                              "Check hasTrackDetails() before!\n");
1097       std::call_once(
1098           covariance_load_flag, [](int v) { covarianceParameterization_.load(v); }, covarianceVersion_);
1099       if (covarianceParameterization_.loadedVersion() != covarianceVersion_) {
1100         throw edm::Exception(edm::errors::UnimplementedFeature)
1101             << "Attempting to load multiple covariance version in same process. "
1102                "This is not supported.";
1103       }
1104       return covarianceParameterization_;
1105     }
1106 
1107     /// check overlap with another Candidate
1108     bool overlap(const reco::Candidate &) const override;
1109     template <typename, typename, typename>
1110     friend struct component;
1111     friend class ::OverlapChecker;
1112     friend class ShallowCloneCandidate;
1113     friend class ShallowClonePtrCandidate;
1114 
1115     enum qualityFlagsShiftsAndMasks {
1116       assignmentQualityMask = 0x7,
1117       assignmentQualityShift = 0,
1118       trackHighPurityMask = 0x8,
1119       trackHighPurityShift = 3,
1120       lostInnerHitsMask = 0x30,
1121       lostInnerHitsShift = 4,
1122       muonFlagsMask = 0x0600,
1123       muonFlagsShift = 9,
1124       egammaFlagsMask = 0x0800,
1125       egammaFlagsShift = 11
1126     };
1127 
1128     /// static to allow unit testing
1129     static uint8_t packTimeError(float timeError);
1130     static float unpackTimeError(uint8_t timeError);
1131     static float unpackTimeNoError(int16_t time);
1132     static int16_t packTimeNoError(float time);
1133     static float unpackTimeWithError(int16_t time, uint8_t timeError);
1134     static int16_t packTimeWithError(float time, float timeError);
1135     static constexpr float MIN_TIMEERROR = 0.002f;      // 2 ps, smallest storable non-zero uncertainty
1136     static constexpr float MIN_TIME_NOERROR = 0.0002f;  // 0.2 ps, smallest non-zero time that can be stored by
1137                                                         // packTimeNoError
1138     static constexpr int EXPO_TIMEERROR = 5;            // power of 2 used in encoding timeError
1139     static constexpr int EXPO_TIME_NOERROR = 6;         // power of 2 used in encoding time without timeError
1140     static constexpr int EXPO_TIME_WITHERROR = -6;      // power of 2 used in encoding time with timeError
1141   public:
1142     uint16_t firstHit_;
1143   };
1144 
1145   typedef std::vector<pat::PackedCandidate> PackedCandidateCollection;
1146   typedef edm::Ref<pat::PackedCandidateCollection> PackedCandidateRef;
1147   typedef edm::RefVector<pat::PackedCandidateCollection> PackedCandidateRefVector;
1148 }  // namespace pat
1149 
1150 #endif