CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalPreshowerSimHitsValidation.cc
Go to the documentation of this file.
1 /*
2  * \file EcalPreshowerSimHitsValidation.cc
3  *
4  * \author C.Rovelli
5  *
6 */
7 
13 
14 using namespace cms;
15 using namespace edm;
16 using namespace std;
17 
19 
20  HepMCLabel(ps.getParameter<std::string>("moduleLabelMC")),
21  g4InfoLabel(ps.getParameter<std::string>("moduleLabelG4")),
22  EEHitsCollection(ps.getParameter<std::string>("EEHitsCollection")),
23  ESHitsCollection(ps.getParameter<std::string>("ESHitsCollection")){
24 
25  HepMCToken = consumes <edm::HepMCProduct> (HepMCLabel);
26  EEHitsToken = consumes <edm::PCaloHitContainer> (edm::InputTag(std::string(g4InfoLabel),std::string(EEHitsCollection)));
27  ESHitsToken = consumes <edm::PCaloHitContainer> (edm::InputTag(std::string(g4InfoLabel),std::string(ESHitsCollection)));
28  // verbosity switch
29  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
30 
31  // get hold of back-end interface
32  dbe_ = 0;
33  dbe_ = edm::Service<DQMStore>().operator->();
34  if ( dbe_ ) {
35  if ( verbose_ ) { dbe_->setVerbose(1); }
36  else { dbe_->setVerbose(0); }
37  }
38 
39  if ( dbe_ ) {
40  if ( verbose_ ) dbe_->showDirStructure();
41  }
42 
43 
44  menESHits1zp_ = 0;
45  menESHits2zp_ = 0;
46  menESHits1zm_ = 0;
47  menESHits2zm_ = 0;
48 
53 
56 
57  meE1alphaE2zp_ = 0;
58  meE1alphaE2zm_ = 0;
59  meEEoverESzp_ = 0;
60  meEEoverESzm_ = 0;
61 
64 
65 
66  Char_t histo[200];
67 
68  if ( dbe_ ) {
69  dbe_->setCurrentFolder("EcalHitsV/EcalSimHitsValidation");
70 
71  sprintf (histo, "ES hits layer 1 multiplicity z+" ) ;
72  menESHits1zp_ = dbe_->book1D(histo, histo, 50, 0., 50. ) ;
73 
74  sprintf (histo, "ES hits layer 2 multiplicity z+" ) ;
75  menESHits2zp_ = dbe_->book1D(histo, histo, 50, 0., 50. ) ;
76 
77  sprintf (histo, "ES hits layer 1 multiplicity z-" ) ;
78  menESHits1zm_ = dbe_->book1D(histo, histo, 50, 0., 50. ) ;
79 
80  sprintf (histo, "ES hits layer 2 multiplicity z-" ) ;
81  menESHits2zm_ = dbe_->book1D(histo, histo, 50, 0., 50. ) ;
82 
83  sprintf (histo, "ES hits energy layer 1 z+" ) ;
84  meESEnergyHits1zp_ = dbe_->book1D(histo, histo, 100, 0., 0.05 ) ;
85 
86  sprintf (histo, "ES hits energy layer 2 z+" ) ;
87  meESEnergyHits2zp_ = dbe_->book1D(histo, histo, 100, 0., 0.05 ) ;
88 
89  sprintf (histo, "ES hits energy layer 1 z-" ) ;
90  meESEnergyHits1zm_ = dbe_->book1D(histo, histo, 100, 0., 0.05 ) ;
91 
92  sprintf (histo, "ES hits energy layer 2 z-" ) ;
93  meESEnergyHits2zm_ = dbe_->book1D(histo, histo, 100, 0., 0.05 ) ;
94 
95  sprintf (histo, "ES hits log10energy spectrum" );
96  meEShitLog10Energy_ = dbe_->book1D(histo, histo, 140, -10., 4.);
97 
98  sprintf (histo, "ES hits log10energy spectrum vs normalized energy" );
99  meEShitLog10EnergyNorm_ = dbe_->bookProfile(histo, histo, 140, -10., 4., 100, 0., 1.);
100 
101  sprintf (histo, "ES E1+07E2 z+" ) ;
102  meE1alphaE2zp_ = dbe_->book1D(histo, histo, 100, 0., 0.05);
103 
104  sprintf (histo, "ES E1+07E2 z-" ) ;
105  meE1alphaE2zm_ = dbe_->book1D(histo, histo, 100, 0., 0.05);
106 
107  sprintf (histo, "EE vs ES z+" ) ;
108  meEEoverESzp_ = dbe_->bookProfile(histo, histo, 250, 0., 500., 200, 0., 200.);
109 
110  sprintf (histo, "EE vs ES z-" ) ;
111  meEEoverESzm_ = dbe_->bookProfile(histo, histo, 250, 0., 500., 200, 0., 200.);
112 
113  sprintf (histo, "ES ene2oEne1 z+" ) ;
114  me2eszpOver1eszp_ = dbe_->book1D(histo, histo, 50, 0., 10.);
115 
116  sprintf (histo, "ES ene2oEne1 z-" ) ;
117  me2eszmOver1eszm_ = dbe_->book1D(histo, histo, 50, 0., 10.);
118  }
119 }
120 
122 
123 }
124 
126 
127 }
128 
130 
131 }
132 
134 
135  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
136 
138  e.getByToken(HepMCToken, MCEvt);
139 
141  e.getByToken(EEHitsToken,EcalHitsEE);
142 
144  e.getByToken(ESHitsToken,EcalHitsES);
145 
146  std::vector<PCaloHit> theEECaloHits;
147  if( EcalHitsEE.isValid() ) {
148  theEECaloHits.insert(theEECaloHits.end(), EcalHitsEE->begin(), EcalHitsEE->end());
149  }
150 
151  std::vector<PCaloHit> theESCaloHits;
152  if( EcalHitsES.isValid() ) {
153  theESCaloHits.insert(theESCaloHits.end(), EcalHitsES->begin(), EcalHitsES->end());
154  }
155 
156  double ESEnergy_ = 0.;
157  //std::map<unsigned int, std::vector<PCaloHit>,std::less<unsigned int> > CaloHitMap;
158 
159  // endcap
160  double EEetzp_ = 0.;
161  double EEetzm_ = 0.;
162  for (std::vector<PCaloHit>::iterator isim = theEECaloHits.begin(); isim != theEECaloHits.end(); ++isim){
163  EEDetId eeid (isim->id()) ;
164  if (eeid.zside() > 0 ) EEetzp_ += isim->energy();
165  if (eeid.zside() < 0 ) EEetzm_ += isim->energy();
166  }
167 
168 
169  uint32_t nESHits1zp = 0;
170  uint32_t nESHits1zm = 0;
171  uint32_t nESHits2zp = 0;
172  uint32_t nESHits2zm = 0;
173  double ESet1zp_ = 0.;
174  double ESet2zp_ = 0.;
175  double ESet1zm_ = 0.;
176  double ESet2zm_ = 0.;
177  std::vector<double> econtr(140, 0. );
178 
179  for (std::vector<PCaloHit>::iterator isim = theESCaloHits.begin();
180  isim != theESCaloHits.end(); ++isim){
181  //CaloHitMap[(*isim).id()].push_back((*isim));
182 
183  ESDetId esid (isim->id()) ;
184 
185  LogDebug("HitInfo")
186  << " CaloHit " << isim->getName() << "\n"
187  << " DetID = "<<isim->id()<< " ESDetId: z side " << esid.zside() << " plane " << esid.plane() << esid.six() << ',' << esid.siy() << ':' << esid.strip() << "\n"
188  << " Time = " << isim->time() << "\n"
189  << " Track Id = " << isim->geantTrackId() << "\n"
190  << " Energy = " << isim->energy();
191 
192  ESEnergy_ += isim->energy();
193  if( isim->energy() > 0 ) {
194  meEShitLog10Energy_->Fill(log10(isim->energy()));
195  int log10i = int( ( log10(isim->energy()) + 10. ) * 10. );
196  if( log10i >=0 && log10i < 140 ) econtr[log10i] += isim->energy();
197  }
198 
199  if (esid.plane() == 1 ) {
200  if (esid.zside() > 0 ) {
201  nESHits1zp++ ;
202  ESet1zp_ += isim->energy();
203  if (meESEnergyHits1zp_) meESEnergyHits1zp_->Fill(isim->energy()) ;
204  }
205  else if (esid.zside() < 0 ) {
206  nESHits1zm++ ;
207  ESet1zm_ += isim->energy();
208  if (meESEnergyHits1zm_) meESEnergyHits1zm_->Fill(isim->energy()) ;
209  }
210  }
211  else if (esid.plane() == 2 ) {
212  if (esid.zside() > 0 ) {
213  nESHits2zp++ ;
214  ESet2zp_ += isim->energy();
215  if (meESEnergyHits2zp_) meESEnergyHits2zp_->Fill(isim->energy()) ;
216  }
217  else if (esid.zside() < 0 ) {
218  nESHits2zm++ ;
219  ESet2zm_ += isim->energy();
220  if (meESEnergyHits2zm_) meESEnergyHits2zm_->Fill(isim->energy()) ;
221  }
222  }
223 
224  }
225 
226  if (menESHits1zp_) menESHits1zp_->Fill(nESHits1zp);
227  if (menESHits1zm_) menESHits1zm_->Fill(nESHits1zm);
228 
229  if (menESHits2zp_) menESHits2zp_->Fill(nESHits2zp);
230  if (menESHits2zm_) menESHits2zm_->Fill(nESHits2zm);
231 
232  if( meEShitLog10EnergyNorm_ && ESEnergy_ != 0 ) {
233  for( int i=0; i<140; i++ ) {
234  meEShitLog10EnergyNorm_->Fill( -10.+(float(i)+0.5)/10., econtr[i]/ESEnergy_ );
235  }
236  }
237 
238 
239  for ( HepMC::GenEvent::particle_const_iterator p = MCEvt->GetEvent()->particles_begin();
240  p != MCEvt->GetEvent()->particles_end(); ++p ) {
241 
242  double htheta = (*p)->momentum().theta();
243  double heta = -99999.;
244  if( tan(htheta * 0.5) > 0 ) {
245  heta = -log(tan(htheta * 0.5));
246  }
247 
248  if ( heta > 1.653 && heta < 2.6 ) {
249 
250  if (meE1alphaE2zp_) meE1alphaE2zp_->Fill(ESet1zp_+0.7*ESet2zp_);
251  if (meEEoverESzp_) meEEoverESzp_ ->Fill((ESet1zp_+0.7*ESet2zp_)/0.00009, EEetzp_);
252  if ((me2eszpOver1eszp_) && (ESet1zp_ != 0.)) me2eszpOver1eszp_->Fill(ESet2zp_/ESet1zp_);
253  }
254  if ( heta < -1.653 && heta > -2.6 ) {
255  if (meE1alphaE2zm_) meE1alphaE2zm_->Fill(ESet1zm_+0.7*ESet2zm_);
256  if (meEEoverESzm_) meEEoverESzm_ ->Fill((ESet1zm_+0.7*ESet2zm_)/0.00009, EEetzm_);
257  if ((me2eszmOver1eszm_) && (ESet1zm_ != 0.)) me2eszmOver1eszm_->Fill(ESet2zm_/ESet1zm_);
258  }
259  }
260 }
261 
262 
263 // LocalWords: EcalSimHitsValidation
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
EcalPreshowerSimHitsValidation(const edm::ParameterSet &ps)
Constructor.
void analyze(const edm::Event &e, const edm::EventSetup &c)
Analyze.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
void Fill(long long x)
edm::EDGetTokenT< edm::HepMCProduct > HepMCToken
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
bool isValid() const
Definition: HandleBase.h:75
edm::EventID id() const
Definition: EventBase.h:60
edm::EDGetTokenT< edm::PCaloHitContainer > EEHitsToken
edm::EDGetTokenT< edm::PCaloHitContainer > ESHitsToken
tuple log
Definition: cmsBatch.py:347