CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PixelResolutionHistograms.h
Go to the documentation of this file.
1 #ifndef FastSimulation_TrackingRecHitProducer_PixelResolutionHistograms_h
2 #define FastSimulation_TrackingRecHitProducer_PixelResolutionHistograms_h 1
3 
4 class TFile;
5 class TH1F;
6 class TH2F;
7 class TAxis;
10 
11 #include <memory>
12 #include <string>
13 
17 
18 static constexpr unsigned int COTBETA_HIST_MAX = 30;
19 static constexpr unsigned int COTALPHA_HIST_MAX = 20;
20 static constexpr unsigned int QBIN_HIST_MAX = 4;
21 
23 public:
24  //--- Constructor to use when generating resolution histograms.
25  // We make empty histograms (which we own), but generator pointers
26  // remain null.
27  //
28  PixelResolutionHistograms(std::string filename, // ROOT file for histograms
29  std::string rootdir, // Subdirectory in the file, "" if none
30  std::string descTitle, // Descriptive title
31  unsigned int detType, // Where we are... (&&& do we need this?)
32  double cotbetaBinWidth, // cot(beta) : bin width,
33  double cotbetaLowEdge, // : low endpoint,
34  int cotbetaBins, // : # of bins
35  double cotalphaBinWidth, // cot(alpha): bin width,
36  double cotalphaLowEdge, // : low endpoint,
37  int cotalphaBins); // : # of bins
38  //int qbinWidth,
39  //int qbins )
40 
41  //--- Constructor to use when reading the histograms from a file (e.g. when
42  // inside a running FastSim job). We get the histograms from a
43  // ROOT file, and we do *not* own them. But we do own the
44  // generators.
45  //
46  PixelResolutionHistograms(std::string filename, // ROOT file for histograms
47  std::string rootdir = "", // ROOT dir, "" if none
48  int detType = -1, // default: read from ROOT file.
49  bool ignore_multi = false, // Forward Big is always single
50  bool ignore_single = false, // Edge does not need single pixels
51  bool ignore_qBin = false); // qBin histograms not used right now (future expansion)
52 
53  //--- Destructor (virtual, just in case)
55 
56  //--- Status after construction (esp.loading from file). Non-zero if there
57  // were problems.
58  inline int status() { return status_; }
59 
60  //--- Fill one entry in one resolution histogram. Use when making histograms.
61  int Fill(double dx,
62  double dy, // the difference wrt true hit
63  double cotalpha,
64  double cotbeta, // cotangent of local angles
65  int qbin, // Qbin = category for how much charge we have
66  int nxpix,
67  int nypix); // length of cluster along x,y (only care if ==1 or not)
68 
69  //--- Get generators, for resolution in X and Y. Use in FastSim.
70  const SimpleHistogramGenerator* getGeneratorX(double cotalpha, double cotbeta, int qbin, bool singlex);
71 
72  const SimpleHistogramGenerator* getGeneratorY(double cotalpha, double cotbeta, int qbin, bool singley);
73 
74 private:
75  // Do we own the histograms, or not?
77 
78  // Where we are.
79  unsigned int detType_; // 1 for barrel, 0 for forward /// May not need this?
80 
81  // Resolution binning
89  int qbins_;
90 
91  // The dummy histogram to hold the binning, and the two cached axes.
93  TAxis* cotbetaAxis_;
94  TAxis* cotalphaAxis_;
95 
96  // Resolution histograms. I (Petar) tried to dynamically allocate
97  // these histograms, but all possible implementations were somewhat
98  // complicated, which would make the code harder to understand,
99  // debug, and thus support in the long term. Since we are here only
100  // booking pointers of histograms, we will instead book larger
101  // matrices, and leave them partially empty. But who cares -- the
102  // wasted memory of a few hundred null pointers is negligible.
103  //
104  // The constructor will then fill only the first cotbetaBins_ x
105  // cotalphaBins_ x qbinBins_ histograms in the matrix, and we'll
106  // ignore the rest.
107  //
113 
114  // File with histograms to load.
115  std::unique_ptr<TFile> file_;
116 
117  // Status of loading. Check if there were errors.
118  int status_;
119 
120  // Identical binning and parameterization for FastSim generators.
126 };
127 #endif
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]
static constexpr unsigned int COTBETA_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_
SimpleHistogramGenerator * resMultiPixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resMultiPixelXGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
int Fill(double dx, double dy, double cotalpha, double cotbeta, int qbin, int nxpix, int nypix)
TH1F * resMultiPixelXHist_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX][QBIN_HIST_MAX]
SimpleHistogramGenerator * resSinglePixelYGen_[COTBETA_HIST_MAX][COTALPHA_HIST_MAX]
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]
static constexpr unsigned int COTALPHA_HIST_MAX
tuple filename
Definition: lut2db_cfg.py:20
static constexpr unsigned int QBIN_HIST_MAX