Belle II Software  release-08-01-10
TRGECLQAMModule.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 //framework headers
10 #include <framework/datastore/StoreObjPtr.h>
11 #include <framework/datastore/StoreArray.h>
12 #include <framework/dataobjects/EventMetaData.h>
13 //trg package headers
14 
15 #include <trg/ecl/modules/trgeclQAM/TRGECLQAMModule.h>
16 
17 #include <string.h>
18 #include <fstream>
19 #include <TFile.h>
20 #include <TString.h>
21 
22 #include <TAxis.h>
23 #include <TTree.h>
24 
25 #include <RooRealVar.h>
26 #include <RooDataHist.h>
27 #include <RooGaussian.h>
28 #include <RooAddModel.h>
29 #include <RooPlot.h>
30 #include <RooArgList.h>
31 #include <RooNovosibirsk.h>
32 
33 using namespace std;
34 
35 namespace Belle2 {
40  //
41  //
42 
44  REG_MODULE(TRGECLQAM);
45  //
46  //
47  //
48  string
49  TRGECLQAMModule::version() const
50  {
51  return string("TRGECLQAMModule 1.00");
52  }
53  //
54  //
55  //
56  TRGECLQAMModule::TRGECLQAMModule()
57  : Module::Module(),
58  m_nevent(1), m_outputfile("")
59  {
60 
61  string desc = "TRGECLQAMModule(" + version() + ")";
62  setDescription(desc);
64 
65  addParam("outputFileName", m_outputfile, "TRGECL QAM file", m_outputfile);
66 
67  m_nRun = 0;
68  m_nevent = 1;
69 
70  //Problem TC Check
71  m_FWD = 0;
72  m_BAR = 0;
73  m_BWD = 0;
74  m_ALL = 0;
75  memset(TCID, 0, sizeof(TCID));
76 
77 
78  //Total Energy
79 
80  m_etot = 0;
81  m_etot_mean = 0;
82  m_etot_error = 0;
83  m_etot_sigma = 0;
84 
85  // std::vector<int> etot;
86  etot.clear();
87 
88  //Caltime
89 
90  m_caltime_mean = 0;
91  m_caltime_error = 0;
92  m_caltime_sigma = 0;
93  caltime.clear();
94 
95  //Cluster E
96 
97  clusterE = 0;
98  cluster.clear();
99 
100  B2DEBUG(100, "TRGECLQAMModule ... created");
101  }
102  //
103  //
104  //
106  {
107 
108  B2DEBUG(100, "TRGECLQAMModule ... destructed ");
109  }
110  //
111  //
112  //
113  void
115  {
116 
117  B2DEBUG(100, "TRGECLQAMModule::initialize ... options");
118 
119 
120  m_TRGECLUnpackerStore.registerInDataStore();
121  m_TRGECLUnpackerEvtStore.registerInDataStore();
122 
123  m_TRGSummary.registerInDataStore();
124  // EvtMeta.registterInDataStore();
125 
126 
127  }
128  //
129  //
130  //
131  void
133  {
134 
135  B2DEBUG(200, "TRGECLQAMModule ... beginRun called ");
136 
137  }
138  //
139  //
140  //
141  void
143  {
144 
145  B2DEBUG(200, "TRGECLFAMMoudle ... event called");
146 
147 
148 
149  //QAM ECL Hit Map
150 
151  // StoreArray<TRGECLUnpackerStore> m_TRGECLUnpackerStore;
152  for (int ii = 0; ii < m_TRGECLUnpackerStore.getEntries(); ii++) {
153  TRGECLUnpackerStore* m_TRGECLUnpacker = m_TRGECLUnpackerStore[ii];
154 
155 
156  int TCId = (m_TRGECLUnpacker->getTCId());
157  int HitEnergy = (m_TRGECLUnpacker->getTCEnergy());
158 
159  m_etot += HitEnergy;
160  if (TCId < 1 || TCId > 576 || HitEnergy == 0) {continue;}
161  TCID[TCId - 1]++;
162  int Caltime = m_TRGECLUnpacker->getTCCALTime();
163  if (Caltime != 0) {
164  caltime.push_back(Caltime);
165  }
166  }
167 
168  etot.push_back(m_etot);
169 
170  m_etot = 0;
171 
172  //Unpacker EvtStore
173  //for ( int ii = 0; ii < m_TRGECLUnpackerSumStore.getEntries(); ii++){
174  TRGECLUnpackerEvtStore* m_TRGECLUnpackerEvt = m_TRGECLUnpackerEvtStore[0];
175 
176  if (m_TRGECLUnpackerEvt->get3DBhabhaS() == 1) {
177  //Cluster Energy 2D Check
178  const int* clE = m_TRGECLUnpackerEvt->getCLEnergy();
179 
180  int clER[6] = {0};
181 
182  for (int i = 0; i < 6; i++) {
183  clER[i] = clE[i];
184  }
185 
186  sort(clER, clER + 6);
187 
188  clusterE = clER[5] + clER[4]; // Ecl1+Ecl2
189  cluster.push_back(clusterE);
190  }
191 
193  m_nRun = bevt->getRun();
194  m_nevent++;
195 
196  }
197  //
198  //
199  //
200  void
202  {
203 
204  double mean_FWD = 0;
205  double mean_BAR = 0;
206  double mean_BWD = 0;
207 
208  //QAM TC Hit Map Check
209  for (int TCId = 1; TCId < 577; TCId++) {
210  if (TCId < 81) { //Forward Endcap
211  mean_FWD += TCID[TCId - 1];
212  } else if (TCId < 513) { //Barrel
213  mean_BAR += TCID[TCId - 1];
214  } else { //Backward Endcap
215  mean_BWD += TCID[TCId - 1];
216  }
217  }
218  mean_FWD /= 80;
219  mean_BAR /= (512 - 80);
220  mean_BWD /= (576 - 512);
221 
222  for (int TCId = 1; TCId < 577; TCId++) {
223  if (TCId < 81) { //Forward Endcap
224  if (TCID[TCId - 1] < mean_FWD * 0.1)m_FWD++;
225  } else if (TCId < 513) { //Barrel
226  if (TCID[TCId - 1] < mean_BAR * 0.1)m_BAR++;
227  } else { //Backward Endcap
228  if (TCID[TCId - 1] < mean_BWD * 0.1)m_BWD++;
229  }
230  }
231 
232  m_ALL = m_FWD + m_BAR + m_BWD;
233  const int etot_size = etot.size();
234  for (int i = 0; i < etot_size; i++) {
235  h_etot->Fill(etot[i]);
236  }
237  //Caluate E_total peak and width
238 
239  RooRealVar* E_tot = new RooRealVar("E_tot", "E_tot", 0, 4000);
240  RooDataHist* data_E_tot = new RooDataHist("data_E_tot", "dataset with E_tot", *E_tot, h_etot);
241  RooRealVar* mean_check = new RooRealVar("mean_check", "Mean for checking", 2000, 1000, 3000);
242  RooRealVar* sigma_check = new RooRealVar("sigma_check", "Sigma for checking", 100, 0, 300);
243  RooGaussian* gauss_check = new RooGaussian("gauss_check", "Gaussian for checking", *E_tot, *mean_check, *sigma_check);
244  gauss_check->fitTo(*data_E_tot);
245  double mean_ch = mean_check->getVal();
246 
247  RooRealVar* mean_E_tot = new RooRealVar("mean_E_tot", "Mean of Gaussian", mean_ch, mean_ch - 200, mean_ch + 200);
248  // RooRealVar sigma_E_tot("sigma_E_tot","Width of Gaussian", 10, 0, 50);
249 
250  RooRealVar* sigma1_E_tot = new RooRealVar("sigma1_E_tot", "Width of Gaussian", 30, 0, 60);
251  RooGaussian* gauss1_E_tot = new RooGaussian("gauss1_E_tot", "gauss(x,mean,sigma)", *E_tot, *mean_E_tot, *sigma1_E_tot);
252  RooRealVar* sigma2_E_tot = new RooRealVar("sigma2_E_tot", "Width of Gaussian", 100, 30, 250);
253  RooGaussian* gauss2_E_tot = new RooGaussian("gauss2_E_tot", "gauss(x,mean,sigma)", *E_tot, *mean_E_tot, *sigma2_E_tot);
254 
255  RooRealVar* frac_core_E_tot = new RooRealVar("frac_core", "core fraction", 0.85, 0.0, 1.0);
256  RooAddModel* gaussm_E_tot = new RooAddModel("gaussm", "core+tail gauss", RooArgList(*gauss1_E_tot, *gauss2_E_tot),
257  RooArgList(*frac_core_E_tot));
258 
259  // RooGaussian gauss_E_tot("gauss_E_tot","gauss(x,mean,sigma)",E_tot,mean_E_tot,sigma_E_tot);
260 
261  gaussm_E_tot->fitTo(*data_E_tot, RooFit::Range(mean_ch - 300, mean_ch + 300));
262  //RooPlot* xframe_E_tot = E_tot->frame();
263  RooPlot* xframe_E_tot = E_tot->frame(RooFit::Title("energy_tot"));
264  //data_E_tot->plotOn(xframe_E_tot);
265  data_E_tot->plotOn(xframe_E_tot, RooFit::MarkerStyle(7));
266  // gauss_E_tot->plotOn(xframe_E_tot);
267  gaussm_E_tot->plotOn(xframe_E_tot, RooFit::Name("Cal_energy"), RooFit::LineColor(2));
268  double v_m_E_tot = mean_E_tot->getVal();
269  double e_m_E_tot = mean_E_tot->getError();
270 
271  double v1_s1_E_tot = sigma1_E_tot->getVal();
272  // double e1_s1_E_tot = sigma1_E_tot->getError();
273  double v2_s1_E_tot = sigma2_E_tot->getVal();
274  // double e2_s1_E_tot = sigma2_E_tot->getError();
275 
276  double frac_E_tot = frac_core_E_tot->getVal();
277  double v_s1_E_tot = v1_s1_E_tot * (frac_E_tot) + v2_s1_E_tot * (1 - frac_E_tot);
278 
279 
280 
281  for (long unsigned int ii = 0; ii < caltime.size(); ii++) {
282  h_caltime->Fill(caltime[ii]);
283  }
284 
285 
286  //Calculate Caltime peak and width
287  RooRealVar T_c("T_c", "T_c", -100, 100);
288  RooDataHist data_T_c("data_T_c", "dataset with T_c", T_c, h_caltime);
289 
290  RooRealVar mean_T_c("mean_T_c", "Mean of Gaussian", 0, -100, 100);
291  RooRealVar sigma1_T_c("sigma1_T_c", "Width of Gaussian", 10, 0, 50);
292  RooGaussian gauss1_T_c("gauss1_T_c", "gauss(x,mean,sigma)", T_c, mean_T_c, sigma1_T_c);
293  RooRealVar sigma2_T_c("sigma2_T_c", "Width of Gaussian", 30, 10, 100);
294  RooGaussian gauss2_T_c("gauss2_T_c", "gauss(x,mean,sigma)", T_c, mean_T_c, sigma2_T_c);
295 
296  RooRealVar frac_core_T_c("frac_core", "core fraction", 0.85, 0.0, 1.0);
297  RooAddModel gaussm_T_c("gaussm", "core+tail gauss", RooArgList(gauss1_T_c, gauss2_T_c), RooArgList(frac_core_T_c));
298 
299  gaussm_T_c.fitTo(data_T_c, RooFit::Range(-45, 45));
300  RooPlot* xframe_T_c = T_c.frame(RooFit::Title("Caltime"));
301  data_T_c.plotOn(xframe_T_c, RooFit::MarkerStyle(7));
302  gaussm_T_c.plotOn(xframe_T_c, RooFit::Name("Caltime"), RooFit::LineColor(2));
303  gaussm_T_c.plotOn(xframe_T_c, RooFit::Components("gauss1_T_c"), RooFit::LineColor(kGreen));
304  gaussm_T_c.plotOn(xframe_T_c, RooFit::Components("gauss2_T_c"), RooFit::LineColor(kMagenta));
305 
306  xframe_T_c->GetXaxis()->SetTitle("Caltime[ns]");
307  xframe_T_c->GetYaxis()->SetTitle("Entries");
308 
309  double v_m_T_c = mean_T_c.getVal();
310  double e_m_T_c = mean_T_c.getError();
311 
312  double v1_s1_T_c = sigma1_T_c.getVal();
313  // double e1_s1_T_c = sigma1_T_c.getError();
314  double v2_s1_T_c = sigma2_T_c.getVal();
315  // double e2_s1_T_c = sigma2_T_c.getError();
316 
317  double frac_T_c = frac_core_T_c.getVal();
318  double v_s1_T_c = v1_s1_T_c * (frac_T_c) + v2_s1_T_c * (1 - frac_T_c);
319 
320 
321 
322  //Cluster E sum for Bhabha skim
323  const int cluster_size = cluster.size();
324  for (int ii = 0; ii < cluster_size; ii++) {
325  h_clusterE->Fill(cluster[ii]);
326  }
327 
328  RooRealVar x("x", "esum", 0, 4000);
329  RooDataHist dh("dh", "ecltrg", x, h_clusterE);
330 
331  RooRealVar mean("mean", "mean", 2000, 1000, 2500);
332  RooRealVar sigma("sigma", "sigma", 150, 0, 300);
333  RooRealVar tail("tail", "tail", 0.45, 0, 1);
334  RooNovosibirsk novo("novo", "", x, mean, sigma, tail);
335  novo.fitTo(dh, RooFit::Extended(0), RooFit::Range(1500, 2200));
336 
337 
338  double clusterE_vm, clusterE_em, clusterE_vs;
339 
340  clusterE_vm = mean.getVal();
341  clusterE_em = mean.getError();
342  clusterE_vs = sigma.getVal();
343 
344 
345  TString outputfile = m_outputfile;
346  TFile file(outputfile, "RECREATE");
347  TTree* tree = new TTree("tree", "tree");
348  int nevent;
349  double FWD;
350  double BAR;
351  double BWD;
352  double ALL;
353 
354  nevent = m_nevent;
355  FWD = m_FWD;
356  BAR = m_BAR;
357  BWD = m_BWD;
358  ALL = m_ALL;
359 
360  m_etot_mean = v_m_E_tot;
361  m_etot_error = e_m_E_tot;
362  m_etot_sigma = v_s1_E_tot;
363 
364  m_caltime_mean = v_m_T_c;
365  m_caltime_error = e_m_T_c;
366  m_caltime_sigma = v_s1_T_c;
367 
368  tree->Branch("m_nRun", &m_nRun, "m_nRun/I");
369  tree->Branch("m_nevent", &nevent, "m_nevent/I");
370  tree->Branch("m_FWD", &FWD, "m_FWD/D");
371  tree->Branch("m_BAR", &BAR, "m_BAR/D");
372  tree->Branch("m_BWD", &BWD, "m_BWD/D");
373  tree->Branch("m_ALL", &ALL, "m_ALL/D");
374  tree->Branch("m_etot_mean", &m_etot_mean, "m_etot_mean/D");
375  tree->Branch("m_etot_error", &m_etot_error, "m_etot_error/D");
376  tree->Branch("m_etot_sigma", &m_etot_sigma, "m_etot_sigma/D");
377  tree->Branch("m_caltime_mean", &m_caltime_mean, "m_caltime_mean/D");
378  tree->Branch("m_caltime_error", &m_caltime_error, "m_caltime_error/D");
379  tree->Branch("m_caltime_sigma", &m_caltime_sigma, "m_caltime_sigma/D");
380  tree->Branch("m_clusterE_mean", &clusterE_vm, "m_clusterE_mean/D");
381  tree->Branch("m_clusterE_error", &clusterE_em, "m_clusterE_error/D");
382  tree->Branch("m_clusterE_sigma", &clusterE_vs, "m_clusterE_sigma/D");
383 
384  tree->Fill();
385  tree->Write();
386  file.Close();
387 
388 
389 
390  B2INFO("ECL Trigger QAM result "
391  << LogVar("Run Number", m_nRun)
392  << LogVar("Total Event", nevent)
393  << LogVar("The # of Low Hit TC in Forward Endcap", FWD)
394  << LogVar("The # of Low Hit TC in Barrel", BAR)
395  << LogVar("The # of Low Hit TC in Backward Endcap", BWD)
396  << LogVar("Total Energy Peak", m_etot_mean)
397  << LogVar("Total Energy Peak width", m_etot_sigma)
398  << LogVar("Caltime Peak", m_caltime_mean)
399  << LogVar("Caltime Peak width", m_caltime_sigma)
400  << LogVar("Cluster Energy Peak", clusterE_vm)
401  << LogVar("Cluster Energy Peak width", clusterE_vs));
402 
403 
404  B2DEBUG(100, "TRGECLQAMModule ... endRun called ");
405 
406  }
407  //
408  //
409  //
411  {
412  B2DEBUG(100, "TRGECLQAMModule ... terminate called ");
413 
414  }
415  //
416  //
417  //
419 } // namespace Belle2
Base class for Modules.
Definition: Module.h:72
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
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
TH1F * h_caltime
Caltime Histogram.
double m_ALL
Total proble TC.
double m_etot_sigma
Total Energy sigma.
StoreArray< TRGECLUnpackerEvtStore > m_TRGECLUnpackerEvtStore
ECL Trigger Unpacker Event output.
int TCID[576]
Hit TCID.
std::vector< int > caltime
Caltime check.
double m_BAR
Barrel problem TC.
double m_BWD
BWD problem TC.
TH1F * h_clusterE
Cluster Energy Histogram.
std::string m_outputfile
Output Root File Name.
double m_etot_mean
Total Energy mean.
std::vector< int > etot
tcenergy check
double m_caltime_mean
Caltime mean.
TH1F * h_etot
Total Energy Histogram.
StoreArray< TRGECLUnpackerStore > m_TRGECLUnpackerStore
ECL Trigger Unpacker TC output.
std::vector< int > cluster
Cluster Energy Vector.
double m_etot_error
Total Energy error.
double m_caltime_error
Caltime error.
int m_nevent
The # of Events.
Double_t m_FWD
FWD problem TC.
StoreObjPtr< TRGSummary > m_TRGSummary
Trigger Summary.
int clusterE
Cluster Energy.
double m_caltime_sigma
Caltime sigma.
int m_etot
Total Energy.
int get3DBhabhaS() const
The method to get 3D Bhabha selection bit.
const int * getCLEnergy() const
The mothod to get Cluster energy.
int getTCCALTime() const
The method to get cal timing.
int getTCEnergy() const
The method to get deposited energy.
int getTCId() const
The method to get cell id.
Class to store variables with their name which were sent to the logging service.
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:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
virtual ~TRGECLQAMModule()
Destrunctor.
virtual void initialize() override
initialize
virtual void event() override
Event.
virtual void endRun() override
End Run.
virtual void terminate() override
terminate
virtual void beginRun() override
begin Run
std::string version() const
version
Abstract base class for different kinds of events.