Belle II Software  release-05-01-25
TRGCDCT3DDQMModule.cc
1 //---------------------------------------------------------------
2 // $Id$
3 //---------------------------------------------------------------
4 // Filename : TRGCDCT3DModule.cc
5 // Section : TRG CDCT3D
6 // Owner :
7 // Email :
8 //---------------------------------------------------------------
9 // Description : A trigger module for TRG CDCT3D
10 //---------------------------------------------------------------
11 // 1.00 : 2017/05/08 : First version
12 //---------------------------------------------------------------
13 #include <trg/cdc/modules/trgcdct3dDQM/TRGCDCT3DDQMModule.h>
14 
15 #include <framework/datastore/StoreObjPtr.h>
16 #include <framework/datastore/StoreArray.h>
17 
18 #include <TDirectory.h>
19 #include <TPostScript.h>
20 #include <TCanvas.h>
21 #include <iostream>
22 
23 using namespace std;
24 using namespace Belle2;
25 
26 REG_MODULE(TRGCDCT3DDQM);
27 
28 
29 TRGCDCT3DDQMModule::TRGCDCT3DDQMModule() : HistoModule()
30 {
31 
32  setDescription("DQM for CDCT3D Trigger system");
34 
35  addParam("T3DMOD", m_T3DMOD,
36  "T3D module number",
37  0);
38  addParam("generatePostscript", m_generatePostscript,
39  "Genarete postscript file or not",
40  false);
41  addParam("postScriptName", m_postScriptName,
42  "postscript file name",
43  string("cdct3ddqm.ps"));
44 
45 
46 }
47 
49 {
50  oldDir = gDirectory;
51  dirDQM = NULL;
52  //dirDQM = oldDir->mkdir("TRGCDCT3D");
53  if (!oldDir->Get("TRGCDCT3D"))dirDQM = oldDir->mkdir("TRGCDCT3D");
54  else dirDQM = (TDirectory*)oldDir->Get("TRGCDCT3D");
55  dirDQM->cd();
56  //dz distribution
57  h_dz = new TH1D(Form("hdz_mod%d", m_T3DMOD), Form("hdz_mod%d", m_T3DMOD), 80, -40, 40);
58  h_dz->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
59  h_dz->GetXaxis()->SetTitle("dz");
60  //phi distribution
61  h_phi = new TH1D(Form("hphi_mod%d", m_T3DMOD), Form("hphi_mod%d", m_T3DMOD), 100, -6, 6);
62  h_phi->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
63  h_phi->GetXaxis()->SetTitle("phi");
64  //tanlambda distribution
65  h_tanlambda = new TH1D(Form("htanlambda_mod%d", m_T3DMOD), Form("htanlambda_mod%d", m_T3DMOD), 100, -2.5, 2.5);
66  h_tanlambda->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
67  h_tanlambda->GetXaxis()->SetTitle("tanlambda");
68  //pt distribution
69  h_pt = new TH1D(Form("hpt_mod%d", m_T3DMOD), Form("hpt_mod%d", m_T3DMOD), 100, 0, 3);
70  h_pt->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
71  h_pt->GetXaxis()->SetTitle("pt");
72 
73  //2D phi distribution
74  h_phi_2D = new TH1D(Form("hphi_2D_mod%d", m_T3DMOD), Form("hphi_2D_mod%d", m_T3DMOD), 100, -6, 6);
75  h_phi_2D->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
76  h_phi_2D->GetXaxis()->SetTitle("phi_2D");
77  //2D pt distribution
78  h_pt_2D = new TH1D(Form("hpt_2D_mod%d", m_T3DMOD), Form("hpt_2D_mod%d", m_T3DMOD), 100, 0, 3);
79  h_pt_2D->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
80  h_pt_2D->GetXaxis()->SetTitle("pt_2D");
81  //TSF1 ID
82  h_ID_TSF1 = new TH1D(Form("hID_TSF1_mod%d", m_T3DMOD), Form("hID_TSF1_mod%d", m_T3DMOD), 160, 160, 319);
83  h_ID_TSF1->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
84  h_ID_TSF1->GetXaxis()->SetTitle("ID_TSF1");
85  //TSF3 ID
86  h_ID_TSF3 = new TH1D(Form("hID_TSF3_mod%d", m_T3DMOD), Form("hID_TSF3_mod%d", m_T3DMOD), 224, 512, 735);
87  h_ID_TSF3->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
88  h_ID_TSF3->GetXaxis()->SetTitle("ID_TSF3");
89  //TSF5 ID
90  h_ID_TSF5 = new TH1D(Form("hID_TSF5_mod%d", m_T3DMOD), Form("hID_TSF5_mod%d", m_T3DMOD), 288, 992, 1279);
91  h_ID_TSF5->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
92  h_ID_TSF5->GetXaxis()->SetTitle("ID_TSF5");
93  //TSF7 ID
94  h_ID_TSF7 = new TH1D(Form("hID_TSF7_mod%d", m_T3DMOD), Form("hID_TSF7_mod%d", m_T3DMOD), 352, 1600, 1951);
95  h_ID_TSF7->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
96  h_ID_TSF7->GetXaxis()->SetTitle("ID_TSF7");
97  //TSF1 rt
98  h_rt_TSF1 = new TH1D(Form("hrt_TSF1_mod%d", m_T3DMOD), Form("hrt_TSF1_mod%d", m_T3DMOD), 512, 0, 511);
99  h_rt_TSF1->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
100  h_rt_TSF1->GetXaxis()->SetTitle("rt_TSF1");
101  //TSF3 rt
102  h_rt_TSF3 = new TH1D(Form("hrt_TSF3_mod%d", m_T3DMOD), Form("hrt_TSF3_mod%d", m_T3DMOD), 512, 0, 511);
103  h_rt_TSF3->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
104  h_rt_TSF3->GetXaxis()->SetTitle("rt_TSF3");
105  //TSF5 rt
106  h_rt_TSF5 = new TH1D(Form("hrt_TSF5_mod%d", m_T3DMOD), Form("hrt_TSF5_mod%d", m_T3DMOD), 512, 0, 511);
107  h_rt_TSF5->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
108  h_rt_TSF5->GetXaxis()->SetTitle("rt_TSF5");
109  //TSF7 rt
110  h_rt_TSF7 = new TH1D(Form("hrt_TSF7_mod%d", m_T3DMOD), Form("hrt_TSF7_mod%d", m_T3DMOD), 512, 0, 511);
111  h_rt_TSF7->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
112  h_rt_TSF7->GetXaxis()->SetTitle("rt_TSF7");
113  //TSF1 validity
114  h_validity_TSF1 = new TH1D(Form("hvalidity_TSF1_mod%d", m_T3DMOD), Form("hvalidity_TSF1_mod%d", m_T3DMOD), 16, 0, 15);
115  h_validity_TSF1->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
116  h_validity_TSF1->GetXaxis()->SetTitle("validity_TSF1");
117  //TSF3 validity
118  h_validity_TSF3 = new TH1D(Form("hvalidity_TSF3_mod%d", m_T3DMOD), Form("hvalidity_TSF3_mod%d", m_T3DMOD), 16, 0, 15);
119  h_validity_TSF3->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
120  h_validity_TSF3->GetXaxis()->SetTitle("validity_TSF3");
121  //TSF5 validity
122  h_validity_TSF5 = new TH1D(Form("hvalidity_TSF5_mod%d", m_T3DMOD), Form("hvalidity_TSF5_mod%d", m_T3DMOD), 16, 0, 15);
123  h_validity_TSF5->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
124  h_validity_TSF5->GetXaxis()->SetTitle("validity_TSF5");
125  //TSF7 validity
126  h_validity_TSF7 = new TH1D(Form("hvalidity_TSF7_mod%d", m_T3DMOD), Form("hvalidity_TSF7_mod%d", m_T3DMOD), 16, 0, 15);
127  h_validity_TSF7->SetTitle(Form("Exp%d Run%d 3Dmodule%d", _exp, _run, m_T3DMOD));
128  h_validity_TSF7->GetXaxis()->SetTitle("validity_TSF7");
129 
130  oldDir->cd();
131 }
132 
134 {
135 
136  dirDQM->cd();
137 
138  h_dz->Reset();
139  h_phi->Reset();
140  h_tanlambda->Reset();
141  h_pt->Reset();
142  h_phi_2D->Reset();
143  h_pt_2D->Reset();
144  h_ID_TSF1->Reset();
145 
146  oldDir->cd();
147 }
148 
150 {
151 
153  _exp = bevt->getExperiment();
154  _run = bevt->getRun();
155 
156  // calls back the defineHisto() function, but the HistoManager module has to be in the path
157  REG_HISTOGRAM
158 
159  char c_name_3D[100];
160  sprintf(c_name_3D, "FirmTRGCDC3DFitterTracks%d", m_T3DMOD);
161  char c_name_2D[100];
162  sprintf(c_name_2D, "FirmTRGCDC2DFinderTracks%d", m_T3DMOD);
163  entAry.isRequired(c_name_3D);
164  entAry_2D.isRequired(c_name_2D);
165  char c_name_TSF[100];
166  sprintf(c_name_TSF, "FirmCDCTriggerSegmentHits%d", m_T3DMOD);
167  entAry_TSF.isRequired(c_name_TSF);
168  if (!entAry || !entAry.getEntries()) return;
169  if (!entAry_2D || !entAry_2D.getEntries()) return;
170  if (!entAry_TSF || !entAry_TSF.getEntries()) return;
171 
172 }
173 
175 {
176  dirDQM->cd();
177 
178  //Draw and save histograms
179  if (m_generatePostscript) {
180  //gStyle->SetOptStat(0);
181  TCanvas c1("c1", "", 0, 0, 500, 300);
182  c1.cd();
183 
184  TPostScript* ps_dz = new TPostScript((m_postScriptName + ".dz_T3DModule" + to_string(m_T3DMOD) + ".ps").c_str(), 112);
185  h_dz->Draw();
186  c1.Update();
187  ps_dz->Close();
188 
189  TPostScript* ps_phi = new TPostScript((m_postScriptName + ".phi_T3DModule" + to_string(m_T3DMOD) + ".ps").c_str(), 112);
190  h_phi->Draw();
191  c1.Update();
192  ps_phi->Close();
193 
194  TPostScript* ps_tanlambda = new TPostScript((m_postScriptName + ".tanlambda_T3DModule" + to_string(m_T3DMOD) + ".ps").c_str(), 112);
195  h_tanlambda->Draw();
196  c1.Update();
197  ps_tanlambda->Close();
198 
199  TPostScript* ps_pt = new TPostScript((m_postScriptName + ".pt_T3DModule" + to_string(m_T3DMOD) + ".ps").c_str(), 112);
200  h_pt->Draw();
201  c1.Update();
202  ps_pt->Close();
203 
204  TPostScript* ps_phi_2D = new TPostScript((m_postScriptName + ".phi_2D_T3DModule" + to_string(m_T3DMOD) + ".ps").c_str(), 112);
205  h_phi_2D->Draw();
206  c1.Update();
207  ps_phi_2D->Close();
208 
209  TPostScript* ps_pt_2D = new TPostScript((m_postScriptName + ".pt_2D_T3DModule" + to_string(m_T3DMOD) + ".ps").c_str(), 112);
210  h_pt_2D->Draw();
211  c1.Update();
212  ps_pt_2D->Close();
213  }
214 
215  oldDir->cd();
216 }
217 
219 {
220 
221  dirDQM->cd();
222 
223  //Fill 3D histo
224  if (!(!entAry || !entAry.getEntries())) {
225  for (int i = 0; i < entAry.getEntries(); i++) {
226  h_dz->Fill(entAry[i]->getZ0());
227  h_phi->Fill(entAry[i]->getPhi0());
228  h_tanlambda->Fill(entAry[i]->getTanLambda());
229  if (entAry[i]->getOmega() != 0)
230  h_pt->Fill(fabs(1. / entAry[i]->getOmega() * 0.3 * 1.5 * 0.01));
231  }
232  }
233  //Fill 2D histo
234  if (!(!entAry_2D || !entAry_2D.getEntries())) {
235  for (int i = 0; i < entAry_2D.getEntries(); i++) {
236  h_phi_2D->Fill(entAry_2D[i]->getPhi0());
237  if (entAry_2D[i]->getOmega() != 0)
238  h_pt_2D->Fill(fabs(1. / entAry_2D[i]->getOmega() * 0.3 * 1.5 * 0.01));
239  }
240  }
241  //Fill stereo TSF histo
242  if (!(!entAry_TSF || !entAry_TSF.getEntries())) {
243  for (int i = 0; i < entAry_TSF.getEntries(); i++) {
244  int sl = entAry_TSF[i]->getISuperLayer();
245  int id = entAry_TSF[i]->getSegmentID();
246  int pr = entAry_TSF[i]->getPriorityPosition();
247  int lr = entAry_TSF[i]->getLeftRight();
248  int rt = entAry_TSF[i]->priorityTime();
249  int validity = lr * 4 + pr;
250  if (lr == 0 || pr == 0) continue;
251  if (sl == 1) {
252  h_ID_TSF1->Fill(id);
253  h_rt_TSF1->Fill(rt);
254  h_validity_TSF1->Fill(validity);
255  } else if (sl == 3) {
256  h_ID_TSF3->Fill(id);
257  h_rt_TSF3->Fill(rt);
258  h_validity_TSF3->Fill(validity);
259  } else if (sl == 5) {
260  h_ID_TSF5->Fill(id);
261  h_rt_TSF5->Fill(rt);
262  h_validity_TSF5->Fill(validity);
263  } else if (sl == 7) {
264  h_ID_TSF7->Fill(id);
265  h_rt_TSF7->Fill(rt);
266  h_validity_TSF7->Fill(validity);
267  }
268 
269 
270  }
271  }
272 
273  oldDir->cd();
274 
275 }
276 
277 
Belle2::TRGCDCT3DDQMModule::entAry_2D
StoreArray< CDCTriggerTrack > entAry_2D
2D data store
Definition: TRGCDCT3DDQMModule.h:104
Belle2::TRGCDCT3DDQMModule::h_ID_TSF1
TH1D * h_ID_TSF1
TSF1 ID of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:55
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::Module::c_ParallelProcessingCertified
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:82
Belle2::TRGCDCT3DDQMModule::h_dz
TH1D * h_dz
dz of T3D in each module
Definition: TRGCDCT3DDQMModule.h:43
Belle2::TRGCDCT3DDQMModule::h_pt_2D
TH1D * h_pt_2D
pt (from 2D) of T3D in each module
Definition: TRGCDCT3DDQMModule.h:53
Belle2::TRGCDCT3DDQMModule::h_ID_TSF3
TH1D * h_ID_TSF3
TSF3 ID of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:57
Belle2::TRGCDCT3DDQMModule::h_ID_TSF7
TH1D * h_ID_TSF7
TSF7 ID of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:61
Belle2::TRGCDCT3DDQMModule::h_rt_TSF7
TH1D * h_rt_TSF7
TSF7 priority time of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:69
Belle2::TRGCDCT3DDQMModule::m_generatePostscript
bool m_generatePostscript
flag to save ps file
Definition: TRGCDCT3DDQMModule.h:86
Belle2::TRGCDCT3DDQMModule::h_validity_TSF3
TH1D * h_validity_TSF3
TSF3 validity of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:73
Belle2::TRGCDCT3DDQMModule::h_validity_TSF5
TH1D * h_validity_TSF5
TSF5 validity of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:75
Belle2::TRGCDCT3DDQMModule::h_rt_TSF1
TH1D * h_rt_TSF1
TSF1 priority time of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:63
Belle2::TRGCDCT3DDQMModule::_exp
unsigned _exp
experiment number
Definition: TRGCDCT3DDQMModule.h:92
Belle2::Module::setPropertyFlags
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:210
Belle2::TRGCDCT3DDQMModule::dirDQM
TDirectory * dirDQM
TDirectories for DQM histograms.
Definition: TRGCDCT3DDQMModule.h:83
Belle2::TRGCDCT3DDQMModule::m_T3DMOD
int m_T3DMOD
T3D module number.
Definition: TRGCDCT3DDQMModule.h:98
Belle2::TRGCDCT3DDQMModule::oldDir
TDirectory * oldDir
TDirectories for DQM histograms.
Definition: TRGCDCT3DDQMModule.h:81
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::TRGCDCT3DDQMModule::h_phi_2D
TH1D * h_phi_2D
phi (from 2D) of T3D in each module
Definition: TRGCDCT3DDQMModule.h:51
Belle2::TRGCDCT3DDQMModule::h_validity_TSF7
TH1D * h_validity_TSF7
TSF7 validity of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:77
Belle2::TRGCDCT3DDQMModule::entAry_TSF
StoreArray< CDCTriggerSegmentHit > entAry_TSF
Stereo TSF data store.
Definition: TRGCDCT3DDQMModule.h:107
Belle2::TRGCDCT3DDQMModule::m_postScriptName
std::string m_postScriptName
name of ps file
Definition: TRGCDCT3DDQMModule.h:89
Belle2::TRGCDCT3DDQMModule::endRun
virtual void endRun() override
End Run.
Definition: TRGCDCT3DDQMModule.cc:174
Belle2::TRGCDCT3DDQMModule::initialize
virtual void initialize() override
initialize
Definition: TRGCDCT3DDQMModule.cc:149
Belle2::Module::addParam
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:562
Belle2::TRGCDCT3DDQMModule::h_ID_TSF5
TH1D * h_ID_TSF5
TSF5 ID of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:59
Belle2::TRGCDCT3DDQMModule::entAry
StoreArray< CDCTriggerTrack > entAry
3D data store
Definition: TRGCDCT3DDQMModule.h:101
Belle2::TRGCDCT3DDQMModule::h_rt_TSF5
TH1D * h_rt_TSF5
TSF5 priority time of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:67
Belle2::TRGCDCT3DDQMModule::h_pt
TH1D * h_pt
pt of T3D in each module
Definition: TRGCDCT3DDQMModule.h:49
Belle2::TRGCDCT3DDQMModule::beginRun
virtual void beginRun() override
begin Run
Definition: TRGCDCT3DDQMModule.cc:133
Belle2::TRGCDCT3DDQMModule::_run
unsigned _run
run number
Definition: TRGCDCT3DDQMModule.h:95
Belle2::TRGCDCT3DDQMModule::h_phi
TH1D * h_phi
phi of T3D in each module
Definition: TRGCDCT3DDQMModule.h:45
Belle2::TRGCDCT3DDQMModule::event
virtual void event() override
Event.
Definition: TRGCDCT3DDQMModule.cc:218
Belle2::TRGCDCT3DDQMModule::h_rt_TSF3
TH1D * h_rt_TSF3
TSF3 priority time of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:65
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Belle2::TRGCDCT3DDQMModule::defineHisto
virtual void defineHisto() override
Define Histogram.
Definition: TRGCDCT3DDQMModule.cc:48
Belle2::TRGCDCT3DDQMModule::h_tanlambda
TH1D * h_tanlambda
tanlambda of T3D in each module
Definition: TRGCDCT3DDQMModule.h:47
Belle2::TRGCDCT3DDQMModule::h_validity_TSF1
TH1D * h_validity_TSF1
TSF1 validity of T3D in each module.
Definition: TRGCDCT3DDQMModule.h:71