Belle II Software  release-08-01-10
SVDDQMClustersOnTrackModule.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 <map>
12 #include <framework/core/HistoModule.h>
13 #include <mdst/dataobjects/SoftwareTriggerResult.h>
14 #include <framework/datastore/StoreObjPtr.h>
15 #include <framework/datastore/StoreArray.h>
16 #include <mdst/dataobjects/Track.h>
17 #include <svd/dataobjects/SVDEventInfo.h>
18 #include <framework/dataobjects/EventT0.h>
19 #include "TList.h"
20 #include "TH1F.h"
21 #include "TH2F.h"
22 
23 namespace Belle2 {
33  class SVDDQMClustersOnTrackModule : public HistoModule { // <- derived from HistoModule class
34 
35  public:
36 
41  /* Destructor */
42  virtual ~SVDDQMClustersOnTrackModule();
45 
47  void initialize() override final;
49  void terminate() override final;
51  void beginRun() override final;
53  void event() override final;
54 
56  void defineHisto() override final;
57 
58  private:
59 
61  bool m_desynchSVDTime = false;
62 
63  std::string m_svdShaperDigitsName;
64  std::string m_svdRecoDigitsName;
65  std::string m_svdClustersName;
66  std::string m_svdEventInfoName;
72 
75 
77  bool m_skipRejectedEvents = true;
78 
79  int m_tb = -1;
82  TList* m_histoList = nullptr;
83 
85  int m_expNumber = 0;
87  int m_runNumber = 0;
88 
91 
93  TH1F** m_clsTrkCharge = nullptr;
94 
96  TH1F** m_clsTrkSNR = nullptr;
97 
99  TH1F** m_clsTrkChargeL3 = nullptr;
100 
102  TH1F** m_clsTrkSNRL3 = nullptr;
103 
105  TH1F* m_clsTrkChargeU3 = nullptr;
107  TH1F* m_clsTrkChargeV3 = nullptr;
109  TH1F* m_clsTrkChargeU456 = nullptr;
111  TH1F* m_clsTrkChargeV456 = nullptr;
112 
114  TH1F* m_clsTrkSNRU3 = nullptr;
116  TH1F* m_clsTrkSNRV3 = nullptr;
118  TH1F* m_clsTrkSNRU456 = nullptr;
120  TH1F* m_clsTrkSNRV456 = nullptr;
121 
123  TH1F* m_stripMaxBinUAll = nullptr;
125  TH1F* m_stripMaxBinVAll = nullptr;
126 
128  TH2F* m_clsTrkTimeUEvtT0 = nullptr;
130  TH2F* m_clsTrkTimeVEvtT0 = nullptr;
132  TH1F* m_clsTrkTimeU3 = nullptr;
134  TH1F* m_clsTrkTimeV3 = nullptr;
135 
137  TH1F* m_cls3TrkTimeU3 = nullptr;
139  TH1F* m_cls3TrkTimeV3 = nullptr;
140 
142  TH1F* m_cls6TrkTimeU3 = nullptr;
144  TH1F* m_cls6TrkTimeV3 = nullptr;
145 
147  TH1F* m_clsTrkTimeU456 = nullptr;
149  TH1F* m_clsTrkTimeV456 = nullptr;
150 
152  TH1F* m_cls3TrkTimeU456 = nullptr;
154  TH1F* m_cls3TrkTimeV456 = nullptr;
155 
157  TH1F* m_cls6TrkTimeU456 = nullptr;
159  TH1F* m_cls6TrkTimeV456 = nullptr;
160 
162  std::map<std::pair<int, int>, int> m_ladderMap;
163 
164  };
165 
167 }
Storage element for the eventwise T0 estimation.
Definition: EventT0.h:29
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
SVD DQM Module for Clusters related to Tracks.
TH1F * m_cls6TrkTimeV456
v Time of clusters related to tracks for layer 4,5,6 sensors for 6 samples
TH1F ** m_clsTrkSNRL3
SNR of clusters related to tracks per layer 3.
TH1F * m_cls6TrkTimeU3
u Time of clusters related to tracks for layer 3 sensors for 3 samples
std::string m_svdClustersName
SVDClusters data object name.
void initialize() override final
Module function initialize.
TH1F ** m_clsTrkChargeL3
charge of clusters related to tracks per layer 3
TH1F ** m_clsTrkCharge
charge of clusters related to tracks per ladder
TH1F * m_clsTrkSNRV456
v SNR of clusters related to tracks for layer 4,5,6 sensors
std::string m_svdShaperDigitsName
SVDShaperDigits data object name.
StoreObjPtr< EventT0 > m_eventT0
EventT0 data object.
bool m_skipRejectedEvents
if true skip events rejected by HLT (default)
TH1F * m_cls3TrkTimeU456
u Time of clusters related to tracks for layer 4,5,6 sensors for 3 samples
TH1F * m_clsTrkTimeU3
u Time of clusters related to tracks for layer 3 sensors
TH1F * m_cls3TrkTimeV456
v Time of clusters related to tracks for layer 4,5,6 sensors for 3 samples
StoreObjPtr< SVDEventInfo > m_svdEventInfo
SVDEventInfo data object.
TH1F * m_clsTrkChargeV3
v charge of clusters related to tracks for layer 3 sensors
void defineHisto() override final
Contains the Histogram definitions
TH1F * m_clsTrkSNRU456
u SNR of clusters related to tracks for layer 4,5,6 sensors
TH1F * m_clsTrkSNRV3
v SNR of clusters related to tracks for layer 3 sensors
TH1F * m_cls6TrkTimeU456
u Time of clusters related to tracks for layer 4,5,6 sensors for 6 samples
TH2F * m_clsTrkTimeUEvtT0
u Time of clusters related to tracks vs EventT0
void terminate() override final
Module function terminate.
TH1F * m_clsTrkChargeU456
u charge of clusters related to tracks for layer 4,5,6 sensors
TH1F * m_cls3TrkTimeV3
v Time of clusters related to tracks for layer 3 sensors for 3 sampes
void event() override final
Module function event.
TH1F * m_clsTrkTimeV3
v Time of clusters related to tracks for layer 3 sensors
SVDDQMClustersOnTrackModule & operator=(const SVDDQMClustersOnTrackModule &)=delete
Operator = (disabled)
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
StoreArray< Track > m_tracks
StoreArray of the Tracks.
int m_tb
choose one trigger bin, or none if the value is -1
TH1F * m_stripMaxBinUAll
u MaxBin of strips related to tracks for all sensors
TList * m_histoList
list of cumulative histograms
std::string m_svdRecoDigitsName
SVDRecoDigits data object name.
std::string m_svdEventInfoName
SVDEventInfo data object name.
TH1F * m_cls6TrkTimeV3
v Time of clusters related to tracks for layer 3 sensors for 3 sampes
TH1F * m_cls3TrkTimeU3
u Time of clusters related to tracks for layer 3 sensors for 3 samples
void beginRun() override final
Module function beginRun.
TH1F * m_clsTrkSNRU3
u SNR of clusters related to tracks for layer 3 sensors
bool m_desynchSVDTime
if TRUE: svdTime back in SVD time reference
TH1F * m_clsTrkChargeU3
u charge of clusters related to tracks for layer 3 sensors
TH1F * m_stripMaxBinVAll
v MaxBin of strips related to tracks for all sensors
SVDDQMClustersOnTrackModule(const SVDDQMClustersOnTrackModule &)=delete
Copy constructor (disabled)
StoreObjPtr< SoftwareTriggerResult > m_resultStoreObjectPointer
Store Object for reading the trigger decision.
TH1F * m_clsTrkTimeU456
u Time of clusters related to tracks for layer 4,5,6 sensors
TH1F ** m_clsTrkSNR
SNR of clusters related to tracks per ladder.
TH2F * m_clsTrkTimeVEvtT0
v Time of clusters related to tracks vs EventT0
TH1F * m_clsTrkChargeV456
v charge of clusters related to tracks for layer 4,5,6 sensors
TH1F * m_clsTrkTimeV456
v Time of clusters related to tracks for layer 4,5,6 sensors
std::map< std::pair< int, int >, int > m_ladderMap
map of ladder index
Stores SVDModeByte object with Trigger time, DAQ mode, Run type & Event type! Also - the information ...
Definition: SVDEventInfo.h:31
Dataobject to store the results of the cut calculations performed by the SoftwareTriggerModule.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Class that bundles various TrackFitResults.
Definition: Track.h:25
Abstract base class for different kinds of events.