CMS 3D CMS Logo

Photon.cc
Go to the documentation of this file.
3 
4 using namespace reco;
5 
7  const Point& caloPos,
8  const PhotonCoreRef & core,
9  const Point & vtx) :
10  RecoCandidate( 0, p4, vtx, 22 ),
11  caloPosition_( caloPos ),
12  photonCore_(core),
13  pixelSeed_(false)
14  {}
15 
16 
17 Photon::Photon( const Photon& rhs ) :
18  RecoCandidate(rhs),
19  caloPosition_(rhs.caloPosition_),
20  photonCore_ ( rhs.photonCore_),
21  pixelSeed_ ( rhs.pixelSeed_ ),
22  fiducialFlagBlock_ ( rhs.fiducialFlagBlock_ ),
23  isolationR04_ ( rhs.isolationR04_),
24  isolationR03_ ( rhs.isolationR03_),
25  showerShapeBlock_ ( rhs.showerShapeBlock_),
26  full5x5_showerShapeBlock_ ( rhs.full5x5_showerShapeBlock_),
27  saturationInfo_ ( rhs.saturationInfo_ ),
28  eCorrections_(rhs.eCorrections_),
29  mipVariableBlock_ (rhs.mipVariableBlock_),
30  pfIsolation_ ( rhs.pfIsolation_ )
31  {}
32 
33 
34 
35 
36 Photon::~Photon() { }
37 
38 Photon * Photon::clone() const {
39  return new Photon( * this );
40 }
41 
42 
43 bool Photon::overlap( const Candidate & c ) const {
44  const RecoCandidate * o = dynamic_cast<const RecoCandidate *>( & c );
45  return ( o != 0 &&
47  );
48  return false;
49 }
50 
51 void Photon::setVertex(const Point & vertex) {
52  math::XYZVectorF direction = caloPosition() - vertex;
53  double energy = this->energy();
54  math::XYZVectorF momentum = direction.unit() * energy;
55  math::XYZTLorentzVector lv(momentum.x(), momentum.y(), momentum.z(), energy );
56  setP4(lv);
58 }
59 
61  return this->photonCore()->superCluster();
62 }
63 
65 
66  const reco::ConversionRefVector & conv2leg = this->photonCore()->conversions();
67  const reco::ConversionRefVector & conv1leg = this->photonCore()->conversionsOneLeg();
68 
69  int origin = -1;
70  bool isEg=false, isPf=false;
71 
72  for (unsigned iConv=0; iConv<conv2leg.size(); iConv++){
73  std::vector<edm::RefToBase<reco::Track> > convtracks = conv2leg[iConv]->tracks();
74  for (unsigned itk=0; itk<convtracks.size(); itk++){
75  if (convTrack==convtracks[itk]) isEg=true;
76  }
77  }
78 
79  for (unsigned iConv=0; iConv<conv1leg.size(); iConv++){
80  std::vector<edm::RefToBase<reco::Track> > convtracks = conv1leg[iConv]->tracks();
81  for (unsigned itk=0; itk<convtracks.size(); itk++){
82  if (convTrack==convtracks[itk]) isPf=true;
83  }
84  }
85 
86  if (isEg) origin=egamma;
87  if (isPf) origin=pflow;
88  if (isEg && isPf) origin=both;
89 
90  return origin;
91 }
92 
93 
94 void Photon::setCorrectedEnergy( P4type type, float newEnergy, float delta_e, bool setToRecoCandidate ) {
95 
96  math::XYZTLorentzVectorD newP4 = p4() ;
97  newP4 *= newEnergy/newP4.e() ;
98  switch(type)
99  {
100  case ecal_standard:
101  eCorrections_.scEcalEnergy = newEnergy;
102  eCorrections_.scEcalEnergyError = delta_e;
103  break ;
104  case ecal_photons:
105  eCorrections_.phoEcalEnergy = newEnergy;
106  eCorrections_.phoEcalEnergyError = delta_e;
107  break ;
108  case regression1:
109  eCorrections_.regression1Energy = newEnergy ;
110  eCorrections_.regression1EnergyError = delta_e;
111  case regression2:
112  eCorrections_.regression2Energy = newEnergy ;
113  eCorrections_.regression2EnergyError = delta_e;
114  break ;
115  default:
116  throw cms::Exception("reco::Photon")<<"unexpected p4 type: "<< type ;
117  }
118  setP4(type, newP4, delta_e, setToRecoCandidate);
119 
120 }
121 
122 
123 
124 
125  float Photon::getCorrectedEnergy( P4type type) const {
126  switch(type)
127  {
128  case ecal_standard:
129  return eCorrections_.scEcalEnergy;
130  break ;
131  case ecal_photons:
132  return eCorrections_.phoEcalEnergy;
133  break ;
134  case regression1:
135  return eCorrections_.regression1Energy;
136  case regression2:
137  return eCorrections_.regression2Energy;
138  break ;
139  default:
140  throw cms::Exception("reco::Photon")<<"unexpected p4 type " << type << " cannot return the energy value: " ;
141  }
142  }
143 
144 
145  float Photon::getCorrectedEnergyError( P4type type) const {
146  switch(type)
147  {
148  case ecal_standard:
149  return eCorrections_.scEcalEnergyError;
150  break ;
151  case ecal_photons:
152  return eCorrections_.phoEcalEnergyError;
153  break ;
154  case regression1:
155  return eCorrections_.regression1EnergyError;
156  case regression2:
157  return eCorrections_.regression2EnergyError;
158  break ;
159  default:
160  throw cms::Exception("reco::Photon")<<"unexpected p4 type " << type << " cannot return the uncertainty on the energy: " ;
161  }
162  }
163 
164 
165 
166 void Photon::setP4(P4type type, const LorentzVector & p4, float error, bool setToRecoCandidate ) {
167 
168 
169  switch(type)
170  {
171  case ecal_standard:
172  eCorrections_.scEcalP4 = p4 ;
173  eCorrections_.scEcalEnergyError = error ;
174  break ;
175  case ecal_photons:
176  eCorrections_.phoEcalP4 = p4 ;
177  eCorrections_.phoEcalEnergyError = error ;
178  break ;
179  case regression1:
180  eCorrections_.regression1P4 = p4 ;
181  eCorrections_.regression1EnergyError = error ;
182  case regression2:
183  eCorrections_.regression2P4 = p4 ;
184  eCorrections_.regression2EnergyError = error ;
185  break ;
186  default:
187  throw cms::Exception("reco::Photon")<<"unexpected p4 type: "<< type ;
188  }
189  if (setToRecoCandidate)
190  {
191  setP4(p4) ;
192  eCorrections_.candidateP4type = type ;
193  }
194 
195 
196 }
197 
198 const Candidate::LorentzVector& Photon::p4( P4type type ) const
199  {
200  switch(type)
201  {
202  case ecal_standard: return eCorrections_.scEcalP4 ;
203  case ecal_photons: return eCorrections_.phoEcalP4 ;
204  case regression1: return eCorrections_.regression1P4 ;
205  case regression2: return eCorrections_.regression2P4 ;
206  default: throw cms::Exception("reco::Photon")<<"unexpected p4 type: "<< type << " cannot return p4 ";
207  }
208  }
type
Definition: HCALResponse.h:21
virtual Photon * clone() const
returns a clone of the candidate
Definition: Photon.py:1
void setVertex(const Point &vertex)
set primary event vertex used to define photon direction
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
reco::SuperClusterRef superCluster() const
Ref to SuperCluster.
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
virtual ~Photon()
destructor
Definition: __init__.py:1
math::XYZTLorentzVector LorentzVector
int conversionTrackProvenance(const edm::RefToBase< reco::Track > &convTrack) const
Photon()
default constructor
Definition: Photon.h:32
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
virtual bool overlap(const Candidate &) const
check overlap with another candidate
double p4[4]
Definition: TauolaWrapper.h:92
math::XYZPoint Point
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:17
virtual void setVertex(const Point &vertex)
set vertex
float getCorrectedEnergyError(P4type type) const
const SuperClusterRef & superCluster() const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
float getCorrectedEnergy(P4type type) const
def checkOverlap(process)
fixed size matrix
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
void setP4(P4type type, const LorentzVector &p4, float p4Error, bool setToRecoCandidate)
virtual reco::SuperClusterRef superCluster() const
reference to a SuperCluster