Belle II Software  release-06-00-14
PXDClusterChargeCollectorModule.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 <pxd/modules/pxdClusterChargeCollector/PXDClusterChargeCollectorModule.h>
10 
11 #include <framework/datastore/RelationArray.h>
12 
13 #include <vxd/geometry/GeoCache.h>
14 #include <pxd/geometry/SensorInfo.h>
15 #include <pxd/reconstruction/PXDGainCalibrator.h>
16 
17 #include <TTree.h>
18 #include <TH2I.h>
19 
20 #include <boost/format.hpp>
21 
22 using namespace std;
23 using boost::format;
24 using namespace Belle2;
25 
26 //-----------------------------------------------------------------
27 // Register the Module
28 //-----------------------------------------------------------------
29 REG_MODULE(PXDClusterChargeCollector)
30 
31 //-----------------------------------------------------------------
32 // Implementation
33 //-----------------------------------------------------------------
34 
36  , m_signal(0), m_run(0), m_exp(0)
37 {
38  // Set module properties
39  setDescription("Calibration collector module for cluster charge related calibrations.");
40  setPropertyFlags(c_ParallelProcessingCertified);
41 
42  addParam("clustersName", m_storeClustersName, "Name of the collection to use for PXDClusters", string(""));
43  addParam("mcParticlesName", m_storeMCParticlesName, "Name of the collection to use for MCParticles", string("MCParticles"));
44  addParam("minClusterCharge", m_minClusterCharge, "Minimum cluster charge cut", int(0));
45  addParam("minClusterSize", m_minClusterSize, "Minimum cluster size cut ", int(2));
46  addParam("maxClusterSize", m_maxClusterSize, "Maximum cluster size cut ", int(6));
47  addParam("nBinsU", m_nBinsU, "Number of gain corrections per sensor along u side", int(4));
48  addParam("nBinsV", m_nBinsV, "Number of gain corrections per sensor along v side", int(6));
49  addParam("chargePayloadName", m_chargeName, "Payload name for Cluster Charge to be read from DB", string(""));
50  addParam("gainPayloadName", m_gainName, "Payload name for Gain to be read from DB", string(""));
51  addParam("fillChargeHistogram", m_fillChargeHistogram, "Flag to fill Charge histograms", bool(false));
52  addParam("matchTrack", m_matchTrack,
53  "Flag to use track matched clusters (=1) and apply theta angle projection to cluster charge (=2)", int(0));
54  addParam("mcSamples", m_mcSamples, "Flag to deal with MC samples", bool(false));
55  addParam("relationCheck", m_relationCheck, "Flag to check relations between PXDClusters and PXDClustersFromTracks", bool(false));
56 
57 }
58 
59 void PXDClusterChargeCollectorModule::prepare() // Do your initialise() stuff here
60 {
61  m_pxdClusters.isRequired(m_storeClustersName);
62  if (m_relationCheck) {
63  std::string storeClustersName2;
64  storeClustersName2 = (m_storeClustersName == "PXDClusters") ?
65  "PXDClustersFromTracks" : "PXDClusters";
66  StoreArray<PXDCluster> pxdClusters2;
67  pxdClusters2.isRequired(storeClustersName2);
68  }
69  m_tracks.isOptional(); //m_storeTracksName);
70  m_recoTracks.isOptional(); //m_storeRecoTracksName);
71 
72  if (m_mcSamples && m_matchTrack > 0) {
73  m_mcParticles.isRequired(m_storeMCParticlesName);
74  } else {
75  m_mcParticles.isOptional();
76  }
77 
78  if (m_nBinsU == 0) {
79  B2WARNING("Number of bins along u side incremented from 0->1");
80  m_nBinsU = 1;
81  }
82 
83  if (m_nBinsV == 0) {
84  B2WARNING("Number of bins along v side incremented from 0->1");
85  m_nBinsV = 1;
86  }
87 
88  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
89  int nPXDSensors = gTools->getNumberOfPXDSensors();
90 
91  //-------------------------------------------------------------------------------------
92  // PXDClusterCounter: Count the number of PXDClusters and store charge for each uBin/vBin pair
93  //-------------------------------------------------------------------------------------
94 
95  auto hPXDClusterCounter = new TH1I("hPXDClusterCounter", "Number of clusters found in data sample",
96  m_nBinsU * m_nBinsV * nPXDSensors, 0,
97  m_nBinsU * m_nBinsV * nPXDSensors);
98  hPXDClusterCounter->GetXaxis()->SetTitle("bin id");
99  hPXDClusterCounter->GetYaxis()->SetTitle("Number of clusters");
100  for (int iSensor = 0; iSensor < nPXDSensors; iSensor++) {
101  for (int uBin = 0; uBin < m_nBinsU; uBin++) {
102  for (int vBin = 0; vBin < m_nBinsV; vBin++) {
103  VxdID id = gTools->getSensorIDFromPXDIndex(iSensor);
104  string sensorDescr = id;
105  hPXDClusterCounter->GetXaxis()->SetBinLabel(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin + 1,
106  str(format("%1%_%2%_%3%") % sensorDescr % uBin % vBin).c_str());
107  }
108  }
109  }
110  registerObject<TH1I>("PXDClusterCounter", hPXDClusterCounter);
111  if (m_fillChargeHistogram) {
112  auto hPXDClusterCharge = new TH2I("hPXDClusterCharge", "Charge of clusters found in data sample",
113  m_nBinsU * m_nBinsV * nPXDSensors, 0,
114  m_nBinsU * m_nBinsV * nPXDSensors,
115  250, 0, 250);
116  hPXDClusterCharge->GetXaxis()->SetTitle("bin id");
117  hPXDClusterCharge->GetYaxis()->SetTitle("Cluster charge [ADU]");
118  registerObject<TH2I>("PXDClusterCharge", hPXDClusterCharge);
119  }
120 
121  //----------------------------------------------------------------------
122  // PXDTrees: One tree to store the calibration data for each grid bin
123  //----------------------------------------------------------------------
124 
125  for (int iSensor = 0; iSensor < nPXDSensors; iSensor++) {
126  for (int uBin = 0; uBin < m_nBinsU; uBin++) {
127  for (int vBin = 0; vBin < m_nBinsV; vBin++) {
128  VxdID id = gTools->getSensorIDFromPXDIndex(iSensor);
129  auto layerNumber = id.getLayerNumber();
130  auto ladderNumber = id.getLadderNumber();
131  auto sensorNumber = id.getSensorNumber();
132  string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
133  auto tree = new TTree(treename.c_str(), treename.c_str());
134  tree->Branch<int>("signal", &m_signal);
135  registerObject<TTree>(treename, tree);
136  }
137  }
138  }
139 
140  auto dbtree = new TTree("dbtree", "dbtree");
141  dbtree->Branch<int>("run", &m_run);
142  dbtree->Branch<int>("exp", &m_exp);
143  dbtree->Branch<PXDClusterChargeMapPar>("chargeMap", &m_chargeMap);
144  dbtree->Branch<PXDGainMapPar>("gainMap", &m_gainMap);
145  registerObject<TTree>("dbtree", dbtree);
146 }
147 
148 void PXDClusterChargeCollectorModule::startRun() // Do your beginRun() stuff here
149 {
150  m_run = m_evtMetaData->getRun();
151  m_exp = m_evtMetaData->getExperiment();
152  if (m_chargeName.length()) {
153  DBObjPtr<PXDClusterChargeMapPar> chargeMap(m_chargeName);
154  m_chargeMap = *chargeMap;
155  } else {
156  m_chargeMap = PXDClusterChargeMapPar();
157  }
158  if (m_gainName.length()) {
159  DBObjPtr<PXDGainMapPar> gainMap(m_gainName);
160  m_gainMap = *gainMap;
161  } else {
162  m_gainMap = PXDGainMapPar();
163  }
164  getObjectPtr<TTree>("dbtree")->Fill();
165 }
166 
167 void PXDClusterChargeCollectorModule::collect() // Do your event() stuff here
168 {
169  // If no input, nothing to do
170  if (!m_pxdClusters) return;
171 
172  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
173 
174  if (!m_matchTrack || m_mcSamples) { // all clusters for data
175 
176  for (auto& cluster : m_pxdClusters) {
177  if (m_relationCheck) {
178  if (m_storeClustersName == "PXDClusters") {
179  RelationVector<PXDCluster> relPXDClusters = cluster.getRelationsTo<PXDCluster>("PXDClustersFromTracks");
180  if (!relPXDClusters.size()) continue;
181  } else if (m_storeClustersName == "PXDClustersFromTracks") {
182  RelationVector<PXDCluster> relPXDClusters = DataStore::getRelationsWithObj<PXDCluster>(&cluster, "PXDClusters");
183  if (!relPXDClusters.size()) continue;
184  }
185  }
186  double correction = 1.0;
187  if (m_matchTrack > 0) { // mc samples
188  RelationVector<MCParticle> mcParticlesPXD = DataStore::getRelationsWithObj<MCParticle>(&cluster, m_storeMCParticlesName);
189  if (!mcParticlesPXD.size()) continue;
190  if (m_matchTrack == 2) {
191  correction = sin(mcParticlesPXD[0]->getMomentum().Theta());
192  }
193  }
194  // Apply cluster selection cuts
195  if (cluster.getCharge() >= m_minClusterCharge && cluster.getSize() >= m_minClusterSize && cluster.getSize() <= m_maxClusterSize) {
196 
197  VxdID sensorID = cluster.getSensorID();
198  const PXD::SensorInfo& Info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
199  auto uID = Info.getUCellID(cluster.getU());
200  auto vID = Info.getVCellID(cluster.getV());
201  auto iSensor = gTools->getPXDSensorIndex(sensorID);
202  auto layerNumber = sensorID.getLayerNumber();
203  auto ladderNumber = sensorID.getLadderNumber();
204  auto sensorNumber = sensorID.getSensorNumber();
205  auto uBin = PXD::PXDGainCalibrator::getInstance().getBinU(sensorID, uID, vID, m_nBinsU);
206  auto vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, vID, m_nBinsV);
207  string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
208 
209  // Compute variables from cluster needed for gain estimation
210  m_signal = cluster.getCharge();
211  // Fill variabels into tree
212  getObjectPtr<TTree>(treename)->Fill();
213  // Increment the counter & store charge (optional)
214  getObjectPtr<TH1I>("PXDClusterCounter")->Fill(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin);
215  if (m_fillChargeHistogram)
216  getObjectPtr<TH2I>("PXDClusterCharge")->Fill(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin,
217  cluster.getCharge()*correction);
218  }
219 
220  }
221 
222  } else { // store only clusters matched with tracks
223 
224  // If no input, nothing to do
225  if (!m_tracks) return;
226 
227  for (const Track& track : m_tracks) {
228 
229  RelationVector<RecoTrack> recoTrack = track.getRelationsTo<RecoTrack>(m_storeRecoTracksName);
230  if (!recoTrack.size()) continue;
231  RelationVector<PXDCluster> pxdClustersTrack = DataStore::getRelationsWithObj<PXDCluster>(recoTrack[0], m_storeClustersName);
232 
233  const TrackFitResult* tfr = track.getTrackFitResultWithClosestMass(Const::pion);
234  double correction = 1.0;
235  if (tfr && m_matchTrack == 2) correction = sin(tfr->getMomentum().Theta());
236 
237  for (auto& cluster : pxdClustersTrack) {
238  if (m_relationCheck) {
239  if (m_storeClustersName == "PXDClustersFromTracks") {
240  RelationVector<PXDCluster> relPXDClusters = DataStore::getRelationsWithObj<PXDCluster>(&cluster, "PXDClusters");
241  if (!relPXDClusters.size()) continue;
242  } else {
243  B2WARNING("Relation check for data only works when clustersName is set to PXDClustersFromTracks");
244  }
245  }
246 
247  // Apply cluster selection cuts
248  if (cluster.getCharge() >= m_minClusterCharge && cluster.getSize() >= m_minClusterSize && cluster.getSize() <= m_maxClusterSize) {
249 
250  VxdID sensorID = cluster.getSensorID();
251  const PXD::SensorInfo& Info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
252  auto uID = Info.getUCellID(cluster.getU());
253  auto vID = Info.getVCellID(cluster.getV());
254  auto iSensor = gTools->getPXDSensorIndex(sensorID);
255  auto layerNumber = sensorID.getLayerNumber();
256  auto ladderNumber = sensorID.getLadderNumber();
257  auto sensorNumber = sensorID.getSensorNumber();
258  auto uBin = PXD::PXDGainCalibrator::getInstance().getBinU(sensorID, uID, vID, m_nBinsU);
259  auto vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, vID, m_nBinsV);
260  string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
261 
262  // Compute variables from cluster needed for gain estimation
263  m_signal = cluster.getCharge();
264  // Fill variabels into tree
265  getObjectPtr<TTree>(treename)->Fill();
266  // Increment the counter & store charge (optional)
267  getObjectPtr<TH1I>("PXDClusterCounter")->Fill(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin);
268  if (m_fillChargeHistogram)
269  getObjectPtr<TH2I>("PXDClusterCharge")->Fill(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin,
270  cluster.getCharge()*correction);
271 
272  }
273 
274  }
275 
276  }
277 
278  }
279 
280 }
Calibration collector module base class.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Collector module for PXD gain calibration.
The payload class for PXD cluster charge calibrations.
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
The payload class for PXD gain corrections.
Definition: PXDGainMapPar.h:43
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:76
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Values of the result of a track fit with a given particle hypothesis.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Class that bundles various TrackFitResults.
Definition: Track.h:25
int getVCellID(double v, bool clamp=false) const
Return the corresponding pixel/strip ID of a given v coordinate.
int getUCellID(double u, double v=0, bool clamp=false) const
Return the corresponding pixel/strip ID of a given u coordinate.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:100
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
#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.