CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
PixelResolutionHistograms Class Reference

#include <PixelResolutionHistograms.h>

Public Member Functions

int Fill (double dx, double dy, double cotalpha, double cotbeta, int qbin, int nxpix, int nypix)
 
const SimpleHistogramGeneratorgetGeneratorX (double cotalpha, double cotbeta, int qbin, bool singlex)
 
const SimpleHistogramGeneratorgetGeneratorY (double cotalpha, double cotbeta, int qbin, bool singley)
 
 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)
 
 PixelResolutionHistograms (std::string filename, std::string rootdir="", int detType=-1, bool ignore_multi=false, bool ignore_single=false, bool ignore_qBin=false)
 
int status ()
 
virtual ~PixelResolutionHistograms ()
 

Private Attributes

TH2F * binningHisto_
 
TAxis * cotalphaAxis_
 
int cotalphaBins_
 
double cotalphaBinWidth_
 
double cotalphaLowEdge_
 
TAxis * cotbetaAxis_
 
int cotbetaBins_
 
double cotbetaBinWidth_
 
double cotbetaLowEdge_
 
unsigned int detType_
 
std::unique_ptr< TFile > file_
 
SimpleHistogramGeneratorqbinGen_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
 
TH1F * qbinHist_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
 
int qbins_
 
int qbinWidth_
 
SimpleHistogramGeneratorresMultiPixelXGen_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
 
TH1F * resMultiPixelXHist_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
 
SimpleHistogramGeneratorresMultiPixelYGen_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
 
TH1F * resMultiPixelYHist_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
 
SimpleHistogramGeneratorresSinglePixelXGen_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
 
TH1F * resSinglePixelXHist_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
 
SimpleHistogramGeneratorresSinglePixelYGen_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
 
TH1F * resSinglePixelYHist_ [COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
 
int status_
 
bool weOwnHistograms_
 

Detailed Description

Definition at line 21 of file PixelResolutionHistograms.h.

Constructor & Destructor Documentation

PixelResolutionHistograms::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 
)

Definition at line 67 of file PixelResolutionHistograms.cc.

References binningHisto_, cotalphaBins_, cotalphaBinWidth_, cotalphaLowEdge_, cotbetaBins_, cotbetaBinWidth_, cotbetaLowEdge_, detType_, file_, trackerHits::histo, cuy::ii, findQualityFiles::jj, GetRecoTauVFromDQM_MC_cff::kk, qbinHist_, qbins_, resMultiPixelXHist_, resMultiPixelYHist_, resSinglePixelXHist_, resSinglePixelYHist_, and fftjetcommon_cfi::title.

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),
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 }
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]
TH1F * qbinHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
std::unique_ptr< TFile > file_
SimpleHistogramGenerator * resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
TH1F * resMultiPixelXHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resSinglePixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
ii
Definition: cuy.py:590
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
PixelResolutionHistograms::PixelResolutionHistograms ( std::string  filename,
std::string  rootdir = "",
int  detType = -1,
bool  ignore_multi = false,
bool  ignore_single = false,
bool  ignore_qBin = false 
)

Definition at line 194 of file PixelResolutionHistograms.cc.

References binningHisto_, cotalphaAxis_, cotalphaBins_, cotalphaBinWidth_, cotalphaLowEdge_, cotbetaAxis_, cotbetaBins_, cotbetaBinWidth_, cotbetaLowEdge_, detType_, file_, trackerHits::histo, cuy::ii, findQualityFiles::jj, GetRecoTauVFromDQM_MC_cff::kk, LOGDEBUG, LOGERROR, LOGINFO, qbinGen_, qbinHist_, qbins_, resMultiPixelXGen_, resMultiPixelXHist_, resMultiPixelYGen_, resMultiPixelYHist_, resSinglePixelXGen_, resSinglePixelXHist_, resSinglePixelYGen_, resSinglePixelYHist_, status_, AlCaHLTBitMon_QueryRunRegistry::string, and fftjetcommon_cfi::title.

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),
210  binningHisto_(nullptr),
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 }
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]
TH1F * qbinHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
std::unique_ptr< TFile > file_
SimpleHistogramGenerator * resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
#define LOGDEBUG
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:590
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
PixelResolutionHistograms::~PixelResolutionHistograms ( )
virtual

delete file_ ; // no need to delete if unique_ptr<> file_ = 0;

Definition at line 418 of file PixelResolutionHistograms.cc.

