CMS 3D CMS Logo

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