CMS 3D CMS Logo

PixelResolutionHistograms.cc
Go to the documentation of this file.
1 //
2 // The implementation of the PixelResolutinHistograms.cc class. Please
3 // look at PixelResolutionHistograms.h header file for the interface.
4 //
5 //------------------------------------------------------------------------------
6 
7 // The switch, undefined in CMSSW release, and defined by standalone compilation script:
8 
9 
10 #ifdef SI_PIXEL_TEMPLATE_STANDALONE
11 //
12 //--- Stand-alone: Include a the header file from the local directory, as well as
13 // dummy implementations of SimpleHistogramGenerator, LogInfo, LogError and LogDebug...
14 //
16 //
17 class TH1F;
18 class TH2F;
20 public:
21  SimpleHistogramGenerator(TH1F * hist) : hist_(hist) {};
22 private:
23  TH1F * hist_; // we don't own it
24 };
25 #define LOGDEBUG std::cout
26 #define LOGERROR std::cout
27 #define LOGINFO std::cout
28 //
29 #else
30 //--- We're inside a CMSSW release: Include the real thing.
31 //
36 //
37 #define LOGDEBUG LogDebug("")
38 #define LOGERROR edm::LogError("Error")
39 #define LOGINFO edm::LogInfo("Info")
40 //
41 #endif
42 
43 
44 
45 
46 // Generic C stuff
47 #include <cmath>
48 #include <iostream>
49 #include <string>
50 
51 // ROOT
52 #include <TFile.h>
53 #include <TH1F.h>
54 #include <TH2F.h>
55 
56 // Global definitions
57 const float cmtomicron = 10000.0;
58 
59 
60 //------------------------------------------------------------------------------
61 // Constructor: Books the FastSim Histograms, given the input parameters
62 // which are provided as arguments. These variables are then const inside
63 // the class. (That is, once we make the histograms, we can't change the
64 // definition of the binning.)
65 //------------------------------------------------------------------------------
67 PixelResolutionHistograms( std::string filename, // ROOT file for histograms
68  std::string rootdir, // Subdirectory in the file, "" if none
69  std::string descTitle, // Descriptive title
70  unsigned int detType, // Where we are... (&&& do we need this?)
71  double cotbetaBinWidth,
72  double cotbetaLowEdge,
73  int cotbetaBins,
74  double cotalphaBinWidth,
75  double cotalphaLowEdge,
76  int cotalphaBins) //,
77  //int qbinWidth,
78  //int qbins )
79  : weOwnHistograms_(true), // we'll be making some histos
80  detType_ ( detType ),
81  cotbetaBinWidth_ ( cotbetaBinWidth ),
82  cotbetaLowEdge_ ( cotbetaLowEdge ),
83  cotbetaBins_ ( cotbetaBins ),
84  cotalphaBinWidth_ ( cotalphaBinWidth ),
85  cotalphaLowEdge_ ( cotalphaLowEdge ),
86  cotalphaBins_ ( cotalphaBins ),
87  qbinWidth_ ( 1 ),
88  qbins_ ( 4 ),
89  binningHisto_(nullptr),
90  resMultiPixelXHist_(), resSinglePixelXHist_(), // all to nullptr
91  resMultiPixelYHist_(), resSinglePixelYHist_(), // all to nullptr
92  qbinHist_(), // all to nullptr
93  file_(nullptr),
94  status_(0),
95  resMultiPixelXGen_(), resSinglePixelXGen_(),
96  resMultiPixelYGen_(), resSinglePixelYGen_(),
97  qbinGen_()
98 {
99 
100  file_ = std::make_unique<TFile>( filename.c_str(), "RECREATE" );
101  //Resolution binning
102  // const double cotbetaBinWidth = 1.0;
103  // const double cotbetaLowEdge = -11.5 ;
104  // const int cotbetaBins = 23;
105  // const double cotalphaBinWidth = 0.08 ;
106  // const double cotalphaLowEdge = -0.36 ;
107  // const int cotalphaBins = 9;
108  // const int qbinWidth = 1;
109  // const int qbins = 4;
110 
111  // Dummy 2D histogram to store binning:
112  binningHisto_ = new TH2F( "ResHistoBinning", descTitle.c_str(),
113  cotbetaBins, cotbetaLowEdge, cotbetaLowEdge + cotbetaBins*cotbetaBinWidth,
114  cotalphaBins, cotalphaLowEdge, cotalphaLowEdge + cotalphaBins*cotalphaBinWidth );
115 
116  // Store detType in the underflow bin
117  binningHisto_->SetBinContent(0, 0, detType_);
118 
119  // All other histograms:
120  Char_t histo[200];
121  Char_t title[200];
122  //
123  //--- Histograms for clusters with multiple pixels hit in a given direction.
124  //
125  for( int ii=0; ii<cotbetaBins_; ii++ ) {
126  for( int jj=0; jj<cotalphaBins_; jj++ ) {
127  for( int kk=0; kk<qbins_; kk++ ) {
128  //
129  sprintf( histo, "hx%d1%02d%d%d", detType_, ii+1, jj+1, kk+1 ); //information of bits of histogram names
130  //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha, qbins
131  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 X",
132  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
133  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_,
134  kk+1 );
135  //
136  resMultiPixelXHist_ [ ii ][ jj ][ kk ] = new TH1F(histo, title, 1000, -0.05, 0.05);
137 
138  sprintf( histo, "hy%d1%02d%d%d", detType_, ii+1, jj+1, kk+1 );
139  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 Y",
140  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
141  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_,
142  kk+1 );
143  //
144  resMultiPixelYHist_ [ ii ][ jj ][ kk ] = new TH1F(histo, title, 1000, -0.05, 0.05);
145  }
146  }
147  }
148 
149 
150  //
151  //--- Histograms for clusters where only a single pixel was hit in a given direction.
152  //
153  for( int ii=0; ii<cotbetaBins_; ii++) {
154  for( int jj=0; jj<cotalphaBins_; jj++) {
155 
156  sprintf( histo, "hx%d0%02d%d", detType_, ii+1, jj+1 ); //information of bits of histogram names
157  //first bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha
158  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 X",
159  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
160  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_ );
161  //
162  resSinglePixelXHist_ [ ii ][ jj ] = new TH1F(histo, title, 1000, -0.05, 0.05);
163 
164  sprintf( histo, "hy%d0%02d%d", detType_, ii+1, jj+1 );
165  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 Y",
166  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
167  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_ );
168  //
169  resSinglePixelYHist_ [ ii ][ jj ] = new TH1F(histo, title, 1000, -0.05, 0.05);
170 
171  sprintf( histo, "hqbin%d%02d%d", detType_, ii+1, jj+1 );
172  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin",
173  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
174  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_ );
175  //
176  qbinHist_ [ ii ][ jj ] = new TH1F(histo, title, 4, -0.49, 3.51);
177 
178  }
179  }
180 
181 }
182 
183 
184 
185 //------------------------------------------------------------------------------
186 // Another constructor: load the histograms from one file.
187 // filename = full path to filename
188 // rootdir = ROOT directory inside the file
189 //
190 // The other parameters are the same (needed later) and must correspond
191 // to the histograms we are loading from the file.
192 //------------------------------------------------------------------------------
195  std::string rootdir,
196  int detType,
197  bool ignore_multi,
198  bool ignore_single,
199  bool ignore_qBin )
200  : weOwnHistograms_(false), // resolution histograms are owned by the ROOT file
201  detType_ (-1),
202  cotbetaBinWidth_ (0),
203  cotbetaLowEdge_ (0),
204  cotbetaBins_ (0),
205  cotalphaBinWidth_ (0),
206  cotalphaLowEdge_ (0),
207  cotalphaBins_ (0),
208  qbinWidth_ (1),
209  qbins_ (4),
211  resMultiPixelXHist_(), resSinglePixelXHist_(), // all to nullptr
212  resMultiPixelYHist_(), resSinglePixelYHist_(), // all to nullptr
213  qbinHist_(), // all to nullptr
214  file_(nullptr),
215  status_(0),
218  qbinGen_()
219 {
220  Char_t histo[200]; // the name of the histogram
221  Char_t title[200]; // histo title, for debugging and sanity checking (compare inside file)
222  TH1F * tmphist = nullptr; // cache for histo pointer
223 
224  //--- Open the file for reading.
225  file_ = std::make_unique<TFile>( filename.c_str() ,"READ");
226  if ( !file_ ) {
227  status_ = 1;
228  LOGERROR << "PixelResolutionHistograms:: Error, file " << filename << " not found.";
229  return; // PixelTemplateSmearerBase will throw an exception upon our return.
230  }
231 
232  //--- The dummy 2D histogram with the binning of cot\beta and cot\alpha:
233  binningHisto_ = (TH2F*) file_->Get( Form( "%s%s" , rootdir.c_str(), "ResHistoBinning" ) );
234  if ( !binningHisto_ ) {
235  status_ = 11;
236  LOGERROR << "PixelResolutionHistograms:: Error, binning histogrram ResHistoBinning not found.";
237  return; // PixelTemplateSmearerBase will throw an exception upon our return.
238  }
239 
240  if ( detType == -1 ) {
241  //--- Fish out detType from the underflow bin:
242  detType_ = binningHisto_->GetBinContent(0, 0);
243  }
244  else {
245  detType_ = detType; // constructor's argument overrides what's in ResHistoBinning histogram.
246  }
247 
248  //--- Now we fill the binning variables:
249  cotbetaAxis_ = binningHisto_->GetXaxis();
250  cotbetaBinWidth_ = binningHisto_->GetXaxis()->GetBinWidth(1); // assume all same width
251  cotbetaLowEdge_ = binningHisto_->GetXaxis()->GetXmin(); // low edge of the first bin
252  cotbetaBins_ = binningHisto_->GetXaxis()->GetNbins();
253  cotalphaAxis_ = binningHisto_->GetYaxis();
254  cotalphaBinWidth_ = binningHisto_->GetYaxis()->GetBinWidth(1); // assume all same width;
255  cotalphaLowEdge_ = binningHisto_->GetYaxis()->GetXmin(); // low edge of the first bin;
256  cotalphaBins_ = binningHisto_->GetYaxis()->GetNbins();
257 
258 
259  if ( !ignore_multi ) {
260  //
261  //--- Histograms for clusters with multiple pixels hit in a given direction.
262  //
263  for( int ii=0; ii<cotbetaBins_; ii++ ) {
264  for( int jj=0; jj<cotalphaBins_; jj++ ) {
265  for( int kk=0; kk<qbins_; kk++ ) {
266  //
267  sprintf( histo, "hx%d1%02d%d%d", detType_, ii+1, jj+1, kk+1 ); //information of bits of histogram names
268  //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha, qbins
269  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 X",
270  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
271  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_,
272  kk+1 );
273  //
274  tmphist = (TH1F*) file_->Get( Form( "%s%s" , rootdir.c_str(), histo ) );
275  if ( !tmphist ) {
276  status_ = 2;
277  LOGERROR << "Failed to find histogram=" << std::string( histo );
278  return;
279  }
280  LOGDEBUG << "Found histo " << std::string(histo)
281  << " with title = " << std::string( tmphist->GetTitle() ) << std::endl;
282  if ( tmphist->GetEntries() < 5 ) {
283  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
284  << " entries. Trouble ahead." << std::endl;
285  }
286  resMultiPixelXHist_ [ ii ][ jj ][ kk ] = tmphist;
287  resMultiPixelXGen_ [ ii ][ jj ][ kk ] = new SimpleHistogramGenerator( tmphist );
288 
289 
290  sprintf( histo, "hy%d1%02d%d%d", detType_, ii+1, jj+1, kk+1 );
291  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 Y",
292  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
293  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_,
294  kk+1 );
295  //
296  tmphist = (TH1F*) file_->Get( Form( "%s%s" , rootdir.c_str(), histo ) );
297  if ( !tmphist ) {
298  status_ = 3;
299  LOGERROR << "Failed to find histogram=" << std::string( histo );
300  return;
301  }
302  LOGDEBUG << "Found histo " << std::string(histo)
303  << " with title = " << std::string( tmphist->GetTitle() ) << std::endl;
304  if ( tmphist->GetEntries() < 5 ) {
305  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
306  << " entries. Trouble ahead." << std::endl;
307  }
308  resMultiPixelYHist_ [ ii ][ jj ][ kk ] = tmphist;
309  resMultiPixelYGen_ [ ii ][ jj ][ kk ] = new SimpleHistogramGenerator( tmphist );
310  }
311  }
312  }
313  //
314  } // if (not ignore multi)
315 
316 
317  //
318  //--- Histograms for clusters where only a single pixel was hit in a given direction.
319  //
320  for( int ii=0; ii<cotbetaBins_; ii++) {
321  for( int jj=0; jj<cotalphaBins_; jj++) {
322 
323  //--- Single pixel, along X.
324  //
325  sprintf( histo, "hx%d0%02d%d", detType_, ii+1, jj+1 ); //information of bits of histogram names
326  //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha
327  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 X",
328  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
329  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_ );
330  //
331  tmphist = (TH1F*) file_->Get( Form( "%s%s" , rootdir.c_str(), histo ) );
332  if ( !tmphist ) {
333  if ( !ignore_single ) {
334  LOGERROR << "Failed to find histogram=" << std::string( histo );
335  status_ = 4;
336  return;
337  }
338  }
339  else {
340  LOGDEBUG << "Found histo " << std::string(histo)
341  << " with title = " << std::string( tmphist->GetTitle() ) << std::endl;
342  LOGDEBUG << "Found histo with title = " << std::string( tmphist->GetTitle() ) << std::endl;
343  if ( tmphist->GetEntries() < 5 ) {
344  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
345  << " entries. Trouble ahead." << std::endl;
346  }
347  resSinglePixelXHist_ [ ii ][ jj ] = tmphist;
348  resSinglePixelXGen_ [ ii ][ jj ] = new SimpleHistogramGenerator( tmphist );
349  }
350 
351 
352  //--- Single pixel, along Y.
353  //
354  sprintf( histo, "hy%d0%02d%d", detType_, ii+1, jj+1 );
355  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 Y",
356  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
357  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_ );
358  //
359  tmphist = (TH1F*) file_->Get( Form( "%s%s" , rootdir.c_str(), histo ) );
360  if ( !tmphist ) {
361  if ( !ignore_single ) {
362  LOGERROR << "Failed to find histogram=" << std::string( histo );
363  status_ = 5;
364  return;
365  }
366  }
367  else {
368  LOGDEBUG << "Found histo " << std::string(histo)
369  << " with title = " << std::string( tmphist->GetTitle() ) << std::endl;
370  if ( tmphist->GetEntries() < 5 ) {
371  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
372  << " entries. Trouble ahead." << std::endl;
373  }
374  resSinglePixelYHist_ [ ii ][ jj ] = tmphist;
375  resSinglePixelYGen_ [ ii ][ jj ] = new SimpleHistogramGenerator( tmphist );
376  }
377 
378 
379  //--- qBin distribution, for this (cotbeta, cotalpha) bin.
380  //
381  sprintf( histo, "hqbin%d%02d%d", detType_, ii+1, jj+1 );
382  sprintf( title, "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin",
383  cotbetaLowEdge_ + ii*cotbetaBinWidth_ , cotbetaLowEdge_ + (ii+1)*cotbetaBinWidth_,
384  cotalphaLowEdge_ +jj*cotalphaBinWidth_, cotalphaLowEdge_ +(jj+1)*cotalphaBinWidth_ );
385  //
386  tmphist = (TH1F*) file_->Get( Form( "%s%s" , rootdir.c_str(), histo ) );
387  if ( !tmphist ) {
388  if ( !ignore_qBin ) {
389  LOGERROR << "Failed to find histogram=" << std::string( histo );
390  status_ = 6;
391  return;
392  }
393  }
394  else {
395  LOGDEBUG << "Found histo " << std::string(histo)
396  << " with title = " << std::string( tmphist->GetTitle() ) << std::endl;
397  if ( tmphist->GetEntries() < 5 ) {
398  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
399  << " entries. Trouble ahead." << std::endl;
400  }
401  qbinHist_ [ ii ][ jj ] = tmphist;
402  qbinGen_ [ ii ][ jj ] = new SimpleHistogramGenerator( tmphist );
403  }
404 
405 
406  }
407  }
408 }
409 
410 
411 
412 
413 //------------------------------------------------------------------------------
414 // Destructor. Use file_ pointer to tell whether we loaded the histograms
415 // from a file (and do not own them), or we built them ourselves and thus need
416 // to delete them.
417 //------------------------------------------------------------------------------
419 {
420  //--- Delete histograms, but only if we own them. If
421  //--- they came from a file, let them be.
422  //
423  if ( ! weOwnHistograms_ ) {
424  //--- Read the histograms from the TFile, the file will take care of them.
425  file_->Close();
428  }
429  else {
430  //--- We made the histograms, so first write them inthe output ROOT file and close it.
431  LOGINFO
432  << "PixelResHistoStore: Writing the histograms to the output file. " // << filename
433  << std::endl;
434  file_->Write();
435  file_->Close();
436 
437  // ROOT file has the ownership, and once the file is closed,
438  // all of these guys are deleted. So, we don't need to do anything.
439  } // else
440 
441 
442 
443  //--- Delete FastSim generators. (It's safe to delete a nullptr.)
444  for( int ii=0; ii<cotbetaBins_; ii++ ) {
445  for( int jj=0; jj<cotalphaBins_; jj++ ) {
446  for( int kk=0; kk<qbins_; kk++ ) {
447  delete resMultiPixelXGen_ [ ii ][ jj ][ kk ];
448  delete resMultiPixelYGen_ [ ii ][ jj ][ kk ];
449  }
450  }
451  }
452  for( int ii=0; ii<cotbetaBins_; ii++ ) {
453  for( int jj=0; jj<cotalphaBins_; jj++ ) {
454  delete resSinglePixelXGen_ [ ii ][ jj ];
455  delete resSinglePixelYGen_ [ ii ][ jj ];
456  delete qbinGen_ [ ii ][ jj ];
457  }
458  }
459 
460 }
461 
462 
463 
464 
465 
466 
467 //------------------------------------------------------------------------------
468 // Fills the appropriate FastSim histograms.
469 // Returns 0 if the relevant histogram(s) were found and filled, 1 if not.
470 //------------------------------------------------------------------------------
471 int
473 Fill( double dx, double dy, double cotalpha, double cotbeta,
474  int qbin, int nxpix, int nypix )
475 {
476  int icotalpha, icotbeta, iqbin ;
477  icotalpha = (int)floor( (cotalpha - cotalphaLowEdge_)/cotalphaBinWidth_ ) ;
478  icotbeta = (int)floor( (cotbeta - cotbetaLowEdge_) /cotbetaBinWidth_ ) ;
479  iqbin = qbin > 2 ? 3 : qbin;
480  if( icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_ ) {
481  qbinHist_[icotbeta][icotalpha]->Fill((double)iqbin);
482  if( nxpix == 1 )
483  resSinglePixelXHist_[icotbeta][icotalpha]->Fill(dx/cmtomicron);
484  else
485  resMultiPixelXHist_[icotbeta][icotalpha][iqbin]->Fill(dx/cmtomicron);
486  if( nypix == 1 )
487  resSinglePixelYHist_[icotbeta][icotalpha]->Fill(dy/cmtomicron);
488  else
489  resMultiPixelYHist_[icotbeta][icotalpha][iqbin]->Fill(dy/cmtomicron);
490  }
491 
492  return 0;
493 }
494 
495 
496 //------------------------------------------------------------------------------
497 // Return the histogram generator for resolution in X. A generator contains
498 // both the histogram and knows how to throw a random number off it. It is
499 // called from FastSim (from PixelTemplateSmearerBase).
500 // If cotalpha or cotbeta are outside of the range, return the end of the range.
501 //------------------------------------------------------------------------------
504 getGeneratorX( double cotalpha, double cotbeta, int qbin, bool single )
505 {
506  int icotalpha, icotbeta, iqbin ;
507  icotalpha = (int)floor( (cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_ ) ;
508  icotbeta = (int)floor( (cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_ ) ;
509  iqbin = qbin > 2 ? 3 : qbin; // if (qbin>2) then = 3, else return qbin
510  //
511  //if( icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_ ) {
512 
513  if (icotalpha < 0)
514  icotalpha = 0;
515  if (icotalpha >= cotalphaBins_)
516  icotalpha = cotalphaBins_ - 1;
517 
518  if (icotbeta < 0)
519  icotbeta = 0;
520  if (icotbeta >= cotbetaBins_)
521  icotbeta = cotbetaBins_ - 1;
522 
523  // At this point we are sure to return *some bin* from the 3D histogram
524 
525  if( single )
526  return resSinglePixelXGen_[icotbeta][icotalpha];
527  else
528  return resMultiPixelXGen_[icotbeta][icotalpha][iqbin];
529 
530  // }
531  //else
532  //return nullptr;
533 }
534 
535 //------------------------------------------------------------------------------
536 // Return the histogram generator for resolution in Y. A generator contains
537 // both the histogram and knows how to throw a random number off it. It is
538 // called from FastSim (from PixelTemplateSmearerBase).
539 // If cotalpha or cotbeta are outside of the range, return the end of the range.
540 //------------------------------------------------------------------------------
543 getGeneratorY( double cotalpha, double cotbeta, int qbin, bool single )
544 {
545  int icotalpha, icotbeta, iqbin ;
546  icotalpha = (int)floor( (cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_ ) ;
547  icotbeta = (int)floor( (cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_ ) ;
548  iqbin = qbin > 2 ? 3 : qbin; // if (qbin>2) then = 3, else return qbin
549  //
550  //if( icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_ ) {
551 
552  if (icotalpha < 0)
553  icotalpha = 0;
554  if (icotalpha >= cotalphaBins_)
555  icotalpha = cotalphaBins_ - 1;
556 
557  if (icotbeta < 0)
558  icotbeta = 0;
559  if (icotbeta >= cotbetaBins_)
560  icotbeta = cotbetaBins_ - 1;
561 
562  // At this point we are sure to return *some bin* from the 3D histogram
563 
564  if( single )
565  return resSinglePixelYGen_[icotbeta][icotalpha];
566  else
567  return resMultiPixelYGen_[icotbeta][icotalpha][iqbin];
568 
569  //}
570  //else
571  //return nullptr;
572 }
TH1F * resSinglePixelYHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
TH1F * resMultiPixelYHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
TH1F * resSinglePixelXHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
SimpleHistogramGenerator * qbinGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
const SimpleHistogramGenerator * getGeneratorX(double cotalpha, double cotbeta, int qbin, bool singlex)
const SimpleHistogramGenerator * getGeneratorY(double cotalpha, double cotbeta, int qbin, bool singley)
TH1F * qbinHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
std::unique_ptr< TFile > file_
const float cmtomicron
SimpleHistogramGenerator * resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
#define LOGDEBUG
int Fill(double dx, double dy, double cotalpha, double cotbeta, int qbin, int nxpix, int nypix)
#define nullptr
TH1F * resMultiPixelXHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resSinglePixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
#define LOGINFO
#define LOGERROR
ii
Definition: cuy.py:588
PixelResolutionHistograms(std::string filename, std::string rootdir, std::string descTitle, unsigned int detType, double cotbetaBinWidth, double cotbetaLowEdge, int cotbetaBins, double cotalphaBinWidth, double cotalphaLowEdge, int cotalphaBins)
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]