9 #include <svd/modules/svdReconstruction/SVDClusterizerDirectModule.h>
11 #include <framework/datastore/StoreArray.h>
12 #include <framework/logging/Logger.h>
14 #include <vxd/geometry/GeoCache.h>
15 #include <svd/geometry/SensorInfo.h>
17 #include <mdst/dataobjects/MCParticle.h>
18 #include <svd/dataobjects/SVDTrueHit.h>
19 #include <svd/dataobjects/SVDShaperDigit.h>
20 #include <svd/dataobjects/SVDCluster.h>
21 #include <mva/dataobjects/DatabaseRepresentationOfWeightfile.h>
23 #include <svd/reconstruction/NNWaveFitTool.h>
25 #include <unordered_map>
46 B2DEBUG(200,
"SVDClusterizerDirectModule ctor");
48 setDescription(
"Clusterize SVDShaperDigits and reconstruct hits");
49 setPropertyFlags(c_ParallelProcessingCertified);
52 addParam(
"Digits", m_storeShaperDigitsName,
53 "6-digits collection name",
string(
""));
54 addParam(
"Clusters", m_storeClustersName,
55 "Cluster collection name",
string(
""));
56 addParam(
"TrueHits", m_storeTrueHitsName,
57 "TrueHit collection name",
string(
""));
58 addParam(
"MCParticles", m_storeMCParticlesName,
59 "MCParticles collection name",
string(
""));
60 addParam(
"SVDEventInfo", m_svdEventInfoName,
61 "SVDEventInfo name",
string(
""));
64 addParam(
"TimeFitterName", m_timeFitterName,
65 "Name of time fitter data file",
string(
"SVDTimeNet_6samples"));
66 addParam(
"CalibratePeak", m_calibratePeak,
"Use calibrattion (vs. default) for peak widths and positions",
bool(
false));
70 addParam(
"NoiseSN", m_cutAdjacent,
71 "SN for digits to be considered for clustering", m_cutAdjacent);
72 addParam(
"SeedSN", m_cutSeed,
73 "SN for digits to be considered as seed", m_cutSeed);
74 addParam(
"ClusterSN", m_cutCluster,
75 "Minimum SN for clusters", m_cutCluster);
76 addParam(
"HeadTailSize", m_sizeHeadTail,
77 "Cluster size at which to switch to Analog head tail algorithm", m_sizeHeadTail);
81 void SVDClusterizerDirectModule::initialize()
93 m_storeSVDEvtInfo.isRequired();
95 if (!m_storeSVDEvtInfo.isOptional(m_svdEventInfoName)) m_svdEventInfoName =
"SVDEventInfoSim";
96 m_storeSVDEvtInfo.isRequired(m_svdEventInfoName);
98 RelationArray relClusterShaperDigits(storeClusters, storeShaperDigits);
99 RelationArray relClusterTrueHits(storeClusters, storeTrueHits);
100 RelationArray relClusterMCParticles(storeClusters, storeMCParticles);
101 RelationArray relShaperDigitTrueHits(storeShaperDigits, storeTrueHits);
102 RelationArray relShaperDigitMCParticles(storeShaperDigits, storeMCParticles);
112 m_storeClustersName = storeClusters.
getName();
113 m_storeShaperDigitsName = storeShaperDigits.
getName();
114 m_storeTrueHitsName = storeTrueHits.
getName();
115 m_storeMCParticlesName = storeMCParticles.
getName();
117 m_relClusterShaperDigitName = relClusterShaperDigits.
getName();
118 m_relClusterTrueHitName = relClusterTrueHits.
getName();
119 m_relClusterMCParticleName = relClusterMCParticles.
getName();
120 m_relShaperDigitTrueHitName = relShaperDigitTrueHits.
getName();
121 m_relShaperDigitMCParticleName = relShaperDigitMCParticles.
getName();
123 B2INFO(
" 1. COLLECTIONS:");
124 B2INFO(
" --> MCParticles: " << m_storeMCParticlesName);
125 B2INFO(
" --> Digits: " << m_storeShaperDigitsName);
126 B2INFO(
" --> Clusters: " << m_storeClustersName);
127 B2INFO(
" --> TrueHits: " << m_storeTrueHitsName);
128 B2INFO(
" --> DigitMCRel: " << m_relShaperDigitMCParticleName);
129 B2INFO(
" --> ClusterMCRel: " << m_relClusterMCParticleName);
130 B2INFO(
" --> ClusterDigitRel: " << m_relClusterShaperDigitName);
131 B2INFO(
" --> DigitTrueRel: " << m_relShaperDigitTrueHitName);
132 B2INFO(
" --> ClusterTrueRel: " << m_relClusterTrueHitName);
133 B2INFO(
" 2. CALIBRATION DATA:");
134 B2INFO(
" --> Time NN: " << m_timeFitterName);
135 B2INFO(
" 4. CLUSTERING:");
136 B2INFO(
" --> Neighbour cut: " << m_cutAdjacent);
137 B2INFO(
" --> Seed cut: " << m_cutSeed);
138 B2INFO(
" --> Cluster charge cut: " << m_cutCluster);
139 B2INFO(
" --> HT for clusters >: " << m_sizeHeadTail);
145 m_fitter.setNetwrok(dbXml->m_data);
148 void SVDClusterizerDirectModule::createRelationLookup(
const RelationArray& relation,
153 if (!relation)
return;
155 lookup.resize(digits);
156 for (
const auto& element : relation) {
157 lookup[element.getFromIndex()] = &element;
162 std::map<unsigned int, float>& relation,
unsigned int index)
165 if (!lookup.empty() && lookup[index]) {
167 const unsigned int size = element.getSize();
169 for (
unsigned int i = 0; i < size; ++i) {
172 if (element.getWeight(i) < 0)
continue;
173 relation[element.getToIndex(i)] += element.getWeight(i);
178 void SVDClusterizerDirectModule::event()
183 if (!storeShaperDigits || !storeShaperDigits.
getEntries() || !m_storeSVDEvtInfo.isValid())
return;
185 SVDModeByte modeByte = m_storeSVDEvtInfo->getModeByte();
187 size_t nDigits = storeShaperDigits.
getEntries();
188 B2DEBUG(90,
"Initial size of StoreDigits array: " << nDigits);
193 RelationArray relShaperDigitMCParticle(storeShaperDigits, storeMCParticles, m_relShaperDigitMCParticleName);
194 RelationArray relShaperDigitTrueHit(storeShaperDigits, storeTrueHits, m_relShaperDigitTrueHitName);
197 storeClusters.
clear();
199 RelationArray relClusterMCParticle(storeClusters, storeMCParticles,
200 m_relClusterMCParticleName);
201 if (relClusterMCParticle) relClusterMCParticle.
clear();
203 RelationArray relClusterShaperDigit(storeClusters, storeShaperDigits,
204 m_relClusterShaperDigitName);
205 if (relClusterShaperDigit) relClusterShaperDigit.
clear();
207 RelationArray relClusterTrueHit(storeClusters, storeTrueHits,
208 m_relClusterTrueHitName);
209 if (relClusterTrueHit) relClusterTrueHit.
clear();
212 createRelationLookup(relShaperDigitMCParticle, m_mcRelation, nDigits);
213 createRelationLookup(relShaperDigitTrueHit, m_trueRelation, nDigits);
219 vector<pair<unsigned short, unsigned short> > sensorDigits;
220 VxdID lastSensorID(0);
221 size_t firstSensorDigit = 0;
222 for (
size_t iDigit = 0; iDigit < nDigits; ++iDigit) {
226 if (sensorID != lastSensorID) {
227 sensorDigits.push_back(make_pair(firstSensorDigit, iDigit));
228 firstSensorDigit = iDigit;
229 lastSensorID = sensorID;
233 sensorDigits.push_back(make_pair(firstSensorDigit, nDigits));
236 for (
auto id_indices : sensorDigits) {
238 unsigned int firstDigit = id_indices.first;
239 unsigned int lastDigit = id_indices.second;
241 const SVDShaperDigit& sampleDigit = *storeShaperDigits[firstDigit];
247 double pitch = isU ? info.getUPitch() : info.getVPitch();
259 vector<pair<size_t, size_t> > stripGroups;
260 unordered_map<size_t, apvSamples> storedNormedSamples;
261 size_t firstClusterDigit = firstDigit;
262 size_t lastClusterDigit = firstDigit;
263 short lastStrip = -2;
265 B2DEBUG(300,
"Clustering digits " << firstDigit <<
" to " << lastDigit);
266 for (
size_t iDigit = firstDigit; iDigit < lastDigit; ++iDigit) {
269 unsigned short currentStrip = digit.
getCellID();
270 B2DEBUG(300,
"Digit " << iDigit <<
", strip: " << currentStrip <<
", lastStrip: " << lastStrip);
271 B2DEBUG(300,
"First CD: " << firstClusterDigit <<
" Last CD: " << lastClusterDigit);
273 float stripNoiseADU = m_noiseCal.getNoise(sensorID, isU, currentStrip);
279 transform(samples.begin(), samples.end(), normedSamples.begin(),
280 bind2nd(divides<float>(), stripNoiseADU));
281 bool validDigit =
pass3Samples(normedSamples, m_cutAdjacent);
286 zeroSuppress(normedSamples, m_cutAdjacent);
287 storedNormedSamples.insert(make_pair(iDigit, normedSamples));
291 bool consecutive = ((currentStrip - lastStrip) == 1);
292 lastStrip = currentStrip;
294 B2DEBUG(300, (validDigit ?
"Valid " :
"Invalid ") << (consecutive ?
"consecutive" :
"gap"));
297 if ((!validDigit || !consecutive) && (firstClusterDigit < lastClusterDigit)) {
298 B2DEBUG(300,
"Saving (" << firstClusterDigit <<
", " << lastClusterDigit <<
")");
299 stripGroups.emplace_back(firstClusterDigit, lastClusterDigit);
305 lastClusterDigit = iDigit + 1;
307 firstClusterDigit = iDigit;
308 lastClusterDigit = iDigit + 1;
311 firstClusterDigit = iDigit + 1;
312 lastClusterDigit = iDigit + 1;
316 if (firstClusterDigit < lastClusterDigit) {
317 B2DEBUG(300,
"Saving (" << firstClusterDigit <<
", " << lastDigit <<
")");
318 stripGroups.emplace_back(firstClusterDigit, lastDigit);
322 os << sensorID <<
"NormedSamples: " << endl;
323 for (
auto item : storedNormedSamples) {
324 os << item.first <<
": ";
325 copy(item.second.begin(), item.second.end(), ostream_iterator<double>(os,
" "));
328 os <<
"StripGroups: " << endl;
329 for (
auto item : stripGroups) {
330 os <<
"(" << item.first <<
", " << item.second <<
"), ";
332 B2DEBUG(300, os.str());
340 vector<unsigned short> stripNumbers;
341 vector<float> stripPositions;
342 vector<float> stripNoises;
343 vector<float> timeShifts;
344 vector<float> waveWidths;
346 for (
auto clusterBounds : stripGroups) {
348 unsigned short clusterSize = clusterBounds.second - clusterBounds.first;
349 assert(clusterSize > 0);
351 stripNumbers.clear();
352 stripPositions.clear();
357 for (
size_t iDigit = clusterBounds.first; iDigit < clusterBounds.second; ++iDigit) {
363 unsigned short stripNo = digit.
getCellID();
364 stripNumbers.push_back(stripNo);
366 double stripNoiseADU = m_noiseCal.getNoise(sensorID, isU, stripNo);
367 stripNoises.push_back(
368 m_pulseShapeCal.getChargeFromADC(sensorID, isU, stripNo, stripNoiseADU)
373 double peakWidth = 270;
374 double timeShift = isU ? 2.5 : -2.2;
375 if (m_calibratePeak) {
376 peakWidth = 1.988 * m_pulseShapeCal.getWidth(sensorID, isU, stripNo);
377 timeShift = m_pulseShapeCal.getPeakTime(sensorID, isU, stripNo)
381 const double triggerBinSep = 4 * 1.96516;
382 double apvPhase = triggerBinSep * (0.5 +
static_cast<int>(modeByte.
getTriggerBin()));
383 timeShift = timeShift + apvPhase;
384 waveWidths.push_back(peakWidth);
385 timeShifts.push_back(timeShift);
386 stripPositions.push_back(
387 isU ? info.getUCellPosition(stripNo) : info.getVCellPosition(stripNo)
389 stripPositions.push_back(
390 isU ? info.getUCellPosition(stripNo) : info.getVCellPosition(stripNo)
396 float clusterNoise = sqrt(
398 * inner_product(stripNoises.begin(), stripNoises.end(), stripNoises.begin(), 0.0)
400 B2DEBUG(200,
"RMS cluster noise: " << clusterNoise);
404 shared_ptr<nnFitterBinData> pStrip;
407 fill(pCluster.begin(), pCluster.end(),
double(1.0));
408 for (
size_t iClusterStrip = 0; iClusterStrip < clusterSize; ++iClusterStrip) {
409 size_t iDigit = clusterBounds.first + iClusterStrip;
411 os1 <<
"Input to NNFitter: iDigit = " << iDigit << endl <<
"Samples: ";
412 copy(storedNormedSamples[iDigit].begin(), storedNormedSamples[iDigit].end(),
413 ostream_iterator<double>(os1,
" "));
415 pStrip = m_fitter.getFit(storedNormedSamples[iDigit], waveWidths[iClusterStrip]);
416 os1 <<
"Output from NNWaveFitter: " << endl;
417 copy(pStrip->begin(), pStrip->end(), ostream_iterator<double>(os1,
" "));
420 fitTool.
shiftInTime(*pStrip, -timeShifts[iClusterStrip]);
421 fitTool.
multiply(pCluster, *pStrip);
422 os1 <<
"Accummulated: " << endl;
423 copy(pCluster.begin(), pCluster.end(), ostream_iterator<double>(os1,
" "));
424 B2DEBUG(200, os1.str());
427 double clusterTime, clusterTimeErr;
428 tie(clusterTime, clusterTimeErr) = fitTool.
getTimeShift(pCluster);
429 B2DEBUG(200,
"Time: " << clusterTime <<
" +/- " << clusterTimeErr);
433 vector<double> stripAmplitudes(stripNoises.size());
434 vector<double> stripAmplitudeErrors(stripNoises.size());
435 double clusterChi2 = 0.0;
436 for (
size_t iClusterStrip = 0; iClusterStrip < clusterSize; ++iClusterStrip) {
437 size_t iDigit = clusterBounds.first + iClusterStrip;
438 double snAmp, snAmpError, chi2;
439 tie(snAmp, snAmpError, chi2) =
440 fitTool.
getAmplitudeChi2(storedNormedSamples[iDigit], clusterTime, waveWidths[iClusterStrip]);
442 stripAmplitudes[iClusterStrip] =
443 stripNoises[iClusterStrip] * snAmp;
444 stripAmplitudeErrors[iClusterStrip] =
445 stripNoises[iClusterStrip] * snAmpError;
447 B2DEBUG(200,
"Digit " << iDigit <<
" Noise: " << stripNoises[iClusterStrip]
448 <<
" Amplitude: " << stripAmplitudes[iClusterStrip]
449 <<
" +/- " << stripAmplitudeErrors[iClusterStrip]
454 float clusterCharge = accumulate(stripAmplitudes.begin(), stripAmplitudes.end(), 0.0);
455 float clusterChargeError = sqrt(
456 inner_product(stripAmplitudeErrors.begin(), stripAmplitudeErrors.end(),
457 stripAmplitudeErrors.begin(), 0.0)
459 float clusterSN = (clusterChargeError > 0) ? clusterCharge / clusterChargeError : clusterCharge;
460 clusterChi2 /= clusterSize;
462 size_t seedIndex = distance(stripAmplitudes.begin(), max_element(
463 stripAmplitudes.begin(), stripAmplitudes.end()));
465 float clusterSeedCharge = stripAmplitudes[seedIndex];
466 B2DEBUG(200,
"Cluster parameters:");
467 B2DEBUG(200,
"Charge: " << clusterCharge <<
" +/- " << clusterChargeError);
468 B2DEBUG(200,
"Seed: " << clusterSeedCharge <<
" +/- " << stripAmplitudeErrors[seedIndex]);
469 B2DEBUG(200,
"S/N: " << clusterSN);
473 float clusterPosition, clusterPositionError;
478 if ((clusterCharge < clusterNoise * m_cutCluster) || (clusterSeedCharge < clusterNoise * m_cutSeed))
481 if (clusterSize < m_sizeHeadTail) {
482 clusterPosition = 1.0 / clusterCharge * inner_product(
483 stripAmplitudes.begin(), stripAmplitudes.end(), stripPositions.begin(), 0.0
486 float phantomCharge = m_cutAdjacent * clusterNoise;
487 if (clusterSize == 1) {
488 clusterPositionError = pitch * phantomCharge / (clusterCharge + phantomCharge);
490 clusterPositionError = pitch * phantomCharge / clusterCharge;
493 float leftStripCharge = stripAmplitudes.front();
494 float leftPos = stripPositions.front();
495 float rightStripCharge = stripAmplitudes.back();
496 float rightPos = stripPositions.back();
497 float centreCharge = (clusterCharge - leftStripCharge - rightStripCharge) / (clusterSize - 2);
498 leftStripCharge = (leftStripCharge < centreCharge) ? leftStripCharge : centreCharge;
499 rightStripCharge = (rightStripCharge < centreCharge) ? rightStripCharge : centreCharge;
500 clusterPosition = 0.5 * (leftPos + rightPos)
501 + 0.5 * (rightStripCharge - leftStripCharge) / centreCharge * pitch;
502 float sn = centreCharge / m_cutAdjacent / clusterNoise;
504 float landauHead = leftStripCharge / centreCharge;
505 double landauTail = rightStripCharge / centreCharge;
506 clusterPositionError = 0.5 * pitch * sqrt(1.0 / sn / sn
507 + 0.5 * landauHead * landauHead + 0.5 * landauTail * landauTail);
510 clusterPosition -= info.getLorentzShift(isU, clusterPosition);
513 map<unsigned int, float> mc_relations;
514 map<unsigned int, float> truehit_relations;
515 vector<pair<unsigned int, float> > digit_weights;
516 digit_weights.reserve(clusterSize);
518 for (
size_t iDigit = clusterBounds.first; iDigit < clusterBounds.second; ++iDigit) {
520 fillRelationMap(m_mcRelation, mc_relations, iDigit);
521 fillRelationMap(m_trueRelation, truehit_relations, iDigit);
523 digit_weights.emplace_back(iDigit, stripAmplitudes[iDigit - clusterBounds.first]);
528 VxdID clusterSensorID = sensorID;
531 SVDCluster(sensorID, isU, clusterPosition, clusterPositionError, clusterTime,
532 clusterTimeErr, clusterCharge, clusterSeedCharge, clusterSize, clusterSN, clusterChi2)
536 if (!mc_relations.empty()) {
537 relClusterMCParticle.
add(clsIndex, mc_relations.begin(), mc_relations.end());
539 if (!truehit_relations.empty()) {
540 relClusterTrueHit.
add(clsIndex, truehit_relations.begin(), truehit_relations.end());
542 relClusterShaperDigit.
add(clsIndex, digit_weights.begin(), digit_weights.end());
547 B2DEBUG(100,
"Number of clusters: " << storeClusters.
getEntries());
Class for accessing objects in the database.
Low-level class to create/modify relations between StoreArrays.
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.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Class to store SVD mode information.
baseType getTriggerBin() const
Get the triggerBin id.
The SVD ShaperDigit class.
VxdID getSensorID() const
Get the sensor ID.
APVFloatSamples getSamples() const
Get array of samples.
short int getCellID() const
Get strip ID.
bool isUStrip() const
Get strip direction.
std::vector< const RelationElement * > RelationLookup
Container for a RelationArray Lookup table.
Specific implementation of SensorInfo for SVD Sensors which provides additional sensor specific infor...
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.
const std::string & getName() const
Return name under which the object is saved in the DataStore.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
T * appendNew()
Construct a new T object at the end of the array.
int getEntries() const
Get the number of objects in the array.
void clear() override
Delete all entries in this array.
Class to uniquely identify a any structure of the PXD and SVD.
void setSegmentNumber(baseType segment)
Set the sensor segment.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Namespace to encapsulate code needed for simulation and reconstrucion of the SVD.
std::array< apvSampleBaseType, nAPVSamples > apvSamples
vector od apvSample BaseType objects
std::vector< double > nnFitterBinData
Vector of values defined for bins, such as bin times or bin probabilities.
bool pass3Samples(const T &a, double thr)
pass 3-samples
Abstract base class for different kinds of events.