44 : fOutputFileName_(
pset.getUntrackedParameter<
string>(
"HistOutFile",
std::
string(
"TestConversions.root"))),
45 fOutputFile_(nullptr) {
63 h_ErecoEMC_ =
new TH1F(
"deltaE",
" delta(reco-mc) energy", 100, 0., 2.);
64 h_deltaPhi_ =
new TH1F(
"deltaPhi",
" delta(reco-mc) phi", 100, -0.1, 0.1);
65 h_deltaEta_ =
new TH1F(
"deltaEta",
" delta(reco-mc) eta", 100, -0.05, 0.05);
68 h_MCphoE_ =
new TH1F(
"MCphoE",
"MC photon energy", 100, 0., 100.);
69 h_MCphoPhi_ =
new TH1F(
"MCphoPhi",
"MC photon phi", 40, -3.14, 3.14);
70 h_MCphoEta_ =
new TH1F(
"MCphoEta",
"MC photon eta", 40, -3., 3.);
73 h_MCConvE_ =
new TH1F(
"MCConvE",
"MC converted photon energy", 100, 0., 100.);
74 h_MCConvPt_ =
new TH1F(
"MCConvPt",
"MC converted photon pt", 100, 0., 100.);
75 h_MCConvEta_ =
new TH1F(
"MCConvEta",
"MC converted photon eta", 50, 0., 2.5);
78 h_scE_ =
new TH1F(
"scE",
"Uncorrected converted photons : SC Energy ", 100, 0., 200.);
79 h_scEta_ =
new TH1F(
"scEta",
"Uncorrected converted photons: SC Eta ", 40, -3., 3.);
80 h_scPhi_ =
new TH1F(
"scPhi",
"Uncorrected converted photons: SC Phi ", 40, -3.14, 3.14);
82 h_phoE_ =
new TH1F(
"phoE",
"Uncorrected converted photons : Energy ", 100, 0., 200.);
83 h_phoEta_ =
new TH1F(
"phoEta",
"Uncorrected converted photons: Eta ", 40, -3., 3.);
84 h_phoPhi_ =
new TH1F(
"phoPhi",
"Uncorrected converted photons: Phi ", 40, -3.14, 3.14);
87 h2_tk_nHitsVsR_ =
new TH2F(
"tknHitsVsR",
"Tracks Hits vs R ", 12, 0., 120., 20, 0.5, 20.5);
88 h2_tk_inPtVsR_ =
new TH2F(
"tkInPtvsR",
"Tracks inner Pt vs R ", 12, 0., 120., 100, 0., 100.);
95 const float PI = 3.1415927;
99 const float R_ECAL = 136.5;
106 float ZEcal =
R_ECAL * sinh(EtaParticle) + Zvertex;
116 if (EtaParticle < 0.0)
118 float Zlen = Zend - Zvertex;
119 float RR = Zlen / sinh(EtaParticle);
120 Theta = atan(RR / Zend);
132 const float etaPhiDistance = 0.01;
135 const float BARL = 1.4442;
136 const float END_LO = 1.566;
137 const float END_HI = 2.5;
142 LogInfo(
"ConvertedPhotonAnalyzer") <<
"ConvertedPhotonAnalyzer Analyzing event number: " <<
e.id()
143 <<
" Global Counter " <<
nEvt_ <<
"\n";
145 std::cout <<
"ConvertedPhotonAnalyzer Analyzing event number: " <<
e.id() <<
" Global Counter " <<
nEvt_ <<
"\n";
151 std::cout <<
"ConvertedPhotonAnalyzer Converted photon collection size " << phoCollection.size() <<
"\n";
154 std::cout <<
" ConvertedPhotonAnalyzer Looking for MC truth "
158 std::vector<SimTrack> theSimTracks;
159 std::vector<SimVertex> theSimVertices;
163 e.getByLabel(
"g4SimHits", SimTk);
164 e.getByLabel(
"g4SimHits", SimVtx);
166 theSimTracks.insert(theSimTracks.end(), SimTk->begin(), SimTk->end());
167 theSimVertices.insert(theSimVertices.end(), SimVtx->begin(), SimVtx->end());
170 std::cout <<
" ConvertedPhotonAnalyzer mcPhotons size " << mcPhotons.size() << std::endl;
177 for (std::vector<PhotonMCTruth>::const_iterator mcPho = mcPhotons.begin(); mcPho != mcPhotons.end(); mcPho++) {
178 float mcPhi = (*mcPho).fourMomentum().phi();
179 float mcEta = (*mcPho).fourMomentum().pseudoRapidity();
182 if (!(fabs(mcEta) <= BARL || (fabs(mcEta) >= END_LO && fabs(mcEta) <= END_HI))) {
186 std::cout <<
" ConvertedPhotonAnalyzer MC Photons before matching " << std::endl;
187 std::cout <<
" ConvertedPhotonAnalyzer Photons isAconversion " << (*mcPho).isAConversion()
188 <<
" mcMatchingPhoton energy " << (*mcPho).fourMomentum().e() <<
" conversion vertex R "
189 << (*mcPho).vertex().perp() <<
" Z " << (*mcPho).vertex().z() <<
" x " << (*mcPho).vertex().x() <<
" y "
190 << (*mcPho).vertex().y() <<
" z " << (*mcPho).vertex().z() << std::endl;
191 std::cout <<
" ConvertedPhotonAnalyzer mcEta " << mcEta <<
" mcPhi " << mcPhi << std::endl;
193 h_MCphoE_->Fill((*mcPho).fourMomentum().e());
198 if ((*mcPho).isAConversion() == 0)
203 h_MCConvEta_->Fill(fabs((*mcPho).fourMomentum().pseudoRapidity()) - 0.001);
208 std::cout <<
" ConvertedPhotonAnalyzer Starting loop over photon candidates "
210 for (reco::ConversionCollection::const_iterator iPho = phoCollection.begin(); iPho != phoCollection.end(); iPho++) {
213 std::cout <<
" ConvertedPhotonAnalyzer Reco SC energy " << (*iPho).caloCluster()[0]->energy() <<
"\n";
215 float phiClu = (*iPho).caloCluster()[0]->phi();
216 float etaClu = (*iPho).caloCluster()[0]->eta();
227 if (
delta >= etaPhiDistance)
237 std::cout <<
" ConvertedPhotonAnalyzer Matching candidate " << std::endl;
239 std::cout <<
" ConvertedPhotonAnalyzer Photons isAconversion " << (*mcPho).isAConversion()
240 <<
" mcMatchingPhoton energy " << (*mcPho).fourMomentum().e()
241 <<
" ConvertedPhotonAnalyzer conversion vertex R " << (*mcPho).vertex().perp() <<
" Z "
242 << (*mcPho).vertex().z() << std::endl;
244 h_ErecoEMC_->Fill((*iPho).caloCluster()[0]->energy() / (*mcPho).fourMomentum().e());
245 h_deltaPhi_->Fill((*iPho).caloCluster()[0]->position().phi() - mcPhi);
246 h_deltaEta_->Fill((*iPho).caloCluster()[0]->position().eta() - mcEta);
248 h_scE_->Fill((*iPho).caloCluster()[0]->energy());
249 h_scEta_->Fill((*iPho).caloCluster()[0]->position().eta());
250 h_scPhi_->Fill((*iPho).caloCluster()[0]->position().phi());
252 for (
unsigned int i = 0;
i < (*iPho).tracks().size();
i++) {
253 std::cout <<
" ConvertedPhotonAnalyzer Reco Track charge " << (*iPho).tracks()[
i]->charge()
254 <<
" Num of RecHits " << (*iPho).tracks()[
i]->recHitsSize() <<
" inner momentum "
255 <<
sqrt((*iPho).tracks()[
i]->innerMomentum().Mag2()) <<
"\n";
257 h2_tk_nHitsVsR_->Fill((*mcPho).vertex().perp(), (*iPho).tracks()[
i]->recHitsSize());
258 h2_tk_inPtVsR_->Fill((*mcPho).vertex().perp(),
sqrt((*iPho).tracks()[
i]->innerMomentum().Mag2()));
272 std::cout <<
"ConvertedPhotonAnalyzer::endJob Analyzed " <<
nEvt_ <<
" events "