CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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")))
23 {
24  EBHitsCollectionToken = consumes<edm::PCaloHitContainer>(edm::InputTag(g4InfoLabel, ps.getParameter<std::string>("EBHitsCollection")));
25  EEHitsCollectionToken = consumes<edm::PCaloHitContainer>(edm::InputTag(g4InfoLabel, ps.getParameter<std::string>("EEHitsCollection")));
26  ESHitsCollectionToken = consumes<edm::PCaloHitContainer>(edm::InputTag(g4InfoLabel, ps.getParameter<std::string>("ESHitsCollection")));
27 
28  // DQM ROOT output
29  outputFile_ = ps.getUntrackedParameter<std::string>("outputFile", "");
30 
31  if ( outputFile_.size() != 0 ) {
32  edm::LogInfo("OutputInfo") << " Ecal SimHits Task histograms will be saved to " << outputFile_.c_str();
33  } else {
34  edm::LogInfo("OutputInfo") << " Ecal SimHits Task histograms will NOT be saved";
35  }
36 
37  // verbosity switch
38  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
39 
40  // DQMServices
41  dbe_ = 0;
42 
43  // get hold of back-end interface
44  dbe_ = edm::Service<DQMStore>().operator->();
45  if ( dbe_ ) {
46  if ( verbose_ ) { dbe_->setVerbose(1); }
47  else { dbe_->setVerbose(0); }
48  }
49 
50  if ( dbe_ ) {
51  if ( verbose_ ) dbe_->showDirStructure();
52  }
53 
54  meGunEnergy_ = 0;
55  meGunEta_ = 0;
56  meGunPhi_ = 0;
60 
61  Char_t histo[200];
62 
63 
64  if ( dbe_ ) {
65  dbe_->setCurrentFolder("EcalHitsV/EcalSimHitsValidation");
66 
67  sprintf (histo, "EcalSimHitsValidation Gun Momentum" ) ;
68  meGunEnergy_ = dbe_->book1D(histo, histo, 100, 0., 1000.);
69 
70  sprintf (histo, "EcalSimHitsValidation Gun Eta" ) ;
71  meGunEta_ = dbe_->book1D(histo, histo, 700, -3.5, 3.5);
72 
73  sprintf (histo, "EcalSimHitsValidation Gun Phi" ) ;
74  meGunPhi_ = dbe_->book1D(histo, histo, 360, 0., 360.);
75 
76  sprintf (histo, "EcalSimHitsValidation Barrel fraction of energy" ) ;
77  meEBEnergyFraction_ = dbe_->book1D(histo, histo, 100 , 0. , 1.1);
78 
79  sprintf (histo, "EcalSimHitsValidation Endcap fraction of energy" ) ;
80  meEEEnergyFraction_ = dbe_->book1D(histo, histo, 100 , 0. , 1.1);
81 
82  sprintf (histo, "EcalSimHitsValidation Preshower fraction of energy" ) ;
83  meESEnergyFraction_ = dbe_->book1D(histo, histo, 60 , 0. , 0.001);
84  }
85 
86 }
87 
89 
90  if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
91 
92 }
93 
95 
96 }
97 
99 
100 }
101 
103 
104  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
105 
106  std::vector<PCaloHit> theEBCaloHits;
107  std::vector<PCaloHit> theEECaloHits;
108  std::vector<PCaloHit> theESCaloHits;
109 
114 
115  e.getByToken(HepMCToken, MCEvt);
116  e.getByToken(EBHitsCollectionToken, EcalHitsEB);
117  e.getByToken(EEHitsCollectionToken, EcalHitsEE);
118  e.getByToken(ESHitsCollectionToken, EcalHitsES);
119 
120  for ( HepMC::GenEvent::particle_const_iterator p = MCEvt->GetEvent()->particles_begin();
121  p != MCEvt->GetEvent()->particles_end(); ++p ) {
122 
123  double htheta = (*p)->momentum().theta();
124  double heta = -99999.;
125  if( tan(htheta * 0.5) > 0 ) {
126  heta = -log(tan(htheta * 0.5));
127  }
128  double hphi = (*p)->momentum().phi();
129  hphi = (hphi>=0) ? hphi : hphi+2*M_PI;
130  hphi = hphi / M_PI * 180.;
131 
132  LogDebug("EventInfo") << "Particle gun type form MC = " << abs((*p)->pdg_id()) << "\n" << "Energy = "<< (*p)->momentum().e() << " Eta = " << heta << " Phi = " << hphi;
133 
134  if (meGunEnergy_) meGunEnergy_->Fill((*p)->momentum().e());
135  if (meGunEta_) meGunEta_ ->Fill(heta);
136  if (meGunPhi_) meGunPhi_ ->Fill(hphi);
137 
138  }
139 
140  double EBEnergy_ = 0.;
141  if ( EcalHitsEB.isValid() ) {
142  theEBCaloHits.insert(theEBCaloHits.end(), EcalHitsEB->begin(), EcalHitsEB->end());
143  for (std::vector<PCaloHit>::iterator isim = theEBCaloHits.begin();
144  isim != theEBCaloHits.end(); ++isim){
145  EBEnergy_ += isim->energy();
146  }
147  }
148 
149  double EEEnergy_ = 0.;
150  if ( EcalHitsEE.isValid() ) {
151  theEECaloHits.insert(theEECaloHits.end(), EcalHitsEE->begin(), EcalHitsEE->end());
152  for (std::vector<PCaloHit>::iterator isim = theEECaloHits.begin();
153  isim != theEECaloHits.end(); ++isim){
154  EEEnergy_ += isim->energy();
155  }
156  }
157 
158  double ESEnergy_ = 0.;
159  if ( EcalHitsES.isValid() ) {
160  theESCaloHits.insert(theESCaloHits.end(), EcalHitsES->begin(), EcalHitsES->end());
161  for (std::vector<PCaloHit>::iterator isim = theESCaloHits.begin();
162  isim != theESCaloHits.end(); ++isim){
163  ESEnergy_ += isim->energy();
164  }
165  }
166 
167  double etot = EBEnergy_ + EEEnergy_ + ESEnergy_ ;
168  double fracEB = 0.0;
169  double fracEE = 0.0;
170  double fracES = 0.0;
171 
172  if (etot>0.0) {
173  fracEB = EBEnergy_/etot;
174  fracEE = EEEnergy_/etot;
175  fracES = ESEnergy_/etot;
176  }
177 
179 
181 
183 }
184 
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:42
MonitorElement * meEBEnergyFraction_
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
T getUntrackedParameter(std::string const &, T const &) const
static std::vector< std::string > checklist log
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void analyze(const edm::Event &e, const edm::EventSetup &c)
Analyze.
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:76
edm::EDGetTokenT< edm::PCaloHitContainer > EEHitsCollectionToken
EcalSimHitsValidation(const edm::ParameterSet &ps)
Constructor.
#define M_PI
edm::EDGetTokenT< edm::PCaloHitContainer > ESHitsCollectionToken
MonitorElement * meEEEnergyFraction_
edm::EventID id() const
Definition: EventBase.h:56
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, const uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE", const bool resetMEsAfterWriting=false)
Definition: DQMStore.cc:2540
MonitorElement * meESEnergyFraction_
~EcalSimHitsValidation()
Destructor.