Belle II Software  release-08-02-00
TRGEFFDQMModule.cc
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 // Own include
10 #include <trg/gdl/modules/trggdlDQM/TRGEFFDQMModule.h>
11 
12 // Dataobject classes
13 #include <TF1.h>
14 #include <Math/Vector3D.h>
15 #include <TVector3.h>
16 #include <TDirectory.h>
17 #include <fstream>
18 #include <math.h>
19 #include <mdst/dataobjects/HitPatternCDC.h>
20 
21 using namespace std;
22 using namespace Belle2;
23 
24 //-----------------------------------------------------------------
25 // Register module
26 //-----------------------------------------------------------------
27 
28 REG_MODULE(TRGEFFDQM);
29 
30 TRGEFFDQMModule::TRGEFFDQMModule() : HistoModule()
31 {
32  // set module description (e.g. insert text)
33  setDescription("Make kinematics dependent efficiency plot");
35 }
36 
38 {
39 }
40 
42 {
43  TDirectory* oldDir = gDirectory;
44  oldDir->mkdir("TRGEFF");
45  oldDir->cd("TRGEFF");
46 
47  m_hPhi = new TH1F("hPhi", "", 36, -180.0, 180.0);
48  m_hPhi_psnecl = new TH1F("hPhi_psnecl", "", 36, -180.0, 180.0);
49  m_hPhi_psnecl_ftdf = new TH1F("hPhi_psnecl_ftdf", "", 36, -180.0, 180.0);
50 
51  m_hPt = new TH1F("hPt", "", 50, 0, 5);
52  m_hPt_psnecl = new TH1F("hPt_psnecl", "", 50, 0, 5);
53  m_hPt_psnecl_ftdf = new TH1F("hPt_psnecl_ftdf", "", 50, 0, 5);
54 
55  m_nobha_hPt = new TH1F("nobha_hPt", "", 50, 0, 5);
56  m_nobha_hPt_psnecl = new TH1F("nobha_hPt_psnecl", "", 50, 0, 5);
57  m_nobha_hPt_psnecl_ftdf = new TH1F("nobha_hPt_psnecl_ftdf", "", 50, 0, 5);
58 
59  m_hP3_y = new TH1F("hP3_y", "", 50, 0, 5);
60  m_hP3_y_psnecl = new TH1F("hP3_y_psnecl", "", 50, 0, 5);
61  m_hP3_y_psnecl_ftdf = new TH1F("hP3_y_psnecl_ftdf", "", 50, 0, 5);
62 
63  m_hP3_z = new TH1F("hP3_z", "", 50, 0, 5);
64  m_hP3_z_psnecl = new TH1F("hP3_z_psnecl", "", 50, 0, 5);
65  m_hP3_z_psnecl_ftdf = new TH1F("hP3_z_psnecl_ftdf", "", 50, 0, 5);
66 
67  m_nobha_hP3_y = new TH1F("nobha_hP3_y", "", 50, 0, 5);
68  m_nobha_hP3_y_psnecl = new TH1F("nobha_hP3_y_psnecl", "", 50, 0, 5);
69  m_nobha_hP3_y_psnecl_ftdf = new TH1F("nobha_hP3_y_psnecl_ftdf", "", 50, 0, 5);
70 
71  m_nobha_hP3_z = new TH1F("nobha_hP3_z", "", 50, 0, 5);
72  m_nobha_hP3_z_psnecl = new TH1F("nobha_hP3_z_psnecl", "", 50, 0, 5);
73  m_nobha_hP3_z_psnecl_ftdf = new TH1F("nobha_hP3_z_psnecl_ftdf", "", 50, 0, 5);
74 
75  m_fyo_dphi = new TH1F("fyo_dphi", "", 18, 0., 180.);
76  m_fyo_dphi_psnecl = new TH1F("fyo_dphi_psnecl", "", 18, 0., 180.);
77  m_fyo_dphi_psnecl_ftdf = new TH1F("fyo_dphi_psnecl_ftdf", "", 18, 0., 180.);
78 
79  m_nobha_fyo_dphi = new TH1F("nobha_fyo_dphi", "", 18, 0., 180.);
80  m_nobha_fyo_dphi_psnecl = new TH1F("nobha_fyo_dphi_psnecl", "", 18, 0., 180.);
81  m_nobha_fyo_dphi_psnecl_ftdf = new TH1F("nobha_fyo_dphi_psnecl_ftdf", "", 18, 0., 180.);
82 
83  m_stt_phi = new TH1F("stt_phi", "", 36, -180., 180.);
84  m_stt_phi_psnecl = new TH1F("stt_phi_psnecl", "", 36, -180., 180.);
85  m_stt_phi_psnecl_ftdf = new TH1F("stt_phi_psnecl_ftdf", "", 36, -180., 180.);
86 
87  m_stt_P3 = new TH1F("stt_P3", "", 50, 0, 5);
88  m_stt_P3_psnecl = new TH1F("stt_P3_psnecl", "", 50, 0, 5);
89  m_stt_P3_psnecl_ftdf = new TH1F("stt_P3_psnecl_ftdf", "", 50, 0, 5);
90 
91  m_stt_theta = new TH1F("stt_theta", "", 18, 0, 180);
92  m_stt_theta_psnecl = new TH1F("stt_theta_psnecl", "", 18, 0, 180);
93  m_stt_theta_psnecl_ftdf = new TH1F("stt_theta_psnecl_ftdf", "", 18, 0, 180);
94 
95  m_nobha_stt_phi = new TH1F("nobha_stt_phi", "", 36, -180., 180.);
96  m_nobha_stt_phi_psnecl = new TH1F("nobha_stt_phi_psnecl", "", 36, -180., 180.);
97  m_nobha_stt_phi_psnecl_ftdf = new TH1F("nobha_stt_phi_psnecl_ftdf", "", 36, -180., 180.);
98 
99  m_nobha_stt_P3 = new TH1F("nobha_stt_P3", "", 50, 0, 5);
100  m_nobha_stt_P3_psnecl = new TH1F("nobha_stt_P3_psnecl", "", 50, 0, 5);
101  m_nobha_stt_P3_psnecl_ftdf = new TH1F("nobha_stt_P3_psnecl_ftdf", "", 50, 0, 5);
102 
103  m_nobha_stt_theta = new TH1F("nobha_stt_theta", "", 18, 0, 180);
104  m_nobha_stt_theta_psnecl = new TH1F("nobha_stt_theta_psnecl", "", 18, 0, 180);
105  m_nobha_stt_theta_psnecl_ftdf = new TH1F("nobha_stt_theta_psnecl_ftdf", "", 18, 0, 180);
106 
107  m_hie_E = new TH1F("hie_E", "", 60, 0, 12);
108  m_hie_E_psnecl = new TH1F("hie_E_psnecl", "", 60, 0, 12);
109  m_hie_E_psnecl_ftdf = new TH1F("hie_E_psnecl_ftdf", "", 60, 0, 12);
110 
111  m_nobha_hie_E = new TH1F("nobha_hie_E", "", 60, 0, 12);
112  m_nobha_hie_E_psnecl = new TH1F("nobha_hie_E_psnecl", "", 60, 0, 12);
113  m_nobha_hie_E_psnecl_ftdf = new TH1F("nobha_hie_E_psnecl_ftdf", "", 60, 0, 12);
114 
115  m_ecltiming_E = new TH1F("ecltiming_E", "", 60, 0, 12);
116  m_ecltiming_E_psnecl = new TH1F("ecltiming_E_psnecl", "", 60, 0, 12);
117  m_ecltiming_E_psnecl_ftdf = new TH1F("ecltiming_E_psnecl_ftdf", "", 60, 0, 12);
118 
119  m_ecltiming_theta = new TH1F("ecltiming_theta", "", 18, 0, 180);
120  m_ecltiming_theta_psnecl = new TH1F("ecltiming_theta_psnecl", "", 18, 0, 180);
121  m_ecltiming_theta_psnecl_ftdf = new TH1F("ecltiming_theta_psnecl_ftdf", "", 18, 0, 180);
122 
123  m_ecltiming_phi = new TH1F("ecltiming_phi", "", 36, -180., 180.);
124  m_ecltiming_phi_psnecl = new TH1F("ecltiming_phi_psnecl", "", 36, -180., 180.);
125  m_ecltiming_phi_psnecl_ftdf = new TH1F("ecltiming_phi_psnecl_ftdf", "", 36, -180., 180.);
126 
127  m_klmhit_phi = new TH1F("klmhit_phi", "", 18, -180., 180.);
128  m_klmhit_phi_psnecl = new TH1F("klmhit_phi_psnecl", "", 18, -180., 180.);
129  m_klmhit_phi_psnecl_ftdf = new TH1F("klmhit_phi_psnecl_ftdf", "", 18, -180., 180.);
130 
131  m_klmhit_theta = new TH1F("klmhit_theta", "", 18, 0, 180);
132  m_klmhit_theta_psnecl = new TH1F("klmhit_theta_psnecl", "", 18, 0, 180);
133  m_klmhit_theta_psnecl_ftdf = new TH1F("klmhit_theta_psnecl_ftdf", "", 18, 0, 180);
134 
135  m_eklmhit_phi = new TH1F("eklmhit_phi", "", 18, -180., 180.);
136  m_eklmhit_phi_psnecl = new TH1F("eklmhit_phi_psnecl", "", 18, -180., 180.);
137  m_eklmhit_phi_psnecl_ftdf = new TH1F("eklmhit_phi_psnecl_ftdf", "", 18, -180., 180.);
138 
139  m_eklmhit_theta = new TH1F("eklmhit_theta", "", 18, 0, 180);
140  m_eklmhit_theta_psnecl = new TH1F("eklmhit_theta_psnecl", "", 18, 0, 180);
141  m_eklmhit_theta_psnecl_ftdf = new TH1F("eklmhit_theta_psnecl_ftdf", "", 18, 0, 180);
142 
143 
144 
145 
146 
147 
148 
149  oldDir->cd();
150 }
151 
153 {
154  REG_HISTOGRAM
155 
156  if (!m_Tracks.isOptional()) {
157  B2WARNING("Missing Tracks array");
158  return;
159  }
160  if (!m_ECLClusters.isOptional()) {
161  B2WARNING("Missing ECLCLuster array");
162  return;
163  }
164  if (!m_KLMClusters.isOptional()) {
165  B2WARNING("Missing KLMCLuster array");
166  return;
167  }
168  if (!m_trgSummary.isOptional()) {
169  B2WARNING("Missing TRGSummary");
170  return;
171  }
172 }
173 
175 {
176  if (!m_RecoTracks.isOptional()) {
177  B2DEBUG(22, "Missing recoTracks array in beginRun() ");
178  return;
179  }
180 
181  m_hPhi->Reset();
182 }
183 
185 {
186  // bool debug = 1; // if(debug)cout<<"line "<<__LINE__<<endl;
187 
188  if (!m_trgSummary.isValid()) {
189  B2WARNING("TRGSummary object not available but require to estimate trg efficiency");
190  return;
191  }
192 
193  if (!m_TrgResult.isValid()) {
194  B2WARNING("SoftwareTriggerResult object not available but require to select bhabha/mumu/hadron events skim");
195  return;
196  }
197 
198  const std::map<std::string, int>& fresults = m_TrgResult->getResults();
199  if ((fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end())
200  || (fresults.find("software_trigger_cut&skim&accept_hadron") == fresults.end())) {
201  B2WARNING("TRGEFFDQMModule: Can't find required bhabha or mumu or hadron trigger identifier");
202  return;
203  }
204 
205  const bool IsBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") == SoftwareTriggerCutResult::c_accept);
206  const bool IsHadron = (m_TrgResult->getResult("software_trigger_cut&skim&accept_hadron") == SoftwareTriggerCutResult::c_accept);
207 
208 
209  /*///////////////////////////////////////////////////////////////////
212 
213  // calculate the total energy
214  double E_ecl_all = 0; // the ECL total energy
215  double E_ecl_hie =
216  0; // the ECL total energy in the thetaID range 2<=ThetaID<=15 (corresponds to 22.49<=theta<=126.80) for ehigh bit
217  for (const auto& test_b2eclcluster : m_ECLClusters) {
218  if (!(test_b2eclcluster.hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons))) continue;
219  double energy = test_b2eclcluster.getEnergyRaw();
220  double theta = test_b2eclcluster.getTheta() / Unit::deg;
221 
222  if (energy < 0.1) continue;
223 
224  E_ecl_all = E_ecl_all + energy;
225  if (theta >= 22.49 && theta <= 126.8) {
226  E_ecl_hie = E_ecl_hie + energy;
227  }
228  }
229 
230  bool trg_hie_psncdc = 0; // for ECL energy trigger, for hie
231  bool trg_hie_Eecl = 0; // for ECL energy trigger, for hie
232  bool trg_ecltiming_psncdc = 0; // for ECL energy trigger, for ecltiming
233  bool trg_ecltiming_Eecl = 0; // for ECL energy trigger. for ecltiming
234  bool trg_nobha_hie_Eecl = 0;
235 
236  for (const auto& b2eclcluster : m_ECLClusters) {
237  if (!(b2eclcluster.hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons))) continue;
238  double phi = b2eclcluster.getPhi() / Unit::deg;
239  double theta = b2eclcluster.getTheta() / Unit::deg;
240  double energy = b2eclcluster.getEnergyRaw();
241 
242  if (energy < 0.1) continue;
243 
244  bool trg_psncdc = m_trgSummary->testPsnm("ffy") || m_trgSummary->testPsnm("fyo") || m_trgSummary->testPsnm("stt");
245  bool trg_hie = m_trgSummary->testFtdl("hie");
246  bool trg_hie_nobha = m_trgSummary->testInput("ehigh"); // remove the bha_veto for hie bit
247  bool trg_ecltiming = m_trgSummary->testFtdl("ecltiming");
248 
249  m_ecltiming_theta->Fill(theta);
250  m_ecltiming_phi->Fill(phi);
251 
252  if (trg_psncdc) {
253  trg_hie_psncdc = 1;
254  trg_ecltiming_psncdc = 1;
255 
256  m_ecltiming_theta_psnecl->Fill(theta);
257  m_ecltiming_phi_psnecl->Fill(phi);
258  }
259 
260  if (trg_psncdc && trg_hie) {
261  trg_hie_Eecl = 1;
262  }
263  if (trg_psncdc && trg_hie_nobha) {
264  trg_nobha_hie_Eecl = 1;
265  }
266  if (trg_psncdc && trg_ecltiming) {
267  trg_ecltiming_Eecl = 1;
268 
269  m_ecltiming_theta_psnecl_ftdf->Fill(theta);
270  m_ecltiming_phi_psnecl_ftdf->Fill(phi);
271  }
272  }
273 
274  // sum of the ECL total energy
275  m_hie_E->Fill(E_ecl_hie);
276  m_nobha_hie_E->Fill(E_ecl_hie);
277  if (trg_hie_psncdc) {
278  m_hie_E_psnecl->Fill(E_ecl_hie);
279  m_nobha_hie_E_psnecl->Fill(E_ecl_hie);
280  }
281  if (trg_hie_Eecl) {
282  m_hie_E_psnecl_ftdf->Fill(E_ecl_hie);
283  }
284  if (trg_nobha_hie_Eecl) {
285  m_nobha_hie_E_psnecl_ftdf->Fill(E_ecl_hie);
286  }
287 
288  m_ecltiming_E->Fill(E_ecl_all);
289  if (trg_ecltiming_psncdc) {
290  m_ecltiming_E_psnecl->Fill(E_ecl_all);
291  }
292  if (trg_ecltiming_Eecl) {
293  m_ecltiming_E_psnecl_ftdf->Fill(E_ecl_all);
294  }
295 
296 
297 
298 
299 
300 
301  /*/
304  for (const auto& b2klmcluster : m_KLMClusters) {
305  int nlayer = b2klmcluster.getLayers();
306  // ROOT::Math::XYZVector position = b2klmcluster.getClusterPosition();
307 
308  if (nlayer <= 6)
309  continue;
310 
311  // double p3 = b2klmcluster.getMomentum().R(); // 3-momentum
312  double theta = b2klmcluster.getMomentum().Theta() / Unit::deg;
313  double phiDegree = b2klmcluster.getMomentum().Phi() / Unit::deg;
314 
315  bool trg_KLMecl = m_trgSummary->testPsnm("hie") || m_trgSummary->testPsnm("c4") || m_trgSummary->testPsnm("eclmumu") ||
316  m_trgSummary->testPsnm("ecltaub2b3") || m_trgSummary->testPsnm("hie4") || m_trgSummary->testPsnm("lml1") ||
317  m_trgSummary->testPsnm("lml2") || m_trgSummary->testPsnm("lml6") || m_trgSummary->testPsnm("lml7") ||
318  m_trgSummary->testPsnm("lml8") || m_trgSummary->testPsnm("lml9") || m_trgSummary->testPsnm("lml10");
319 
320  bool trg_klmhit = m_trgSummary->testFtdl("klmhit");
321  bool trg_eklmhit = m_trgSummary->testFtdl("eklmhit");
322 
323  m_klmhit_phi->Fill(phiDegree);
324  m_klmhit_theta->Fill(theta);
325  m_eklmhit_phi->Fill(phiDegree);
326  m_eklmhit_theta->Fill(theta);
327 
328  if (trg_KLMecl) {
329  m_klmhit_phi_psnecl->Fill(phiDegree);
330  m_klmhit_theta_psnecl->Fill(theta);
331  m_eklmhit_phi_psnecl->Fill(phiDegree);
332  m_eklmhit_theta_psnecl->Fill(theta);
333  }
334  if (trg_KLMecl && trg_klmhit) {
335  m_klmhit_phi_psnecl_ftdf->Fill(phiDegree);
336  m_klmhit_theta_psnecl_ftdf->Fill(theta);
337  }
338  if (trg_KLMecl && trg_eklmhit) {
339  m_eklmhit_phi_psnecl_ftdf->Fill(phiDegree);
340  m_eklmhit_theta_psnecl_ftdf->Fill(theta);
341  }
342  }
343 
344 
345 
346  /*///////////////////////////////////////////////////////////////////
349  // int N_tracks = m_Tracks.getEntries(); // number of tracks
350 
351  vector<double> p_stt_P3_psnecl_ftdf, p_stt_P3_psnecl, p_stt_P3, phi_fyo_dphi, phi_fyo_dphi_psnecl, phi_fyo_dphi_psnecl_ftdf ;
352  p_stt_P3_psnecl_ftdf.clear();
353  p_stt_P3_psnecl.clear();
354  p_stt_P3.clear();
355 
356  phi_fyo_dphi.clear();
357  phi_fyo_dphi_psnecl.clear();
358  phi_fyo_dphi_psnecl_ftdf.clear();
359 
360  vector<double> p_nobha_stt_P3_psnecl_ftdf, p_nobha_stt_P3_psnecl, p_nobha_stt_P3, phi_nobha_fyo_dphi, phi_nobha_fyo_dphi_psnecl,
361  phi_nobha_fyo_dphi_psnecl_ftdf ;
362  p_nobha_stt_P3_psnecl_ftdf.clear();
363  p_nobha_stt_P3_psnecl.clear();
364  p_nobha_stt_P3.clear();
365 
366  phi_nobha_fyo_dphi.clear();
367  phi_nobha_fyo_dphi_psnecl.clear();
368  phi_nobha_fyo_dphi_psnecl_ftdf.clear();
369 
370  int nitrack = 0; // the i track of the m_Tracks
371  for (const auto& b2track : m_Tracks) {
372  const Belle2::TrackFitResult* fitresult = b2track.getTrackFitResultWithClosestMass(Const::pion);
373  if (!fitresult) {
374  B2WARNING("No track fit result found.");
375  nitrack++;
376  continue;
377  }
378 
379  // require high NDF track
380  int ndf = fitresult->getNDF();
381  if (ndf < 20) {
382  nitrack++;
383  continue;
384  }
385 
386  // IP tracks at barrel
387  if (fabs(fitresult->getD0()) < 1.0 && fabs(fitresult->getZ0()) < 1.0 && fitresult->getHitPatternCDC().getLastLayer() > 50
388  && fitresult->getHitPatternCDC().getFirstLayer() < 5) {
389  double phiDegree = fitresult->getPhi() / Unit::deg;
390  double pt = fitresult->getTransverseMomentum();
391  double p3 = fitresult->getMomentum().R(); // 3-momentum
392  double theta = fitresult->getMomentum().Theta() / Unit::deg;
393 
394  bool trg_psnecl = m_trgSummary->testPsnm("hie") || m_trgSummary->testPsnm("c4") || m_trgSummary->testPsnm("eclmumu") ||
395  m_trgSummary->testPsnm("ecltaub2b3") || m_trgSummary->testPsnm("hie4") || m_trgSummary->testPsnm("lml1") ||
396  m_trgSummary->testPsnm("lml2") || m_trgSummary->testPsnm("lml6") || m_trgSummary->testPsnm("lml7") ||
397  m_trgSummary->testPsnm("lml8") || m_trgSummary->testPsnm("lml9") || m_trgSummary->testPsnm("lml10");
398 
399  bool trg_ftdf = m_trgSummary->testFtdl("f");
400 
401  // for f bit, reomove the Bhabha_veto
402  bool trg_itdt2 = (m_trgSummary->testInput("t2_0") || m_trgSummary->testInput("t2_1") || m_trgSummary->testInput("t2_2")
403  || m_trgSummary->testInput("t2_3"));
404 
405  // for z, reomove the Bhabha_veto
406  bool trg_itdt3 = (m_trgSummary->testInput("t3_0") || m_trgSummary->testInput("t3_1") || m_trgSummary->testInput("t3_2")
407  || m_trgSummary->testInput("t3_3"));
408 
409  // for y, reomove the Bhabha_veto
410  bool trg_itdt4 = (m_trgSummary->testInput("ty_0") || m_trgSummary->testInput("ty_1") || m_trgSummary->testInput("ty_2")
411  || m_trgSummary->testInput("ty_3"));
412 
413  // (t3>0 and !bhaveto and !veto) for z
414  bool trg_ftdz = m_trgSummary->testFtdl("z");
415 
416  // (ty>0 and !bhaveto and !veto) for y
417  bool trg_ftdy = m_trgSummary->testFtdl("y");
418 
419  // typ and !bha veto and !veto
420  bool trg_stt = m_trgSummary->testFtdl("stt");
421 
422  // remove bha_veto
423  bool trg_stt_nobha = m_trgSummary->testInput("typ") ;
424 
425 
426  // require pt > 0.3 GeV
427  if (pt > 0.3) {
428  m_hPhi->Fill(phiDegree);
429  if (trg_psnecl) {
430  m_hPhi_psnecl->Fill(phiDegree);
431  }
432  if (trg_psnecl && trg_ftdf) {
433  m_hPhi_psnecl_ftdf->Fill(phiDegree);
434  }
435  }
436 
437  m_hPt->Fill(pt);
438  m_nobha_hPt->Fill(pt);
439 
440  m_hP3_z->Fill(p3);
441  m_hP3_y->Fill(p3);
442 
443  m_nobha_hP3_z->Fill(p3);
444  m_nobha_hP3_y->Fill(p3);
445 
446  m_stt_theta->Fill(theta);
447  m_stt_phi->Fill(phiDegree);
448  p_stt_P3.push_back(p3);
449 
450  m_nobha_stt_theta->Fill(theta);
451  m_nobha_stt_phi->Fill(phiDegree);
452  p_nobha_stt_P3.push_back(p3);
453 
454 
455  if (trg_psnecl) {
456  m_hPt_psnecl->Fill(pt);
457  m_nobha_hPt_psnecl->Fill(pt);
458 
459  m_hP3_z_psnecl->Fill(p3); // for z bit
460  m_hP3_y_psnecl->Fill(p3); // for y bit
461 
462  m_nobha_hP3_z_psnecl->Fill(p3); // remove bhabha veto for z bit
463  m_nobha_hP3_y_psnecl->Fill(p3); // remove bhabha veto for y bit
464 
465  m_stt_phi_psnecl->Fill(phiDegree);
466  p_stt_P3_psnecl.push_back(p3);
467  m_stt_theta_psnecl->Fill(theta);
468 
469  m_nobha_stt_phi_psnecl->Fill(phiDegree);
470  p_nobha_stt_P3_psnecl.push_back(p3);
471  m_nobha_stt_theta_psnecl->Fill(theta);
472  }
473 
474  if (trg_psnecl && trg_ftdf) {
475  m_hPt_psnecl_ftdf->Fill(pt);
476  }
477  if (trg_psnecl && trg_itdt2) {
478  m_nobha_hPt_psnecl_ftdf->Fill(pt);
479  }
480 
481  if (trg_psnecl && trg_ftdz) {
482  m_hP3_z_psnecl_ftdf->Fill(p3);
483  }
484  if (trg_psnecl && trg_ftdy) {
485  m_hP3_y_psnecl_ftdf->Fill(p3);
486  }
487  if (trg_psnecl && trg_itdt3) {
488  m_nobha_hP3_z_psnecl_ftdf->Fill(p3);
489  }
490  if (trg_psnecl && trg_itdt4) {
491  m_nobha_hP3_y_psnecl_ftdf->Fill(p3);
492  }
493 
494  if (trg_psnecl && trg_stt) {
495  m_stt_phi_psnecl_ftdf->Fill(phiDegree);
496  p_stt_P3_psnecl_ftdf.push_back(p3);
497  m_stt_theta_psnecl_ftdf->Fill(theta);
498  }
499  if (trg_psnecl && trg_stt_nobha) {
500  m_nobha_stt_phi_psnecl_ftdf->Fill(phiDegree);
501  p_nobha_stt_P3_psnecl_ftdf.push_back(p3);
502  m_nobha_stt_theta_psnecl_ftdf->Fill(theta);
503  }
504 
505 
506  // the following is for fyo \deleta_phi
507  int njtrack = 0; // the j track of the m_Tracks
508  for (const auto& j_b2track : m_Tracks) {
509  if (nitrack >= njtrack) {
510  njtrack++;
511  continue;
512  }
513 
514  const Belle2::TrackFitResult* jfitresult = j_b2track.getTrackFitResultWithClosestMass(Const::pion);
515  if (!jfitresult) {
516  B2WARNING("No track fit result found.");
517  njtrack++;
518  continue;
519  }
520 
521  // require high NDF track
522  int jndf = jfitresult->getNDF();
523  if (jndf < 20) {
524  njtrack++;
525  continue;
526  }
527 
528  // IP tracks at barrel
529  if (fabs(jfitresult->getD0()) < 1.0 && fabs(jfitresult->getZ0()) < 1.0 && jfitresult->getHitPatternCDC().getLastLayer() > 50
530  && jfitresult->getHitPatternCDC().getFirstLayer() < 5) {
531  double jrk_phiDegree = jfitresult->getPhi() / Unit::deg;
532  double deltea_phi = fabs(phiDegree - jrk_phiDegree);
533  double dphi = deltea_phi;
534 
535  if (deltea_phi > 180)
536  dphi = 360 - deltea_phi;
537 
538  bool trg_fyo = m_trgSummary->testFtdl("fyo");
539  bool trg_fyo_nobha = (m_trgSummary->testInput("t2_1") || m_trgSummary->testInput("t2_2") || m_trgSummary->testInput("t2_3"))
540  &&
541  (m_trgSummary->testInput("ty_0") || m_trgSummary->testInput("ty_1") || m_trgSummary->testInput("ty_2")
542  || m_trgSummary->testInput("ty_3")) &&
543  m_trgSummary->testInput("cdc_open90");
544 
545  // cout<<"i = "<<nitrack <<" j= "<<njtrack <<" phiDegree = "<<phiDegree<<" jrk_phiDegree = "<<jrk_phiDegree<<" dphi "<<dphi<<endl;
546 
547  phi_fyo_dphi.push_back(dphi);
548  phi_nobha_fyo_dphi.push_back(dphi);
549 
550  if (trg_psnecl) {
551  phi_fyo_dphi_psnecl.push_back(dphi);
552  phi_nobha_fyo_dphi_psnecl.push_back(dphi);
553  }
554  if (trg_psnecl && trg_fyo) {
555  phi_fyo_dphi_psnecl_ftdf.push_back(dphi);
556  }
557  if (trg_psnecl && trg_fyo_nobha) {
558  phi_nobha_fyo_dphi_psnecl_ftdf.push_back(dphi);
559  }
560  }
561  njtrack++;
562  }
563  }
564  nitrack++;
565  }
566 
567 
568  // the largest cdc_open angle in an event for fyo bit
569  if (phi_fyo_dphi_psnecl_ftdf.size() != 0) {
570  auto max_it = std::max_element(phi_fyo_dphi_psnecl_ftdf.begin(), phi_fyo_dphi_psnecl_ftdf.end());
571  double max_value = *max_it;
572  m_fyo_dphi_psnecl_ftdf->Fill(max_value);
573  }
574  if (phi_fyo_dphi_psnecl.size() != 0) {
575  auto max_it = std::max_element(phi_fyo_dphi_psnecl.begin(), phi_fyo_dphi_psnecl.end());
576  double max_value = *max_it;
577  m_fyo_dphi_psnecl->Fill(max_value);
578  }
579  if (phi_fyo_dphi.size() != 0) {
580  auto max_it = std::max_element(phi_fyo_dphi.begin(), phi_fyo_dphi.end());
581  double max_value = *max_it;
582  m_fyo_dphi->Fill(max_value);
583  }
584 
585  //
586  if (phi_nobha_fyo_dphi_psnecl_ftdf.size() != 0) {
587  auto max_it = std::max_element(phi_nobha_fyo_dphi_psnecl_ftdf.begin(), phi_nobha_fyo_dphi_psnecl_ftdf.end());
588  double max_value = *max_it;
589  m_nobha_fyo_dphi_psnecl_ftdf->Fill(max_value);
590  }
591  if (phi_nobha_fyo_dphi_psnecl.size() != 0) {
592  auto max_it = std::max_element(phi_nobha_fyo_dphi_psnecl.begin(), phi_nobha_fyo_dphi_psnecl.end());
593  double max_value = *max_it;
594  m_nobha_fyo_dphi_psnecl->Fill(max_value);
595  }
596  if (phi_nobha_fyo_dphi.size() != 0) {
597  auto max_it = std::max_element(phi_nobha_fyo_dphi.begin(), phi_nobha_fyo_dphi.end());
598  double max_value = *max_it;
599  m_nobha_fyo_dphi->Fill(max_value);
600  }
601 
602 
603  // the largest momentum track p in an event for stt bit
604  if (p_stt_P3_psnecl_ftdf.size() != 0) {
605  auto max_it = std::max_element(p_stt_P3_psnecl_ftdf.begin(), p_stt_P3_psnecl_ftdf.end());
606  double max_value = *max_it;
607  m_stt_P3_psnecl_ftdf->Fill(max_value);
608  }
609  if (p_stt_P3_psnecl.size() != 0) {
610  auto max_it = std::max_element(p_stt_P3_psnecl.begin(), p_stt_P3_psnecl.end());
611  double max_value = *max_it;
612  m_stt_P3_psnecl->Fill(max_value);
613  }
614  if (p_stt_P3.size() != 0) {
615  auto max_it = std::max_element(p_stt_P3.begin(), p_stt_P3.end());
616  double max_value = *max_it;
617  m_stt_P3->Fill(max_value);
618  }
619 
620  //
621  if (p_nobha_stt_P3_psnecl_ftdf.size() != 0) {
622  auto max_it = std::max_element(p_nobha_stt_P3_psnecl_ftdf.begin(), p_nobha_stt_P3_psnecl_ftdf.end());
623  double max_value = *max_it;
624  m_nobha_stt_P3_psnecl_ftdf->Fill(max_value);
625  }
626  if (p_nobha_stt_P3_psnecl.size() != 0) {
627  auto max_it = std::max_element(p_nobha_stt_P3_psnecl.begin(), p_nobha_stt_P3_psnecl.end());
628  double max_value = *max_it;
629  m_nobha_stt_P3_psnecl->Fill(max_value);
630  }
631  if (p_nobha_stt_P3.size() != 0) {
632  auto max_it = std::max_element(p_nobha_stt_P3.begin(), p_nobha_stt_P3.end());
633  double max_value = *max_it;
634  m_nobha_stt_P3->Fill(max_value);
635  }
636 
637 
638 
639 
640 }
641 
642 
643 
644 void TRGEFFDQMModule::endRun()
645 {
646 }
647 
648 void TRGEFFDQMModule::terminate()
649 {
650 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ 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:80
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
void initialize() override
Initialize the Module.
TH1F * m_hPhi_psnecl_ftdf
Histogram of cdc phi of IP tracks with ecl and f bit.
virtual ~TRGEFFDQMModule()
Destructor.
void event() override
Event processor.
StoreObjPtr< TRGSummary > m_trgSummary
Trigger summary.
TH1F * m_hPhi_psnecl
Histogram of cdc phi of IP tracks with ecl bit.
void beginRun() override
Called when entering a new run.
StoreArray< ECLCluster > m_ECLClusters
ECL Clusters.
TH1F * m_hPhi
Histogram of cdc phi of IP tracks.
StoreArray< KLMCluster > m_KLMClusters
KLM Clusters.
StoreArray< RecoTrack > m_RecoTracks
RecoTracks.
StoreArray< Track > m_Tracks
Tracks.
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
Store array for Trigger selection.
void defineHisto() override
Histogram definitions.
static const double deg
degree to radians
Definition: Unit.h:109
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
@ c_accept
Accept this event.
Abstract base class for different kinds of events.