Belle II Software development
SVDB4CommissioningPlotsModule.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#include <svd/modules/svdPerformance/SVDB4CommissioningPlotsModule.h>
9#include <framework/datastore/RelationVector.h>
10#include <time.h>
11#include <mdst/dataobjects/HitPatternVXD.h>
12#include <vxd/geometry/GeoCache.h>
13
14using namespace std;
15using namespace Belle2;
16
17REG_MODULE(SVDB4CommissioningPlots);
18
20 , m_nTracks(0), m_Pvalue(), m_mom(0), m_nSVDhits(0)
21{
22
23 setDescription("This module check performances of SVD reconstruction on Commissioning data");
24
25 addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("SVDB4CommissioningPlots_output.root"));
26
27 addParam("RecoDigitsName", m_RecoDigitName, "Name of RecoDigit Store Array.", std::string("SVDRecoDigits"));
28 addParam("ClustersName", m_ClusterName, "Name of Cluster Store Array.", std::string("SVDClusters"));
29 addParam("TrackListName", m_TrackName, "Name of Track Store Array.", std::string("Tracks"));
30 addParam("TrackFitResultListName", m_TrackFitResultName, "Name of TracksFitResult Store Array.", std::string("TrackFitResults"));
31 addParam("plotRecoDigits", m_plotRecoDigits,
32 "Set true to produce the plots for RecoDigits (false by default)", m_plotRecoDigits);
33}
34
36{
37
38 m_svdRecos.isOptional(m_RecoDigitName);
40 m_Tracks.isOptional(m_TrackName);
41 m_recoTracks.isOptional();
42 m_tfr.isOptional(m_TrackFitResultName);
43
44 B2INFO(" RecoDigits: " << m_RecoDigitName);
45 B2INFO(" Clusters: " << m_ClusterName);
46 B2INFO(" Tracks: " << m_TrackName);
47 B2INFO(" TrackFitResults: " << m_TrackFitResultName);
48
49
50
51 m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
52
53 m_ntracks = 0;
54}
55
57{
58 m_nEvents = 0;
59
60 //RECO DIGITS
61 TH1F hRecoCharge("reco_charge_L@layerL@ladderS@sensor@view",
62 "Charge of RecoDigits in @layer.@ladder.@sensor @view/@side side",
63 100, 0, 1000);
64 hRecoCharge.GetXaxis()->SetTitle("charge (ke-)");
65 h_recoCharge = new SVDHistograms<TH1F>(hRecoCharge);
66
67 TH1F hRecoEnergy("reco_energy_L@layerL@ladderS@sensor@view",
68 "Energy of RecoDigits in @layer.@ladder.@sensor @view/@side side",
69 100, 0, 360);
70 hRecoEnergy.GetXaxis()->SetTitle("energy (keV)");
71 h_recoEnergy = new SVDHistograms<TH1F>(hRecoEnergy);
72
73 TH1F hRecoTime("reco_time_L@layerL@ladderS@sensor@view",
74 "Time of RecoDigits in @layer.@ladder.@sensor @view/@side side",
75 200, -100, 100);
76 hRecoTime.GetXaxis()->SetTitle("time (ns)");
77 h_recoTime = new SVDHistograms<TH1F>(hRecoTime);
78
79 TH1F hRecoNoise("reco_noise_L@layerL@ladderS@sensor@view",
80 "Noise of RecoDigits in @layer.@ladder.@sensor @view/@side side",
81 200, 300, 1800);
82 hRecoNoise.GetXaxis()->SetTitle("strip noise (e-)");
83 h_recoNoise = new SVDHistograms<TH1F>(hRecoNoise);
84
85 //CLUSTER NOT RELATED TO TRACKS
86 TH1F hClusterCharge("cluster_charge_L@layerL@ladderS@sensor@view",
87 "Charge of Clusters in @layer.@ladder.@sensor @view/@side side",
88 100, 0, 100);
89 hClusterCharge.GetXaxis()->SetTitle("charge (ke-)");
90 h_clusterCharge = new SVDHistograms<TH1F>(hClusterCharge);
91
92 TH1F hClusterSize("cluster_size_L@layerL@ladderS@sensor@view",
93 "Clusters Size for @layer.@ladder.@sensor @view/@side side",
94 20, 0, 20);
95 hClusterSize.GetXaxis()->SetTitle("cluster size");
96 h_clusterSize = new SVDHistograms<TH1F>(hClusterSize);
97
98 TH1F hClusterSNR("cluster_SNR_L@layerL@ladderS@sensor@view",
99 "SNR of Clusters in @layer.@ladder.@sensor @view/@side side",
100 100, 0, 140);
101 hClusterSNR.GetXaxis()->SetTitle("SNR");
102 h_clusterSNR = new SVDHistograms<TH1F>(hClusterSNR);
103
104 TH1F hClusterEnergy("cluster_energy_L@layerL@ladderS@sensor@view",
105 "Energy of Clusters in @layer.@ladder.@sensor @view/@side side",
106 100, 0, 360);
107 hClusterEnergy.GetXaxis()->SetTitle("energy (keV)");
108 h_clusterEnergy = new SVDHistograms<TH1F>(hClusterEnergy);
109
110 TH1F hClusterTime("cluster_time_L@layerL@ladderS@sensor@view",
111 "Time of Clusters in @layer.@ladder.@sensor @view/@side side",
112 200, -100, 100);
113 hClusterTime.GetXaxis()->SetTitle("time (ns)");
114 h_clusterTime = new SVDHistograms<TH1F>(hClusterTime);
115
116 //CLUSTER RELATED TO TRACKS
117 TH1F hClusterTrkCharge("clusterTrk_charge_L@layerL@ladderS@sensor@view",
118 "Charge of Clusters Related to Tracks in @layer.@ladder.@sensor @view/@side side",
119 100, 0, 100);
120 hClusterTrkCharge.GetXaxis()->SetTitle("charge (ke-)");
121 h_clusterTrkCharge = new SVDHistograms<TH1F>(hClusterTrkCharge);
122
123 TH1F hClusterTrkSize("clusterTrk_size_L@layerL@ladderS@sensor@view",
124 "Cluster Size for @layer.@ladder.@sensor @view/@side side",
125 20, 0, 20);
126 hClusterTrkSize.GetXaxis()->SetTitle("cluster size");
127 h_clusterTrkSize = new SVDHistograms<TH1F>(hClusterTrkSize);
128
129 TH1F hClusterTrkSNR("clusterTrk_SNR_L@layerL@ladderS@sensor@view",
130 "SNR of Clusters Related to Tracks in @layer.@ladder.@sensor @view/@side side",
131 100, 0, 140);
132 hClusterTrkSNR.GetXaxis()->SetTitle("SNR");
133 h_clusterTrkSNR = new SVDHistograms<TH1F>(hClusterTrkSNR);
134
135 TH1F hClusterTrkEnergy("clusterTrk_energy_L@layerL@ladderS@sensor@view",
136 "Energy of Clusters Related to Tracks in @layer.@ladder.@sensor @view/@side side",
137 100, 0, 360);
138 hClusterTrkEnergy.GetXaxis()->SetTitle("energy (keV)");
139 h_clusterTrkEnergy = new SVDHistograms<TH1F>(hClusterTrkEnergy);
140
141 TH1F hClusterTrkTime("clusterTrk_time_L@layerL@ladderS@sensor@view",
142 "Time of Clusters Related to Tracks in @layer.@ladder.@sensor @view/@side side",
143 200, -100, 100);
144 hClusterTrkTime.GetXaxis()->SetTitle("time (ns)");
145 h_clusterTrkTime = new SVDHistograms<TH1F>(hClusterTrkTime);
146
147 TH1F hClusterTrkInterstripPos("clusterTrk_interstripPos_L@layerL@ladderS@sensor@view",
148 "Interstrip Position of Clusters Related to Tracks in @layer.@ladder.@sensor @view/@side side",
149 400, 0, 1);
150 hClusterTrkInterstripPos.GetXaxis()->SetTitle("interstrip position");
151 h_clusterTrkInterstripPos = new SVDHistograms<TH1F>(hClusterTrkInterstripPos);
152
153 //tracks
154 m_nTracks = new TH1F("h1nTracks", "number of Tracks per event", 50, 0, 50);
155 m_nTracks->GetXaxis()->SetTitle("n Tracks");
156 m_Pvalue = new TH1F("h1pValue", "Tracks p value", 100, 0, 1);
157 m_Pvalue->GetXaxis()->SetTitle("p value");
158 m_mom = new TH1F("h1momentum", " Tracks Momentum", 200, 0, 10);
159 m_mom->GetXaxis()->SetTitle("p (GeV/c)");
160 m_nSVDhits = new TH1F("h1nSVDhits", "# SVD hits per track", 20, 0, 20);
161 m_nSVDhits->GetXaxis()->SetTitle("# SVD hits");
162
163
164}
165
167{
168
169 // StoreObjPtr<EventMetaData> eventMD;
170
171 m_nEvents++;
172 float c_eTOkeV = 3.6 / 1000; //keV = e * c_eTOkeV
173
174 //tracks
175 if (m_Tracks) {
176 m_nTracks->Fill(m_Tracks.getEntries());
177 m_ntracks += m_Tracks.getEntries();
178 }
179 for (Track& track : m_Tracks) {
180
181 const TrackFitResult* tfr = track.getTrackFitResultWithClosestMass(Const::pion);
182 if (tfr) {
183 m_Pvalue->Fill(tfr->getPValue());
184 m_mom->Fill(tfr->getMomentum().R());
185 m_nSVDhits->Fill((tfr->getHitPatternVXD()).getNSVDHits());
186 }
187
188 RelationVector<RecoTrack> theRC = DataStore::getRelationsWithObj<RecoTrack>(&track);
189 RelationVector<SVDCluster> svdClustersTrack = DataStore::getRelationsWithObj<SVDCluster>(theRC[0]);
190
191 for (int cl = 0 ; cl < (int)svdClustersTrack.size(); cl++) {
192
193 float clCharge = svdClustersTrack[cl]->getCharge();
194 float clEnergy = svdClustersTrack[cl]->getCharge() * c_eTOkeV;
195 int clSize = svdClustersTrack[cl]->getSize();
196 float clSN = svdClustersTrack[cl]->getSNR();
197 float clTime = svdClustersTrack[cl]->getClsTime();
198 float clPosition = svdClustersTrack[cl]->getPosition();
199 VxdID::baseType theVxdID = (VxdID::baseType)svdClustersTrack[cl]->getSensorID();
200 int side = svdClustersTrack[cl]->isUCluster();
201
202 const VXD::SensorInfoBase* aSensorInfo = &VXD::GeoCache::getInstance().getSensorInfo(theVxdID);
203 float pitch = aSensorInfo->getVPitch();
204 if (side == 1)
205 pitch = aSensorInfo->getUPitch();
206 float clInterstripPos = fmod(clPosition, pitch) / pitch;
207
208 h_clusterTrkCharge->fill(theVxdID, side, clCharge / 1000.);
209 h_clusterTrkSize->fill(theVxdID, side, clSize);
210 h_clusterTrkSNR->fill(theVxdID, side, clSN);
211 h_clusterTrkEnergy->fill(theVxdID, side, clEnergy);
212 h_clusterTrkTime->fill(theVxdID, side, clTime);
213 h_clusterTrkInterstripPos->fill(theVxdID, side, clInterstripPos);
214
215 // if(svdClustersTrack.size()%2 != 0)
216 // B2INFO(" Event Number = "<<eventMD->getEvent()<<" "<<theVxdID<<"."<<side<<" position = "<<svdClustersTrack[cl]->getPosition()<<" charge = "<<svdClustersTrack[cl]->getCharge());
217
218 }
219 }
220
221 if (m_Tracks)
222 B2DEBUG(29, "%%%%%%%% NEW EVENT, number of Tracks = " << m_Tracks.getEntries());
223
224 //reco digits
225 if (m_plotRecoDigits) {
226 if (m_svdRecos.isValid()) {
227 for (int digi = 0 ; digi < m_svdRecos.getEntries(); digi++) {
228
229 VxdID::baseType theVxdID = (VxdID::baseType)m_svdRecos[digi]->getSensorID();
230 int side = m_svdRecos[digi]->isUStrip();
231 int cellID = m_svdRecos[digi]->getCellID();
232
233 float thisNoise = m_NoiseCal.getNoiseInElectrons(theVxdID, side, cellID);
234
235 h_recoNoise->fill(theVxdID, side, thisNoise);
236 h_recoCharge->fill(theVxdID, side, m_svdRecos[digi]->getCharge() / 1000.);
237 h_recoEnergy->fill(theVxdID, side, m_svdRecos[digi]->getCharge()*c_eTOkeV);
238 h_recoTime->fill(theVxdID, side, m_svdRecos[digi]->getTime());
239 }
240 }
241 }
242
243
244 //clusters NOT related to tracks
245 for (int cl = 0 ; cl < m_svdClusters.getEntries(); cl++) {
246
247 float clCharge = m_svdClusters[cl]->getCharge();
248 float clEnergy = m_svdClusters[cl]->getCharge() * c_eTOkeV;
249 int clSize = m_svdClusters[cl]->getSize();
250 float clTime = m_svdClusters[cl]->getClsTime();
251 float clSN = m_svdClusters[cl]->getSNR();
252
253 RelationVector<RecoTrack> theRC = DataStore::getRelationsWithObj<RecoTrack>(m_svdClusters[cl]);
254
255 bool isAssigned = false;
256 // for(int r = 0; r<(int)theRC.size(); r++){
257 if (theRC.size() > 0) {
258 RelationVector<Track> theT = DataStore::getRelationsWithObj<Track>(theRC[0]);
259 if (theT.size() > 0)
260 isAssigned = true;
261 // }
262 }
263
264 if (isAssigned)
265 continue;
266
267 VxdID::baseType theVxdID = (VxdID::baseType)m_svdClusters[cl]->getSensorID();
268 int side = m_svdClusters[cl]->isUCluster();
269
270 h_clusterCharge->fill(theVxdID, side, clCharge / 1000.);
271 h_clusterSize->fill(theVxdID, side, clSize);
272 h_clusterSNR->fill(theVxdID, side, clSN);
273 h_clusterEnergy->fill(theVxdID, side, clEnergy);
274 h_clusterTime->fill(theVxdID, side, clTime);
275 }
276}
277
278
280{
281 B2INFO("SVDB4CommissioningPlotsModule::endRun(), writing the histograms");
282
283 if (m_rootFilePtr != nullptr) {
284 m_rootFilePtr->cd();
285
286 TDirectory* oldDir = gDirectory;
287
288 TDirectory* dir_track = oldDir->mkdir("tracks");
289 dir_track->cd();
290 m_nTracks->Write();
291 m_Pvalue->Write();
292 m_mom->Write();
293 m_nSVDhits->Write();
294
295 TDirectory* dir_reco = oldDir->mkdir("recoDigits");
296 TDirectory* dir_clusterAssigned = oldDir->mkdir("clusters_assigned");
297 TDirectory* dir_clusterNotAssigned = oldDir->mkdir("clusters_not_assigned");
298
300
301 for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
302 TString layerName = Form("layer%d", layer.getLayerNumber());
303 TDirectory* dir_reco_layer = dir_reco->mkdir(layerName.Data());
304 TDirectory* dir_clusterAssigned_layer = dir_clusterAssigned->mkdir(layerName.Data());
305 TDirectory* dir_clusterNotAssigned_layer = dir_clusterNotAssigned->mkdir(layerName.Data());
306 for (auto ladder : geoCache.getLadders(layer))
307 for (Belle2::VxdID sensor : geoCache.getSensors(ladder))
308 for (int view = SVDHistograms<TH1F>::VIndex ; view < SVDHistograms<TH1F>::UIndex + 1; view++) {
309
310 if (m_plotRecoDigits) {
311 dir_reco_layer->cd();
312 (h_recoCharge->getHistogram(sensor, view))->Write();
313 (h_recoEnergy->getHistogram(sensor, view))->Write();
314 (h_recoTime->getHistogram(sensor, view))->Write();
315 (h_recoNoise->getHistogram(sensor, view))->Write();
316 }
317 dir_clusterAssigned_layer->cd();
318 (h_clusterTrkCharge->getHistogram(sensor, view))->Write();
319 (h_clusterTrkSNR->getHistogram(sensor, view))->Write();
320 (h_clusterTrkSize->getHistogram(sensor, view))->Write();
321 (h_clusterTrkEnergy->getHistogram(sensor, view))->Write();
322 (h_clusterTrkTime->getHistogram(sensor, view))->Write();
323 (h_clusterTrkInterstripPos->getHistogram(sensor, view))->Write();
324 dir_clusterNotAssigned_layer->cd();
325 (h_clusterCharge->getHistogram(sensor, view))->Write();
326 (h_clusterSNR->getHistogram(sensor, view))->Write();
327 (h_clusterSize->getHistogram(sensor, view))->Write();
328 (h_clusterEnergy->getHistogram(sensor, view))->Write();
329 (h_clusterTime->getHistogram(sensor, view))->Write();
330 }
331 }
332
333 m_rootFilePtr->Close();
334
335 }
336
337}
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
SVDHistograms< TH1F > * h_clusterTrkSNR
SVDClusters SNR.
bool m_plotRecoDigits
Produce plots for SVDRecoDigits when True.
SVDHistograms< TH1F > * h_clusterTrkCharge
SVDClusters Charge.
virtual void initialize() override
check StoreArrays & create rootfile
StoreArray< SVDCluster > m_svdClusters
SVDCluster StoreArray.
virtual void event() override
fill histograms
StoreArray< TrackFitResult > m_tfr
TrackFitResult StoreArray.
SVDNoiseCalibrations m_NoiseCal
SVDNoise Calibrations db object.
virtual void endRun() override
write histogrmas
SVDHistograms< TH1F > * h_recoEnergy
SVDRecoDigits Energy.
SVDHistograms< TH1F > * h_clusterTrkEnergy
SVDClusters Energy.
SVDHistograms< TH1F > * h_clusterEnergy
SVDClusters Energy.
SVDHistograms< TH1F > * h_clusterSNR
SVDClusters SNR.
std::string m_TrackName
Track StoreArray name.
std::string m_RecoDigitName
SVDRecoDigit StoreArray name.
SVDHistograms< TH1F > * h_recoNoise
SVDRecoDigits Noise.
virtual void beginRun() override
create histograms
std::string m_TrackFitResultName
TrackFitResult name.
SVDHistograms< TH1F > * h_clusterSize
SVDClusters Size.
SVDHistograms< TH1F > * h_recoCharge
SVDRecoDigits Charge.
SVDHistograms< TH1F > * h_clusterTime
SVDClusters Time.
SVDHistograms< TH1F > * h_recoTime
SVDRecoDigits Time.
StoreArray< RecoTrack > m_recoTracks
RecoTrack StoreArray.
std::string m_ClusterName
SVDCluster StoreArray name.
StoreArray< SVDRecoDigit > m_svdRecos
SVDRecoDigit StoreArray.
SVDHistograms< TH1F > * h_clusterTrkTime
SVDClusters Time.
SVDHistograms< TH1F > * h_clusterCharge
SVDClusters Charge.
StoreArray< Track > m_Tracks
Track StoreArray.
TFile * m_rootFilePtr
pointer at root file used for storing histograms
SVDHistograms< TH1F > * h_clusterTrkSize
SVDClusters Size.
SVDHistograms< TH1F > * h_clusterTrkInterstripPos
SVDClusters InterstripPosition.
template class for SVd histograms
Definition: SVDHistograms.h:24
void fill(const VxdID &vxdID, int view, Types ... args)
fill the histogram for
Definition: SVDHistograms.h:77
H * getHistogram(const VxdID &vxdID, int view)
get a reference to the histogram for
Definition: SVDHistograms.h:56
float getNoiseInElectrons(const VxdID &sensorID, const bool &isU, const unsigned short &strip) const
This method provides the correct noise conversion into electrons, taking into account that the noise ...
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
Values of the result of a track fit with a given particle hypothesis.
double getPValue() const
Getter for Chi2 Probability of the track fit.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
Class to facilitate easy access to sensor information of the VXD like coordinate transformations or p...
Definition: GeoCache.h:39
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
Base class to provide Sensor Information for PXD and SVD.
double getUPitch(double v=0) const
Return the pitch of the sensor.
double getVPitch(double v=0) const
Return the pitch of the sensor.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
unsigned short baseType
The base integer type for VxdID.
Definition: VxdID.h:36
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
Abstract base class for different kinds of events.
STL namespace.