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 
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 
42  //--- Constructor to use when reading the histograms from a file (e.g. when
43  // inside a running FastSim job). We get the histograms from a
44  // ROOT file, and we do *not* own them. But we do own the
45  // generators.
46  //
47  PixelResolutionHistograms( std::string filename, // ROOT file for histograms
48  std::string rootdir = "", // ROOT dir, "" if none
49  int detType = -1, // default: read from ROOT file.
50  bool ignore_multi = false, // Forward Big is always single
51  bool ignore_single = false, // Edge does not need single pixels
52  bool ignore_qBin = false ); // qBin histograms not used right now (future expansion)
53 
54 
55  //--- Destructor (virtual, just in case)
57 
58  //--- Status after construction (esp.loading from file). Non-zero if there
59  // were problems.
60  inline int status() { return status_ ; }
61 
62  //--- Fill one entry in one resolution histogram. Use when making histograms.
63  int Fill( double dx, double dy, // the difference wrt true hit
64  double cotalpha, double cotbeta, // cotangent of local angles
65  int qbin, // Qbin = category for how much charge we have
66  int nxpix, int nypix ); // length of cluster along x,y (only care if ==1 or not)
67 
68 
69  //--- Get generators, for resolution in X and Y. Use in FastSim.
70  const SimpleHistogramGenerator * getGeneratorX( double cotalpha,
71  double cotbeta,
72  int qbin,
73  bool singlex );
74 
75  const SimpleHistogramGenerator * getGeneratorY( double cotalpha,
76  double cotbeta,
77  int qbin,
78  bool singley );
79 
80 
81  private:
82  // Do we own the histograms, or not?
84 
85  // Where we are.
86  unsigned int detType_ ; // 1 for barrel, 0 for forward /// May not need this?
87 
88  // Resolution binning
90  double cotbetaLowEdge_ ;
95  int qbinWidth_ ;
96  int qbins_ ;
97 
98  // The dummy histogram to hold the binning, and the two cached axes.
99  TH2F * binningHisto_ ;
100  TAxis * cotbetaAxis_ ;
101  TAxis * cotalphaAxis_ ;
102 
103  // Resolution histograms. I (Petar) tried to dynamically allocate
104  // these histograms, but all possible implementations were somewhat
105  // complicated, which would make the code harder to understand,
106  // debug, and thus support in the long term. Since we are here only
107  // booking pointers of histograms, we will instead book larger
108  // matrices, and leave them partially empty. But who cares -- the
109  // wasted memory of a few hundred null pointers is negligible.
110  //
111  // The constructor will then fill only the first cotbetaBins_ x
112  // cotalphaBins_ x qbinBins_ histograms in the matrix, and we'll
113  // ignore the rest.
114  //
120 
121  // File with histograms to load.
122  std::unique_ptr<TFile> file_ ;
123 
124  // Status of loading. Check if there were errors.
125  int status_ ;
126 
127  // Identical binning and parameterization for FastSim generators.
133 
134 };
135 #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)
#define constexpr
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
static constexpr unsigned int QBIN_HIST_MAX