References cotalphaBins_, cotbetaBins_, file_, Fill(), cuy::ii, findQualityFiles::jj, GetRecoTauVFromDQM_MC_cff::kk, LOGINFO, qbinGen_, qbins_, resMultiPixelXGen_, resMultiPixelYGen_, resSinglePixelXGen_, resSinglePixelYGen_, and weOwnHistograms_.

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 }
SimpleHistogramGenerator * qbinGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
std::unique_ptr< TFile > file_
SimpleHistogramGenerator * resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resSinglePixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
#define LOGINFO
ii
Definition: cuy.py:590
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]

Member Function Documentation

int PixelResolutionHistograms::Fill ( double  dx,
double  dy,
double  cotalpha,
double  cotbeta,
int  qbin,
int  nxpix,
int  nypix 
)

Definition at line 473 of file PixelResolutionHistograms.cc.

References cmtomicron, cotalphaBinWidth_, cotalphaLowEdge_, cotbetaBins_, cotbetaBinWidth_, cotbetaLowEdge_, getGeneratorX(), createfilelist::int, qbinHist_, resMultiPixelXHist_, resMultiPixelYHist_, resSinglePixelXHist_, and resSinglePixelYHist_.

Referenced by status(), and ~PixelResolutionHistograms().

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 }
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]
TH1F * qbinHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
const float cmtomicron
TH1F * resMultiPixelXHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
const SimpleHistogramGenerator * PixelResolutionHistograms::getGeneratorX ( double  cotalpha,
double  cotbeta,
int  qbin,
bool  singlex 
)

Definition at line 504 of file PixelResolutionHistograms.cc.

References cotalphaBins_, cotalphaBinWidth_, cotalphaLowEdge_, cotbetaBins_, cotbetaBinWidth_, cotbetaLowEdge_, getGeneratorY(), createfilelist::int, resMultiPixelXGen_, and resSinglePixelXGen_.

Referenced by Fill(), and status().

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 }
SimpleHistogramGenerator * resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
const SimpleHistogramGenerator * PixelResolutionHistograms::getGeneratorY ( double  cotalpha,
double  cotbeta,
int  qbin,
bool  singley 
)

Definition at line 543 of file PixelResolutionHistograms.cc.

References cotalphaBins_, cotalphaBinWidth_, cotalphaLowEdge_, cotbetaBins_, cotbetaBinWidth_, cotbetaLowEdge_, createfilelist::int, resMultiPixelYGen_, and resSinglePixelYGen_.

Referenced by getGeneratorX(), and status().

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 }
SimpleHistogramGenerator * resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resSinglePixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
int PixelResolutionHistograms::status ( void  )
inline

Member Data Documentation

TH2F* PixelResolutionHistograms::binningHisto_
private

Definition at line 99 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

TAxis* PixelResolutionHistograms::cotalphaAxis_
private

Definition at line 101 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

int PixelResolutionHistograms::cotalphaBins_
private
double PixelResolutionHistograms::cotalphaBinWidth_
private
double PixelResolutionHistograms::cotalphaLowEdge_
private
TAxis* PixelResolutionHistograms::cotbetaAxis_
private

Definition at line 100 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

int PixelResolutionHistograms::cotbetaBins_
private
double PixelResolutionHistograms::cotbetaBinWidth_
private
double PixelResolutionHistograms::cotbetaLowEdge_
private
unsigned int PixelResolutionHistograms::detType_
private

Definition at line 86 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

std::unique_ptr<TFile> PixelResolutionHistograms::file_
private
SimpleHistogramGenerator* PixelResolutionHistograms::qbinGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private
TH1F* PixelResolutionHistograms::qbinHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private

Definition at line 119 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

int PixelResolutionHistograms::qbins_
private
int PixelResolutionHistograms::qbinWidth_
private

Definition at line 95 of file PixelResolutionHistograms.h.

SimpleHistogramGenerator* PixelResolutionHistograms::resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
private
TH1F* PixelResolutionHistograms::resMultiPixelXHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
private

Definition at line 115 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

SimpleHistogramGenerator* PixelResolutionHistograms::resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
private
TH1F* PixelResolutionHistograms::resMultiPixelYHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
private

Definition at line 117 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

SimpleHistogramGenerator* PixelResolutionHistograms::resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private
TH1F* PixelResolutionHistograms::resSinglePixelXHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private

Definition at line 116 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

SimpleHistogramGenerator* PixelResolutionHistograms::resSinglePixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private
TH1F* PixelResolutionHistograms::resSinglePixelYHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private

Definition at line 118 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

int PixelResolutionHistograms::status_
private

Definition at line 125 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms(), and status().

bool PixelResolutionHistograms::weOwnHistograms_
private

Definition at line 83 of file PixelResolutionHistograms.h.

Referenced by ~PixelResolutionHistograms().