Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:06

0001 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0002 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
0003 
0004 using namespace reco;
0005 
0006 Photon::Photon(const LorentzVector& p4, const Point& caloPos, const PhotonCoreRef& core, const Point& vtx)
0007     : RecoCandidate(0, p4, vtx, 22),
0008       caloPosition_(caloPos),
0009       photonCore_(core),
0010       pixelSeed_(false),
0011       haloTaggerMVAVal_(99) {}
0012 
0013 Photon::Photon(const Photon& rhs)
0014     : RecoCandidate(rhs),
0015       caloPosition_(rhs.caloPosition_),
0016       photonCore_(rhs.photonCore_),
0017       pixelSeed_(rhs.pixelSeed_),
0018       fiducialFlagBlock_(rhs.fiducialFlagBlock_),
0019       isolationR04_(rhs.isolationR04_),
0020       isolationR03_(rhs.isolationR03_),
0021       showerShapeBlock_(rhs.showerShapeBlock_),
0022       full5x5_showerShapeBlock_(rhs.full5x5_showerShapeBlock_),
0023       saturationInfo_(rhs.saturationInfo_),
0024       eCorrections_(rhs.eCorrections_),
0025       mipVariableBlock_(rhs.mipVariableBlock_),
0026       pfIsolation_(rhs.pfIsolation_),
0027       pfID_(rhs.pfID_),
0028       haloTaggerMVAVal_(rhs.haloTaggerMVAVal_) {}
0029 
0030 Photon::~Photon() {}
0031 
0032 Photon* Photon::clone() const { return new Photon(*this); }
0033 
0034 bool Photon::overlap(const Candidate& c) const {
0035   const RecoCandidate* o = dynamic_cast<const RecoCandidate*>(&c);
0036   return (o != nullptr && (checkOverlap(superCluster(), o->superCluster())));
0037   return false;
0038 }
0039 
0040 void Photon::setVertex(const Point& vertex) {
0041   math::XYZVectorF direction = caloPosition() - vertex;
0042   double energy = this->energy();
0043   math::XYZVectorF momentum = direction.unit() * energy;
0044   math::XYZTLorentzVector lv(momentum.x(), momentum.y(), momentum.z(), energy);
0045   setP4(lv);
0046   LeafCandidate::setVertex(vertex);
0047 }
0048 
0049 reco::SuperClusterRef Photon::superCluster() const { return this->photonCore()->superCluster(); }
0050 
0051 int Photon::conversionTrackProvenance(const edm::RefToBase<reco::Track>& convTrack) const {
0052   const reco::ConversionRefVector& conv2leg = this->photonCore()->conversions();
0053   const reco::ConversionRefVector& conv1leg = this->photonCore()->conversionsOneLeg();
0054 
0055   int origin = -1;
0056   bool isEg = false, isPf = false;
0057 
0058   for (unsigned iConv = 0; iConv < conv2leg.size(); iConv++) {
0059     std::vector<edm::RefToBase<reco::Track> > convtracks = conv2leg[iConv]->tracks();
0060     for (unsigned itk = 0; itk < convtracks.size(); itk++) {
0061       if (convTrack == convtracks[itk])
0062         isEg = true;
0063     }
0064   }
0065 
0066   for (unsigned iConv = 0; iConv < conv1leg.size(); iConv++) {
0067     std::vector<edm::RefToBase<reco::Track> > convtracks = conv1leg[iConv]->tracks();
0068     for (unsigned itk = 0; itk < convtracks.size(); itk++) {
0069       if (convTrack == convtracks[itk])
0070         isPf = true;
0071     }
0072   }
0073 
0074   if (isEg)
0075     origin = egamma;
0076   if (isPf)
0077     origin = pflow;
0078   if (isEg && isPf)
0079     origin = both;
0080 
0081   return origin;
0082 }
0083 
0084 void Photon::setCorrectedEnergy(P4type type, float newEnergy, float delta_e, bool setToRecoCandidate) {
0085   math::XYZTLorentzVectorD newP4 = p4();
0086   newP4 *= newEnergy / newP4.e();
0087   switch (type) {
0088     case ecal_standard:
0089       eCorrections_.scEcalEnergy = newEnergy;
0090       eCorrections_.scEcalEnergyError = delta_e;
0091       break;
0092     case ecal_photons:
0093       eCorrections_.phoEcalEnergy = newEnergy;
0094       eCorrections_.phoEcalEnergyError = delta_e;
0095       break;
0096     case regression1:
0097       eCorrections_.regression1Energy = newEnergy;
0098       eCorrections_.regression1EnergyError = delta_e;
0099       [[fallthrough]];
0100     case regression2:
0101       eCorrections_.regression2Energy = newEnergy;
0102       eCorrections_.regression2EnergyError = delta_e;
0103       break;
0104     default:
0105       throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type;
0106   }
0107   setP4(type, newP4, delta_e, setToRecoCandidate);
0108 }
0109 
0110 float Photon::getCorrectedEnergy(P4type type) const {
0111   switch (type) {
0112     case ecal_standard:
0113       return eCorrections_.scEcalEnergy;
0114       break;
0115     case ecal_photons:
0116       return eCorrections_.phoEcalEnergy;
0117       break;
0118     case regression1:
0119       return eCorrections_.regression1Energy;
0120     case regression2:
0121       return eCorrections_.regression2Energy;
0122       break;
0123     default:
0124       throw cms::Exception("reco::Photon") << "unexpected p4 type " << type << " cannot return the energy value: ";
0125   }
0126 }
0127 
0128 float Photon::getCorrectedEnergyError(P4type type) const {
0129   switch (type) {
0130     case ecal_standard:
0131       return eCorrections_.scEcalEnergyError;
0132       break;
0133     case ecal_photons:
0134       return eCorrections_.phoEcalEnergyError;
0135       break;
0136     case regression1:
0137       return eCorrections_.regression1EnergyError;
0138     case regression2:
0139       return eCorrections_.regression2EnergyError;
0140       break;
0141     default:
0142       throw cms::Exception("reco::Photon")
0143           << "unexpected p4 type " << type << " cannot return the uncertainty on the energy: ";
0144   }
0145 }
0146 
0147 void Photon::setP4(P4type type, const LorentzVector& p4, float error, bool setToRecoCandidate) {
0148   switch (type) {
0149     case ecal_standard:
0150       eCorrections_.scEcalP4 = p4;
0151       eCorrections_.scEcalEnergyError = error;
0152       break;
0153     case ecal_photons:
0154       eCorrections_.phoEcalP4 = p4;
0155       eCorrections_.phoEcalEnergyError = error;
0156       break;
0157     case regression1:
0158       eCorrections_.regression1P4 = p4;
0159       eCorrections_.regression1EnergyError = error;
0160       [[fallthrough]];
0161     case regression2:
0162       eCorrections_.regression2P4 = p4;
0163       eCorrections_.regression2EnergyError = error;
0164       break;
0165     default:
0166       throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type;
0167   }
0168   if (setToRecoCandidate) {
0169     setP4(p4);
0170     eCorrections_.candidateP4type = type;
0171   }
0172 }
0173 
0174 const Candidate::LorentzVector& Photon::p4(P4type type) const {
0175   switch (type) {
0176     case ecal_standard:
0177       return eCorrections_.scEcalP4;
0178     case ecal_photons:
0179       return eCorrections_.phoEcalP4;
0180     case regression1:
0181       return eCorrections_.regression1P4;
0182     case regression2:
0183       return eCorrections_.regression2P4;
0184     default:
0185       throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type << " cannot return p4 ";
0186   }
0187 }
0188 
0189 void Photon::hcalToRun2EffDepth() {
0190   auto& ss1 = showerShapeBlock_;
0191   auto& ss2 = full5x5_showerShapeBlock_;
0192   auto& iv1 = isolationR03_;
0193   auto& iv2 = isolationR04_;
0194 
0195   for (uint id = 2u; id < ss1.hcalOverEcal.size(); ++id) {
0196     ss1.hcalOverEcal[1] += ss1.hcalOverEcal[id];
0197     ss1.hcalOverEcalBc[1] += ss1.hcalOverEcalBc[id];
0198 
0199     ss1.hcalOverEcal[id] = 0.f;
0200     ss1.hcalOverEcalBc[id] = 0.f;
0201 
0202     ss2.hcalOverEcal[1] += ss2.hcalOverEcal[id];
0203     ss2.hcalOverEcalBc[1] += ss2.hcalOverEcalBc[id];
0204 
0205     ss2.hcalOverEcal[id] = 0.f;
0206     ss2.hcalOverEcalBc[id] = 0.f;
0207 
0208     iv1.hcalRecHitSumEt[1] += iv1.hcalRecHitSumEt[id];
0209     iv1.hcalRecHitSumEtBc[1] += iv1.hcalRecHitSumEtBc[id];
0210 
0211     iv1.hcalRecHitSumEt[id] = 0.f;
0212     iv1.hcalRecHitSumEtBc[id] = 0.f;
0213 
0214     iv2.hcalRecHitSumEt[1] += iv2.hcalRecHitSumEt[id];
0215     iv2.hcalRecHitSumEtBc[1] += iv2.hcalRecHitSumEtBc[id];
0216 
0217     iv2.hcalRecHitSumEt[id] = 0.f;
0218     iv2.hcalRecHitSumEtBc[id] = 0.f;
0219   }
0220 }