CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
ContainmentCorrectionAnalyzer Class Reference

#include <ContainmentCorrectionAnalyzer.h>

Inheritance diagram for ContainmentCorrectionAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
 ContainmentCorrectionAnalyzer (const edm::ParameterSet &)
 
virtual void endJob ()
 
 ~ContainmentCorrectionAnalyzer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

float ecalEta (float EtaParticle, float Zvertex, float plane_Radius)
 
void fillMcTruth (std::vector< SimTrack > &theSimTracks, std::vector< SimVertex > &theSimVertices)
 
std::vector< EcalSimPhotonMCTruthfindMcTruth (std::vector< SimTrack > &theSimTracks, std::vector< SimVertex > &theSimVertices)
 

Private Attributes

edm::EDGetTokenT
< SuperClusterCollection
BarrelSuperClusterCollection_
 
std::vector< float > e1
 
std::vector< float > e25
 
std::vector< float > e9
 
edm::EDGetTokenT
< SuperClusterCollection
EndcapSuperClusterCollection_
 
std::map< unsigned, unsigned > geantToIndex_
 
TH1F * h_EB_converted
 
TH1F * h_EB_e25EtrueReference
 
TH1F * h_EB_e9EtrueReference
 
TH1F * h_EB_eRecoEtrueReference
 
TH1F * h_EB_eTrue
 
TH1F * h_EE_converted
 
TH1F * h_EE_e25EtrueReference
 
TH1F * h_EE_e9EtrueReference
 
TH1F * h_EE_eRecoEtrueReference
 
TH1F * h_EE_eTrue
 
std::vector< int > iMC
 
std::vector< int > isConverted
 
std::vector< float > mcEnergy
 
std::vector< float > mcEta
 
std::vector< float > mcPhi
 
std::vector< float > mcPt
 
int nMCphotons
 
int nRECOphotons
 
edm::EDGetTokenT
< EcalRecHitCollection
reducedBarrelRecHitCollection_
 
edm::EDGetTokenT
< EcalRecHitCollection
reducedEndcapRecHitCollection_
 
std::vector< int > seedXtal
 
edm::EDGetTokenT
< SimTrackContainer
SimTrackCollection_
 
edm::EDGetTokenT
< SimVertexContainer
SimVertexCollection_
 
std::vector< float > superClusterEnergy
 
std::vector< float > superClusterEt
 
std::vector< float > superClusterEta
 
std::vector< float > superClusterPhi
 
std::vector< float > x_vtx
 
std::vector< float > y_vtx
 
std::vector< float > z_vtx
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 53 of file ContainmentCorrectionAnalyzer.h.

Constructor & Destructor Documentation

ContainmentCorrectionAnalyzer::ContainmentCorrectionAnalyzer ( const edm::ParameterSet ps)
explicit

Definition at line 9 of file ContainmentCorrectionAnalyzer.cc.

References edm::ParameterSet::getParameter().

9  {
10 
11  BarrelSuperClusterCollection_ = consumes<reco::SuperClusterCollection>(ps.getParameter<InputTag>("BarrelSuperClusterCollection"));
12  EndcapSuperClusterCollection_ = consumes<reco::SuperClusterCollection>(ps.getParameter<InputTag>("EndcapSuperClusterCollection"));
13  reducedBarrelRecHitCollection_ = consumes<EcalRecHitCollection>(ps.getParameter<InputTag>("reducedBarrelRecHitCollection"));
14  reducedEndcapRecHitCollection_ = consumes<EcalRecHitCollection>(ps.getParameter<InputTag>("reducedEndcapRecHitCollection"));
15  SimTrackCollection_ = consumes<SimTrackContainer>(ps.getParameter<InputTag>("simTrackCollection"));
16  SimVertexCollection_ = consumes<SimVertexContainer>(ps.getParameter<InputTag>("simVertexCollection"));
17 }
T getParameter(std::string const &) const
edm::EDGetTokenT< EcalRecHitCollection > reducedEndcapRecHitCollection_
edm::EDGetTokenT< SuperClusterCollection > EndcapSuperClusterCollection_
edm::EDGetTokenT< SimTrackContainer > SimTrackCollection_
edm::EDGetTokenT< SimVertexContainer > SimVertexCollection_
edm::EDGetTokenT< SuperClusterCollection > BarrelSuperClusterCollection_
edm::EDGetTokenT< EcalRecHitCollection > reducedBarrelRecHitCollection_
ContainmentCorrectionAnalyzer::~ContainmentCorrectionAnalyzer ( )

