Belle II Software development
TOPChannelT0MCModule.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 header.
10#include <top/modules/TOPChannelT0MC/TOPChannelT0MCModule.h>
11
12// framework - DataStore
13#include <framework/datastore/StoreArray.h>
14
15// database classes
16#include <framework/database/DBStore.h>
17
18#include <TFile.h>
19#include <TH1F.h>
20#include <TTree.h>
21#include <sstream>
22
23using namespace std;
24
25namespace Belle2 {
30 //-----------------------------------------------------------------
32 //-----------------------------------------------------------------
33
34 REG_MODULE(TOPChannelT0MC);
35
36 //-----------------------------------------------------------------
37 // Implementation
38 //-----------------------------------------------------------------
39
41 {
42 std::vector<double> frange = {100, 0., 1.};
43 // set module description (e.g. insert text)
44 setDescription("TOP channel T0 Calibration of MC extraction");
45 // Add parameters
46 addParam("outputFile", m_outputFile, "Output root file name",
47 string(""));
48 addParam("fitRange", m_fitRange, "fit range[nbins, xmin, xmax]", frange);
49
50 }
51
53 {
54 m_digits.isRequired();
55 }
56
58 {
59
60 for (auto& digit : m_digits) {
61 unsigned channel = digit.getChannel();
62 if (channel < c_NumChannels) {
63 auto histo = m_histo[channel];
64 if (!histo) {
65 stringstream ss;
66 ss << "chan" << channel ;
67 string name;
68 ss >> name;
69 string title = "Times " + name;
70 histo = new TH1F(name.c_str(), title.c_str(), (int)m_fitRange[0], m_fitRange[1], m_fitRange[2]);
71 m_histo[channel] = histo;
72 }
73 if (digit.getTime() < m_fitRange[1] || digit.getTime() > m_fitRange[2]) continue;
74 histo->Fill(digit.getTime());
75 }
76 }
77 }
78
80 {
81 auto file = new TFile(m_outputFile.c_str(), "RECREATE");
82 auto otree = new TTree("t0MC", "extract channel t0 info. from MC");
83
84 unsigned channel = 0;
85 double maxpos = 0;
86
87 otree->Branch("maxpos", &maxpos, "maxpos/D");
88 otree->Branch("channel", &channel, "channel/I");
89
90 for (int i = 0; i < c_NumChannels; i++) {
91 channel = i;
92 maxpos = m_histo[i]->GetXaxis()->GetBinCenter(m_histo[i]->GetMaximumBin());
93 otree->Fill();
94 }
95 otree->Write();
96 delete otree;
97 file->Close();
98 delete file;
99 }
100
101} // end Belle2 namespace
102
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
Module()
Constructor.
Definition Module.cc:30
std::vector< double > m_fitRange
fit range [nbins, xmin, xmax]
std::string m_outputFile
output root file name
StoreArray< TOPDigit > m_digits
collection of digits
TH1F * m_histo[c_NumChannels]
profile histograms
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void terminate() override
Termination action.
Abstract base class for different kinds of events.
STL namespace.