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 22 of file PixelResolutionHistograms.h.

Constructor & Destructor Documentation

◆ PixelResolutionHistograms() [1/2]

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 62 of file PixelResolutionHistograms.cc.

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

74  : weOwnHistograms_(true), // we'll be making some histos
75  detType_(detType),
76  cotbetaBinWidth_(cotbetaBinWidth),
77  cotbetaLowEdge_(cotbetaLowEdge),
78  cotbetaBins_(cotbetaBins),
79  cotalphaBinWidth_(cotalphaBinWidth),
80  cotalphaLowEdge_(cotalphaLowEdge),
81  cotalphaBins_(cotalphaBins),
82  qbinWidth_(1),
83  qbins_(4),
84  binningHisto_(nullptr),
86  resSinglePixelXHist_(), // all to nullptr
88  resSinglePixelYHist_(), // all to nullptr
89  qbinHist_(), // all to nullptr
90  file_(nullptr),
91  status_(0),
96  qbinGen_() {
97  file_ = std::make_unique<TFile>(filename.c_str(), "RECREATE");
98  //Resolution binning
99  // const double cotbetaBinWidth = 1.0;
100  // const double cotbetaLowEdge = -11.5 ;
101  // const int cotbetaBins = 23;
102  // const double cotalphaBinWidth = 0.08 ;
103  // const double cotalphaLowEdge = -0.36 ;
104  // const int cotalphaBins = 9;
105  // const int qbinWidth = 1;
106  // const int qbins = 4;
107 
108  // Dummy 2D histogram to store binning:
109  binningHisto_ = new TH2F("ResHistoBinning",
110  descTitle.c_str(),
111  cotbetaBins,
112  cotbetaLowEdge,
113  cotbetaLowEdge + cotbetaBins * cotbetaBinWidth,
114  cotalphaBins,
115  cotalphaLowEdge,
116  cotalphaLowEdge + cotalphaBins * cotalphaBinWidth);
117 
118  // Store detType in the underflow bin
119  binningHisto_->SetBinContent(0, 0, detType_);
120 
121  // All other histograms:
122  Char_t histo[200];
123  Char_t title[200];
124  //
125  //--- Histograms for clusters with multiple pixels hit in a given direction.
126  //
127  for (int ii = 0; ii < cotbetaBins_; ii++) {
128  for (int jj = 0; jj < cotalphaBins_; jj++) {
129  for (int kk = 0; kk < qbins_; kk++) {
130  //
131  sprintf(histo, "hx%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1); //information of bits of histogram names
132  //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha, qbins
133  sprintf(title,
134  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 X",
139  kk + 1);
140  //
141  resMultiPixelXHist_[ii][jj][kk] = new TH1F(histo, title, 1000, -0.05, 0.05);
142 
143  sprintf(histo, "hy%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1);
144  sprintf(title,
145  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 Y",
150  kk + 1);
151  //
152  resMultiPixelYHist_[ii][jj][kk] = new TH1F(histo, title, 1000, -0.05, 0.05);
153  }
154  }
155  }
156 
157  //
158  //--- Histograms for clusters where only a single pixel was hit in a given direction.
159  //
160  for (int ii = 0; ii < cotbetaBins_; ii++) {
161  for (int jj = 0; jj < cotalphaBins_; jj++) {
162  sprintf(histo, "hx%d0%02d%d", detType_, ii + 1, jj + 1); //information of bits of histogram names
163  //first bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha
164  sprintf(title,
165  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 X",
170  //
171  resSinglePixelXHist_[ii][jj] = new TH1F(histo, title, 1000, -0.05, 0.05);
172 
173  sprintf(histo, "hy%d0%02d%d", detType_, ii + 1, jj + 1);
174  sprintf(title,
175  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 Y",
180  //
181  resSinglePixelYHist_[ii][jj] = new TH1F(histo, title, 1000, -0.05, 0.05);
182 
183  sprintf(histo, "hqbin%d%02d%d", detType_, ii + 1, jj + 1);
184  sprintf(title,
185  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin",
190  //
191  qbinHist_[ii][jj] = new TH1F(histo, title, 4, -0.49, 3.51);
192  }
193  }
194 }
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:589
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]

◆ PixelResolutionHistograms() [2/2]

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 204 of file PixelResolutionHistograms.cc.

References binningHisto_, cotalphaAxis_, cotalphaBins_, cotalphaBinWidth_, cotalphaLowEdge_, cotbetaAxis_, cotbetaBins_, cotbetaBinWidth_, cotbetaLowEdge_, detType_, file_, corrVsCorr::filename, timingPdfMaker::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 runGCPTkAlMap::title.

210  : weOwnHistograms_(false), // resolution histograms are owned by the ROOT file
211  detType_(-1),
212  cotbetaBinWidth_(0),
213  cotbetaLowEdge_(0),
214  cotbetaBins_(0),
216  cotalphaLowEdge_(0),
217  cotalphaBins_(0),
218  qbinWidth_(1),
219  qbins_(4),
220  binningHisto_(nullptr),
222  resSinglePixelXHist_(), // all to nullptr
224  resSinglePixelYHist_(), // all to nullptr
225  qbinHist_(), // all to nullptr
226  file_(nullptr),
227  status_(0),
232  qbinGen_() {
233  Char_t histo[200]; // the name of the histogram
234  Char_t title[200]; // histo title, for debugging and sanity checking (compare inside file)
235  TH1F* tmphist = nullptr; // cache for histo pointer
236 
237  //--- Open the file for reading.
238  file_ = std::make_unique<TFile>(filename.c_str(), "READ");
239  if (!file_) {
240  status_ = 1;
241  LOGERROR << "PixelResolutionHistograms:: Error, file " << filename << " not found.";
242  return; // PixelTemplateSmearerBase will throw an exception upon our return.
243  }
244 
245  //--- The dummy 2D histogram with the binning of cot\beta and cot\alpha:
246  binningHisto_ = (TH2F*)file_->Get(Form("%s%s", rootdir.c_str(), "ResHistoBinning"));
247  if (!binningHisto_) {
248  status_ = 11;
249  LOGERROR << "PixelResolutionHistograms:: Error, binning histogrram ResHistoBinning not found.";
250  return; // PixelTemplateSmearerBase will throw an exception upon our return.
251  }
252 
253  if (detType == -1) {
254  //--- Fish out detType from the underflow bin:
255  detType_ = binningHisto_->GetBinContent(0, 0);
256  } else {
257  detType_ = detType; // constructor's argument overrides what's in ResHistoBinning histogram.
258  }
259 
260  //--- Now we fill the binning variables:
261  cotbetaAxis_ = binningHisto_->GetXaxis();
262  cotbetaBinWidth_ = binningHisto_->GetXaxis()->GetBinWidth(1); // assume all same width
263  cotbetaLowEdge_ = binningHisto_->GetXaxis()->GetXmin(); // low edge of the first bin
264  cotbetaBins_ = binningHisto_->GetXaxis()->GetNbins();
265  cotalphaAxis_ = binningHisto_->GetYaxis();
266  cotalphaBinWidth_ = binningHisto_->GetYaxis()->GetBinWidth(1); // assume all same width;
267  cotalphaLowEdge_ = binningHisto_->GetYaxis()->GetXmin(); // low edge of the first bin;
268  cotalphaBins_ = binningHisto_->GetYaxis()->GetNbins();
269 
270  if (!ignore_multi) {
271  //
272  //--- Histograms for clusters with multiple pixels hit in a given direction.
273  //
274  for (int ii = 0; ii < cotbetaBins_; ii++) {
275  for (int jj = 0; jj < cotalphaBins_; jj++) {
276  for (int kk = 0; kk < qbins_; kk++) {
277  //
278  sprintf(histo, "hx%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1); //information of bits of histogram names
279  //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha, qbins
280  sprintf(title,
281  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 X",
286  kk + 1);
287  //
288  tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
289  if (!tmphist) {
290  status_ = 2;
291  LOGERROR << "Failed to find histogram=" << std::string(histo);
292  return;
293  }
294  LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
295  << std::endl;
296  if (tmphist->GetEntries() < 5) {
297  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
298  << " entries. Trouble ahead." << std::endl;
299  }
300  resMultiPixelXHist_[ii][jj][kk] = tmphist;
302 
303  sprintf(histo, "hy%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1);
304  sprintf(title,
305  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 Y",
310  kk + 1);
311  //
312  tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
313  if (!tmphist) {
314  status_ = 3;
315  LOGERROR << "Failed to find histogram=" << std::string(histo);
316  return;
317  }
318  LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
319  << std::endl;
320  if (tmphist->GetEntries() < 5) {
321  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
322  << " entries. Trouble ahead." << std::endl;
323  }
324  resMultiPixelYHist_[ii][jj][kk] = tmphist;
326  }
327  }
328  }
329  //
330  } // if (not ignore multi)
331 
332  //
333  //--- Histograms for clusters where only a single pixel was hit in a given direction.
334  //
335  for (int ii = 0; ii < cotbetaBins_; ii++) {
336  for (int jj = 0; jj < cotalphaBins_; jj++) {
337  //--- Single pixel, along X.
338  //
339  sprintf(histo, "hx%d0%02d%d", detType_, ii + 1, jj + 1); //information of bits of histogram names
340  //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha
341  sprintf(title,
342  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 X",
347  //
348  tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
349  if (!tmphist) {
350  if (!ignore_single) {
351  LOGERROR << "Failed to find histogram=" << std::string(histo);
352  status_ = 4;
353  return;
354  }
355  } else {
356  LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
357  << std::endl;
358  LOGDEBUG << "Found histo with title = " << std::string(tmphist->GetTitle()) << std::endl;
359  if (tmphist->GetEntries() < 5) {
360  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
361  << " entries. Trouble ahead." << std::endl;
362  }
363  resSinglePixelXHist_[ii][jj] = tmphist;
365  }
366 
367  //--- Single pixel, along Y.
368  //
369  sprintf(histo, "hy%d0%02d%d", detType_, ii + 1, jj + 1);
370  sprintf(title,
371  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 Y",
376  //
377  tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
378  if (!tmphist) {
379  if (!ignore_single) {
380  LOGERROR << "Failed to find histogram=" << std::string(histo);
381  status_ = 5;
382  return;
383  }
384  } else {
385  LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
386  << std::endl;
387  if (tmphist->GetEntries() < 5) {
388  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
389  << " entries. Trouble ahead." << std::endl;
390  }
391  resSinglePixelYHist_[ii][jj] = tmphist;
393  }
394 
395  //--- qBin distribution, for this (cotbeta, cotalpha) bin.
396  //
397  sprintf(histo, "hqbin%d%02d%d", detType_, ii + 1, jj + 1);
398  sprintf(title,
399  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin",
404  //
405  tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
406  if (!tmphist) {
407  if (!ignore_qBin) {
408  LOGERROR << "Failed to find histogram=" << std::string(histo);
409  status_ = 6;
410  return;
411  }
412  } else {
413  LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
414  << std::endl;
415  if (tmphist->GetEntries() < 5) {
416  LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
417  << " entries. Trouble ahead." << std::endl;
418  }
419  qbinHist_[ii][jj] = tmphist;
420  qbinGen_[ii][jj] = new SimpleHistogramGenerator(tmphist);
421  }
422  }
423  }
424 }
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:589
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]

◆ ~PixelResolutionHistograms()

PixelResolutionHistograms::~PixelResolutionHistograms ( )
virtual

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

Definition at line 431 of file PixelResolutionHistograms.cc.

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

431  {
432  //--- Delete histograms, but only if we own them. If
433  //--- they came from a file, let them be.
434  //
435  if (!weOwnHistograms_) {
436  //--- Read the histograms from the TFile, the file will take care of them.
437  file_->Close();
440  } else {
441  //--- We made the histograms, so first write them inthe output ROOT file and close it.
442  LOGINFO << "PixelResHistoStore: Writing the histograms to the output file. " // << filename
443  << std::endl;
444  file_->Write();
445  file_->Close();
446 
447  // ROOT file has the ownership, and once the file is closed,
448  // all of these guys are deleted. So, we don't need to do anything.
449  } // else
450 
451  //--- Delete FastSim generators. (It's safe to delete a nullptr.)
452  for (int ii = 0; ii < cotbetaBins_; ii++) {
453  for (int jj = 0; jj < cotalphaBins_; jj++) {
454  for (int kk = 0; kk < qbins_; kk++) {
455  delete resMultiPixelXGen_[ii][jj][kk];
456  delete resMultiPixelYGen_[ii][jj][kk];
457  }
458  }
459  }
460  for (int ii = 0; ii < cotbetaBins_; ii++) {
461  for (int jj = 0; jj < cotalphaBins_; jj++) {
462  delete resSinglePixelXGen_[ii][jj];
463  delete resSinglePixelYGen_[ii][jj];
464  delete qbinGen_[ii][jj];
465  }
466  }
467 }
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:589
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]

Member Function Documentation

◆ Fill()

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_, PVValHelper::dx, PVValHelper::dy, createfilelist::int, qbinHist_, resMultiPixelXHist_, resMultiPixelYHist_, resSinglePixelXHist_, and resSinglePixelYHist_.

474  {
475  int icotalpha, icotbeta, iqbin;
476  icotalpha = (int)floor((cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_);
477  icotbeta = (int)floor((cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_);
478  iqbin = qbin > 2 ? 3 : qbin;
479  if (icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_) {
480  qbinHist_[icotbeta][icotalpha]->Fill((double)iqbin);
481  if (nxpix == 1)
482  resSinglePixelXHist_[icotbeta][icotalpha]->Fill(dx / cmtomicron);
483  else
484  resMultiPixelXHist_[icotbeta][icotalpha][iqbin]->Fill(dx / cmtomicron);
485  if (nypix == 1)
486  resSinglePixelYHist_[icotbeta][icotalpha]->Fill(dy / cmtomicron);
487  else
488  resMultiPixelYHist_[icotbeta][icotalpha][iqbin]->Fill(dy / cmtomicron);
489  }
490 
491  return 0;
492 }
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]

◆ getGeneratorX()

const SimpleHistogramGenerator * PixelResolutionHistograms::getGeneratorX ( double  cotalpha,
double  cotbeta,
int  qbin,
bool  singlex 
)

Definition at line 500 of file PixelResolutionHistograms.cc.

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

503  {
504  int icotalpha, icotbeta, iqbin;
505  icotalpha = (int)floor((cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_);
506  icotbeta = (int)floor((cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_);
507  iqbin = qbin > 2 ? 3 : qbin; // if (qbin>2) then = 3, else return qbin
508  //
509  //if( icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_ ) {
510 
511  if (icotalpha < 0)
512  icotalpha = 0;
513  if (icotalpha >= cotalphaBins_)
514  icotalpha = cotalphaBins_ - 1;
515 
516  if (icotbeta < 0)
517  icotbeta = 0;
518  if (icotbeta >= cotbetaBins_)
519  icotbeta = cotbetaBins_ - 1;
520 
521  // At this point we are sure to return *some bin* from the 3D histogram
522 
523  if (single)
524  return resSinglePixelXGen_[icotbeta][icotalpha];
525  else
526  return resMultiPixelXGen_[icotbeta][icotalpha][iqbin];
527 
528  // }
529  //else
530  //return nullptr;
531 }
SimpleHistogramGenerator * resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
int single(int argc, char *argv[])
Definition: DMRsingle.cc:18
SimpleHistogramGenerator * resSinglePixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]

◆ getGeneratorY()

const SimpleHistogramGenerator * PixelResolutionHistograms::getGeneratorY ( double  cotalpha,
double  cotbeta,
int  qbin,
bool  singley 
)

Definition at line 539 of file PixelResolutionHistograms.cc.

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

542  {
543  int icotalpha, icotbeta, iqbin;
544  icotalpha = (int)floor((cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_);
545  icotbeta = (int)floor((cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_);
546  iqbin = qbin > 2 ? 3 : qbin; // if (qbin>2) then = 3, else return qbin
547  //
548  //if( icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_ ) {
549 
550  if (icotalpha < 0)
551  icotalpha = 0;
552  if (icotalpha >= cotalphaBins_)
553  icotalpha = cotalphaBins_ - 1;
554 
555  if (icotbeta < 0)
556  icotbeta = 0;
557  if (icotbeta >= cotbetaBins_)
558  icotbeta = cotbetaBins_ - 1;
559 
560  // At this point we are sure to return *some bin* from the 3D histogram
561 
562  if (single)
563  return resSinglePixelYGen_[icotbeta][icotalpha];
564  else
565  return resMultiPixelYGen_[icotbeta][icotalpha][iqbin];
566 
567  //}
568  //else
569  //return nullptr;
570 }
SimpleHistogramGenerator * resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
int single(int argc, char *argv[])
Definition: DMRsingle.cc:18
SimpleHistogramGenerator * resSinglePixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]

◆ status()

int PixelResolutionHistograms::status ( void  )
inline

Definition at line 58 of file PixelResolutionHistograms.h.

References status_.

Member Data Documentation

◆ binningHisto_

TH2F* PixelResolutionHistograms::binningHisto_
private

Definition at line 92 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

◆ cotalphaAxis_

TAxis* PixelResolutionHistograms::cotalphaAxis_
private

Definition at line 94 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

◆ cotalphaBins_

int PixelResolutionHistograms::cotalphaBins_
private

◆ cotalphaBinWidth_

double PixelResolutionHistograms::cotalphaBinWidth_
private

◆ cotalphaLowEdge_

double PixelResolutionHistograms::cotalphaLowEdge_
private

◆ cotbetaAxis_

TAxis* PixelResolutionHistograms::cotbetaAxis_
private

Definition at line 93 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

◆ cotbetaBins_

int PixelResolutionHistograms::cotbetaBins_
private

◆ cotbetaBinWidth_

double PixelResolutionHistograms::cotbetaBinWidth_
private

◆ cotbetaLowEdge_

double PixelResolutionHistograms::cotbetaLowEdge_
private

◆ detType_

unsigned int PixelResolutionHistograms::detType_
private

Definition at line 79 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms().

◆ file_

std::unique_ptr<TFile> PixelResolutionHistograms::file_
private

◆ qbinGen_

SimpleHistogramGenerator* PixelResolutionHistograms::qbinGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private

◆ qbinHist_

TH1F* PixelResolutionHistograms::qbinHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
private

Definition at line 112 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

◆ qbins_

int PixelResolutionHistograms::qbins_
private

◆ qbinWidth_

int PixelResolutionHistograms::qbinWidth_
private

Definition at line 88 of file PixelResolutionHistograms.h.

◆ resMultiPixelXGen_

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

◆ resMultiPixelXHist_

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

Definition at line 108 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

◆ resMultiPixelYGen_

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

◆ resMultiPixelYHist_

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

Definition at line 110 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

◆ resSinglePixelXGen_

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

◆ resSinglePixelXHist_

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

Definition at line 109 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

◆ resSinglePixelYGen_

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

◆ resSinglePixelYHist_

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

Definition at line 111 of file PixelResolutionHistograms.h.

Referenced by Fill(), and PixelResolutionHistograms().

◆ status_

int PixelResolutionHistograms::status_
private

Definition at line 118 of file PixelResolutionHistograms.h.

Referenced by PixelResolutionHistograms(), and status().

◆ weOwnHistograms_

bool PixelResolutionHistograms::weOwnHistograms_
private

Definition at line 76 of file PixelResolutionHistograms.h.

Referenced by ~PixelResolutionHistograms().