Definition at line 19 of file ContainmentCorrectionAnalyzer.cc.

19 { }

Member Function Documentation

void ContainmentCorrectionAnalyzer::analyze ( const edm::Event evt,
const edm::EventSetup es 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 38 of file ContainmentCorrectionAnalyzer.cc.

References deltaR(), alignCSCRings::e, eta(), HcalObjRepresent::Fill(), edm::EventSetup::get(), edm::Event::getByToken(), i, edm::EventBase::id(), edm::ESHandleBase::isValid(), edm::HandleBase::isValid(), ConfigFiles::l, edm::EDConsumerBase::Labels::module, interactiveExample::photons, pi, reco::CaloCluster::position(), funct::pow(), reco::SuperCluster::preshowerEnergy(), edm::Handle< T >::product(), edm::ESHandle< class >::product(), reco::SuperCluster::rawEnergy(), reco::SuperCluster::seed(), mathSSE::sqrt(), and ecaldqm::topology().

38  {
39 
40  LogInfo("ContainmentCorrectionAnalyzer") << "Analyzing event " << evt.id() << "\n";
41 
42  // taking the needed collections
43  std::vector<SimTrack> theSimTracks;
46  Labels l;
48 
49  if (SimTk.isValid() ) theSimTracks.insert(theSimTracks.end(),SimTk->begin(),SimTk->end());
50  else { LogError("ContainmentCorrectionAnalyzer") << "Error! can't get collection with label " << l.module; }
51 
52  std::vector<SimVertex> theSimVertexes;
54  evt.getByToken(SimVertexCollection_,SimVtx);
56 
57 
58  if (SimVtx.isValid()) theSimVertexes.insert(theSimVertexes.end(),SimVtx->begin(),SimVtx->end());
59  else {LogError("ContainmentCorrectionAnalyzer") << "Error! can't get collection with label " << l.module; }
60 
61  const reco::SuperClusterCollection* BarrelSuperClusters = 0;
62  Handle<reco::SuperClusterCollection> pHybridBarrelSuperClusters;
63  evt.getByToken(BarrelSuperClusterCollection_, pHybridBarrelSuperClusters);
65 
66  if (pHybridBarrelSuperClusters.isValid()) { BarrelSuperClusters = pHybridBarrelSuperClusters.product(); }
67  else {LogError("ContainmentCorrectionAnalyzer") << "Error! can't get collection with label " << l.module; }
68 
69  const reco::SuperClusterCollection* EndcapSuperClusters = 0;
70  Handle<reco::SuperClusterCollection> pMulti5x5EndcapSuperClusters;
71  evt.getByToken(EndcapSuperClusterCollection_, pMulti5x5EndcapSuperClusters);
73 
74  if (pMulti5x5EndcapSuperClusters.isValid()) EndcapSuperClusters = pMulti5x5EndcapSuperClusters.product();
75  else {LogError("ContainmentCorrectionAnalyzer") << "Error! can't get collection with label " << l.module; }
76 
77  const EcalRecHitCollection *ebRecHits = 0;
79  evt.getByToken( reducedBarrelRecHitCollection_, pEBRecHits );
81 
82  if (pEBRecHits.isValid()) ebRecHits = pEBRecHits.product();
83  else {LogError("ContainmentCorrectionAnalyzer") << "Error! can't get collection with label " << l.module; }
84 
85  const EcalRecHitCollection *eeRecHits = 0;
87  evt.getByToken( reducedEndcapRecHitCollection_, pEERecHits );
89 
90  if (pEERecHits.isValid()) eeRecHits = pEERecHits.product();
91  else {LogError("ContainmentCorrectionAnalyzer") << "Error! can't get collection with label " << l.module; }
92 
93  const CaloTopology *topology = 0;
94  ESHandle<CaloTopology> pTopology;
95  es.get<CaloTopologyRecord>().get(pTopology);
96  if(pTopology.isValid()) topology = pTopology.product();
97 
98  std::vector<EcalSimPhotonMCTruth> photons=findMcTruth(theSimTracks,theSimVertexes);
99 
100  nMCphotons = 0;
101  nRECOphotons = 0;
102 
103  int mc_size = photons.size();
104  mcEnergy.resize(mc_size);
105  mcEta.resize(mc_size);
106  mcPhi.resize(mc_size);
107  mcPt.resize(mc_size);
108  isConverted.resize(mc_size);
109  x_vtx.resize(mc_size);
110  y_vtx.resize(mc_size);
111  z_vtx.resize(mc_size);
112 
113  superClusterEnergy.resize(mc_size);
114  superClusterEta.resize(mc_size);
115  superClusterPhi.resize(mc_size);
116  superClusterEt.resize(mc_size);
117  e1.resize(mc_size);
118  e9.resize(mc_size);
119  e25.resize(mc_size);
120  seedXtal.resize(mc_size);
121  iMC.resize(mc_size);
122 
123 
124  // loop over MC truth photons
125  for (unsigned int ipho=0;ipho<photons.size();ipho++) {
126 
127  math::XYZTLorentzVectorD vtx = photons[ipho].primaryVertex();
128  double phiTrue = photons[ipho].fourMomentum().phi();
129  double vtxPerp = sqrt(vtx.x()*vtx.x() + vtx.y()*vtx.y());
130  double etaTrue = ecalEta(photons[ipho].fourMomentum().eta(), vtx.z(), vtxPerp);
131  double etTrue = photons[ipho].fourMomentum().e()/cosh(etaTrue);
132  nMCphotons++;
133  mcEnergy[nMCphotons-1]=photons[ipho].fourMomentum().e();
134  mcEta[nMCphotons-1]=etaTrue;
135  mcPhi[nMCphotons-1]=phiTrue;
136  mcPt[nMCphotons-1]=etTrue;
137  isConverted[nMCphotons-1]=photons[ipho].isAConversion();
138  x_vtx[nMCphotons-1]=vtx.x();
139  y_vtx[nMCphotons-1]=vtx.y();
140  z_vtx[nMCphotons-1]=vtx.z();
141 
142  // check histos for MC truth
143  if(std::fabs(etaTrue) < 1.479) {
144  h_EB_eTrue ->Fill(photons[ipho].fourMomentum().e());
145  h_EB_converted ->Fill(photons[ipho].isAConversion());
146  }
147  if(std::fabs(etaTrue) >= 1.6) {
148  h_EE_eTrue ->Fill(photons[ipho].fourMomentum().e());
149  h_EE_converted ->Fill(photons[ipho].isAConversion());
150  }
151 
152  // barrel
153  if(std::fabs(etaTrue) < 1.479) {
154  double etaCurrent; // , etaFound = 0; // UNUSED
155  double phiCurrent;
156  // double etCurrent, etFound = 0; // UNUSED
157  const reco::SuperCluster* nearSC = 0;
158 
159  double closestParticleDistance = 999;
160  for(reco::SuperClusterCollection::const_iterator aClus = BarrelSuperClusters->begin();
161  aClus != BarrelSuperClusters->end(); aClus++) {
162  etaCurrent = aClus->position().eta();
163  phiCurrent = aClus->position().phi();
164  // etCurrent = aClus->energy()/std::cosh(etaCurrent); // UNUSED
165  double deltaR = std::sqrt(std::pow(etaCurrent-etaTrue,2)+std::pow(phiCurrent-phiTrue,2));
166  if(deltaR < closestParticleDistance) {
167  // etFound = etCurrent; // UNUSED
168  // etaFound = etaCurrent; // UNUSED
169  closestParticleDistance = deltaR;
170  nearSC=&(*aClus);
171  }
172  }
173 
174  if(closestParticleDistance < 0.3) {
175  // Is a matched particle dumping informations
176  nRECOphotons++;
178  superClusterEta[nRECOphotons-1]=nearSC->position().eta();
179  superClusterPhi[nRECOphotons-1]=nearSC->position().phi();
180  superClusterEt[nRECOphotons-1]=nearSC->rawEnergy()/std::cosh(nearSC->position().eta());
182 
183  reco::CaloClusterPtr theSeed = nearSC->seed();
184  seedXtal[nRECOphotons-1] = EcalClusterTools::getMaximum(*theSeed, ebRecHits).first;
185  e1[nRECOphotons-1] = EcalClusterTools::eMax(*theSeed, ebRecHits);
186  e9[nRECOphotons-1] = EcalClusterTools::e3x3(*theSeed, ebRecHits, topology );
187  e25[nRECOphotons-1] = EcalClusterTools::e5x5(*theSeed, ebRecHits, topology );
188  }
189  }
190 
191  // endcap
192  if(std::fabs(etaTrue) >= 1.6) {
193  double etaCurrent; // , etaFound = 0; // UNUSED
194  double phiCurrent;
195  // double etCurrent, etFound = 0; // UNUSED
196  const reco::SuperCluster* nearSC = 0;
197 
198  double closestParticleDistance = 999;
199  for(reco::SuperClusterCollection::const_iterator aClus = EndcapSuperClusters->begin();
200  aClus != EndcapSuperClusters->end(); aClus++) {
201  etaCurrent = aClus->position().eta();
202  phiCurrent = aClus->position().phi();
203  // etCurrent = aClus->energy()/std::cosh(etaCurrent);
204  double deltaR = std::sqrt(std::pow(etaCurrent-etaTrue,2)+std::pow(phiCurrent-phiTrue,2));
205  if(deltaR < closestParticleDistance) {
206  // etFound = etCurrent; // UNUSED
207  // etaFound = etaCurrent; // UNUSED
208  closestParticleDistance = deltaR;
209  nearSC=&(*aClus);
210  }
211  }
212 
213  if(closestParticleDistance < 0.3) {
214  //Is a matched particle dumping informations
215  nRECOphotons++;
216  float psEnergy = nearSC->preshowerEnergy();
217  superClusterEnergy[nRECOphotons-1]=nearSC->rawEnergy() + psEnergy;
218  superClusterEta[nRECOphotons-1]=nearSC->position().eta();
219  superClusterPhi[nRECOphotons-1]=nearSC->position().phi();
220  superClusterEt[nRECOphotons-1]=(nearSC->rawEnergy() + psEnergy)/std::cosh(nearSC->position().eta());
222 
223  reco::CaloClusterPtr theSeed = nearSC->seed();
224  seedXtal[nRECOphotons-1] = EcalClusterTools::getMaximum(*theSeed, eeRecHits).first;
225  e1[nRECOphotons-1] = EcalClusterTools::eMax(*theSeed, eeRecHits) + psEnergy;
226  e9[nRECOphotons-1] = EcalClusterTools::e3x3(*theSeed, eeRecHits, topology ) + psEnergy;
227  e25[nRECOphotons-1] = EcalClusterTools::e5x5(*theSeed, eeRecHits, topology ) + psEnergy;
228  }
229  }
230  }
231 
232 
233  // containment analysis for unconverted photons in the reference region only
234  for (int i=0;i<nRECOphotons;i++) {
235 
236  // barrel
237  if (fabs(superClusterEta[i]) < 1.479 ) {
238  if (isConverted[iMC[i]] != 1){
239  int ietaAbs = (seedXtal[i]>>9)&0x7F;
240  int iphi = seedXtal[i]&0x1FF;
241  if (ietaAbs > 5 && ietaAbs < 21 && ((iphi % 20) > 5) && ((iphi % 20) < 16) ) {
245  }
246  }
247  }
248 
249  // endcap
250  if (fabs(superClusterEta[i]) > 1.6 ) {
251  if (isConverted[iMC[i]] != 1) {
252  if ( fabs(superClusterEta[i]) > 1.7 && fabs(superClusterEta[i] < 2.3) &&
253  ((superClusterPhi[i] > -CLHEP::pi/2. + 0.1 && superClusterPhi[i] < CLHEP::pi/2. - 0.1) ||
254  (superClusterPhi[i] > CLHEP::pi/2. + 0.1) ||
255  (superClusterPhi[i] < -CLHEP::pi/2. - 0.1)
256  ))
257  {
261  }
262  }
263  }
264 
265  } // loop over reco photons
266 }
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:124
int i
Definition: DBlmapReader.cc:9
std::vector< EcalSimPhotonMCTruth > findMcTruth(std::vector< SimTrack > &theSimTracks, std::vector< SimVertex > &theSimVertices)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
CaloTopology const * topology(0)
edm::EDGetTokenT< EcalRecHitCollection > reducedEndcapRecHitCollection_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
edm::EDGetTokenT< SuperClusterCollection > EndcapSuperClusterCollection_
edm::EDGetTokenT< SimTrackContainer > SimTrackCollection_
T eta() const
const Double_t pi
float ecalEta(float EtaParticle, float Zvertex, float plane_Radius)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:48
edm::EDGetTokenT< SimVertexContainer > SimVertexCollection_
bool isValid() const
Definition: HandleBase.h:76
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
edm::EDGetTokenT< SuperClusterCollection > BarrelSuperClusterCollection_
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:47
edm::EDGetTokenT< EcalRecHitCollection > reducedBarrelRecHitCollection_
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:81
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
edm::EventID id() const
Definition: EventBase.h:56
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:66
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:50
bool isValid() const
Definition: ESHandle.h:37
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void ContainmentCorrectionAnalyzer::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 21 of file ContainmentCorrectionAnalyzer.cc.

References TFileService::make().

21  {
22 
24 
25  // Define reference histograms
26  h_EB_eRecoEtrueReference = fs->make<TH1F>("EB_eRecoEtrueReference","EB_eRecoEtrueReference",440,0.,1.1);
27  h_EB_e9EtrueReference = fs->make<TH1F>("EB_e9EtrueReference", "EB_e9EtrueReference", 440,0.,1.1);
28  h_EB_e25EtrueReference = fs->make<TH1F>("EB_e25EtrueReference", "EB_e25EtrueReference", 440,0.,1.1);
29  h_EE_eRecoEtrueReference = fs->make<TH1F>("EE_eRecoEtrueReference","EE_eRecoEtrueReference",440,0.,1.1);
30  h_EE_e9EtrueReference = fs->make<TH1F>("EE_e9EtrueReference", "EE_e9EtrueReference", 440,0.,1.1);
31  h_EE_e25EtrueReference = fs->make<TH1F>("EE_e25EtrueReference", "EE_e25EtrueReference", 440,0.,1.1);
32  h_EB_eTrue = fs->make<TH1F>("EB_eTrue", "EB_eTrue", 41,40.,60.);
33  h_EE_eTrue = fs->make<TH1F>("EE_eTrue", "EE_eTrue", 41,40.,60.);
34  h_EB_converted = fs->make<TH1F>("EB_converted", "EB_converted", 2,0.,2.);
35  h_EE_converted = fs->make<TH1F>("EE_converted", "EE_converted", 2,0.,2.);
36 }
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
float ContainmentCorrectionAnalyzer::ecalEta ( float  EtaParticle,
float  Zvertex,
float  plane_Radius 
)
private

Definition at line 270 of file ContainmentCorrectionAnalyzer.cc.

References ETA, etaBarrelEndcap, fff_deleter::log, Geom::pi(), R_ECAL, funct::tan(), and Z_Endcap.

270  {
271 
272  const float R_ECAL = 136.5;
273  const float Z_Endcap = 328.0;
274  const float etaBarrelEndcap = 1.479;
275 
276  if(EtaParticle != 0.) {
277  float Theta = 0.0 ;
278  float ZEcal = (R_ECAL-plane_Radius)*sinh(EtaParticle)+Zvertex;
279 
280  if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
281  if(Theta<0.0) Theta = Theta+Geom::pi() ;
282 
283  float ETA = - log(tan(0.5*Theta));
284 
285  if( fabs(ETA) > etaBarrelEndcap ) {
286  float Zend = Z_Endcap ;
287  if(EtaParticle<0.0 ) Zend = -Zend ;
288  float Zlen = Zend - Zvertex ;
289  float RR = Zlen/sinh(EtaParticle);
290  Theta = atan((RR+plane_Radius)/Zend);
291  if(Theta<0.0) Theta = Theta+Geom::pi() ;
292  ETA = - log(tan(0.5*Theta));
293  }
294 
295  return ETA;
296  }
297  else {
298  LogWarning("") << "[ContainmentCorrectionAnalyzer::ecalEta] Warning: Eta equals to zero, not correcting" ;
299  return EtaParticle;
300  }
301 }
#define ETA
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
static const float etaBarrelEndcap
static const float Z_Endcap
double pi()
Definition: Pi.h:31
static const float R_ECAL
void ContainmentCorrectionAnalyzer::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 268 of file ContainmentCorrectionAnalyzer.cc.

268 { }
void ContainmentCorrectionAnalyzer::fillMcTruth ( std::vector< SimTrack > &  theSimTracks,
std::vector< SimVertex > &  theSimVertices 
)
private

Definition at line 436 of file ContainmentCorrectionAnalyzer.cc.

436  {
437 
438  unsigned nVtx = simVertices.size();
439  unsigned nTks = simTracks.size();
440  if ( nVtx == 0 ) return;
441  // create a map associating geant particle id and position in the event SimTrack vector
442  for( unsigned it=0; it<nTks; ++it ) {
443  geantToIndex_[ simTracks[it].trackId() ] = it;
444  }
445 }
std::map< unsigned, unsigned > geantToIndex_
std::vector< EcalSimPhotonMCTruth > ContainmentCorrectionAnalyzer::findMcTruth ( std::vector< SimTrack > &  theSimTracks,
std::vector< SimVertex > &  theSimVertices 
)
private

store this electron since it's from a converted photon

Definition at line 304 of file ContainmentCorrectionAnalyzer.cc.

References funct::abs(), if(), SimVertex::parentIndex(), CoreSimVertex::position(), and query::result.

304  {
305 
306  std::vector<EcalSimPhotonMCTruth> result;
307 
308  geantToIndex_.clear();
309  // int idTrk1_[10]; // UNUSED
310  // int idTrk2_[10]; // UNUSED
311 
312  // Local variables
313  // const int SINGLE=1; // UNUSED
314  // const int DOUBLE=2; // UNUSED
315  // const int PYTHIA=3; // UNUSED
316  const int ELECTRON_FLAV=1;
317  const int PIZERO_FLAV=2;
318  const int PHOTON_FLAV=3;
319 
320  // int ievtype=0; // UNUSED
321  int ievflav=0;
322  std::vector<SimTrack*> photonTracks;
323  std::vector<SimTrack*> pizeroTracks;
324  std::vector<const SimTrack *> trkFromConversion;
325  SimVertex primVtx;
326  std::vector<int> convInd;
327 
328  fillMcTruth(theSimTracks, theSimVertices);
329  int iPV=-1;
330  int partType1=0;
331  int partType2=0;
332  std::vector<SimTrack>::iterator iFirstSimTk = theSimTracks.begin();
333  if ( !(*iFirstSimTk).noVertex() ) {
334  iPV = (*iFirstSimTk).vertIndex();
335  int vtxId = (*iFirstSimTk).vertIndex();
336  primVtx = theSimVertices[vtxId];
337  partType1 = (*iFirstSimTk).type();
338  }
339 
340  // Look at a second track
341  iFirstSimTk++;
342  if ( iFirstSimTk!= theSimTracks.end() ) {
343  if ( (*iFirstSimTk).vertIndex() == iPV) {
344  partType2 = (*iFirstSimTk).type();
345  }
346  }
347  int npv=0;
348  int iPho=0;
349  for (std::vector<SimTrack>::iterator iSimTk = theSimTracks.begin(); iSimTk != theSimTracks.end(); ++iSimTk){
350  if ( (*iSimTk).noVertex() ) continue;
351  // int vertexId = (*iSimTk).vertIndex(); // UNUSED
352  // SimVertex vertex = theSimVertices[vertexId]; // UNUSED
353  if ( (*iSimTk).vertIndex() == iPV ) {
354  npv++;
355  if ( (*iSimTk).type() == 22) {
356  convInd.push_back(0);
357  photonTracks.push_back( &(*iSimTk) );
358  // math::XYZTLorentzVectorD momentum = (*iSimTk).momentum(); // UNUSED
359  }
360  }
361  }
362 
363  if(npv > 4) { // ievtype = PYTHIA; // UNUSED
364  } else if(npv == 1) {
365  if( abs(partType1) == 11 ) { /* ievtype = SINGLE; ==UNUSED== */ ievflav = ELECTRON_FLAV; }
366  else if(partType1 == 111) { /* ievtype = SINGLE; ==UNUSED== */ ievflav = PIZERO_FLAV; }
367  else if(partType1 == 22) { /* ievtype = SINGLE; ==UNUSED== */ ievflav = PHOTON_FLAV; }
368  } else if(npv == 2) {
369  if ( abs(partType1) == 11 && abs(partType2) == 11 ) { /* ievtype = DOUBLE; ==UNUSED== */ ievflav = ELECTRON_FLAV; }
370  else if(partType1 == 111 && partType2 == 111) { /* ievtype = DOUBLE; ==UNUSED== */ ievflav = PIZERO_FLAV; }
371  else if(partType1 == 22 && partType2 == 22) { /* ievtype = DOUBLE; ==UNUSED== */ ievflav = PHOTON_FLAV; }
372  }
373 
374  // Look into converted photons
375  int isAconversion=0;
376  if(ievflav == PHOTON_FLAV) {
377 
378  int nConv=0;
379  int iConv=0;
380  iPho=0;
381  for (std::vector<SimTrack*>::iterator iPhoTk = photonTracks.begin(); iPhoTk != photonTracks.end(); ++iPhoTk){
382  trkFromConversion.clear();
383  for (std::vector<SimTrack>::iterator iSimTk = theSimTracks.begin(); iSimTk != theSimTracks.end(); ++iSimTk){
384  if ( (*iSimTk).noVertex() ) continue;
385  if ( (*iSimTk).vertIndex() == iPV ) continue;
386  if ( abs((*iSimTk).type()) != 11 ) continue;
387  int vertexId = (*iSimTk).vertIndex();
388  SimVertex vertex = theSimVertices[vertexId];
389  int motherId=-1;
390  if ( vertex.parentIndex() ) {
391  unsigned motherGeantId = vertex.parentIndex();
392  std::map<unsigned, unsigned >::iterator association = geantToIndex_.find( motherGeantId );
393  if(association != geantToIndex_.end() )
394  motherId = association->second;
395  //int motherType = motherId == -1 ? 0 : theSimTracks[motherId].type();
396 
397  if ( theSimTracks[motherId].trackId() == (*iPhoTk)->trackId() ) {
399  trkFromConversion.push_back(&(*iSimTk ) );
400  }
401  }
402  } // loop over the SimTracks
403 
404  if ( trkFromConversion.size() > 0 ) {
405  isAconversion=1;
406  nConv++;
407  convInd[iPho]=nConv;
408  int convVtxId = trkFromConversion[0]->vertIndex();
409  SimVertex convVtx = theSimVertices[convVtxId];
410  math::XYZTLorentzVectorD vtxPosition = convVtx.position();
411  // math::XYZTLorentzVectorD momentum = (*iPhoTk)->momentum(); // UNUSED
412  if ( nConv <= 10) {
413  if ( trkFromConversion.size() > 1) {
414  // idTrk1_[iConv]= trkFromConversion[0]->trackId(); // UNUSED
415  // idTrk2_[iConv]= trkFromConversion[1]->trackId(); // UNUSED
416  } else {
417  // idTrk1_[iConv]=trkFromConversion[0]->trackId(); // UNUSED
418  // idTrk2_[iConv]=-1; // UNUSED
419  }
420  }
421  iConv++;
422 
423  result.push_back( EcalSimPhotonMCTruth(isAconversion, (*iPhoTk)->momentum(), vtxPosition.pt(), vtxPosition.z() , vtxPosition, primVtx.position(), trkFromConversion ));
424  } else {
425  isAconversion=0;
426  math::XYZTLorentzVectorD vtxPosition(0.,0.,0.,0.);
427  result.push_back( EcalSimPhotonMCTruth(isAconversion, (*iPhoTk)->momentum(), vtxPosition.pt(), vtxPosition.z() , vtxPosition, primVtx.position(), trkFromConversion ));
428  }
429  iPho++;
430  } // loop over the primary photons
431  } // Event with one or two photons
432 
433  return result;
434 }
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
std::map< unsigned, unsigned > geantToIndex_
int parentIndex() const
Definition: SimVertex.h:33
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:26
if(conf.exists("allCellsPositionCalc"))
void fillMcTruth(std::vector< SimTrack > &theSimTracks, std::vector< SimVertex > &theSimVertices)

Member Data Documentation

edm::EDGetTokenT<SuperClusterCollection> ContainmentCorrectionAnalyzer::BarrelSuperClusterCollection_
private

Definition at line 70 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::e1
private

Definition at line 81 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::e25
private

Definition at line 81 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::e9
private

Definition at line 81 of file ContainmentCorrectionAnalyzer.h.

edm::EDGetTokenT<SuperClusterCollection> ContainmentCorrectionAnalyzer::EndcapSuperClusterCollection_
private

Definition at line 71 of file ContainmentCorrectionAnalyzer.h.

std::map<unsigned, unsigned> ContainmentCorrectionAnalyzer::geantToIndex_
private

Definition at line 85 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EB_converted
private

Definition at line 95 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EB_e25EtrueReference
private

Definition at line 89 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EB_e9EtrueReference
private

Definition at line 88 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EB_eRecoEtrueReference
private

Definition at line 87 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EB_eTrue
private

Definition at line 93 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EE_converted
private

Definition at line 96 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EE_e25EtrueReference
private

Definition at line 92 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EE_e9EtrueReference
private

Definition at line 91 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EE_eRecoEtrueReference
private

Definition at line 90 of file ContainmentCorrectionAnalyzer.h.

TH1F* ContainmentCorrectionAnalyzer::h_EE_eTrue
private

Definition at line 94 of file ContainmentCorrectionAnalyzer.h.

std::vector<int> ContainmentCorrectionAnalyzer::iMC
private

Definition at line 83 of file ContainmentCorrectionAnalyzer.h.

std::vector<int> ContainmentCorrectionAnalyzer::isConverted
private

Definition at line 77 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::mcEnergy
private

Definition at line 76 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::mcEta
private

Definition at line 76 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::mcPhi
private

Definition at line 76 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::mcPt
private

Definition at line 76 of file ContainmentCorrectionAnalyzer.h.

int ContainmentCorrectionAnalyzer::nMCphotons
private

Definition at line 75 of file ContainmentCorrectionAnalyzer.h.

int ContainmentCorrectionAnalyzer::nRECOphotons
private

Definition at line 79 of file ContainmentCorrectionAnalyzer.h.

edm::EDGetTokenT<EcalRecHitCollection> ContainmentCorrectionAnalyzer::reducedBarrelRecHitCollection_
private

Definition at line 72 of file ContainmentCorrectionAnalyzer.h.

edm::EDGetTokenT<EcalRecHitCollection> ContainmentCorrectionAnalyzer::reducedEndcapRecHitCollection_
private

Definition at line 73 of file ContainmentCorrectionAnalyzer.h.

std::vector<int> ContainmentCorrectionAnalyzer::seedXtal
private

Definition at line 82 of file ContainmentCorrectionAnalyzer.h.

edm::EDGetTokenT<SimTrackContainer> ContainmentCorrectionAnalyzer::SimTrackCollection_
private

Definition at line 68 of file ContainmentCorrectionAnalyzer.h.

edm::EDGetTokenT<SimVertexContainer> ContainmentCorrectionAnalyzer::SimVertexCollection_
private

Definition at line 69 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::superClusterEnergy
private

Definition at line 80 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::superClusterEt
private

Definition at line 80 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::superClusterEta
private

Definition at line 80 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::superClusterPhi
private

Definition at line 80 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::x_vtx
private

Definition at line 78 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::y_vtx
private

Definition at line 78 of file ContainmentCorrectionAnalyzer.h.

std::vector<float> ContainmentCorrectionAnalyzer::z_vtx
private

Definition at line 78 of file ContainmentCorrectionAnalyzer.h.