Belle II Software development
TOPCommonT0LLCollectorModule.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#include <top/modules/collectors/TOPCommonT0LLCollectorModule.h>
10#include <top/geometry/TOPGeometryPar.h>
11#include <top/reconstruction_cpp/TOPTrack.h>
12#include <top/reconstruction_cpp/PDFConstructor.h>
13#include <top/reconstruction_cpp/TOPRecoManager.h>
14
15// framework aux
16#include <framework/gearbox/Const.h>
17#include <framework/logging/Logger.h>
18
19// root
20#include <TRandom.h>
21#include <TH1F.h>
22#include <TH1D.h>
23#include <TH2F.h>
24
25using namespace std;
26
27namespace Belle2 {
33 using namespace TOP;
34
35 //-----------------------------------------------------------------
37 //-----------------------------------------------------------------
38
39 REG_MODULE(TOPCommonT0LLCollector);
40
41 //-----------------------------------------------------------------
42 // Implementation
43 //-----------------------------------------------------------------
44
46 {
47 // set module description and processing properties
48 setDescription("A collector for common T0 calibration with dimuons or Bhabha's using "
49 "neg. log likelihood minimization (method LL)");
51
52 // module parameters
53 addParam("bunchesPerSSTclk", m_bunchesPerSSTclk,
54 "number of bunches per SST clock period", 24);
55 addParam("numBins", m_numBins, "number of bins of the search region", 200);
56 addParam("timeRange", m_timeRange,
57 "time range in which to search for the minimum [ns]", 10.0);
58 addParam("sigmaSmear", m_sigmaSmear,
59 "sigma in [ns] for additional smearing of PDF", 0.0);
60 addParam("sample", m_sample,
61 "sample type: one of dimuon or bhabha", std::string("dimuon"));
62 addParam("deltaEcms", m_deltaEcms,
63 "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
64 addParam("dr", m_dr, "cut on POCA in r", 2.0);
65 addParam("dz", m_dz, "cut on POCA in abs(z)", 4.0);
66 addParam("minZ", m_minZ,
67 "minimal local z of extrapolated hit", -130.0);
68 addParam("maxZ", m_maxZ,
69 "maximal local z of extrapolated hit", 130.0);
70 addParam("pdfOption", m_pdfOption,
71 "PDF option, one of 'rough', 'fine', 'optimal'", std::string("rough"));
72 }
73
74
76 {
77 // input collections
78
79 m_digits.isRequired();
80 m_tracks.isRequired();
81 m_extHits.isRequired();
82 m_recBunch.isRequired();
83
84 // bunch time separation
85
86 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
87 m_bunchTimeSep = geo->getNominalTDC().getSyncTimeBase() / m_bunchesPerSSTclk;
88
89 // Parse PDF option
90
91 if (m_pdfOption == "rough") {
93 } else if (m_pdfOption == "fine") {
95 } else if (m_pdfOption == "optimal") {
97 } else {
98 B2ERROR("Unknown PDF option '" << m_pdfOption << "'");
99 }
100
101 // set track selector
102
103 if (m_sample == "dimuon" or m_sample == "bhabha") {
108 } else {
109 B2ERROR("Invalid sample type '" << m_sample << "'");
110 }
111
112 // create and register histograms
113
114 double tmin = -m_timeRange / 2;
115 double tmax = m_timeRange / 2;
116 for (unsigned i = 0; i < c_numSets; i++) {
117 string name = "chi2_set" + to_string(i);
118 auto h = new TH1D(name.c_str(), "chi2 scan; t0 [ns]; chi2", m_numBins, tmin, tmax);
119 registerObject<TH1D>(name, h);
120 m_names.push_back(name);
121 }
122
123 auto h1 = new TH1F("tracks_per_set", "tracks per sample; sample number; num tracks",
124 c_numSets, 0, c_numSets);
125 registerObject<TH1F>("tracks_per_set", h1);
126
127 auto h2 = new TH1F("numHits", "Number of photons per slot",
128 c_numModules, 0.5, static_cast<float>(c_numModules) + 0.5);
129 h2->SetXTitle("slot number");
130 h2->SetYTitle("hits per slot");
131 registerObject<TH1F>("numHits", h2);
132
133 auto h3 = new TH2F("timeHits", "Photon times vs. boardstacks",
134 c_numModules * 4, 0.5, static_cast<float>(c_numModules) + 0.5, 200, 0.0, 20.0);
135 h3->SetXTitle("slot number");
136 h3->SetYTitle("time [ns]");
137 registerObject<TH2F>("timeHits", h3);
138
139 // this one is needed primarely to pass bunch time separation to the algorithm,
140 // since DB interface doesn't work there
141 auto h4 = new TH1F("offset", "current offset from input files; offset [ns]",
142 200, -m_bunchTimeSep / 2, m_bunchTimeSep / 2);
143 registerObject<TH1F>("offset", h4);
144
145 }
146
147
149 {
150 // bunch must be reconstructed
151
152 if (not m_recBunch.isValid()) return;
153 if (not m_recBunch->isReconstructed()) return;
154
156 double timeMin = TOPRecoManager::getMinTime();
157 double timeMax = TOPRecoManager::getMaxTime();
158
159 // correct-back digits for common T0
160
161 double T0 = 0;
162 if (m_commonT0->isCalibrated()) T0 = m_commonT0->getT0();
163 for (auto& digit : m_digits) {
164 if (digit.isCommonT0Calibrated()) digit.subtractT0(-T0);
165 if (digit.hasStatus(TOPDigit::c_BunchOffsetSubtracted)) digit.subtractT0(-m_recBunch->getAverageOffset());
166 }
167
168 // loop over reconstructed tracks, make a selection and accumulate log likelihoods
169
170 int ntra = 0;
171 for (const auto& track : m_tracks) {
172
173 // track selection
174 TOPTrack trk(track);
175 if (not trk.isValid()) continue;
176
177 if (not m_selector.isSelected(trk)) continue;
178
179 // construct PDF
181 if (not pdfConstructor.isValid()) continue;
182
183 // minimization procedure: accumulate
184 int sub = gRandom->Integer(c_numSets); // generate sub-sample number
185 auto h = getObjectPtr<TH1D>(m_names[sub]);
186 for (int ibin = 0; ibin < h->GetNbinsX(); ibin++) {
187 double t0 = h->GetBinCenter(ibin + 1);
188 double chi = h->GetBinContent(ibin + 1);
189 chi += -2 * pdfConstructor.getLogL(t0, m_sigmaSmear).logL;
190 h->SetBinContent(ibin + 1, chi);
191 }
192 auto h1 = getObjectPtr<TH1F>("tracks_per_set");
193 h1->Fill(sub);
194
195 // fill histograms of hits
196 auto h2 = getObjectPtr<TH1F>("numHits");
197 auto h3 = getObjectPtr<TH2F>("timeHits");
198 for (const auto& digit : m_digits) {
199 if (digit.getHitQuality() != TOPDigit::c_Good) continue;
200 if (digit.getModuleID() != trk.getModuleID()) continue;
201 if (digit.getTime() < timeMin) continue;
202 if (digit.getTime() > timeMax) continue;
203 h2->Fill(digit.getModuleID());
204 int bs = digit.getBoardstackNumber();
205 h3->Fill((digit.getModuleID() * 4 + bs - 1.5) / 4.0, digit.getTime());
206 }
207 ntra++;
208 }
209
210 // fill just another control histogram
211
212 if (ntra > 0) {
213 auto h4 = getObjectPtr<TH1F>("offset");
214 h4->Fill(m_recBunch->getCurrentOffset());
215 }
216
217 }
218
220}
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
StoreObjPtr< TOPRecBunch > m_recBunch
reconstructed bunch
TOP::TrackSelector m_selector
track selection utility
double m_sigmaSmear
additional smearing of PDF in [ns]
DBObjPtr< TOPCalCommonT0 > m_commonT0
common T0 calibration constants
int m_bunchesPerSSTclk
number of bunches per SST clock
int m_numBins
number of bins to which search region is divided
double m_maxZ
maximal local z of extrapolated hit
TOP::PDFConstructor::EPDFOption m_PDFOption
PDF option.
@ c_numSets
number of statistically independent subsamples
double m_minZ
minimal local z of extrapolated hit
std::vector< std::string > m_names
histogram names of chi2 scans
StoreArray< Track > m_tracks
collection of tracks
double m_timeRange
time range in which to search for the minimum [ns]
StoreArray< TOPDigit > m_digits
collection of digits
StoreArray< ExtHit > m_extHits
collection of extrapolated hits
double m_bunchTimeSep
bunch separation in time [ns]
PDF construction and log likelihood determination for a given track and particle hypothesis.
LogL getLogL() const
Returns extended log likelihood (using the default time window)
bool isValid() const
Checks the object status.
@ c_Optimal
y dependent only where necessary
@ c_Fine
y dependent everywhere
@ c_Rough
no dependence on y
const TOPGeometry * getGeometry() const
Returns pointer to geometry object using basf2 units.
static TOPGeometryPar * Instance()
Static method to obtain the pointer to its instance.
static void setDefaultTimeWindow()
Sets default time window (functions getMinTime(), getMaxTime() will then return default values from D...
static double getMaxTime()
Returns time window upper edge.
static double getMinTime()
Returns time window lower edge.
Reconstructed track at TOP.
Definition: TOPTrack.h:39
bool isValid() const
Checks if track is successfully constructed.
Definition: TOPTrack.h:137
int getModuleID() const
Returns slot ID.
Definition: TOPTrack.h:143
Utility for the track selection - used in various calibration modules.
Definition: TrackSelector.h:27
void setDeltaEcms(double deltaEcms)
Sets cut on c.m.s.
Definition: TrackSelector.h:63
void setCutOnPOCA(double dr, double dz)
Sets cut on point of closest approach to (0, 0, 0)
Definition: TrackSelector.h:70
bool isSelected(const TOPTrack &track) const
Returns selection status.
void setCutOnLocalZ(double minZ, double maxZ)
Sets cut on local z coordinate (module frame) of the track extrapolated to TOP.
Definition: TrackSelector.h:81
const Const::ChargedStable & getChargedStable() const
Returns track hypothesis.
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 void collect() final
Replacement for event().
virtual void prepare() final
Replacement for initialize().
Abstract base class for different kinds of events.
STL namespace.
double logL
extended log likelihood