Belle II Software development
TOPPDFCheckerModule.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/TOPPDFChecker/TOPPDFCheckerModule.h>
11
12// TOP headers.
13#include <top/geometry/TOPGeometryPar.h>
14#include <top/reconstruction_cpp/TOPTrack.h>
15#include <top/reconstruction_cpp/PDFConstructor.h>
16#include <top/reconstruction_cpp/TOPRecoManager.h>
17
18// DataStore classes
19#include <tracking/dataobjects/ExtHit.h>
20#include <top/dataobjects/TOPBarHit.h>
21
22// framework aux
23#include <framework/gearbox/Unit.h>
24#include <framework/gearbox/Const.h>
25#include <framework/logging/Logger.h>
26
27// ROOT
28#include <TRandom.h>
29#include <iostream>
30
31
32using namespace std;
33using namespace ROOT::Math;
34
35namespace Belle2 {
40 using namespace TOP;
41
42 //-----------------------------------------------------------------
44 //-----------------------------------------------------------------
45
46 REG_MODULE(TOPPDFChecker);
47
48 //-----------------------------------------------------------------
49 // Implementation
50 //-----------------------------------------------------------------
51
53
54 {
55 // set module description
56 setDescription("Module for checking analytic PDF used in likelihood calculation");
58
59 // Add parameters
60 addParam("minTime", m_minTime,
61 "histogram lower bound in time [ns]", 0.0);
62 addParam("maxTime", m_maxTime,
63 "histogram upper bound in time [ns]", 50.0);
64 addParam("numBins", m_numBins,
65 "histogram number of bins in time", 1000);
66
67 }
68
70 {
71
72 // time vs. pixel ID
73 m_hits = new TH2F("hits", "photon hits", 512, 0.5, 512.5,
75 m_hits->SetXTitle("pixel ID");
76 m_hits->SetYTitle("time [ns]");
77
78 m_pdf = new TH2F("pdf", "PDF", 512, 0.5, 512.5,
80 m_pdf->SetXTitle("pixel ID");
81 m_pdf->SetYTitle("time [ns]");
82
83 // time vs pixel column
84 m_hitsCol = new TH2F("hitsCol", "photon hits", 64, 0.5, 64.5,
86 m_hitsCol->SetXTitle("pixel column");
87 m_hitsCol->SetYTitle("time [ns]");
88
89 m_pdfCol = new TH2F("pdfCol", "PDF", 64, 0.5, 64.5,
91 m_pdfCol->SetXTitle("pixel column");
92 m_pdfCol->SetYTitle("time [ns]");
93
94 }
95
97 {
98 // Register histograms (calls back defineHisto)
99 REG_HISTOGRAM;
100
101 // input
102
103 m_digits.isRequired();
104 m_tracks.isRequired();
105
106 StoreArray<ExtHit> extHits;
107 extHits.isRequired();
108
109 StoreArray<MCParticle> mcParticles;
110 mcParticles.isRequired();
111
112 StoreArray<TOPBarHit> barHits;
113 barHits.isRequired();
114
115 }
116
117
119 {
121 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
122
123 // loop over reconstructed tracks, call reconstruction and fill histograms
124
125 int numTrk = 0;
126 for (const auto& track : m_tracks) {
127 TOPTrack trk(track);
128 if (not trk.isValid()) continue;
129 if (not trk.getMCParticle()) continue;
130 if (not trk.getBarHit()) continue;
131
132 auto chargedStable = Const::chargedStableSet.find(abs(trk.getPDGCode()));
133 if (chargedStable == Const::invalidParticle) continue;
134
135 PDFConstructor pdfConstructor(trk, chargedStable, PDFConstructor::c_Fine);
136 if (not pdfConstructor.isValid()) continue;
137 numTrk++;
138
139 // average values - to print in terminate()
140 const auto& module = geo->getModule(trk.getModuleID());
141 m_avrgMomentum += module.momentumToLocal(trk.getExtHit()->getMomentum());
142 m_avrgPosition += static_cast<XYZVector>(module.pointToLocal(static_cast<XYZPoint>(trk.getExtHit()->getPosition())));
143 m_numTracks++;
144 m_slotIDs.emplace(trk.getModuleID());
145 m_PDGCodes.emplace(trk.getPDGCode());
146
147 // histogram photons in a slot crossed by the track
148 for (const auto& digit : m_digits) {
149 if (digit.getModuleID() == trk.getModuleID() and digit.getTime() < m_maxTime) {
150 if (not isFromThisParticle(digit, trk.getMCParticle())) continue;
151 m_hits->Fill(digit.getPixelID(), digit.getTime());
152 m_hitsCol->Fill(digit.getPixelCol(), digit.getTime());
153 }
154 }
155
156 // histogram PDF using MC approach
157 for (const auto& signalPDF : pdfConstructor.getSignalPDF()) {
158 int pixelID = signalPDF.getPixelID();
159 for (const auto& peak : signalPDF.getPDFPeaks()) {
160 double numPhot = pdfConstructor.getExpectedSignalPhotons() * peak.nph;
161 double sigma = sqrt(peak.wid);
162 for (int i = 0; i < gRandom->Poisson(numPhot); i++) {
163 double time = gRandom->Gaus(peak.t0, sigma);
164 m_pdf->Fill(pixelID, time);
165 int pixelCol = (pixelID - 1) % 64 + 1;
166 m_pdfCol->Fill(pixelCol, time);
167 }
168 }
169 }
170
171 }
172
173 if (numTrk == 0) {
174 B2WARNING("No track hitting the bars");
175 } else if (numTrk > 1) {
176 B2WARNING("More than one track hits the bars");
177 }
178
179 }
180
181
183 {
186
187 cout << "Average particle parameters at entrance to bar (in local frame):" << endl;
188 cout << " slot ID: ";
189 for (auto slot : m_slotIDs) cout << slot << " ";
190 cout << endl;
191 cout << " PDG code: ";
192 for (auto pdg : m_PDGCodes) cout << pdg << " ";
193 cout << endl;
194 cout << " position: x = " << m_avrgPosition.X()
195 << ", y = " << m_avrgPosition.Y()
196 << ", z = " << m_avrgPosition.Z() << endl;
197 cout << " momentum: p = " << m_avrgMomentum.R()
198 << ", theta = " << m_avrgMomentum.Theta() / Unit::deg
199 << ", phi = " << m_avrgMomentum.Phi() / Unit::deg << endl;
200 cout << "Number of particles: " << m_numTracks << endl;
201 cout << "Photons per particle: simulation = " << m_hits->GetSum() / m_numTracks
202 << ", PDF = " << m_pdf->GetSum() / m_numTracks << endl;
203
204 }
205
206
208 const MCParticle* particle)
209 {
210 const auto particles = digit.getRelationsWith<MCParticle>();
211 for (unsigned i = 0; i < particles.size(); ++i) {
212 if (particles[i] == particle and particles.weight(i) > 0) return true;
213 }
214 return false;
215 }
216
217
219} // end Belle2 namespace
220
const ParticleType & find(int pdg) const
Returns particle in set with given PDG code, or invalidParticle if not found.
Definition: Const.h:571
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:618
static const ParticleType invalidParticle
Invalid particle, used internally.
Definition: Const.h:681
ROOT::Math::XYZVector getMomentum() const
Get momentum at this extrapolation hit.
Definition: ExtHit.h:151
ROOT::Math::XYZVector getPosition() const
Get position of this extrapolation hit.
Definition: ExtHit.h:144
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
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
RelationVector< T > getRelationsWith(const std::string &name="", const std::string &namedRelation="") const
Get the relations between this object and another store array.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Class to store TOP digitized hits (output of TOPDigitizer or raw data unpacker) relations to TOPSimHi...
Definition: TOPDigit.h:24
TH2F * m_pdfCol
histogram of PDF projected to pixel columns
TH2F * m_hits
histogram of photon hits
double m_maxTime
histogram upper bound in time [ns]
double m_minTime
histogram lower bound in time [ns]
TH2F * m_hitsCol
histogram of photon hits projected to pixel columns
ROOT::Math::XYZPoint m_avrgPosition
average particle position at bar entrance (bar frame)
int m_numBins
number of bins in time
int m_numTracks
number of tracks
StoreArray< Track > m_tracks
collection of tracks
std::set< int > m_PDGCodes
particle PDG codes
TH2F * m_pdf
histogram of PDF
StoreArray< TOPDigit > m_digits
collection of digits
std::set< int > m_slotIDs
slot ID's that are hit by particle
ROOT::Math::XYZVector m_avrgMomentum
average particle momentum at bar entrance (bar frame)
PDF construction and log likelihood determination for a given track and particle hypothesis.
bool isValid() const
Checks the object status.
const std::vector< SignalPDF > & getSignalPDF() const
Returns signal PDF.
double getExpectedSignalPhotons() const
Returns the expected number of signal photons within the default time window.
@ c_Fine
y dependent everywhere
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 setTimeWindow(double minTime, double maxTime)
Sets time window.
Reconstructed track at TOP.
Definition: TOPTrack.h:39
int getPDGCode() const
Returns PDG code of associated MCParticle (returns 0 if none)
Definition: TOPTrack.h:228
const MCParticle * getMCParticle() const
Returns MC particle assigned to this track (if any)
Definition: TOPTrack.h:222
bool isValid() const
Checks if track is successfully constructed.
Definition: TOPTrack.h:137
const ExtHit * getExtHit() const
Returns extrapolated hit (track entrance to the bar)
Definition: TOPTrack.h:216
int getModuleID() const
Returns slot ID.
Definition: TOPTrack.h:143
const TOPBarHit * getBarHit() const
Returns bar hit of MC particle assigned to this track (if any)
Definition: TOPTrack.h:238
static const double deg
degree to radians
Definition: Unit.h:109
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void terminate() override
Termination action.
bool isFromThisParticle(const TOPDigit &digit, const MCParticle *particle)
Checks if digit comes from given MC particle.
virtual void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
Abstract base class for different kinds of events.
STL namespace.