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