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  haloTaggerMVAVal_(rhs.haloTaggerMVAVal_) {}
28 
29 Photon::~Photon() {}
30 
31 Photon* Photon::clone() const { return new Photon(*this); }
32 
33 bool Photon::overlap(const Candidate& c) const {
34  const RecoCandidate* o = dynamic_cast<const RecoCandidate*>(&c);
35  return (o != nullptr && (checkOverlap(superCluster(), o->superCluster())));
36  return false;
37 }
38 
39 void Photon::setVertex(const Point& vertex) {
40  math::XYZVectorF direction = caloPosition() - vertex;
41  double energy = this->energy();
42  math::XYZVectorF momentum = direction.unit() * energy;
44  setP4(lv);
46 }
47 
48 reco::SuperClusterRef Photon::superCluster() const { return this->photonCore()->superCluster(); }
49 
51  const reco::ConversionRefVector& conv2leg = this->photonCore()->conversions();
52  const reco::ConversionRefVector& conv1leg = this->photonCore()->conversionsOneLeg();
53 
54  int origin = -1;
55  bool isEg = false, isPf = false;
56 
57  for (unsigned iConv = 0; iConv < conv2leg.size(); iConv++) {
58  std::vector<edm::RefToBase<reco::Track> > convtracks = conv2leg[iConv]->tracks();
59  for (unsigned itk = 0; itk < convtracks.size(); itk++) {
60  if (convTrack == convtracks[itk])
61  isEg = true;
62  }
63  }
64 
65  for (unsigned iConv = 0; iConv < conv1leg.size(); iConv++) {
66  std::vector<edm::RefToBase<reco::Track> > convtracks = conv1leg[iConv]->tracks();
67  for (unsigned itk = 0; itk < convtracks.size(); itk++) {
68  if (convTrack == convtracks[itk])
69  isPf = true;
70  }
71  }
72 
73  if (isEg)
74  origin = egamma;
75  if (isPf)
76  origin = pflow;
77  if (isEg && isPf)
78  origin = both;
79 
80  return origin;
81 }
82 
83 void Photon::setCorrectedEnergy(P4type type, float newEnergy, float delta_e, bool setToRecoCandidate) {
84  math::XYZTLorentzVectorD newP4 = p4();
85  newP4 *= newEnergy / newP4.e();
86  switch (type) {
87  case ecal_standard:
88  eCorrections_.scEcalEnergy = newEnergy;
90  break;
91  case ecal_photons:
92  eCorrections_.phoEcalEnergy = newEnergy;
94  break;
95  case regression1:
96  eCorrections_.regression1Energy = newEnergy;
98  [[fallthrough]];
99  case regression2:
100  eCorrections_.regression2Energy = newEnergy;
102  break;
103  default:
104  throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type;
105  }
106  setP4(type, newP4, delta_e, setToRecoCandidate);
107 }
108 
109 float Photon::getCorrectedEnergy(P4type type) const {
110  switch (type) {
111  case ecal_standard:
113  break;
114  case ecal_photons:
116  break;
117  case regression1:
119  case regression2:
121  break;
122  default:
123  throw cms::Exception("reco::Photon") << "unexpected p4 type " << type << " cannot return the energy value: ";
124  }
125 }
126 
127 float Photon::getCorrectedEnergyError(P4type type) const {
128  switch (type) {
129  case ecal_standard:
131  break;
132  case ecal_photons:
134  break;
135  case regression1:
137  case regression2:
139  break;
140  default:
141  throw cms::Exception("reco::Photon")
142  << "unexpected p4 type " << type << " cannot return the uncertainty on the energy: ";
143  }
144 }
145 
146 void Photon::setP4(P4type type, const LorentzVector& p4, float error, bool setToRecoCandidate) {
147  switch (type) {
148  case ecal_standard:
151  break;
152  case ecal_photons:
155  break;
156  case regression1:
159  [[fallthrough]];
160  case regression2:
163  break;
164  default:
165  throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type;
166  }
167  if (setToRecoCandidate) {
168  setP4(p4);
170  }
171 }
172 
173 const Candidate::LorentzVector& Photon::p4(P4type type) const {
174  switch (type) {
175  case ecal_standard:
176  return eCorrections_.scEcalP4;
177  case ecal_photons:
178  return eCorrections_.phoEcalP4;
179  case regression1:
181  case regression2:
183  default:
184  throw cms::Exception("reco::Photon") << "unexpected p4 type: " << type << " cannot return p4 ";
185  }
186 }
187 
189  auto& ss1 = showerShapeBlock_;
190  auto& ss2 = full5x5_showerShapeBlock_;
191  auto& iv1 = isolationR03_;
192  auto& iv2 = isolationR04_;
193 
194  for (uint id = 2u; id < ss1.hcalOverEcal.size(); ++id) {
195  ss1.hcalOverEcal[1] += ss1.hcalOverEcal[id];
196  ss1.hcalOverEcalBc[1] += ss1.hcalOverEcalBc[id];
197 
198  ss1.hcalOverEcal[id] = 0.f;
199  ss1.hcalOverEcalBc[id] = 0.f;
200 
201  ss2.hcalOverEcal[1] += ss2.hcalOverEcal[id];
202  ss2.hcalOverEcalBc[1] += ss2.hcalOverEcalBc[id];
203 
204  ss2.hcalOverEcal[id] = 0.f;
205  ss2.hcalOverEcalBc[id] = 0.f;
206 
207  iv1.hcalRecHitSumEt[1] += iv1.hcalRecHitSumEt[id];
208  iv1.hcalRecHitSumEtBc[1] += iv1.hcalRecHitSumEtBc[id];
209 
210  iv1.hcalRecHitSumEt[id] = 0.f;
211  iv1.hcalRecHitSumEtBc[id] = 0.f;
212 
213  iv2.hcalRecHitSumEt[1] += iv2.hcalRecHitSumEt[id];
214  iv2.hcalRecHitSumEtBc[1] += iv2.hcalRecHitSumEtBc[id];
215 
216  iv2.hcalRecHitSumEt[id] = 0.f;
217  iv2.hcalRecHitSumEtBc[id] = 0.f;
218  }
219 }
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:605
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:603
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:331
IsolationVariables isolationR04_
Definition: Photon.h:602
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
EnergyCorrections eCorrections_
Definition: Photon.h:607
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:88
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:50
LorentzVector regression1P4
Definition: Photon.h:328
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
std::array< float, 7 > hcalRecHitSumEt
Definition: Photon.h:419
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:604
void setP4(P4type type, const LorentzVector &p4, float p4Error, bool setToRecoCandidate)
double energy() const final
energy