Belle II Software  release-08-01-10
DQMHistoModuleBase.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 #include <framework/core/HistoModule.h>
12 #include <framework/geometry/B2Vector3.h>
13 #include <mdst/dataobjects/Track.h>
14 #include <tracking/dataobjects/RecoTrack.h>
15 #include <framework/core/ModuleParam.templateDetails.h>
16 
17 #include <TH1F.h>
18 #include <TH2F.h>
19 
20 namespace Belle2 {
30 
31  public:
32 
37 
42 
46  virtual void initialize() override;
47 
51  virtual void beginRun() override;
52 
56  virtual void event() override;
57 
60  virtual void defineHisto() override;
61 
63  void runningOnHLT() {m_hltDQM = true;};
64 
67  virtual TH1F* Create(std::string name, std::string title, int nbinsx, double xlow, double xup, std::string xTitle,
68  std::string yTitle);
70  virtual TH2F* Create(std::string name, std::string title, int nbinsx, double xlow, double xup, int nbinsy, double ylow, double yup,
71  std::string xTitle, std::string yTitle, std::string zTitle);
72 
82  virtual TH1F** CreateLayers(boost::format nameTemplate, boost::format titleTemplate, int nbinsx, double xlow, double xup,
83  std::string xTitle, std::string yTitle);
85  virtual TH2F** CreateLayers(boost::format nameTemplate, boost::format titleTemplate, int nbinsx, double xlow, double xup,
86  int nbinsy, double ylow, double yup, std::string xTitle, std::string yTitle, std::string zTitle);
87 
97  virtual TH1F** CreateSensors(boost::format nameTemplate, boost::format titleTemplate, int nbinsx, double xlow, double xup,
98  std::string xTitle, std::string yTitle);
100  virtual TH2F** CreateSensors(boost::format nameTemplate, boost::format titleTemplate, int nbinsx, double xlow, double xup,
101  int nbinsy, double ylow, double yup, std::string xTitle, std::string yTitle, std::string zTitle);
102 
108  virtual void FillTrackIndexes(int iTrack, int iTrackVXD, int iTrackCDC, int iTrackVXDCDC);
110  virtual void FillHitNumbers(int nPXD, int nSVD, int nCDC);
112  virtual void FillMomentumAngles(const TrackFitResult* tfr);
114  virtual void FillMomentumCoordinates(const TrackFitResult* tfr);
116  virtual void FillHelixParametersAndCorrelations(const TrackFitResult* tfr);
118  virtual void FillTrackFitStatus(const genfit::FitStatus* tfs);
120  virtual void FillTRClusterCorrelations(float phi_deg, float phiPrev_deg, float theta_deg, float thetaPrev_deg,
121  int correlationIndex);
123  virtual void FillTRClusterHitmap(float phi_deg, float theta_deg, int layerIndex);
125  virtual void FillUBResidualsPXD(const B2Vector3D& residual_um);
127  virtual void FillUBResidualsSVD(const B2Vector3D& residual_um);
129  virtual void FillHalfShellsPXD(const B2Vector3D& globalResidual_um, bool isNotYang);
131  virtual void FillHalfShellsSVD(const B2Vector3D& globalResidual_um, bool isNotMat);
133  virtual void FillUB1DResidualsSensor(const B2Vector3D& residual_um, int sensorIndex);
135  virtual void FillUB2DResidualsSensor(const B2Vector3D& residual_um, int sensorIndex);
138  protected:
140  static std::string SensorNameDescription(VxdID sensorID);
142  static std::string SensorTitleDescription(VxdID sensorID);
143 
148  static void ComputeMean(TH1F* output, TH2F* input, bool onX = true);
149 
155  void ProcessHistogramParameterChange(const std::string& name, const std::string& parameter, const std::string& value);
157  void EditHistogramParameter(TH1* histogram, const std::string& parameter, std::string value);
158 
164  virtual void DefineTracks();
166  virtual void DefineHits();
168  virtual void DefineMomentumAngles();
170  virtual void DefineMomentumCoordinates();
174  virtual void DefineTrackFitStatus();
176  virtual void DefineTRClusters();
178  virtual void DefineUBResidualsVXD();
180  virtual void DefineHalfShellsVXD();
182  virtual void Define1DSensors();
184  virtual void Define2DSensors();
189  std::vector<TH1*> m_histograms;
191  bool histogramsDefined = false;
193  bool m_hltDQM = false;
194 
196  std::vector<std::tuple<std::string, std::string, std::string>> m_histogramParameterChanges;
201 
203  TH1F* m_PValue = nullptr;
205  TH1F* m_Chi2 = nullptr;
207  TH1F* m_NDF = nullptr;
209  TH1F* m_Chi2NDF = nullptr;
211  TH2F* m_UBResidualsPXD = nullptr;
213  TH2F* m_UBResidualsSVD = nullptr;
215  TH2F** m_UBResidualsSensor = nullptr;
217  TH1F* m_UBResidualsPXDU = nullptr;
219  TH1F* m_UBResidualsSVDU = nullptr;
221  TH1F** m_UBResidualsSensorU = nullptr;
223  TH1F* m_UBResidualsPXDV = nullptr;
225  TH1F* m_UBResidualsSVDV = nullptr;
226 
229  TH1F* m_UBResidualsPXDX_Yin = nullptr;
231  TH1F* m_UBResidualsPXDX_Yang = nullptr;
233  TH1F* m_UBResidualsSVDX_Pat = nullptr;
235  TH1F* m_UBResidualsSVDX_Mat = nullptr;
236 
238  TH1F* m_UBResidualsPXDY_Yin = nullptr;
240  TH1F* m_UBResidualsPXDY_Yang = nullptr;
242  TH1F* m_UBResidualsSVDY_Pat = nullptr;
244  TH1F* m_UBResidualsSVDY_Mat = nullptr;
245 
247  TH1F* m_UBResidualsPXDZ_Yin = nullptr;
249  TH1F* m_UBResidualsPXDZ_Yang = nullptr;
251  TH1F* m_UBResidualsSVDZ_Pat = nullptr;
253  TH1F* m_UBResidualsSVDZ_Mat = nullptr;
254 
256  TH1F** m_UBResidualsSensorV = nullptr;
258  TH2F** m_TRClusterHitmap = nullptr;
260  TH2F** m_TRClusterCorrelationsPhi = nullptr;
263 
265  TH1F* m_MomPhi = nullptr;
267  TH1F* m_MomTheta = nullptr;
269  TH1F* m_MomCosTheta = nullptr;
271  TH1F* m_MomX = nullptr;
273  TH1F* m_MomY = nullptr;
275  TH1F* m_MomZ = nullptr;
277  TH1F* m_MomPt = nullptr;
279  TH1F* m_Mom = nullptr;
281  TH1F* m_D0 = nullptr;
283  TH2F* m_PhiD0 = nullptr;
285  TH1F* m_Z0 = nullptr;
287  TH2F* m_D0Z0 = nullptr;
289  TH1F* m_Phi = nullptr;
291  TH1F* m_TanLambda = nullptr;
293  TH1F* m_Omega = nullptr;
294 
296  TH1F* m_HitsPXD = nullptr;
298  TH1F* m_HitsSVD = nullptr;
300  TH1F* m_HitsCDC = nullptr;
302  TH1F* m_Hits = nullptr;
304  TH1F* m_TracksVXD = nullptr;
306  TH1F* m_TracksCDC = nullptr;
308  TH1F* m_TracksVXDCDC = nullptr;
310  TH1F* m_Tracks = nullptr;
311 
312  private:
314  inline bool checkVariableForNANOrINF(const double var)
315  {
316  return std::isnan(var) or std::isinf(var);
317  }
318  };
320 }
This class serves as a base for the TrackDQMModule and AlignDQMModule (and possibly other DQM histogr...
TH1F * m_Hits
Number of all hits in tracks.
virtual void Define1DSensors()
Define 1D histograms with unbiased residuals for individual sensors.
TH1F * m_UBResidualsSVDV
Unbiased residuals for SVD v.
bool histogramsDefined
True if the defineHisto() was called.
TH1F * m_MomTheta
Track momentum Pt.Theta.
TH1F * m_Tracks
Number of all finding tracks.
bool checkVariableForNANOrINF(const double var)
Check a variable for whether or not it is NAN or INF.
TH1F * m_UBResidualsPXDY_Yang
Unbiased residuals in Y for PXD for Yang.
TH1F * m_MomPt
Track momentum Pt.
TH1F * m_TanLambda
TanLambda - the slope of the track in the r-z plane.
TH1F ** m_UBResidualsSensorU
Unbiased residuals for PXD and SVD u per sensor.
virtual void FillHitNumbers(int nPXD, int nSVD, int nCDC)
Fill histograms with numbers of hits.
TH1F * m_UBResidualsPXDZ_Yang
Unbiased residuals in Z for PXD for Yang.
static std::string SensorTitleDescription(VxdID sensorID)
Creates string description of the sensor from given sensor ID to be used in a histogram title.
TH1F * m_UBResidualsPXDU
Unbiased residuals for PXD u.
std::vector< std::tuple< std::string, std::string, std::string > > m_histogramParameterChanges
Used for changing parameters of histograms via the ProcessHistogramParameterChange function.
virtual void DefineHalfShellsVXD()
Define histograms with unbiased residuals for half-shells for PXD and SVD sensors.
TH1F * m_TracksCDC
Number of tracks only with CDC.
TH1F * m_UBResidualsSVDY_Pat
Unbiased residuals in Y for PXD for Pat.
TH1F * m_HitsSVD
Number of hits on SVD.
virtual void initialize() override
Initializer.
virtual void DefineUBResidualsVXD()
Define histograms with unbiased residuals in PXD and SVD sensors.
std::vector< TH1 * > m_histograms
All histograms created via the Create- functions are automatically added to this set.
virtual void DefineTRClusters()
Define histograms with correlations between neighbor layers and cluster hitmap in IP angle range.
TH2F ** m_TRClusterHitmap
Track related clusters - hitmap in IP angle range.
void ProcessHistogramParameterChange(const std::string &name, const std::string &parameter, const std::string &value)
Process one change in histogram parameters.
TH1F * m_UBResidualsSVDZ_Mat
Unbiased residuals in Z for PXD for Mat.
TH1F * m_Omega
Omega - the curvature of the track.
virtual void event() override
This method is called for each event.
TH2F * m_UBResidualsPXD
Unbiased residuals for PXD u vs v.
static std::string SensorNameDescription(VxdID sensorID)
Creates string description of the sensor from given sensor ID to be used in a histogram name.
TH1F * m_TracksVXD
Number of tracks only with VXD.
TH1F * m_UBResidualsSVDX_Mat
Unbiased residuals in X for PXD for Mat.
virtual TH1F ** CreateSensors(boost::format nameTemplate, boost::format titleTemplate, int nbinsx, double xlow, double xup, std::string xTitle, std::string yTitle)
Function to create array of TH1F histograms, one for each sensor.
TH1F * m_UBResidualsPXDX_Yang
Unbiased residuals in X for PXD for Yang.
virtual void DefineTrackFitStatus()
Define histograms which require FitStatus.
virtual void FillHalfShellsSVD(const B2Vector3D &globalResidual_um, bool isNotMat)
Fill histograms with unbiased residuals for half-shells for SVD sensors.
std::string m_tracksStoreArrayName
StoreArray name where Tracks are written.
TH1F * m_Mom
Track momentum Magnitude.
TH1F * m_UBResidualsPXDX_Yin
half-shells
TH1F * m_UBResidualsSVDU
Unbiased residuals for SVD u.
TH2F * m_PhiD0
d0 vs Phi - the signed distance to the IP in the r-phi plane
virtual void FillUB2DResidualsSensor(const B2Vector3D &residual_um, int sensorIndex)
Fill 2D histograms with unbiased residuals for individual sensors.
virtual void FillUBResidualsPXD(const B2Vector3D &residual_um)
Fill histograms with unbiased residuals in PXD sensors.
virtual TH1F * Create(std::string name, std::string title, int nbinsx, double xlow, double xup, std::string xTitle, std::string yTitle)
Function to create TH1F and add it to the vector of histograms (m_histograms).
TH1F * m_MomX
Track momentum Pt.X.
virtual void FillTrackFitStatus(const genfit::FitStatus *tfs)
Fill histograms which require FitStatus.
virtual void FillTRClusterHitmap(float phi_deg, float theta_deg, int layerIndex)
Fill cluster hitmap in IP angle range.
TH1F * m_UBResidualsSVDX_Pat
Unbiased residuals in X for PXD for Pat.
TH2F ** m_TRClusterCorrelationsPhi
Track related clusters - neighbor correlations in Phi.
TH1F * m_UBResidualsSVDZ_Pat
Unbiased residuals in Z for PXD for Pat.
TH2F * m_UBResidualsSVD
Unbiased residuals for SVD u vs v.
TH1F * m_MomPhi
Track momentum Pt.Phi.
TH1F * m_Phi
Phi - the angle of the transverse momentum in the r-phi plane, with CDF naming convention.
virtual void beginRun() override
Called when entering a new run.
TH1F * m_MomCosTheta
Track momentum Pt.CosTheta.
virtual void DefineHits()
Define histograms with numbers of hits.
virtual void FillHelixParametersAndCorrelations(const TrackFitResult *tfr)
Fill histograms with helix parameters and their correlations.
TH2F ** m_UBResidualsSensor
Unbiased residuals for PXD and SVD u vs v per sensor.
TH1F * m_HitsCDC
Number of hits on CDC.
virtual void DefineMomentumCoordinates()
Define histograms with track momentum Pt.
virtual void FillMomentumCoordinates(const TrackFitResult *tfr)
Fill histograms with track momentum Pt.
TH1F ** m_UBResidualsSensorV
Unbiased residuals for PXD and SVD v per sensor.
static void ComputeMean(TH1F *output, TH2F *input, bool onX=true)
Creates a graph of means by given axis from given TH2F histogram.
virtual void DefineTracks()
All the following Define- functions should be used in the defineHisto() function to define histograms...
bool m_hltDQM
True if the DQM module is run on HLT.
TH2F * m_D0Z0
z0 vs d0 - signed distance to the IP in r-phi vs.
virtual void FillHalfShellsPXD(const B2Vector3D &globalResidual_um, bool isNotYang)
Fill histograms with unbiased residuals for half-shells for PXD sensors.
virtual void Define2DSensors()
Define 2D histograms with unbiased residuals for individual sensors.
TH1F * m_UBResidualsPXDY_Yin
Unbiased residuals in Y for PXD for Yin.
TH1F * m_MomY
Track momentum Pt.Y.
TH1F * m_UBResidualsSVDY_Mat
Unbiased residuals in Y for PXD for Mat.
virtual void DefineMomentumAngles()
Define histograms with track momentum Pt.
virtual TH1F ** CreateLayers(boost::format nameTemplate, boost::format titleTemplate, int nbinsx, double xlow, double xup, std::string xTitle, std::string yTitle)
Function to create array of TH1F histograms, one for each layer.
virtual void FillTrackIndexes(int iTrack, int iTrackVXD, int iTrackCDC, int iTrackVXDCDC)
Fill histograms with track indexes.
TH1F * m_HitsPXD
Number of hits on PXD.
void runningOnHLT()
function called when the module is run on HLT
TH1F * m_UBResidualsPXDZ_Yin
Unbiased residuals in Z for PXD for Yin.
std::string m_recoTracksStoreArrayName
StoreArray name where RecoTracks are written.
TH1F * m_TracksVXDCDC
Number of full tracks with VXD+CDC.
TH1F * m_MomZ
Track momentum Pt.Z.
TH1F * m_D0
d0 - the signed distance to the IP in the r-phi plane
TH1F * m_UBResidualsPXDV
Unbiased residuals for PXD v.
TH2F ** m_TRClusterCorrelationsTheta
Track related clusters - neighbor corelations in Theta.
virtual void FillTRClusterCorrelations(float phi_deg, float phiPrev_deg, float theta_deg, float thetaPrev_deg, int correlationIndex)
Fill histograms with correlations between neighbor layers.
virtual void FillMomentumAngles(const TrackFitResult *tfr)
Fill histograms with track momentum Pt.
virtual void FillUBResidualsSVD(const B2Vector3D &residual_um)
Fill histograms with unbiased residuals in SVD sensors.
virtual void DefineHelixParametersAndCorrelations()
Define histograms with helix parameters and their correlations.
TH1F * m_Z0
z0 - the z0 coordinate of the perigee (beam spot position)
virtual void FillUB1DResidualsSensor(const B2Vector3D &residual_um, int sensorIndex)
Fill 1D histograms with unbiased residuals for individual sensors.
void EditHistogramParameter(TH1 *histogram, const std::string &parameter, std::string value)
On given histogram sets given parameter to given value.
virtual void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Values of the result of a track fit with a given particle hypothesis.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
Abstract base class for different kinds of events.