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>
16#include <framework/gearbox/Const.h>
17#include <framework/logging/Logger.h>
48 setDescription(
"A collector for common T0 calibration with dimuons or Bhabha's using "
49 "neg. log likelihood minimization (method LL)");
54 "number of bunches per SST clock period", 24);
57 "time range in which to search for the minimum [ns]", 10.0);
59 "sigma in [ns] for additional smearing of PDF", 0.0);
61 "sample type: one of dimuon or bhabha", std::string(
"dimuon"));
63 "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
67 "minimal local z of extrapolated hit", -130.0);
69 "maximal local z of extrapolated hit", 130.0);
71 "PDF option, one of 'rough', 'fine', 'optimal'", std::string(
"rough"));
98 B2ERROR(
"Unknown PDF option '" <<
m_pdfOption <<
"'");
109 B2ERROR(
"Invalid sample type '" <<
m_sample <<
"'");
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);
123 auto h1 =
new TH1F(
"tracks_per_set",
"tracks per sample; sample number; num tracks",
125 registerObject<TH1F>(
"tracks_per_set", h1);
127 auto h2 =
new TH1F(
"numHits",
"Number of photons per slot",
129 h2->SetXTitle(
"slot number");
130 h2->SetYTitle(
"hits per slot");
131 registerObject<TH1F>(
"numHits", h2);
133 auto h3 =
new TH2F(
"timeHits",
"Photon times vs. boardstacks",
135 h3->SetXTitle(
"slot number");
136 h3->SetYTitle(
"time [ns]");
137 registerObject<TH2F>(
"timeHits", h3);
141 auto h4 =
new TH1F(
"offset",
"current offset from input files; offset [ns]",
143 registerObject<TH1F>(
"offset", h4);
153 if (not
m_recBunch->isReconstructed())
return;
164 if (digit.isCommonT0Calibrated()) digit.subtractT0(-T0);
165 if (digit.hasStatus(TOPDigit::c_BunchOffsetSubtracted)) digit.subtractT0(-
m_recBunch->getAverageOffset());
171 for (
const auto& track :
m_tracks) {
175 if (not trk.
isValid())
continue;
181 if (not pdfConstructor.
isValid())
continue;
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);
190 h->SetBinContent(ibin + 1, chi);
192 auto h1 = getObjectPtr<TH1F>(
"tracks_per_set");
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());
213 auto h4 = getObjectPtr<TH1F>(
"offset");
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
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
double m_dz
cut on POCA in z
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
@ c_numModules
number of modules
double m_minZ
minimal local z of extrapolated hit
std::vector< std::string > m_names
histogram names of chi2 scans
double m_dr
cut on POCA in r
StoreArray< Track > m_tracks
collection of tracks
double m_deltaEcms
c.m.s energy window
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
std::string m_pdfOption
PDF option name.
double m_bunchTimeSep
bunch separation in time [ns]
std::string m_sample
sample type
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.
bool isValid() const
Checks if track is successfully constructed.
int getModuleID() const
Returns slot ID.
Utility for the track selection - used in various calibration modules.
void setDeltaEcms(double deltaEcms)
Sets cut on c.m.s.
void setCutOnPOCA(double dr, double dz)
Sets cut on point of closest approach to (0, 0, 0)
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.
const Const::ChargedStable & getChargedStable() const
Returns track hypothesis.
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
virtual void collect() final
Replacement for event().
TOPCommonT0LLCollectorModule()
Constructor.
virtual void prepare() final
Replacement for initialize().
Abstract base class for different kinds of events.
double logL
extended log likelihood