Belle II Software  release-06-01-15
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 "TDirectory.h"
11 
12 using namespace Belle2;
13 
14 //-----------------------------------------------------------------
15 // Register the Module
16 //-----------------------------------------------------------------
17 REG_MODULE(SVDDQMEfficiency)
18 
19 //-----------------------------------------------------------------
20 // Implementation
21 //-----------------------------------------------------------------
22 
23 SVDDQMEfficiencyModule::SVDDQMEfficiencyModule() : HistoModule(), m_geoCache(VXD::GeoCache::getInstance())
24 {
25  // Set module properties
26  setDescription("Create basic histograms to compute the average sensor efficiency.");
27 
28  // What exactly is needed for this to be true?
29  setPropertyFlags(c_ParallelProcessingCertified);
30 
31  // Parameter definitions
32  addParam("Clusters", m_svdClustersName, "name of StoreArray with SVD cluster.", std::string(""));
33  addParam("Intercepts", m_interceptsName, "name of StoreArray with SVDIntercepts.", std::string(""));
34 
35  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed.",
36  std::string("SVDEfficiency"));
37 
38  addParam("binsU", m_u_bins, "histogram bins in u direction.", int(4));
39  addParam("binsV", m_v_bins, "histogram bins in v direction.", int(6));
40 
41 
42  addParam("saveExpertHistos", m_saveExpertHistos, "if True, save additional histograms.", bool(false));
43 
44  addParam("minSVDHits", m_minSVDHits, "Number of SVD hits required in a track to be considered.", 1u);
45  addParam("minCDCHits", m_minCDCHits, "Number of CDC hits required in a track to be considered.", 20u);
46 
47  addParam("pValCut", m_pcut, "Set a cut on the track p-value.", double(0));
48  // addParam("d0Cut", m_d0cut, "|d0| smaller than the cut", double(0.5));
49  // addParam("z0Cut", m_z0cut, "|z0| smaller than the cut", double(2));
50 
51  addParam("momCut", m_momCut, "Set a cut on the track momentum.", double(0));
52  addParam("ptCut", m_ptCut, "Set a cut on the track transverse momentum.", double(1));
53 
54  addParam("fiducialU", m_fiducialU, "Fiducial Area, U direction.", float(0.5));
55  addParam("fiducialV", m_fiducialV, "Fiducial Area, V direction.", float(0.5));
56  addParam("maxHalfResidU", m_maxResidU, "half window for cluster search around intercept, U direction.", float(0.05));
57  addParam("maxHalfResidV", m_maxResidV, "half window for cluster search around intercept, V direction.", float(0.05));
58 
59 }
60 
61 
63 {
64  //calls the define histogram function
65  REG_HISTOGRAM;
66 
67  //register the required arrays
68  //Register as optional so validation for cases where they are not available still succeeds, but module will not do any meaningful work without them
70  m_intercepts.isOptional(m_interceptsName);
71  m_recoTracks.isOptional();
72 
73 }
74 
75 
77 {
78  if (!m_svdClusters.isValid()) {
79  B2INFO("SVDClusters array is missing, no SVD efficiencies");
80  return;
81  }
82  if (!m_intercepts.isValid()) {
83  B2INFO("SVDIntercepts array is missing, no SVD efficiencies");
84  return;
85  }
86 
87  if (!m_recoTracks.isValid()) {
88  B2INFO("RecoTracks array is missing, no SVD efficiencies");
89  return;
90  }
91 
92 
93  //intercepts
94  for (int inter = 0 ; inter < m_intercepts.getEntries(); inter++) {
95 
96  if (!isGoodIntercept(m_intercepts[inter]))
97  continue;
98 
99  B2DEBUG(10, "this intercept is related to a good track");
100 
101  VxdID::baseType theVxdID = (VxdID::baseType)m_intercepts[inter]->getSensorID();
102  double coorU = m_intercepts[inter]->getCoorU();
103  double coorV = m_intercepts[inter]->getCoorV();
105  int cellU = info.getUCellID(coorU);
106  int cellV = info.getVCellID(coorV);
107 
108  const VXD::SensorInfoBase& theSensorInfo = m_geoCache.getSensorInfo(theVxdID);
109  if (theSensorInfo.inside(coorU, coorV, -m_fiducialU, -m_fiducialV)) {
110 
111  //This track should be on the sensor
112  if (m_saveExpertHistos)
113  m_h_track_hits[theVxdID]->Fill(cellU, cellV);
114 
115  m_TrackHits->fill(theVxdID, 0, 1);
116  m_TrackHits->fill(theVxdID, 1, 1);
117 
118  bool foundU = false;
119  bool foundV = false;
120 
121  //loop on clusters
122  for (int cls = 0 ; cls < m_svdClusters.getEntries(); cls++) {
123 
124  VxdID::baseType clVxdID = (VxdID::baseType)m_svdClusters[cls]->getSensorID();
125  if (clVxdID != theVxdID)
126  continue;
127 
128  double maxResid = m_maxResidV;
129  double interCoor = coorV;
130  double resid = interCoor - m_svdClusters[cls]->getPosition();
131  if (m_svdClusters[cls]->isUCluster()) {
132  interCoor = coorU;
133  maxResid = m_maxResidU;
134  resid = interCoor - m_svdClusters[cls]->getPosition(coorV);
135  }
136 
137 
138  if (abs(resid) < maxResid) {
139  if (m_svdClusters[cls]->isUCluster()) {
140  foundU = true;
141  } else
142  foundV = true;
143  }
144 
145  if (foundU && foundV)
146  break;
147  }
148 
149  if (foundU) {
150  m_MatchedHits->fill(theVxdID, 1, 1);
151  if (m_saveExpertHistos) m_h_matched_clusterU[theVxdID]->Fill(cellU, cellV);
152  }
153 
154  if (foundV) {
155  m_MatchedHits->fill(theVxdID, 0, 1);
156  if (m_saveExpertHistos)m_h_matched_clusterV[theVxdID]->Fill(cellU, cellV);
157  }
158 
159  }
160  }
161 
162 }
163 
164 
165 
167 {
168  // Create a separate histogram directory and cd into it.
169  TDirectory* oldDir = gDirectory;
170  if (m_histogramDirectoryName != "") {
171  oldDir->mkdir(m_histogramDirectoryName.c_str());
172  oldDir->cd(m_histogramDirectoryName.c_str());
173  }
174  m_TrackHits = new SVDSummaryPlots("TrackHits@view", "Number of Tracks intercepting the @view/@side Side");
175  m_MatchedHits = new SVDSummaryPlots("MatchedHits@view", "Number of Matched Clusters on the @view/@side Side");
176 
177  if (!m_saveExpertHistos) {
178  oldDir->cd();
179  return;
180  }
181 
182  std::vector<VxdID> sensors = m_geoCache.getListOfSensors();
183  for (VxdID& avxdid : sensors) {
185  if (info.getType() != VXD::SensorInfoBase::SVD) continue;
186  //Only interested in SVD sensors
187 
188  TString buff = (std::string)avxdid;
189  buff.ReplaceAll(".", "_");
190 
191  int nu = info.getUCells();
192  int nv = info.getVCells();
193 
194  m_h_track_hits[avxdid] = new TH2D("track_hits_" + buff, "tracks through sensor " + buff,
195  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
196  m_h_matched_clusterU[avxdid] = new TH2D("matched_clusterU_" + buff, "track intersections with a matched U cluster" + buff,
197  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
198  m_h_matched_clusterV[avxdid] = new TH2D("matched_clusterV_" + buff, "track intersections with a matched V cluster" + buff,
199  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
200  }
201  // cd back to root directory
202  oldDir->cd();
203 }
204 
205 
206 
208 {
209 
210  RelationVector<RecoTrack> theRC = DataStore::getRelationsWithObj<RecoTrack>(inter);
211 
212  if (theRC.size() == 0)
213  return false;
214 
215 
216  //If fit failed assume position pointed to is useless anyway
217  if (!theRC[0]->wasFitSuccessful()) return false;
218 
219  if (theRC[0]->getNumberOfSVDHits() < m_minSVDHits) return false;
220 
221  if (theRC[0]->getNumberOfCDCHits() < m_minCDCHits) return false;
222 
223  const genfit::FitStatus* fitstatus = theRC[0]->getTrackFitStatus();
224  if (fitstatus->getPVal() < m_pcut) return false;
225  // if (fabs(fitstatus->getD0() < m_d0cut) return false;
226  // if (fabs(fitstatus->getZ0() < m_z0cut) return false;
227 
228 
229 
230  genfit::MeasuredStateOnPlane trackstate;
231  trackstate = theRC[0]->getMeasuredStateOnPlaneFromFirstHit();
232  if (trackstate.getMom().Mag() < m_momCut) return false;
233 
234  if (trackstate.getMom().Perp() < m_ptCut) return false;
235 
236  return true;
237 
238 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
Creates the basic histograms for SVD Efficiency DQM.
float m_fiducialV
stay away from the U border by m_fiducialU (in cm)
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
StoreArray< SVDIntercept > m_intercepts
SVDIntercept StoreArray.
VXD::GeoCache & m_geoCache
BelleII Geometry.
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
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
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.
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
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:58
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:66
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
virtual double getPVal() const
Get the p value of the fit.
Definition: FitStatus.h:128
#StateOnPlane with additional covariance matrix.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.