CMS 3D CMS Logo

Photon.cc
Go to the documentation of this file.
3 
4 using namespace reco;
5 
6 Photon::Photon(const LorentzVector& p4, const Point& caloPos, const PhotonCoreRef& core, const Point& vtx)
7  : RecoCandidate(0, p4, vtx, 22),
8  caloPosition_(caloPos),
9  photonCore_(core),
10  pixelSeed_(false),
11  haloTaggerMVAVal_(99) {}
12 
13 Photon::Photon(const Photon& rhs)
14  : RecoCandidate(rhs),
15  caloPosition_(rhs.caloPosition_),
16  photonCore_(rhs.photonCore_),
17  pixelSeed_(rhs.pixelSeed_),
18  fiducialFlagBlock_(rhs.fiducialFlagBlock_),
19  isolationR04_(rhs.isolationR04_),
20  isolationR03_(rhs.isolationR03_),
21  showerShapeBlock_(rhs.showerShapeBlock_),
22  full5x5_showerShapeBlock_(rhs.full5x5_showerShapeBlock_),
23  saturationInfo_(rhs.saturationInfo_),
24  eCorrections_(rhs.eCorrections_),
25  mipVariableBlock_(rhs.mipVariableBlock_),
26  pfIsolation_(rhs.pfIsolation_),
27  pfID_(rhs.pfID_),
28  haloTaggerMVAVal_(rhs.haloTaggerMVAVal_) {}
29 
30 Photon::~Photon() {}
31 
32 Photon* Photon::clone() const { return new Photon(*this); }
33 
34 bool Photon::overlap(const Candidate& c) const {
35  const RecoCandidate* o = dynamic_cast<const RecoCandidate*>(&c);
36  return (o != nullptr && (checkOverlap(superCluster(), o->superCluster())));
37  return false;
38 }
39 
40 void Photon::setVertex(const Point& vertex) {
41  math::XYZVectorF direction = caloPosition() - vertex;
42  double energy = this->energy();
43  math::XYZVectorF momentum = direction.unit() * energy;
45  setP4(lv);
47 }
48 
49 reco::SuperClusterRef Photon::superCluster() const { return this->photonCore()->superCluster(); }
50 
52  const reco::ConversionRefVector& conv2leg = this->photonCore()->conversions();
53  const reco::ConversionRefVector& conv1leg = this->photonCore()->conversionsOneLeg();
54 
55  int origin = -1;
56  bool isEg = false, isPf = false;
57 
58  for (unsigned iConv = 0; iConv < conv2leg.size(); iConv++) {
59  std::vector<edm::RefToBase<reco::Track> > convtracks = conv2leg[iConv]->tracks();
60  for (unsigned itk = 0; itk < convtracks.size(); itk++) {
61  if (convTrack == convtracks[itk])
62  isEg = true;
63  }
64  }
65 
66  for (unsigned iConv = 0; iConv < conv1leg.size(); iConv++) {
67  std::vector<edm::RefToBase<reco::Track> > convtracks = conv1leg[iConv]->tracks();
68  for (unsigned itk = 0; itk < convtracks.size(); itk++) {
69  if (convTrack == convtracks[itk])
70  isPf = true;
71  }
72  }
73 
74  if (isEg)
75  origin = egamma;
76  if (isPf)
77  origin = pflow;
78  if (isEg && isPf)
79  origin = both;
80 
81  return origin;
82 }
83 
84 void Photon::setCorrectedEnergy(P4type type, float newEnergy, float delta_e, bool setToRecoCandidate) {
85  math::XYZTLorentzVectorD newP4 = p4();
86  newP4 *= newEnergy / newP4.e();
87  switch (type) {
88  case ecal_standard:
89  eCorrections_.scEcalEnergy = newEnergy;
91  break;
92  case ecal_photons:
93  eCorrections_.phoEcalEnergy = newEnergy;
95  break;
96  case regression1:
97  eCorrections_.regression1Energy = newEnergy;
99  [[fallthrough]];
100  case regression2:
101  eCorrections_.regression2Energy = newEnergy;
103  break;
104  default:
105  throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type;
106  }
107  setP4(type, newP4, delta_e, setToRecoCandidate);
108 }
109 
110 float Photon::getCorrectedEnergy(P4type type) const {
111  switch (type) {
112  case ecal_standard:
114  break;
115  case ecal_photons:
117  break;
118  case regression1:
120  case regression2:
122  break;
123  default:
124  throw cms::Exception("reco::Photon") << "unexpected p4 type " << type << " cannot return the energy value: ";
125  }
126 }
127 
128 float Photon::getCorrectedEnergyError(P4type type) const {
129  switch (type) {
130  case ecal_standard:
132  break;
133  case ecal_photons:
135  break;
136  case regression1:
138  case regression2:
140  break;
141  default:
142  throw cms::Exception("reco::Photon")
143  << "unexpected p4 type " << type << " cannot return the uncertainty on the energy: ";
144  }
145 }
146 
147 void Photon::setP4(P4type type, const LorentzVector& p4, float error, bool setToRecoCandidate) {
148  switch (type) {
149  case ecal_standard:
152  break;
153  case ecal_photons:
156  break;
157  case regression1:
160  [[fallthrough]];
161  case regression2:
164  break;
165  default:
166  throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type;
167  }
168  if (setToRecoCandidate) {
169  setP4(p4);
171  }
172 }
173 
174 const Candidate::LorentzVector& Photon::p4(P4type type) const {
175  switch (type) {
176  case ecal_standard:
177  return eCorrections_.scEcalP4;
178  case ecal_photons:
179  return eCorrections_.phoEcalP4;
180  case regression1:
182  case regression2:
184  default:
185  throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type << " cannot return p4 ";
186  }
187 }
188 
190  auto& ss1 = showerShapeBlock_;
191  auto& ss2 = full5x5_showerShapeBlock_;
192  auto& iv1 = isolationR03_;
193  auto& iv2 = isolationR04_;
194 
195  for (uint id = 2u; id < ss1.hcalOverEcal.size(); ++id) {
196  ss1.hcalOverEcal[1] += ss1.hcalOverEcal[id];
197  ss1.hcalOverEcalBc[1] += ss1.hcalOverEcalBc[id];
198 
199  ss1.hcalOverEcal[id] = 0.f;
200  ss1.hcalOverEcalBc[id] = 0.f;
201 
202  ss2.hcalOverEcal[1] += ss2.hcalOverEcal[id];
203  ss2.hcalOverEcalBc[1] += ss2.hcalOverEcalBc[id];
204 
205  ss2.hcalOverEcal[id] = 0.f;
206  ss2.hcalOverEcalBc[id] = 0.f;
207 
208  iv1.hcalRecHitSumEt[1] += iv1.hcalRecHitSumEt[id];
209  iv1.hcalRecHitSumEtBc[1] += iv1.hcalRecHitSumEtBc[id];
210 
211  iv1.hcalRecHitSumEt[id] = 0.f;
212  iv1.hcalRecHitSumEtBc[id] = 0.f;
213 
214  iv2.hcalRecHitSumEt[1] += iv2.hcalRecHitSumEt[id];
215  iv2.hcalRecHitSumEtBc[1] += iv2.hcalRecHitSumEtBc[id];
216 
217  iv2.hcalRecHitSumEt[id] = 0.f;
218  iv2.hcalRecHitSumEtBc[id] = 0.f;
219  }
220 }
Photon * clone() const override
returns a clone of the candidate
float getCorrectedEnergy(P4type type) const
Vector momentum() const final
spatial momentum vector
ShowerShape full5x5_showerShapeBlock_
Definition: Photon.h:608
Definition: Photon.py:1
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
Definition: __init__.py:1
int conversionTrackProvenance(const edm::RefToBase< reco::Track > &convTrack) const
bool overlap(const Candidate &) const override
check overlap with another candidate
const Point & vertex() const override
vertex position (overwritten by PF...)
IsolationVariables isolationR03_
Definition: Photon.h:606
const LorentzVector & p4() const final
four-momentum Lorentz vector
Photon()
default constructor
Definition: Photon.h:32
void setVertex(const Point &vertex) override
set vertex
void setVertex(const Point &vertex) override
set primary event vertex used to define photon direction
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
LorentzVector regression2P4
Definition: Photon.h:334
IsolationVariables isolationR04_
Definition: Photon.h:605
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
EnergyCorrections eCorrections_
Definition: Photon.h:610
float getCorrectedEnergyError(P4type type) const
void hcalToRun2EffDepth()
math::XYZTLorentzVector LorentzVector
math::XYZPointF caloPosition() const
position in ECAL: this is th SC position if r9<0.93. If r8>0.93 is position of seed BasicCluster taki...
Definition: Photon.h:91
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
~Photon() override
destructor
bool checkOverlap(const R &r1, const R &r2) const
check if two components overlap
Definition: RecoCandidate.h:67
reco::PhotonCoreRef photonCore() const
returns a reference to the core photon object
Definition: Photon.h:53
LorentzVector regression1P4
Definition: Photon.h:331
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
std::array< float, 7 > hcalRecHitSumEt
Definition: Photon.h:422
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
fixed size matrix
Structure Point Contains parameters of Gaussian fits to DMRs.
ShowerShape showerShapeBlock_
Definition: Photon.h:607
void setP4(P4type type, const LorentzVector &p4, float p4Error, bool setToRecoCandidate)
double energy() const final
energy