Belle II Software development
PXDClusterizerModule.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/pxdReconstruction/PXDClusterizerModule.h>
10
11#include <framework/datastore/DataStore.h>
12#include <framework/datastore/StoreArray.h>
13#include <framework/datastore/RelationArray.h>
14#include <framework/logging/Logger.h>
15
16#include <vxd/geometry/GeoCache.h>
17
18#include <mdst/dataobjects/MCParticle.h>
19#include <mdst/dataobjects/EventLevelTrackingInfo.h>
20#include <pxd/dataobjects/PXDDigit.h>
21#include <pxd/dataobjects/PXDCluster.h>
22#include <pxd/dataobjects/PXDTrueHit.h>
23#include <pxd/geometry/SensorInfo.h>
24
25#include <pxd/reconstruction/PXDClusterPositionEstimator.h>
26
27#include <pxd/utilities/PXDUtilities.h>
28
29using namespace std;
30using namespace Belle2;
31using namespace Belle2::PXD;
32
33//-----------------------------------------------------------------
34// Register the Module
35//-----------------------------------------------------------------
36REG_MODULE(PXDClusterizer);
37
38//-----------------------------------------------------------------
39// Implementation
40//-----------------------------------------------------------------
41
43 , m_elNoise(0.7), m_cutSeed(5.0), m_cutAdjacent(3.0), m_cutCluster(8.0)
44 , m_cutAdjacentSignal(0), m_sizeHeadTail(3), m_clusterCacheSize(0)
45{
46 //Set module properties
47 setDescription("Cluster PXDHits");
49
50 addParam("ElectronicNoise", m_elNoise,
51 "Noise added by the electronics, set in ADU", m_elNoise);
52 addParam("NoiseSN", m_cutAdjacent,
53 "SN for digits to be considered for clustering", m_cutAdjacent);
54 addParam("SeedSN", m_cutSeed, "SN for digits to be considered as seed",
55 m_cutSeed);
56 addParam("ClusterSN", m_cutCluster, "Minimum SN for clusters", m_cutCluster);
57 addParam("ClusterCacheSize", m_clusterCacheSize,
58 "Maximum desired number of sensor rows", 0);
59 addParam("HeadTailSize", m_sizeHeadTail,
60 "Minimum cluster size to switch to Analog head tail algorithm for cluster center",
62 addParam("Digits", m_storeDigitsName, "Digits collection name", string(""));
63 addParam("Clusters", m_storeClustersName, "Cluster collection name",
64 string(""));
65 addParam("TrueHits", m_storeTrueHitsName, "TrueHit collection name",
66 string(""));
67 addParam("MCParticles", m_storeMCParticlesName, "MCParticles collection name",
68 string(""));
69 addParam("ErrorFromDB", m_errorFromDB, "Assign cluster position error from DB", true);
70 addParam("PositionErrorUPayloadName", m_positionErrorUName, "Payload name for cluster position error in U",
71 string("PXDClusterPositionErrorUPar"));
72 addParam("PositionErrorVPayloadName", m_positionErrorVName, "Payload name for cluster position error in V",
73 string("PXDClusterPositionErrorVPar"));
74 addParam("createPXDClustersForAbortedTrackingEvents", m_createPXDClustersForAbortedTrackingEvents,
75 "Create PXDClusters for events where either the SVDSpacePointCreator abort flag or the VXDTF2 and SVDCKF abort flags are set.",
77
78}
79
81{
82 //Register collections
87
89 storeDigits.isRequired();
90 storeTrueHits.isOptional();
91 storeMCParticles.isOptional();
92
93 RelationArray relClusterDigits(storeClusters, storeDigits);
94 RelationArray relClusterTrueHits(storeClusters, storeTrueHits);
95 RelationArray relClusterMCParticles(storeClusters, storeMCParticles);
96 RelationArray relDigitMCParticles(storeDigits, storeMCParticles);
97 RelationArray relDigitTrueHits(storeDigits, storeTrueHits);
98
99 relClusterDigits.registerInDataStore();
100 //Relations to Simulation objects are only needed if the parent one already
101 //exists
102 if (relDigitTrueHits.isOptional()) {
103 relClusterTrueHits.registerInDataStore();
104 }
105 if (relDigitMCParticles.isOptional()) {
106 relClusterMCParticles.registerInDataStore();
107 }
108
109 //Now save all names to speed creation later, mainly if they had default names
110 m_storeClustersName = storeClusters.getName();
111 m_storeDigitsName = storeDigits.getName();
112 m_storeTrueHitsName = storeTrueHits.getName();
113 m_storeMCParticlesName = storeMCParticles.getName();
114
115 m_relClusterDigitName = relClusterDigits.getName();
116 m_relClusterTrueHitName = relClusterTrueHits.getName();
117 m_relClusterMCParticleName = relClusterMCParticles.getName();
118 m_relDigitTrueHitName = relDigitTrueHits.getName();
119 m_relDigitMCParticleName = relDigitMCParticles.getName();
120
121 B2DEBUG(20,
122 "PXDClusterizer Parameters (in default system units, *=cannot be set directly):");
123 B2DEBUG(20, " --> ElectronicNoise: " << m_elNoise);
124 B2DEBUG(20, " --> NoiseSN: " << m_cutAdjacent);
125 B2DEBUG(20, " --> SeedSN: " << m_cutSeed);
126 B2DEBUG(20, " --> ClusterSN: " << m_cutCluster);
127 B2DEBUG(20,
128 " --> MCParticles: " << DataStore::arrayName<MCParticle>(m_storeMCParticlesName));
129 B2DEBUG(20,
130 " --> Digits: " << DataStore::arrayName<PXDDigit>(m_storeDigitsName));
131 B2DEBUG(20,
132 " --> Clusters: " << DataStore::arrayName<PXDCluster>(m_storeClustersName));
133 B2DEBUG(20,
134 " --> TrueHits: " << DataStore::arrayName<PXDTrueHit>(m_storeTrueHitsName));
135 B2DEBUG(20, " --> DigitMCRel: " << m_relDigitMCParticleName);
136 B2DEBUG(20, " --> ClusterMCRel: " << m_relClusterMCParticleName);
137 B2DEBUG(20, " --> ClusterDigitRel: " << m_relClusterDigitName);
138 B2DEBUG(20, " --> DigitTrueRel: " << m_relDigitTrueHitName);
139 B2DEBUG(20, " --> ClusterTrueRel: " << m_relClusterTrueHitName);
140
141
144 if (m_clusterCacheSize > 0)
145 m_cache = std::unique_ptr<ClusterCache>(new ClusterCache(m_clusterCacheSize));
146 else
147 m_cache = std::unique_ptr<ClusterCache>(new ClusterCache());
148
149 // Cluster position error from DB
150 if (m_errorFromDB) {
151 if (!m_positionErrorUName.size() || !m_positionErrorVName.size()) {
152 B2WARNING("You chose to use cluster position erros from DB but did not provide PositionErrorUPayloadName ("
153 << m_positionErrorUName << ") and/or PositionErrorVPayloadName (" << m_positionErrorVName
154 << "). Disabling DB option.");
155 m_errorFromDB = false;
156 }
157 m_clusterPositionErrorUPar = std::make_unique<DBObjPtr<PXDClusterPositionErrorPar>>(m_positionErrorUName);
158 m_clusterPositionErrorVPar = std::make_unique<DBObjPtr<PXDClusterPositionErrorPar>>(m_positionErrorVName);
159 if (m_clusterPositionErrorUPar == nullptr || m_clusterPositionErrorVPar == nullptr) {
160 B2FATAL("DB objects for ClusterPositionError not valid");
161 }
162 }
163
164 m_eventLevelTrackingInfo.isOptional();
165
166}
167
169{
170 // Need to check if payload should be re-loaded at run change (currently just one revision for the entire period)
171 if (m_errorFromDB && (!m_clusterPositionErrorUPar->isValid() || !m_clusterPositionErrorVPar.get()->isValid())) {
172 B2FATAL("DB objects for ClusterPositionError not valid for this run");
173 }
174}
175
177{
178 // Abort in case SVDSpacePointCreator was aborted (high occupancy events) or if both VXDTF2 and SVDCKF were aborted
179 // as we do not add PXD hits to CDC standalone tracks, so we don't need to run the PXDClusterizer.
180 // This veto can be overwritten by setting m_createPXDClustersForAbortedTrackingEvents to true to create them regardless.
181 if (m_eventLevelTrackingInfo.isValid()) {
183 (m_eventLevelTrackingInfo->hasSVDSpacePointCreatorAbortionFlag() or
184 (m_eventLevelTrackingInfo->hasSVDCKFAbortionFlag() and m_eventLevelTrackingInfo->hasVXDTF2AbortionFlag()))) {
185 return;
186 }
187 }
188
189 const StoreArray<MCParticle> storeMCParticles(m_storeMCParticlesName);
190 const StoreArray<PXDTrueHit> storeTrueHits(m_storeTrueHitsName);
191 const StoreArray<PXDDigit> storeDigits(m_storeDigitsName);
193
194 storeClusters.clear();
195
196 RelationArray relClusterMCParticle(storeClusters, storeMCParticles,
198 if (relClusterMCParticle) relClusterMCParticle.clear();
199
200 RelationArray relClusterDigit(storeClusters, storeDigits,
202 if (relClusterDigit) relClusterDigit.clear();
203
204 RelationArray relClusterTrueHit(storeClusters, storeTrueHits,
206 if (relClusterTrueHit) relClusterTrueHit.clear();
207
208 int nPixels = storeDigits.getEntries();
209 if (nPixels == 0)
210 return;
211
212 //Build lookup tables for relations
213 RelationArray relDigitMCParticle(storeDigits, storeMCParticles, m_relDigitMCParticleName);
214 RelationArray relDigitTrueHit(storeDigits, storeTrueHits, m_relDigitTrueHitName);
215 createRelationLookup(relDigitMCParticle, m_mcRelation, storeDigits.getEntries());
216 createRelationLookup(relDigitTrueHit, m_trueRelation, storeDigits.getEntries());
217
218 m_cache->clear();
219
220 //We require all pixels are already sorted and directly cluster them. Once
221 //the sensorID changes, we write out all existing clusters and continue.
222 VxdID sensorID(0);
223 //To check sorting
224 Pixel lastPixel;
225 for (int i = 0; i < nPixels; i++) {
226 const PXDDigit* const storeDigit = storeDigits[i];
227 Pixel px(storeDigit, i);
228 //New sensor, write clusters
229 if (sensorID != storeDigit->getSensorID()) {
230 writeClusters(sensorID);
231 sensorID = storeDigit->getSensorID();
232 //Load the correct noise map for the new sensor
233 m_noiseMap.setSensorID(sensorID);
234 } else if (px <= lastPixel) {
235 //Check for sorting as precaution
236 B2FATAL("Pixels are not sorted correctly, please include the "
237 "PXDDigitSorter module before running the Clusterizer or fix "
238 "the input to be ordered by v,u in ascending order");
239 }
240 //Remember last pixel to check sorting
241 lastPixel = px;
242
243 //Ignore digits with not enough signal. Has to be done after the check of
244 //the SensorID to make sure we compare to the correct noise level
245 if (!m_noiseMap(px, m_cutAdjacent)) continue;
246
247 // TODO: If we would like to cluster across dead channels we would need a
248 // sorted list of dead pixels. Assuming we have such a list m_deadChannels
249 // containing Pixel instances with all dead channels of this sensor and
250 // m_currentDeadChannel is an iterator to that list (initialized to
251 // begin(m_deadChannels) when the sensorID changes), we would do
252 //
253 // while(m_currentDeadChannel != end(m_deadChannels) && *m_currentDeadChannel < px){
254 // m_cache->findCluster(m_currentDeadChannel->getU(), m_currentDeadChannel->getV());
255 // m_currentDeadChannel++;
256 // }
257 //
258 // This would have the effect of marking the pixel address as belonging
259 // to a cluster but would not modify the clusters itself (except for
260 // possible merging of clusters) so we would end up with clusters that
261 // contain holes or are disconnected.
262
263 // Find correct cluster and add pixel to cluster
264 try {
265 m_cache->findCluster(px.getU(), px.getV()).add(px);
266 } catch (std::out_of_range& e) {
267 B2WARNING("PXD clustering: Ignoring pixel " << px.getU() << "," << px.getV() << ": " << e.what());
268 }
269 }
270 writeClusters(sensorID);
271}
272
274{
275 lookup.clear();
276 //If we don't have a relation we don't build a lookuptable
277 if (!relation) return;
278 //Resize to number of digits and set all values
279 lookup.resize(digits);
280 for (const RelationElement& element : relation) {
281 lookup[element.getFromIndex()] = &element;
282 }
283}
284
285void PXDClusterizerModule::fillRelationMap(const RelationLookup& lookup, std::map<unsigned int, float>& relation,
286 unsigned int index)
287{
288 //If the lookup table is not empty and the element is set
289 if (!lookup.empty() && lookup[index]) {
290 const RelationElement& element = *lookup[index];
291 const unsigned int size = element.getSize();
292 //Add all Relations to the map
293 for (unsigned int i = 0; i < size; ++i) {
294 //negative weights are from ignored particles, we don't like them and
295 //thus ignore them :D
296 if (element.getWeight(i) < 0) continue;
297 relation[element.getToIndex(i)] += element.getWeight(i);
298 }
299 }
300}
301
303{
304 if (m_cache->empty())
305 return;
306
307 //Get all datastore elements
308 const StoreArray<MCParticle> storeMCParticles(m_storeMCParticlesName);
309 const StoreArray<PXDDigit> storeDigits(m_storeDigitsName);
310 const StoreArray<PXDTrueHit> storeTrueHits(m_storeTrueHitsName);
312 RelationArray relClusterMCParticle(storeClusters, storeMCParticles,
314 RelationArray relClusterDigit(storeClusters, storeDigits,
316 RelationArray relClusterTrueHit(storeClusters, storeTrueHits,
318
319 //Get Geometry information
320 const SensorInfo& info = dynamic_cast<const SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(
321 sensorID));
322
323 map<unsigned int, float> mc_relations;
324 map<unsigned int, float> truehit_relations;
325 vector<pair<unsigned int, float> > digit_weights;
326
327
328 for (ClusterCandidate& cls : *m_cache) {
329 //Check for noise cuts
330 if (!(cls.size() > 0 && m_noiseMap(cls.getCharge(), m_cutCluster) && m_noiseMap(cls.getSeed(), m_cutSeed))) continue;
331
332 double rho(0);
333 ClusterProjection projU, projV;
334 mc_relations.clear();
335 truehit_relations.clear();
336 digit_weights.clear();
337 digit_weights.reserve(cls.size());
338
339 const Pixel& seed = cls.getSeed();
340
341 for (const PXD::Pixel& px : cls.pixels()) {
342 //Add the Pixel information to the two projections
343 projU.add(px.getU(), info.getUCellPosition(px.getU()), px.getCharge());
344 projV.add(px.getV(), info.getVCellPosition(px.getV()), px.getCharge());
345
346 //Obtain relations from MCParticles
347 fillRelationMap(m_mcRelation, mc_relations, px.getIndex());
348 //Obtain relations from PXDTrueHits
349 fillRelationMap(m_trueRelation, truehit_relations, px.getIndex());
350 //Save the weight of the digits for the Cluster->Digit relation
351 digit_weights.emplace_back(px.getIndex(), px.getCharge());
352 }
353 projU.finalize();
354 projV.finalize();
355
356 const double pitchU = info.getUPitch();
357 const double pitchV = info.getVPitch(projV.getPos());
358 // Calculate shape correlation coefficient: only for non-trivial shapes
359 if (projU.getSize() > 1 && projV.getSize() > 1) {
360 // Add in-pixel position noise to smear the correlation
361 double posUU = cls.getCharge() * pitchU * pitchU / 12.0;
362 double posVV = cls.getCharge() * pitchV * pitchV / 12.0;
363 double posUV(0);
364 for (const Pixel& px : cls.pixels()) {
365 const double du = info.getUCellPosition(px.getU()) - projU.getPos();
366 const double dv = info.getVCellPosition(px.getV()) - projV.getPos();
367 posUU += px.getCharge() * du * du;
368 posVV += px.getCharge() * dv * dv;
369 posUV += px.getCharge() * du * dv;
370 }
371 rho = posUV / sqrt(posUU * posVV);
372 }
373
374 //Calculate position and error with u as primary axis, fixed pitch size
375 calculatePositionError(cls, projU, projV, pitchU, pitchU, pitchU);
376 //Calculate position and error with v as primary axis, possibly different pitch sizes
377 calculatePositionError(cls, projV, projU, info.getVPitch(projV.getMinPos()), pitchV, info.getVPitch(projV.getMaxPos()));
378
379 if (m_errorFromDB) { // Overwrite cluster position error with value from DB (keep the above calculation untouched for now)
380 unsigned int uID = info.getUCellID(projU.getPos());
381 unsigned int vID = info.getVCellID(projV.getPos());
382 bool isUedge = PXD::isClusterAtUEdge(sensorID, projU.getMinCell(), projU.getMaxCell());
383 bool isVedge = (PXD::isClusterAtVEdge(sensorID, projV.getMinCell(), projV.getMaxCell())
384 || PXD::isClusterAtLadderJoint(sensorID, projV.getMinCell(), projV.getMaxCell()));
385 assignPositionErrorFromDB(projU, **m_clusterPositionErrorUPar, sensorID, uID, vID, pitchU, isUedge, isVedge);
386 assignPositionErrorFromDB(projV, **m_clusterPositionErrorVPar, sensorID, uID, vID, pitchV, isUedge, isVedge);
387 }
388
389 ROOT::Math::XYZVector lorentzShift = info.getLorentzShift(projU.getPos(), projV.getPos());
390 projU.setPos(projU.getPos() - lorentzShift.X());
391 projV.setPos(projV.getPos() - lorentzShift.Y());
392 B2DEBUG(20, "Lorentz shift: " << lorentzShift.X() << " " << lorentzShift.Y());
393
394 // Pre classification of cluster looking at pitch type of pixels and if they touch sensor edges
395 int clusterkind = PXDClusterPositionEstimator::getInstance().getClusterkind(cls.pixels(), sensorID);
396
397 // Compute sorted set of pixel
398 // FIXME: I am not 100% sure if cls.pixels() are sorted
399 set<Pixel> pixelSet(cls.pixels().begin(), cls.pixels().end());
400
401 // Compute classifier variables needed for later retrival of position correction in PXD CKF
402 vector<float> sectorEtaValues = {0, 0, 0, 0};
403 sectorEtaValues[0] = PXDClusterPositionEstimator::getInstance().computeEta(pixelSet, projV.getMinCell(), projV.getSize(), +1.0,
404 +1.0);
405 sectorEtaValues[1] = PXDClusterPositionEstimator::getInstance().computeEta(pixelSet, projV.getMinCell(), projV.getSize(), -1.0,
406 +1.0);
407 sectorEtaValues[2] = PXDClusterPositionEstimator::getInstance().computeEta(pixelSet, projV.getMinCell(), projV.getSize(), -1.0,
408 -1.0);
409 sectorEtaValues[3] = PXDClusterPositionEstimator::getInstance().computeEta(pixelSet, projV.getMinCell(), projV.getSize(), +1.0,
410 -1.0);
411
412 vector<int> sectorShapeIndices = { -1, -1, -1, -1};
413 sectorShapeIndices[0] = PXDClusterPositionEstimator::getInstance().computeShapeIndex(pixelSet, projU.getMinCell(),
414 projV.getMinCell(), projV.getSize(), +1.0, +1.0);
415 sectorShapeIndices[1] = PXDClusterPositionEstimator::getInstance().computeShapeIndex(pixelSet, projU.getMinCell(),
416 projV.getMinCell(), projV.getSize(), -1.0, +1.0);
417 sectorShapeIndices[2] = PXDClusterPositionEstimator::getInstance().computeShapeIndex(pixelSet, projU.getMinCell(),
418 projV.getMinCell(), projV.getSize(), -1.0, -1.0);
419 sectorShapeIndices[3] = PXDClusterPositionEstimator::getInstance().computeShapeIndex(pixelSet, projU.getMinCell(),
420 projV.getMinCell(), projV.getSize(), +1.0, -1.0);
421
422 //Store Cluster into Datastore ...
423 int clsIndex = storeClusters.getEntries();
424 storeClusters.appendNew(sensorID, projU.getPos(), projV.getPos(), projU.getError(), projV.getError(),
425 rho, cls.getCharge(), seed.getCharge(),
426 cls.size(), projU.getSize(), projV.getSize(), projU.getMinCell(), projV.getMinCell(), clusterkind,
427 sectorEtaValues, sectorShapeIndices
428 );
429
430 //Create Relations to this Digit
431 if (!mc_relations.empty()) relClusterMCParticle.add(clsIndex, mc_relations.begin(), mc_relations.end());
432 if (!truehit_relations.empty()) relClusterTrueHit.add(clsIndex, truehit_relations.begin(), truehit_relations.end());
433 relClusterDigit.add(clsIndex, digit_weights.begin(), digit_weights.end());
434 }
435
436 m_cache->clear();
437}
438
440 const ClusterProjection& secondary, double minPitch, double centerPitch, double maxPitch)
441{
442 if (primary.getSize() >= static_cast<unsigned int>(m_sizeHeadTail)) { //Analog head tail
443 //Average charge in the central area
444 const double centerCharge = primary.getCenterCharge() / (primary.getSize() - 2);
445 const double minCharge = (primary.getMinCharge() < centerCharge) ? primary.getMinCharge() : centerCharge;
446 const double maxCharge = (primary.getMaxCharge() < centerCharge) ? primary.getMaxCharge() : centerCharge;
447 primary.setPos(0.5 * (primary.getMinPos() + primary.getMaxPos()
448 + (maxCharge * maxPitch - minCharge * minPitch) / centerCharge));
449 const double snHead = centerCharge / m_cutAdjacentSignal / minPitch;
450 const double snTail = centerCharge / m_cutAdjacentSignal / maxPitch;
451 const double landauHead = minCharge / centerCharge * minPitch;
452 const double landauTail = maxCharge / centerCharge * maxPitch;
453 primary.setError(0.5 * sqrt(1.0 / snHead / snHead + 1.0 / snTail / snTail
454 + 0.5 * landauHead * landauHead + 0.5 * landauTail * landauTail));
455 } else if (primary.getSize() <= 2) { // Add a phantom charge to second strip
456 primary.setError(centerPitch * (secondary.getSize() + 2) * m_cutAdjacentSignal / (primary.getCharge() +
457 (secondary.getSize() + 3) * m_cutAdjacentSignal));
458 } else {
459 const double sn = cls.getSeedCharge() / m_elNoise;
460 primary.setError(2.0 * centerPitch / sn);
461 }
462}
463
465 VxdID sensorID, unsigned int uCell, unsigned int vCell, double centerPitch,
466 bool isAtUEdge, bool isAtVEdge, bool isAdjacentDead)
467{
468 // Get bins from cell ID
469 unsigned short uBin = PXD::getBinU(sensorID, uCell, vCell, errorPar.getBinsU());
470 unsigned short vBin = PXD::getBinV(sensorID, vCell, errorPar.getBinsV());
471 // Get error from DB [in units of pix]
472 double error = errorPar.getContent(sensorID.getID(), uBin, vBin, primary.getSize());
473 double sf = 1.;
474 // Apply additional factor if at sensor edges or adjacent to daed rows/colums
475 if (isAtUEdge) sf *= errorPar.getSensorUEdgeFactor(sensorID.getID(), uBin, vBin, primary.getSize());
476 if (isAtVEdge) sf *= errorPar.getSensorVEdgeFactor(sensorID.getID(), uBin, vBin, primary.getSize());
477 if (isAdjacentDead) sf *= errorPar.getDeadNeighbourFactor(sensorID.getID(), uBin, vBin, primary.getSize());
478 // Set error (convert to [um])
479 if (error) primary.setError(sf * error * centerPitch); // zero means default values to use analytic error
480
481}
@ c_ErrorIfAlreadyRegistered
If the object/array was already registered, produce an error (aborting initialisation).
Definition: DataStore.h:72
Base class for Modules.
Definition: Module.h:72
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
The payload class for PXD cluster position error.
float getContent(unsigned short sensorID, unsigned short globalID) const
Get content.
unsigned short getBinsV() const
Get number of bins along sensor v side.
float getDeadNeighbourFactor(unsigned short sensorID, unsigned short globalID) const
Get scaling factor when neighbouring dead rows/column.
unsigned short getBinsU() const
Get number of bins along sensor u side.
float getSensorVEdgeFactor(unsigned short sensorID, unsigned short globalID) const
Get scaling factor at sensor edge in V.
float getSensorUEdgeFactor(unsigned short sensorID, unsigned short globalID) const
Get scaling factor at sensor edge in U.
The PXD digit class.
Definition: PXDDigit.h:27
VxdID getSensorID() const
Get the sensor ID.
Definition: PXDDigit.h:64
Class to remember recently assigned clusters This class will remember the current and the last pixel ...
Definition: ClusterCache.h:77
Class representing a possible cluster during clustering of the PXD It supports merging of different c...
float getSeedCharge() const
get the seed charge of the cluster
Helper struct to collect information about the 1D projection of a Pixel cluster.
double getMinPos() const
Return the position of the minimum cell of the cluster.
void add(unsigned int cell, float position, float charge)
Add Pixel information to the projection.
double getMaxCharge() const
Return the charge in the maximum cell of the cluster.
void finalize()
Finish calculation of center of gravity and set correct cluster size.
unsigned int getMinCell() const
Return the minimum cell part of the cluster.
double getCharge() const
Return the total charge of the cluster.
unsigned int getMaxCell() const
Return the maximum cell part of the cluster.
void setPos(double pos)
Set the position of the cluster.
double getMinCharge() const
Return the charge in the minimum cell of the cluster.
void setError(double error)
Set the error of the cluster.
double getError() const
Return the error of the cluster.
double getMaxPos() const
Return the position of the maximum cell of the cluster.
double getPos() const
Return the projected position of the cluster.
double getCenterCharge() const
Return the center charge of the cluster, that is total charge minus minimum and maximum cell charge.
unsigned int getSize() const
Return the projected size of the cluster.
void setNoiseLevel(float noise)
Set the noise level.
Definition: NoiseMap.h:33
virtual void setSensorID(VxdID)
Set the sensorID currently used.
Definition: NoiseMap.h:38
static PXDClusterPositionEstimator & getInstance()
Main (and only) way to access the PXDClusterPositionEstimator.
float computeEta(const std::set< Pixel > &pixels, int vStart, int vSize, double thetaU, double thetaV) const
Return the normed charge ratio between head and tail pixels (size>=2) or the charge of the seed (size...
int computeShapeIndex(const std::set< Pixel > &pixels, int uStart, int vStart, int vSize, double thetaU, double thetaV) const
Return the shape index of the pixels.
int getClusterkind(const PXDCluster &cluster) const
Return kind of cluster needed to find cluster position correction.
NoiseMap m_noiseMap
Noise map for the currently active sensor.
int m_clusterCacheSize
Size of cluster Cache (0 = default)
std::unique_ptr< DBObjPtr< PXDClusterPositionErrorPar > > m_clusterPositionErrorVPar
DB object for cluster posotion errors in V.
double m_cutAdjacentSignal
Signal in ADU for Adjacent cut, basically m_elNoise*m_cutAdjacent.
virtual void initialize() override
Initialize the module.
std::string m_positionErrorVName
Name of the DB payload containing cluster posotion errors in V.
StoreObjPtr< EventLevelTrackingInfo > m_eventLevelTrackingInfo
StoreObject to access the event level tracking information.
std::string m_relDigitMCParticleName
Name of the relation between PXDDigits and MCParticles.
virtual void event() override
do the clustering
void calculatePositionError(const ClusterCandidate &cls, ClusterProjection &primary, const ClusterProjection &secondary, double minPitch, double centerPitch, double maxPitch)
Calculate position and error for a given cluster.
std::vector< const RelationElement * > RelationLookup
Container for a RelationArray Lookup table.
std::string m_relClusterDigitName
Name of the relation between PXDClusters and PXDDigits.
bool m_createPXDClustersForAbortedTrackingEvents
bool to override the EventLevelTrackingInfo abort flag decision
double m_cutCluster
Cluster cut in sigma.
std::string m_storeTrueHitsName
Name of the collection to use for the PXDTrueHits.
PXDClusterizerModule()
Constructor defining the parameters.
void fillRelationMap(const RelationLookup &lookup, std::map< unsigned int, float > &relation, unsigned int index)
Add the relation from a given PXDDigit index to a map.
std::string m_storeMCParticlesName
Name of the collection to use for the MCParticles.
RelationLookup m_trueRelation
Lookup table for PXDDigit->PXDTrueHit relation.
std::string m_positionErrorUName
Name of the DB payload containing cluster posotion errors in U.
virtual void beginRun() override
do at every run change
std::unique_ptr< DBObjPtr< PXDClusterPositionErrorPar > > m_clusterPositionErrorUPar
DB object for cluster posotion errors in U.
void assignPositionErrorFromDB(ClusterProjection &primary, PXDClusterPositionErrorPar errorPar, VxdID sensorID, unsigned int uCell, unsigned int vCell, double centerPitch, bool isAtUEdge=false, bool isAtVEdge=false, bool isAdjacentDead=false)
Assign position error for a given cluster from DB.
int m_sizeHeadTail
Size of the cluster at which we switch from Center of Gravity to Analog Head Tail.
std::string m_storeDigitsName
Name of the collection to use for the PXDDigits.
void createRelationLookup(const RelationArray &relation, RelationLookup &lookup, size_t digits)
Create lookup maps for Relations We do not use the RelationIndex as we know much more about the relat...
std::string m_storeClustersName
Name of the collection to use for the PXDClusters.
void writeClusters(VxdID sensorID)
Write clusters to collection.
std::string m_relDigitTrueHitName
Name of the relation between PXDDigits and PXDTrueHits.
std::string m_relClusterMCParticleName
Name of the relation between PXDClusters and MCParticles.
RelationLookup m_mcRelation
Lookup table for PXDDigit->MCParticle relation.
std::unique_ptr< ClusterCache > m_cache
cache of the last seen clusters to speed up clustering
bool m_errorFromDB
Flag to set cluster position error from DB (default = true)
double m_cutAdjacent
Noise cut in sigma.
std::string m_relClusterTrueHitName
Name of the relation between PXDClusters and PXDTrueHits.
Class to represent one pixel, used in clustering for fast access.
Definition: Pixel.h:36
unsigned short getU() const
Return the CellID in u.
Definition: Pixel.h:66
unsigned short getV() const
Return the CellID in v.
Definition: Pixel.h:68
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
void clear() override
Clear all elements from the relation.
Class to store a single element of a relation.
const std::string & getName() const
Return name under which the object is saved in the DataStore.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
void clear() override
Delete all entries in this array.
Definition: StoreArray.h:207
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne 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
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getID() const
Get the unique id.
Definition: VxdID.h:94
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
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
unsigned short getBinU(VxdID id, unsigned int uid, unsigned int vid, unsigned short nBinsU)
Function to return a bin number for equal sized binning in U.
Definition: PXDUtilities.h:161
bool isClusterAtUEdge(VxdID id, unsigned int umin, unsigned int umax)
Helper function to check if one of the end pixels are at the edge of the sensor.
Definition: PXDUtilities.h:125
bool isClusterAtVEdge(VxdID id, unsigned int vmin, unsigned int vmax)
Helper function to check if one of the end pixels are at the edge of the sensor.
Definition: PXDUtilities.h:136
bool isClusterAtLadderJoint(VxdID id, unsigned int vmin, unsigned int vmax)
Helper function to check if one of the end pixels are at the ladder joint.
Definition: PXDUtilities.h:148
unsigned short getBinV(VxdID id, unsigned int vid, unsigned short nBinsV)
Function to return a bin number for equal sized binning in V.
Definition: PXDUtilities.h:171
Abstract base class for different kinds of events.
STL namespace.