Belle II Software development
SVDDQMEfficiencyModule.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 <svd/modules/svdDQM/SVDDQMEfficiencyModule.h>
10#include <svd/dataobjects/SVDEventInfo.h>
11
12#include "TDirectory.h"
13
14#include <cmath>
15
16using namespace Belle2;
17
18//-----------------------------------------------------------------
19// Register the Module
20//-----------------------------------------------------------------
21REG_MODULE(SVDDQMEfficiency);
22
23//-----------------------------------------------------------------
24// Implementation
25//-----------------------------------------------------------------
26
28{
29 // Set module properties
30 setDescription("Create basic histograms to compute the average sensor efficiency.");
31
32 // What exactly is needed for this to be true?
34
35 // Parameter definitions
36 addParam("Clusters", m_svdClustersName, "name of StoreArray with SVD cluster.", std::string(""));
37 addParam("Intercepts", m_interceptsName, "name of StoreArray with SVDIntercepts.", std::string(""));
38
39 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed.",
40 std::string("SVDEfficiency"));
41
42 addParam("binsU", m_u_bins, "histogram bins in u direction.", int(4));
43 addParam("binsV", m_v_bins, "histogram bins in v direction.", int(6));
44
45
46 addParam("saveExpertHistos", m_saveExpertHistos, "if True, save additional histograms.", bool(false));
47
48 addParam("minSVDHits", m_minSVDHits, "Number of SVD hits required in a track to be considered.", 1u);
49 addParam("minCDCHits", m_minCDCHits, "Number of CDC hits required in a track to be considered.", 20u);
50
51 addParam("pValCut", m_pcut, "Set a cut on the track p-value.", double(0));
52 // addParam("d0Cut", m_d0cut, "|d0| smaller than the cut", double(0.5));
53 // addParam("z0Cut", m_z0cut, "|z0| smaller than the cut", double(2));
54
55 addParam("momCut", m_momCut, "Set a cut on the track momentum.", double(0));
56 addParam("ptCut", m_ptCut, "Set a cut on the track transverse momentum.", double(1));
57
58 addParam("fiducialU", m_fiducialU, "Fiducial Area, U direction.", float(0.5));
59 addParam("fiducialV", m_fiducialV, "Fiducial Area, V direction.", float(0.5));
60 addParam("maxHalfResidU", m_maxResidU, "half window for cluster search around intercept, U direction.", float(0.05));
61 addParam("maxHalfResidV", m_maxResidV, "half window for cluster search around intercept, V direction.", float(0.05));
62 addParam("samples3", m_3Samples, "if True 3 samples histograms analysis is performed", bool(false));
63}
64
65
67{
68 //calls the define histogram function
69 REG_HISTOGRAM;
70
71 //register the required arrays
72 //Register as optional so validation for cases where they are not available still succeeds, but module will not do any meaningful work without them
75 m_recoTracks.isOptional();
76
77}
78
79
81{
82 if (!m_svdClusters.isValid()) {
83 B2INFO("SVDClusters array is missing, no SVD efficiencies");
84 return;
85 }
86 if (!m_intercepts.isValid()) {
87 B2INFO("SVDIntercepts array is missing, no SVD efficiencies");
88 return;
89 }
90
91 if (!m_recoTracks.isValid()) {
92 B2INFO("RecoTracks array is missing, no SVD efficiencies");
93 return;
94 }
95
96 std::string m_svdEventInfoName = "SVDEventInfo";
97 StoreObjPtr<SVDEventInfo> eventinfo(m_svdEventInfoName);
98
99 int nSamples = 0;
100 if (eventinfo.isValid())
101 nSamples = eventinfo->getNSamples();
102 else
103 return;
104
105 //intercepts
106 for (int inter = 0 ; inter < m_intercepts.getEntries(); inter++) {
107
108 if (!isGoodIntercept(m_intercepts[inter]))
109 continue;
110
111 B2DEBUG(10, "this intercept is related to a good track");
112
113 VxdID::baseType theVxdID = (VxdID::baseType)m_intercepts[inter]->getSensorID();
114 double coorU = m_intercepts[inter]->getCoorU();
115 double coorV = m_intercepts[inter]->getCoorV();
116 VXD::SensorInfoBase info = m_geoCache.getSensorInfo(theVxdID);
117 int cellU = info.getUCellID(coorU);
118 int cellV = info.getVCellID(coorV);
119
120 const VXD::SensorInfoBase& theSensorInfo = m_geoCache.getSensorInfo(theVxdID);
121 if (theSensorInfo.inside(coorU, coorV, -m_fiducialU, -m_fiducialV)) {
122
123 //This track should be on the sensor
125 m_h_track_hits[theVxdID]->Fill(cellU, cellV);
126
127 m_TrackHits->fill(theVxdID, 0, 1);
128 m_TrackHits->fill(theVxdID, 1, 1);
129
130 if (m_3Samples) {
131 if (nSamples == 3) {
132 m_TrackHits3Sample->fill(theVxdID, 0, 1);
133 m_TrackHits3Sample->fill(theVxdID, 1, 1);
134 } else {
135 m_TrackHits6Sample->fill(theVxdID, 0, 1);
136 m_TrackHits6Sample->fill(theVxdID, 1, 1);
137 }
138 }
139
140 bool foundU = false;
141 bool foundV = false;
142
143 //loop on clusters
144 for (int cls = 0 ; cls < m_svdClusters.getEntries(); cls++) {
145
146 VxdID::baseType clVxdID = (VxdID::baseType)m_svdClusters[cls]->getSensorID();
147 if (clVxdID != theVxdID)
148 continue;
149
150 double maxResid = m_maxResidV;
151 double interCoor = coorV;
152 double resid = interCoor - m_svdClusters[cls]->getPosition();
153 if (m_svdClusters[cls]->isUCluster()) {
154 interCoor = coorU;
155 maxResid = m_maxResidU;
156 resid = interCoor - m_svdClusters[cls]->getPosition(coorV);
157 }
158
159
160 if (std::abs(resid) < maxResid) {
161 if (m_svdClusters[cls]->isUCluster()) {
162 foundU = true;
163 } else
164 foundV = true;
165 }
166
167 if (foundU && foundV)
168 break;
169 }
170
171 if (foundU) {
172 m_MatchedHits->fill(theVxdID, 1, 1);
173 if (m_3Samples) {
174 if (nSamples == 3)
175 m_MatchedHits3Sample->fill(theVxdID, 1, 1);
176 else
177 m_MatchedHits6Sample->fill(theVxdID, 1, 1);
178 }
179
180 if (m_saveExpertHistos) {
181 m_h_matched_clusterU[theVxdID]->Fill(cellU, cellV);
182 if (m_3Samples) {
183 if (nSamples == 3)
184 m_h_matched3_clusterU[theVxdID]->Fill(cellU, cellV);
185 else
186 m_h_matched6_clusterU[theVxdID]->Fill(cellU, cellV);
187 }
188 }
189 }
190
191 if (foundV) {
192 m_MatchedHits->fill(theVxdID, 0, 1);
193 if (m_3Samples) {
194 if (nSamples == 3)
195 m_MatchedHits3Sample->fill(theVxdID, 0, 1);
196 else
197 m_MatchedHits6Sample->fill(theVxdID, 0, 1);
198 }
199
200 if (m_saveExpertHistos) {
201 m_h_matched_clusterV[theVxdID]->Fill(cellU, cellV);
202 if (m_3Samples) {
203 if (nSamples == 3)
204 m_h_matched3_clusterV[theVxdID]->Fill(cellU, cellV);
205 else
206 m_h_matched6_clusterV[theVxdID]->Fill(cellU, cellV);
207 }
208 }
209 }
210
211 }
212 }
213}
214
216{
217 // Create a separate histogram directory and cd into it.
218 TDirectory* oldDir = gDirectory;
219 if (m_histogramDirectoryName != "") {
220 oldDir->mkdir(m_histogramDirectoryName.c_str());
221 oldDir->cd(m_histogramDirectoryName.c_str());
222 }
223 m_TrackHits = new SVDSummaryPlots("TrackHits@view", "Number of Tracks intercepting the @view/@side Side");
224 m_MatchedHits = new SVDSummaryPlots("MatchedHits@view", "Number of Matched Clusters on the @view/@side Side");
225
226 if (m_3Samples) {
227 m_TrackHits3Sample = new SVDSummaryPlots("TrackHits3@view", "Number of Tracks intercepting the @view/@side Side for 3 samples");
228 m_TrackHits6Sample = new SVDSummaryPlots("TrackHits6@view", "Number of Tracks intercepting the @view/@side Side for 6 samples");
229 m_MatchedHits3Sample = new SVDSummaryPlots("MatchedHits3@view", "Number of Matched Clusters on the @view/@side Side for 3 samples");
230 m_MatchedHits6Sample = new SVDSummaryPlots("MatchedHits6@view", "Number of Matched Clusters on the @view/@side Side for 6 samples");
231 }
232
233 if (!m_saveExpertHistos) {
234 oldDir->cd();
235 return;
236 }
237
238 std::vector<VxdID> sensors = m_geoCache.getListOfSensors();
239 for (VxdID& avxdid : sensors) {
240 VXD::SensorInfoBase info = m_geoCache.getSensorInfo(avxdid);
241 if (info.getType() != VXD::SensorInfoBase::SVD) continue;
242 //Only interested in SVD sensors
243
244 TString buff = (std::string)avxdid;
245 buff.ReplaceAll(".", "_");
246
247 int nu = info.getUCells();
248 int nv = info.getVCells();
249
250 m_h_track_hits[avxdid] = new TH2D("track_hits_" + buff, "tracks through sensor " + buff,
251 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
252 m_h_matched_clusterU[avxdid] = new TH2D("matched_clusterU_" + buff, "track intersections with a matched U cluster" + buff,
253 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
254 m_h_matched_clusterV[avxdid] = new TH2D("matched_clusterV_" + buff, "track intersections with a matched V cluster" + buff,
255 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
256
257 if (m_3Samples) {
258 m_h_matched3_clusterU[avxdid] = new TH2D("matched3_clusterU_" + buff,
259 "track intersections with a matched U cluster for 3 samples" + buff,
260 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
261 m_h_matched3_clusterV[avxdid] = new TH2D("matched3_clusterV_" + buff,
262 "track intersections with a matched V cluster for 3 samples" + buff,
263 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
264
265 m_h_matched6_clusterU[avxdid] = new TH2D("matched6_clusterU_" + buff,
266 "track intersections with a matched U cluster for 6 samples" + buff,
267 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
268 m_h_matched6_clusterV[avxdid] = new TH2D("matched6_clusterV_" + buff,
269 "track intersections with a matched V cluster for 6 samples" + buff,
270 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
271 }
272 }
273 // cd back to root directory
274 oldDir->cd();
275}
276
277
278
280{
281
283
284 if (theRC.size() == 0)
285 return false;
286
287 //If fit failed assume position pointed to is useless anyway
288 if (!theRC[0]->wasFitSuccessful()) return false;
289
290 if (theRC[0]->getNumberOfSVDHits() < m_minSVDHits) return false;
291
292 if (theRC[0]->getNumberOfCDCHits() < m_minCDCHits) return false;
293
294 const genfit::FitStatus* fitstatus = theRC[0]->getTrackFitStatus();
295 if (fitstatus->getPVal() < m_pcut) return false;
296
297 genfit::MeasuredStateOnPlane trackstate;
298 trackstate = theRC[0]->getMeasuredStateOnPlaneFromFirstHit();
299 if (trackstate.getMom().Mag() < m_momCut) return false;
300
301 if (trackstate.getMom().Perp() < m_ptCut) return false;
302
303 return true;
304
305}
static RelationVector< T > getRelationsWithObj(const TObject *object, const std::string &name="", const std::string &namedRelation="")
Get the relations between an object and other objects in a store array.
Definition DataStore.h:412
HistoModule()
Constructor.
Definition HistoModule.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
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
float m_fiducialV
stay away from the U border by m_fiducialU (in cm)
SVDDQMEfficiencyModule()
Constructor: Sets the description, the properties and the parameters of the module.
std::string m_svdClustersName
SVDClusters StoreArray name.
std::map< VxdID, TH2D * > m_h_track_hits
track hits histogram map to sensorID
void initialize() override final
initializes the need store arrays, trees and histograms
unsigned int m_minCDCHits
Required hits in CDC for tracks.
bool m_saveExpertHistos
save additional histograms id set True
std::map< VxdID, TH2D * > m_h_matched6_clusterV
matched V-hits histogram map to sensorID for 6 samples
StoreArray< SVDIntercept > m_intercepts
SVDIntercept StoreArray.
SVDSummaryPlots * m_MatchedHits3Sample
matched hits summary plot for 3 samples
VXD::GeoCache & m_geoCache
BelleII Geometry.
std::map< VxdID, TH2D * > m_h_matched3_clusterV
matched V-hits histogram map to sensorID for 3 samples
unsigned int m_minSVDHits
Required hits in SVD strips for tracks.
StoreArray< SVDCluster > m_svdClusters
SVDCluster StoreArray.
double m_momCut
Cut on fitted track momentum.
int m_u_bins
number of U-bins for expert histogram
std::map< VxdID, TH2D * > m_h_matched_clusterU
matched U-hits histogram map to sensorID
std::map< VxdID, TH2D * > m_h_matched6_clusterU
matched U-hits histogram map to sensorID for 6 samples
float m_fiducialU
stay away from the U border by m_fiducialU (in cm)
void defineHisto() override final
actually defines the trees and histograms
double m_pcut
pValue-Cut for tracks
bool isGoodIntercept(SVDIntercept *inter)
returns true if the track related to the intercept passes the selection cuts
void event() override final
main function which fills trees and histograms
std::string m_histogramDirectoryName
name of the directory where to store the histograms
float m_maxResidV
max distance cut in cm V side
bool m_3Samples
if true enable 3 samples histograms analysis
SVDSummaryPlots * m_MatchedHits6Sample
matched hits summary plot for 6 samples
std::map< VxdID, TH2D * > m_h_matched3_clusterU
matched U-hits histogram map to sensorID for 3 samples
int m_v_bins
number of V-bins for expert histogram
SVDSummaryPlots * m_MatchedHits
matched hits summary plot
SVDSummaryPlots * m_TrackHits
track hits summary plot
StoreArray< RecoTrack > m_recoTracks
RecoTrack StoreArray.
double m_ptCut
Cut on fitted track pt.
std::string m_interceptsName
SVDIntercepts StoreArray name.
SVDSummaryPlots * m_TrackHits3Sample
track hits summary plot for 3 samples
std::map< VxdID, TH2D * > m_h_matched_clusterV
matched V-hits histogram map to sensorID
SVDSummaryPlots * m_TrackHits6Sample
track hits summary plot for 6 samples
float m_maxResidU
max distance cut in cm U side
SVDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an SVD ...
class to summarize SVD quantities per sensor and side
Type-safe access to single objects in the data store.
Definition StoreObjPtr.h:96
Base class to provide Sensor Information for PXD and SVD.
bool inside(double u, double v, double uTolerance=DBL_EPSILON, double vTolerance=DBL_EPSILON) const
Check whether a given point is inside the active area.
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
Namespace to provide code needed by both Vertex Detectors, PXD and SVD, and also testbeam telescopes.
Definition GeoCache.h:34
Abstract base class for different kinds of events.