CMS 3D CMS Logo

MCPizeroAnalyzer.cc
Go to the documentation of this file.
19 
20 #include "TFile.h"
21 #include "TH1.h"
22 #include "TH2.h"
23 #include "TTree.h"
24 #include "TVector3.h"
25 #include "TProfile.h"
26 
27 #include <iostream>
28 #include <map>
29 #include <vector>
30 
32 public:
33  //
34  explicit MCPizeroAnalyzer(const edm::ParameterSet&);
35  ~MCPizeroAnalyzer() override;
36 
37  void analyze(const edm::Event&, const edm::EventSetup&) override;
38  void beginJob() override;
39  void endJob() override;
40 
41 private:
42  float etaTransformation(float a, float b);
43  float phiNormalization(float& a);
44 
45  //
47 
49  TFile* fOutputFile_;
50 
51  int nEvt_;
52  int nMatched_;
53 
55  double mcPhi_;
56  double mcEta_;
57 
62 
63  TH1F* h_MCPizE_;
64  TH1F* h_MCPizEta_;
68  TH1F* h_MCPizPhi_;
71 
72  TH1F* h_MCEleE_;
73  TH1F* h_MCEleEta_;
74  TH1F* h_MCElePhi_;
75  TH1F* h_BremFrac_;
77 
79 
80  TH1F* h_MCPhoE_;
81  TH1F* h_MCPhoEta_;
82  TH1F* h_MCPhoPhi_;
87 };
88 
91 
92 using namespace std;
93 
95  : fOutputFileName_(pset.getUntrackedParameter<string>("HistOutFile", std::string("TestConversions.root"))),
96  fOutputFile_(nullptr) {}
97 
99 
101  nEvt_ = 0;
102 
104 
105  fOutputFile_ = new TFile(fOutputFileName_.c_str(), "RECREATE");
106 
108  h_MCPizE_ = new TH1F("MCPizE", "MC piz energy", 100, 0., 200.);
109  h_MCPizPhi_ = new TH1F("MCPizPhi", "MC piz phi", 40, -3.14, 3.14);
110  h_MCPizEta_ = new TH1F("MCPizEta", "MC piz eta", 40, -3., 3.);
111  h_MCPizUnEta_ = new TH1F("MCPizUnEta", "MC un piz eta", 40, -3., 3.);
112  h_MCPiz1ConEta_ = new TH1F("MCPiz1ConEta", "MC con piz eta: at least one converted photon", 40, -3., 3.);
113  h_MCPiz2ConEta_ = new TH1F("MCPiz2ConEta", "MC con piz eta: two converted photons", 40, -3., 3.);
114  h_MCPizMass1_ = new TH1F("MCPizMass1", "Piz mass unconverted ", 100, 0., 200);
115  h_MCPizMass2_ = new TH1F("MCPizMass2", "Piz mass converted ", 100, 0., 200);
116 
117  // All Photons from Pizeros
118  h_MCPhoE_ = new TH1F("MCPhoE", "MC photon energy", 100, 0., 200.);
119  h_MCPhoPhi_ = new TH1F("MCPhoPhi", "MC photon phi", 40, -3.14, 3.14);
120  h_MCPhoEta_ = new TH1F("MCPhoEta", "MC photon eta", 40, -3., 3.);
121 
122  // Converted photons
123  h_MCConvPhoE_ = new TH1F("MCConvPhoE", "MC converted photon energy", 100, 0., 200.);
124  h_MCConvPhoPhi_ = new TH1F("MCConvPhoPhi", "MC converted photon phi", 40, -3.14, 3.14);
125  h_MCConvPhoEta_ = new TH1F("MCConvPhoEta", "MC converted photon eta", 40, -3., 3.);
126  h_MCConvPhoR_ = new TH1F("MCConvPhoR", "MC converted photon R", 120, 0., 120.);
127  // Electrons from converted photons
128  h_MCEleE_ = new TH1F("MCEleE", "MC ele energy", 100, 0., 200.);
129  h_MCElePhi_ = new TH1F("MCElePhi", "MC ele phi", 40, -3.14, 3.14);
130  h_MCEleEta_ = new TH1F("MCEleEta", "MC ele eta", 40, -3., 3.);
131  h_BremFrac_ = new TH1F("bremFrac", "brem frac ", 50, 0., 1.);
132  h_BremEnergy_ = new TH1F("bremE", "Brem energy", 100, 0., 200.);
133 
134  h_EleEvsPhoE_ = new TH2F("eleEvsPhoE", "eleEvsPhoE", 100, 0., 200., 100, 0., 200.);
135 
136  return;
137 }
138 
139 float MCPizeroAnalyzer::etaTransformation(float EtaParticle, float Zvertex) {
140  //---Definitions
141  const float PI = 3.1415927;
142  //UNUSED const float TWOPI = 2.0*PI;
143 
144  //---Definitions for ECAL
145  const float R_ECAL = 136.5;
146  const float Z_Endcap = 328.0;
147  const float etaBarrelEndcap = 1.479;
148 
149  //---ETA correction
150 
151  float Theta = 0.0;
152  float ZEcal = R_ECAL * sinh(EtaParticle) + Zvertex;
153 
154  if (ZEcal != 0.0)
155  Theta = atan(R_ECAL / ZEcal);
156  if (Theta < 0.0)
157  Theta = Theta + PI;
158  float ETA = -log(tan(0.5 * Theta));
159 
160  if (fabs(ETA) > etaBarrelEndcap) {
161  float Zend = Z_Endcap;
162  if (EtaParticle < 0.0)
163  Zend = -Zend;
164  float Zlen = Zend - Zvertex;
165  float RR = Zlen / sinh(EtaParticle);
166  Theta = atan(RR / Zend);
167  if (Theta < 0.0)
168  Theta = Theta + PI;
169  ETA = -log(tan(0.5 * Theta));
170  }
171  //---Return the result
172  return ETA;
173  //---end
174 }
175 
177  //---Definitions
178  const float PI = 3.1415927;
179  const float TWOPI = 2.0 * PI;
180 
181  if (phi > PI) {
182  phi = phi - TWOPI;
183  }
184  if (phi < -PI) {
185  phi = phi + TWOPI;
186  }
187 
188  // cout << " Float_t PHInormalization out " << PHI << endl;
189  return phi;
190 }
191 
193  using namespace edm;
194  //UNUSED const float etaPhiDistance=0.01;
195  // Fiducial region
196  //UNUSED const float TRK_BARL =0.9;
197  //UNUSED const float BARL = 1.4442; // DAQ TDR p.290
198  //UNUSED const float END_LO = 1.566;
199  //UNUSED const float END_HI = 2.5;
200  // Electron mass
201  //UNUSED const Float_t mElec= 0.000511;
202 
203  nEvt_++;
204  LogInfo("MCPizeroAnalyzer") << "MCPizeroAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_
205  << "\n";
206  // LogDebug("MCPizeroAnalyzer") << "MCPizeroAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
207  std::cout << "MCPizeroAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ << "\n";
208 
210  std::cout << " MCPizeroAnalyzer Looking for MC truth "
211  << "\n";
212 
213  //get simtrack info
214  std::vector<SimTrack> theSimTracks;
215  std::vector<SimVertex> theSimVertices;
216 
219  e.getByLabel("g4SimHits", SimTk);
220  e.getByLabel("g4SimHits", SimVtx);
221 
222  theSimTracks.insert(theSimTracks.end(), SimTk->begin(), SimTk->end());
223  theSimVertices.insert(theSimVertices.end(), SimVtx->begin(), SimVtx->end());
224  std::cout << " MCPizeroAnalyzer This Event has " << theSimTracks.size() << " sim tracks " << std::endl;
225  std::cout << " MCPizeroAnalyzer This Event has " << theSimVertices.size() << " sim vertices " << std::endl;
226  if (theSimTracks.empty())
227  std::cout << " Event number " << e.id() << " has NO sim tracks " << std::endl;
228 
229  std::vector<PizeroMCTruth> MCPizeroeros = thePizeroMCTruthFinder_->find(theSimTracks, theSimVertices);
230  std::cout << " MCPizeroAnalyzer MCPizeroeros size " << MCPizeroeros.size() << std::endl;
231 
232  for (std::vector<PizeroMCTruth>::const_iterator iPiz = MCPizeroeros.begin(); iPiz != MCPizeroeros.end(); ++iPiz) {
233  h_MCPizE_->Fill((*iPiz).fourMomentum().e());
234  h_MCPizEta_->Fill((*iPiz).fourMomentum().pseudoRapidity());
235  h_MCPizPhi_->Fill((*iPiz).fourMomentum().phi());
236 
237  std::vector<PhotonMCTruth> mcPhotons = (*iPiz).photons();
238  std::cout << " MCPizeroAnalyzer mcPhotons size " << mcPhotons.size() << std::endl;
239 
240  float px = mcPhotons[0].fourMomentum().x() + mcPhotons[1].fourMomentum().x();
241  float py = mcPhotons[0].fourMomentum().y() + mcPhotons[1].fourMomentum().y();
242  float pz = mcPhotons[0].fourMomentum().z() + mcPhotons[1].fourMomentum().z();
243  float e = mcPhotons[0].fourMomentum().e() + mcPhotons[1].fourMomentum().e();
244  float invM = sqrt(e * e - px * px - py * py - pz * pz) * 1000;
245  h_MCPizMass1_->Fill(invM);
246 
247  int converted = 0;
248  for (std::vector<PhotonMCTruth>::const_iterator iPho = mcPhotons.begin(); iPho != mcPhotons.end(); ++iPho) {
249  h_MCPhoE_->Fill((*iPho).fourMomentum().e());
250  h_MCPhoEta_->Fill((*iPho).fourMomentum().pseudoRapidity());
251  h_MCPhoPhi_->Fill((*iPho).fourMomentum().phi());
252  if ((*iPho).isAConversion()) {
253  converted++;
254 
255  h_MCConvPhoE_->Fill((*iPho).fourMomentum().e());
256  h_MCConvPhoEta_->Fill((*iPho).fourMomentum().pseudoRapidity());
257  h_MCConvPhoPhi_->Fill((*iPho).fourMomentum().phi());
258  h_MCConvPhoR_->Fill((*iPho).vertex().perp());
259 
260  std::vector<ElectronMCTruth> mcElectrons = (*iPho).electrons();
261  std::cout << " MCPizeroAnalyzer mcElectrons size " << mcElectrons.size() << std::endl;
262 
263  for (std::vector<ElectronMCTruth>::const_iterator iEl = mcElectrons.begin(); iEl != mcElectrons.end(); ++iEl) {
264  if ((*iEl).fourMomentum().e() < 30)
265  continue;
266  h_MCEleE_->Fill((*iEl).fourMomentum().e());
267  h_MCEleEta_->Fill((*iEl).fourMomentum().pseudoRapidity());
268  h_MCElePhi_->Fill((*iEl).fourMomentum().phi());
269 
270  h_EleEvsPhoE_->Fill((*iPho).fourMomentum().e(), (*iEl).fourMomentum().e());
271 
272  float totBrem = 0;
273  for (unsigned int iBrem = 0; iBrem < (*iEl).bremVertices().size(); ++iBrem)
274  totBrem += (*iEl).bremMomentum()[iBrem].e();
275 
276  h_BremFrac_->Fill(totBrem / (*iEl).fourMomentum().e());
277  h_BremEnergy_->Fill(totBrem);
278  }
279  }
280  }
281 
282  if (converted > 0) {
283  h_MCPiz1ConEta_->Fill((*iPiz).fourMomentum().pseudoRapidity());
284  if (converted == 2)
285  h_MCPiz2ConEta_->Fill((*iPiz).fourMomentum().pseudoRapidity());
286  } else {
287  h_MCPizUnEta_->Fill((*iPiz).fourMomentum().pseudoRapidity());
288  }
289  }
290 }
291 
293  fOutputFile_->Write();
294  fOutputFile_->Close();
295 
296  edm::LogInfo("MCPizeroAnalyzer") << "Analyzed " << nEvt_ << "\n";
297  std::cout << "MCPizeroAnalyzer::endJob Analyzed " << nEvt_ << " events "
298  << "\n";
299 
300  return;
301 }
MCPizeroAnalyzer::mcPhi_
double mcPhi_
global variable for the MC photon
Definition: MCPizeroAnalyzer.cc:55
PI
Definition: PayloadInspector.h:21
Handle.h
EDAnalyzer.h
MessageLogger.h
MCPizeroAnalyzer::etaTransformation
float etaTransformation(float a, float b)
Definition: MCPizeroAnalyzer.cc:139
MCPizeroAnalyzer::h_MCEleEta_
TH1F * h_MCEleEta_
Definition: MCPizeroAnalyzer.cc:73
TrackerGeometry.h
PizeroMCTruthFinder.h
MCPizeroAnalyzer::h_MCPizPhi_
TH1F * h_MCPizPhi_
Definition: MCPizeroAnalyzer.cc:68
ESHandle.h
PI
#define PI
Definition: QcdUeDQM.h:37
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition: AlignableModifier.h:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
MCPizeroAnalyzer
Definition: MCPizeroAnalyzer.cc:31
MCPizeroAnalyzer::nEvt_
int nEvt_
Definition: MCPizeroAnalyzer.cc:51
MCPizeroAnalyzer::h_MCConvPhoR_
TH1F * h_MCConvPhoR_
Definition: MCPizeroAnalyzer.cc:86
MCPizeroAnalyzer::h_MCConvPhoPhi_
TH1F * h_MCConvPhoPhi_
Definition: MCPizeroAnalyzer.cc:85
MCPizeroAnalyzer::thePizeroMCTruthFinder_
PizeroMCTruthFinder * thePizeroMCTruthFinder_
Definition: MCPizeroAnalyzer.cc:46
MCPizeroAnalyzer::h_MCEleE_
TH1F * h_MCEleE_
Definition: MCPizeroAnalyzer.cc:72
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
MCPizeroAnalyzer::h_MCPizMass1_
TH1F * h_MCPizMass1_
Definition: MCPizeroAnalyzer.cc:69
edm::Handle
Definition: AssociativeIterator.h:50
MCPizeroAnalyzer::h_MCPizEta_
TH1F * h_MCPizEta_
Definition: MCPizeroAnalyzer.cc:64
PizeroMCTruthFinder
Definition: PizeroMCTruthFinder.h:15
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ElectronMCTruth.h
MCPizeroAnalyzer::fOutputFile_
TFile * fOutputFile_
Definition: MCPizeroAnalyzer.cc:49
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
SimVertex.h
MCPizeroAnalyzer::~MCPizeroAnalyzer
~MCPizeroAnalyzer() override
Definition: MCPizeroAnalyzer.cc:98
IdealMagneticFieldRecord.h
MCPizeroAnalyzer::h_MCPiz2ConEta_
TH1F * h_MCPiz2ConEta_
Definition: MCPizeroAnalyzer.cc:67
MCPizeroAnalyzer::h_BremFrac_
TH1F * h_BremFrac_
Definition: MCPizeroAnalyzer.cc:75
b
double b
Definition: hdecay.h:118
MCPizeroAnalyzer::fOutputFileName_
std::string fOutputFileName_
Definition: MCPizeroAnalyzer.cc:48
MCPizeroAnalyzer::h_MCElePhi_
TH1F * h_MCElePhi_
Definition: MCPizeroAnalyzer.cc:74
MCPizeroAnalyzer::endJob
void endJob() override
Definition: MCPizeroAnalyzer.cc:292
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PizeroMCTruth.h
edm::ParameterSet
Definition: ParameterSet.h:47
a
double a
Definition: hdecay.h:119
Event.h
PVValHelper::phi
Definition: PVValidationHelpers.h:69
TWOPI
#define TWOPI
Definition: DQMSourcePi0.cc:36
MCPizeroAnalyzer::h_MCPizMass2_
TH1F * h_MCPizMass2_
Definition: MCPizeroAnalyzer.cc:70
funct::tan
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
MCPizeroAnalyzer::SimHitLabel
std::string SimHitLabel
Definition: MCPizeroAnalyzer.cc:61
MCPizeroAnalyzer::h_MCPhoEta_
TH1F * h_MCPhoEta_
Definition: MCPizeroAnalyzer.cc:81
MCPizeroAnalyzer::h_MCPizE_
TH1F * h_MCPizE_
Definition: MCPizeroAnalyzer.cc:63
MCPizeroAnalyzer::HepMCLabel
std::string HepMCLabel
Definition: MCPizeroAnalyzer.cc:58
MagneticField.h
edm::EventSetup
Definition: EventSetup.h:58
MCPizeroAnalyzer::SimVtxLabel
std::string SimVtxLabel
Definition: MCPizeroAnalyzer.cc:60
MCPizeroAnalyzer::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: MCPizeroAnalyzer.cc:192
MCPizeroAnalyzer::h_MCPhoPhi_
TH1F * h_MCPhoPhi_
Definition: MCPizeroAnalyzer.cc:82
MCPizeroAnalyzer::nMatched_
int nMatched_
Definition: MCPizeroAnalyzer.cc:52
fileCollector2.converted
converted
Definition: fileCollector2.py:211
PizeroMCTruthFinder::find
std::vector< PizeroMCTruth > find(const std::vector< SimTrack > &simTracks, const std::vector< SimVertex > &simVertices)
Definition: PizeroMCTruthFinder.cc:18
DDAxes::phi
MCPizeroAnalyzer::h_BremEnergy_
TH1F * h_BremEnergy_
Definition: MCPizeroAnalyzer.cc:76
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
MCPizeroAnalyzer::h_EleEvsPhoE_
TH2F * h_EleEvsPhoE_
Definition: MCPizeroAnalyzer.cc:78
std
Definition: JetResolutionObject.h:76
MCPizeroAnalyzer::phiNormalization
float phiNormalization(float &a)
Definition: MCPizeroAnalyzer.cc:176
MCPizeroAnalyzer::MCPizeroAnalyzer
MCPizeroAnalyzer(const edm::ParameterSet &)
Definition: MCPizeroAnalyzer.cc:94
ETA
#define ETA
Definition: GenericBenchmark.cc:28
EventSetup.h
MCPizeroAnalyzer::h_MCConvPhoE_
TH1F * h_MCConvPhoE_
Definition: MCPizeroAnalyzer.cc:83
Exception.h
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
MCPizeroAnalyzer::SimTkLabel
std::string SimTkLabel
Definition: MCPizeroAnalyzer.cc:59
SimTrack.h
MCPizeroAnalyzer::h_MCConvPhoEta_
TH1F * h_MCConvPhoEta_
Definition: MCPizeroAnalyzer.cc:84
MCPizeroAnalyzer::h_MCPiz1ConEta_
TH1F * h_MCPiz1ConEta_
Definition: MCPizeroAnalyzer.cc:66
HepMCProduct.h
MCPizeroAnalyzer::beginJob
void beginJob() override
Definition: MCPizeroAnalyzer.cc:100
ZEcal
static constexpr float ZEcal
Definition: L1TkEmParticleProducer.cc:40
edm::Event
Definition: Event.h:73
MCPizeroAnalyzer::mcEta_
double mcEta_
Definition: MCPizeroAnalyzer.cc:56
SimTrackContainer.h
MCPizeroAnalyzer::h_MCPizUnEta_
TH1F * h_MCPizUnEta_
Definition: MCPizeroAnalyzer.cc:65
Z_Endcap
static constexpr float Z_Endcap
Definition: ECALPositionCalculator.cc:11
SimVertexContainer.h
R_ECAL
static constexpr float R_ECAL
Definition: ECALPositionCalculator.cc:10
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
MCPizeroAnalyzer::h_MCPhoE_
TH1F * h_MCPhoE_
Definition: MCPizeroAnalyzer.cc:80
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
etaBarrelEndcap
static constexpr float etaBarrelEndcap
Definition: ECALPositionCalculator.cc:12