CMS 3D CMS Logo

EcalSimHitsValidation.cc
Go to the documentation of this file.
1 /*
2  * \file EcalSimHitsValidation.cc
3  *
4  * \author C.Rovelli
5  *
6  */
7 
15 
16 using namespace cms;
17 using namespace edm;
18 using namespace std;
19 
21  : g4InfoLabel(ps.getParameter<std::string>("moduleLabelG4")),
22  HepMCToken(consumes<edm::HepMCProduct>(ps.getParameter<std::string>("moduleLabelMC"))) {
24  consumes<edm::PCaloHitContainer>(edm::InputTag(g4InfoLabel, ps.getParameter<std::string>("EBHitsCollection")));
26  consumes<edm::PCaloHitContainer>(edm::InputTag(g4InfoLabel, ps.getParameter<std::string>("EEHitsCollection")));
28  consumes<edm::PCaloHitContainer>(edm::InputTag(g4InfoLabel, ps.getParameter<std::string>("ESHitsCollection")));
29 
30  // DQM ROOT output
31  outputFile_ = ps.getUntrackedParameter<std::string>("outputFile", "");
32 
33  if (!outputFile_.empty()) {
34  edm::LogInfo("OutputInfo") << " Ecal SimHits Task histograms will be saved to " << outputFile_.c_str();
35  } else {
36  edm::LogInfo("OutputInfo") << " Ecal SimHits Task histograms will NOT be saved";
37  }
38 
39  // verbosity switch
40  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
41 
42  // DQMServices
43  dbe_ = nullptr;
44 
45  // get hold of back-end interface
46  dbe_ = edm::Service<DQMStore>().operator->();
47 
48  if (dbe_) {
49  if (verbose_)
50  dbe_->showDirStructure();
51  }
52 
53  meGunEnergy_ = nullptr;
54  meGunEta_ = nullptr;
55  meGunPhi_ = nullptr;
56  meEBEnergyFraction_ = nullptr;
57  meEEEnergyFraction_ = nullptr;
58  meESEnergyFraction_ = nullptr;
59 
60  Char_t histo[200];
61 
62  if (dbe_) {
63  dbe_->setCurrentFolder("EcalHitsV/EcalSimHitsValidation");
64 
65  sprintf(histo, "EcalSimHitsValidation Gun Momentum");
66  meGunEnergy_ = dbe_->book1D(histo, histo, 100, 0., 1000.);
67 
68  sprintf(histo, "EcalSimHitsValidation Gun Eta");
69  meGunEta_ = dbe_->book1D(histo, histo, 700, -3.5, 3.5);
70 
71  sprintf(histo, "EcalSimHitsValidation Gun Phi");
72  meGunPhi_ = dbe_->book1D(histo, histo, 360, 0., 360.);
73 
74  sprintf(histo, "EcalSimHitsValidation Barrel fraction of energy");
75  meEBEnergyFraction_ = dbe_->book1D(histo, histo, 100, 0., 1.1);
76 
77  sprintf(histo, "EcalSimHitsValidation Endcap fraction of energy");
78  meEEEnergyFraction_ = dbe_->book1D(histo, histo, 100, 0., 1.1);
79 
80  sprintf(histo, "EcalSimHitsValidation Preshower fraction of energy");
81  meESEnergyFraction_ = dbe_->book1D(histo, histo, 60, 0., 0.001);
82  }
83 }
84 
86  if (!outputFile_.empty() && dbe_)
88 }
89 
91 
93 
95  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
96 
97  std::vector<PCaloHit> theEBCaloHits;
98  std::vector<PCaloHit> theEECaloHits;
99  std::vector<PCaloHit> theESCaloHits;
100 
105 
106  e.getByToken(HepMCToken, MCEvt);
107  e.getByToken(EBHitsCollectionToken, EcalHitsEB);
108  e.getByToken(EEHitsCollectionToken, EcalHitsEE);
109  e.getByToken(ESHitsCollectionToken, EcalHitsES);
110 
111  for (HepMC::GenEvent::particle_const_iterator p = MCEvt->GetEvent()->particles_begin();
112  p != MCEvt->GetEvent()->particles_end();
113  ++p) {
114  double htheta = (*p)->momentum().theta();
115  double heta = -99999.;
116  if (tan(htheta * 0.5) > 0) {
117  heta = -log(tan(htheta * 0.5));
118  }
119  double hphi = (*p)->momentum().phi();
120  hphi = (hphi >= 0) ? hphi : hphi + 2 * M_PI;
121  hphi = hphi / M_PI * 180.;
122 
123  LogDebug("EventInfo") << "Particle gun type form MC = " << abs((*p)->pdg_id()) << "\n"
124  << "Energy = " << (*p)->momentum().e() << " Eta = " << heta << " Phi = " << hphi;
125 
126  if (meGunEnergy_)
127  meGunEnergy_->Fill((*p)->momentum().e());
128  if (meGunEta_)
129  meGunEta_->Fill(heta);
130  if (meGunPhi_)
131  meGunPhi_->Fill(hphi);
132  }
133 
134  double EBEnergy_ = 0.;
135  if (EcalHitsEB.isValid()) {
136  theEBCaloHits.insert(theEBCaloHits.end(), EcalHitsEB->begin(), EcalHitsEB->end());
137  for (std::vector<PCaloHit>::iterator isim = theEBCaloHits.begin(); isim != theEBCaloHits.end(); ++isim) {
138  EBEnergy_ += isim->energy();
139  }
140  }
141 
142  double EEEnergy_ = 0.;
143  if (EcalHitsEE.isValid()) {
144  theEECaloHits.insert(theEECaloHits.end(), EcalHitsEE->begin(), EcalHitsEE->end());
145  for (std::vector<PCaloHit>::iterator isim = theEECaloHits.begin(); isim != theEECaloHits.end(); ++isim) {
146  EEEnergy_ += isim->energy();
147  }
148  }
149 
150  double ESEnergy_ = 0.;
151  if (EcalHitsES.isValid()) {
152  theESCaloHits.insert(theESCaloHits.end(), EcalHitsES->begin(), EcalHitsES->end());
153  for (std::vector<PCaloHit>::iterator isim = theESCaloHits.begin(); isim != theESCaloHits.end(); ++isim) {
154  ESEnergy_ += isim->energy();
155  }
156  }
157 
158  double etot = EBEnergy_ + EEEnergy_ + ESEnergy_;
159  double fracEB = 0.0;
160  double fracEE = 0.0;
161  double fracES = 0.0;
162 
163  if (etot > 0.0) {
164  fracEB = EBEnergy_ / etot;
165  fracEE = EEEnergy_ / etot;
166  fracES = ESEnergy_ / etot;
167  }
168 
170  meEBEnergyFraction_->Fill(fracEB);
171 
173  meEEEnergyFraction_->Fill(fracEE);
174 
176  meESEnergyFraction_->Fill(fracES);
177 }
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:38
MonitorElement * meEBEnergyFraction_
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:40
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
~EcalSimHitsValidation() override
Destructor.
void Fill(long long x)
edm::EDGetTokenT< edm::PCaloHitContainer > EBHitsCollectionToken
edm::EDGetTokenT< edm::HepMCProduct > HepMCToken
MonitorElement * meGunEnergy_
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< edm::PCaloHitContainer > EEHitsCollectionToken
EcalSimHitsValidation(const edm::ParameterSet &ps)
Constructor.
#define M_PI
Namespace of DDCMS conversion namespace.
edm::EDGetTokenT< edm::PCaloHitContainer > ESHitsCollectionToken
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
MonitorElement * meEEEnergyFraction_
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Analyze.
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
MonitorElement * meESEnergyFraction_
void save(std::string const &filename, std::string const &path="", std::string const &pattern="", std::string const &rewrite="", uint32_t run=0, uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, std::string const &fileupdate="RECREATE")
Definition: DQMStore.cc:2244
void endJob(void) override