9 #include <svd/modules/svdReconstruction/SVDNNClusterizerModule.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/SVDRecoDigit.h>
21 #include <svd/dataobjects/SVDCluster.h>
22 #include <mva/dataobjects/DatabaseRepresentationOfWeightfile.h>
24 #include <svd/reconstruction/NNWaveFitTool.h>
46 B2DEBUG(200,
"SVDNNClusterizerModule ctor");
48 setDescription(
"Clusterize SVDRecoDigits and reconstruct hits");
49 setPropertyFlags(c_ParallelProcessingCertified);
52 addParam(
"Digits", m_storeRecoDigitsName,
53 "RecoDigits 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(
""));
62 addParam(
"TimeFitterName", m_timeFitterName,
63 "Name of time fitter data file",
string(
"SVDTimeNet_6samples"));
64 addParam(
"CalibratePeak", m_calibratePeak,
"Use calibrattion (vs. default) for peak widths and positions",
bool(
false));
68 addParam(
"NoiseSN", m_cutAdjacent,
69 "SN for digits to be considered for clustering", m_cutAdjacent);
70 addParam(
"SeedSN", m_cutSeed,
71 "SN for digits to be considered as seed", m_cutSeed);
72 addParam(
"ClusterSN", m_cutCluster,
73 "Minimum SN for clusters", m_cutCluster);
74 addParam(
"HeadTailSize", m_sizeHeadTail,
75 "Cluster size at which to switch to Analog head tail algorithm", m_sizeHeadTail);
79 void SVDNNClusterizerModule::initialize()
92 RelationArray relClusterRecoDigits(storeClusters, storeRecoDigits);
93 RelationArray relClusterTrueHits(storeClusters, storeTrueHits);
94 RelationArray relClusterMCParticles(storeClusters, storeMCParticles);
95 RelationArray relRecoDigitTrueHits(storeRecoDigits, storeTrueHits);
96 RelationArray relRecoDigitMCParticles(storeRecoDigits, storeMCParticles);
106 m_storeClustersName = storeClusters.
getName();
107 m_storeRecoDigitsName = storeRecoDigits.
getName();
108 m_storeTrueHitsName = storeTrueHits.
getName();
109 m_storeMCParticlesName = storeMCParticles.
getName();
111 m_relClusterRecoDigitName = relClusterRecoDigits.
getName();
112 m_relClusterTrueHitName = relClusterTrueHits.
getName();
113 m_relClusterMCParticleName = relClusterMCParticles.
getName();
114 m_relRecoDigitTrueHitName = relRecoDigitTrueHits.
getName();
115 m_relRecoDigitMCParticleName = relRecoDigitMCParticles.
getName();
117 B2INFO(
" 1. COLLECTIONS:");
118 B2INFO(
" --> MCParticles: " << m_storeMCParticlesName);
119 B2INFO(
" --> Digits: " << m_storeRecoDigitsName);
120 B2INFO(
" --> Clusters: " << m_storeClustersName);
121 B2INFO(
" --> TrueHits: " << m_storeTrueHitsName);
122 B2INFO(
" --> DigitMCRel: " << m_relRecoDigitMCParticleName);
123 B2INFO(
" --> ClusterMCRel: " << m_relClusterMCParticleName);
124 B2INFO(
" --> ClusterDigitRel: " << m_relClusterRecoDigitName);
125 B2INFO(
" --> DigitTrueRel: " << m_relRecoDigitTrueHitName);
126 B2INFO(
" --> ClusterTrueRel: " << m_relClusterTrueHitName);
127 B2INFO(
" 2. CALIBRATION DATA:");
128 B2INFO(
" --> Time NN: " << m_timeFitterName);
129 B2INFO(
" 4. CLUSTERING:");
130 B2INFO(
" --> Neighbour cut: " << m_cutAdjacent);
131 B2INFO(
" --> Seed cut: " << m_cutSeed);
132 B2INFO(
" --> Cluster charge cut: " << m_cutCluster);
133 B2INFO(
" --> HT for clusters >: " << m_sizeHeadTail);
139 m_fitter.setNetwrok(dbXml->m_data);
142 void SVDNNClusterizerModule::createRelationLookup(
const RelationArray& relation,
147 if (!relation)
return;
149 lookup.resize(digits);
150 for (
const auto& element : relation) {
151 lookup[element.getFromIndex()] = &element;
156 std::map<unsigned int, float>& relation,
unsigned int index)
159 if (!lookup.empty() && lookup[index]) {
161 const unsigned int size = element.getSize();
163 for (
unsigned int i = 0; i < size; ++i) {
166 if (element.getWeight(i) < 0)
continue;
167 relation[element.getToIndex(i)] += element.getWeight(i);
172 void SVDNNClusterizerModule::event()
176 if (!storeRecoDigits || !storeRecoDigits.
getEntries())
return;
178 size_t nDigits = storeRecoDigits.
getEntries();
179 B2DEBUG(90,
"Initial size of RecoDigits array: " << nDigits);
184 RelationArray relRecoDigitMCParticle(storeRecoDigits, storeMCParticles, m_relRecoDigitMCParticleName);
185 RelationArray relRecoDigitTrueHit(storeRecoDigits, storeTrueHits, m_relRecoDigitTrueHitName);
188 storeClusters.
clear();
190 RelationArray relClusterMCParticle(storeClusters, storeMCParticles,
191 m_relClusterMCParticleName);
192 if (relClusterMCParticle) relClusterMCParticle.
clear();
194 RelationArray relClusterRecoDigit(storeClusters, storeRecoDigits,
195 m_relClusterRecoDigitName);
196 if (relClusterRecoDigit) relClusterRecoDigit.
clear();
198 RelationArray relClusterTrueHit(storeClusters, storeTrueHits,
199 m_relClusterTrueHitName);
200 if (relClusterTrueHit) relClusterTrueHit.
clear();
203 createRelationLookup(relRecoDigitMCParticle, m_mcRelation, nDigits);
204 createRelationLookup(relRecoDigitTrueHit, m_trueRelation, nDigits);
210 vector<pair<unsigned short, unsigned short> > sensorDigits;
211 VxdID lastSensorID(0);
212 size_t firstSensorDigit = 0;
213 for (
size_t iDigit = 0; iDigit < nDigits; ++iDigit) {
217 if (sensorID != lastSensorID) {
218 sensorDigits.push_back(make_pair(firstSensorDigit, iDigit));
219 firstSensorDigit = iDigit;
220 lastSensorID = sensorID;
224 sensorDigits.push_back(make_pair(firstSensorDigit, nDigits));
227 for (
auto id_indices : sensorDigits) {
229 unsigned int firstDigit = id_indices.first;
230 unsigned int lastDigit = id_indices.second;
232 const SVDRecoDigit& sampleRecoDigit = *storeRecoDigits[firstDigit];
234 bool isU = sampleRecoDigit.
isUStrip();
249 vector<pair<size_t, size_t> > stripGroups;
250 size_t firstClusterDigit = firstDigit;
251 size_t lastClusterDigit = firstDigit;
252 short lastStrip = -2;
254 B2DEBUG(300,
"Clustering digits " << firstDigit <<
" to " << lastDigit);
255 for (
size_t iDigit = firstDigit; iDigit < lastDigit; ++iDigit) {
257 const SVDRecoDigit& recoDigit = *storeRecoDigits[iDigit];
258 unsigned short currentStrip = recoDigit.
getCellID();
259 B2DEBUG(300,
"Digit " << iDigit <<
", strip: " << currentStrip <<
", lastStrip: " << lastStrip);
260 B2DEBUG(300,
"First CD: " << firstClusterDigit <<
" Last CD: " << lastClusterDigit);
263 bool consecutive = ((currentStrip - lastStrip) == 1);
264 lastStrip = currentStrip;
266 B2DEBUG(300, (consecutive ?
"consecutive" :
"gap"));
269 if (!consecutive && (firstClusterDigit < lastClusterDigit)) {
270 B2DEBUG(300,
"Saving (" << firstClusterDigit <<
", " << lastClusterDigit <<
")");
271 stripGroups.emplace_back(firstClusterDigit, lastClusterDigit);
276 lastClusterDigit = iDigit + 1;
278 firstClusterDigit = iDigit;
279 lastClusterDigit = iDigit + 1;
283 if (firstClusterDigit < lastClusterDigit) {
284 B2DEBUG(300,
"Saving (" << firstClusterDigit <<
", " << lastDigit <<
")");
285 stripGroups.emplace_back(firstClusterDigit, lastDigit);
289 os <<
"StripGroups: " << endl;
290 for (
auto item : stripGroups) {
291 os <<
"(" << item.first <<
", " << item.second <<
"), ";
293 B2DEBUG(300, os.str());
302 vector<unsigned short> stripNumbers;
303 vector<float> stripPositions;
304 vector<float> stripNoises;
305 vector<float> stripGains;
306 vector<float> timeShifts;
307 vector<float> waveWidths;
308 vector<apvSamples> storedNormedSamples;
309 vector<SVDRecoDigit::OutputProbArray> storedPDFs;
312 for (
auto clusterBounds : stripGroups) {
314 unsigned short clusterSize = clusterBounds.second - clusterBounds.first;
315 assert(clusterSize > 0);
317 stripNumbers.clear();
318 stripPositions.clear();
324 for (
size_t iDigit = clusterBounds.first; iDigit < clusterBounds.second; ++iDigit) {
328 const SVDRecoDigit& recoDigit = *storeRecoDigits[iDigit];
330 unsigned short stripNo = recoDigit.
getCellID();
331 stripNumbers.push_back(stripNo);
333 double stripNoiseADU = m_noiseCal.getNoise(sensorID, isU, stripNo);
334 stripNoises.push_back(
335 m_pulseShapeCal.getChargeFromADC(sensorID, isU, stripNo, stripNoiseADU)
340 double peakWidth = 270;
341 double timeShift = isU ? 4.0 : 0.0;
342 if (m_calibratePeak) {
343 peakWidth = 1.988 * m_pulseShapeCal.getWidth(sensorID, isU, stripNo);
344 timeShift = m_pulseShapeCal.getPeakTime(sensorID, isU, stripNo)
347 waveWidths.push_back(peakWidth);
348 timeShifts.push_back(timeShift);
349 stripPositions.push_back(
356 B2FATAL(
"Missing SVDRecoDigits->SVDShaperDigits relation. This should not happen.");
358 transform(samples.begin(), samples.end(), normedSamples.begin(),
359 bind2nd(divides<float>(), stripNoiseADU));
363 zeroSuppress(normedSamples, m_cutAdjacent);
364 storedNormedSamples.emplace_back(normedSamples);
370 float clusterNoise = sqrt(
372 * inner_product(stripNoises.begin(), stripNoises.end(), stripNoises.begin(), 0.0)
374 B2DEBUG(200,
"RMS cluster noise: " << clusterNoise);
378 shared_ptr<nnFitterBinData> pStrip;
381 fill(pCluster.begin(), pCluster.end(),
double(1.0));
383 for (
size_t iClusterStrip = 0; iClusterStrip < clusterSize; ++iClusterStrip) {
384 size_t iDigit = clusterBounds.first + iClusterStrip;
386 os1 <<
"Input to NNFitter: iDigit = " << iDigit << endl <<
"Samples: ";
387 copy(storedNormedSamples[iClusterStrip].begin(), storedNormedSamples[iClusterStrip].end(),
388 ostream_iterator<double>(os1,
" "));
390 os1 <<
"PDF from RecoDigit: " << endl;
391 copy(storedPDFs[iClusterStrip].begin(), storedPDFs[iClusterStrip].end(), ostream_iterator<double>(os1,
" "));
393 fitTool.
multiply(pCluster, storedPDFs[iClusterStrip]);
394 os1 <<
"Accummulated: " << endl;
395 copy(pCluster.begin(), pCluster.end(), ostream_iterator<double>(os1,
" "));
396 B2DEBUG(200, os1.str());
399 double clusterTime, clusterTimeErr;
400 tie(clusterTime, clusterTimeErr) = fitTool.
getTimeShift(pCluster);
401 B2DEBUG(200,
"Time: " << clusterTime <<
" +/- " << clusterTimeErr);
405 vector<double> stripAmplitudes(stripNoises.size());
406 vector<double> stripAmplitudeErrors(stripNoises.size());
407 double clusterChi2 = 0.0;
408 for (
size_t iClusterStrip = 0; iClusterStrip < clusterSize; ++iClusterStrip) {
409 size_t iDigit = clusterBounds.first + iClusterStrip;
410 double snAmp, snAmpError, chi2;
411 tie(snAmp, snAmpError, chi2) =
412 fitTool.
getAmplitudeChi2(storedNormedSamples[iClusterStrip], clusterTime, waveWidths[iClusterStrip]);
414 stripAmplitudes[iClusterStrip] = stripNoises[iClusterStrip] * snAmp;
415 stripAmplitudeErrors[iClusterStrip] = stripNoises[iClusterStrip] * snAmpError;
417 B2DEBUG(200,
"Digit " << iDigit <<
" Noise: " << stripNoises[iClusterStrip]
418 <<
" Amplitude: " << stripAmplitudes[iClusterStrip]
419 <<
" +/- " << stripAmplitudeErrors[iClusterStrip]
424 float clusterCharge = accumulate(stripAmplitudes.begin(), stripAmplitudes.end(), 0.0);
425 float clusterChargeError = sqrt(
426 inner_product(stripAmplitudeErrors.begin(), stripAmplitudeErrors.end(),
427 stripAmplitudeErrors.begin(), 0.0)
429 float clusterSN = (clusterChargeError > 0) ? clusterCharge / clusterChargeError : clusterCharge;
431 clusterChi2 /= clusterSize;
433 size_t seedIndex = distance(stripAmplitudes.begin(), max_element(
434 stripAmplitudes.begin(), stripAmplitudes.end()));
435 float clusterSeedCharge = stripAmplitudes[seedIndex];
436 B2DEBUG(200,
"Cluster parameters:");
437 B2DEBUG(200,
"Charge: " << clusterCharge <<
" +/- " << clusterChargeError);
438 B2DEBUG(200,
"Seed: " << clusterSeedCharge <<
" +/- " << stripAmplitudeErrors[seedIndex]);
439 B2DEBUG(200,
"S/N: " << clusterSN);
440 B2DEBUG(200,
"chi2: " << clusterChi2);
443 float clusterPosition, clusterPositionError;
448 if ((clusterCharge < clusterNoise * m_cutCluster) || (clusterSeedCharge < clusterNoise * m_cutSeed))
451 if (clusterSize < m_sizeHeadTail) {
452 clusterPosition = 1.0 / clusterCharge * inner_product(
453 stripAmplitudes.begin(), stripAmplitudes.end(), stripPositions.begin(), 0.0
456 float phantomCharge = m_cutAdjacent * clusterNoise;
457 if (clusterSize == 1) {
458 clusterPositionError = pitch * phantomCharge / (clusterCharge + phantomCharge);
460 clusterPositionError = pitch * phantomCharge / clusterCharge;
463 float leftStripCharge = stripAmplitudes.front();
464 float leftPos = stripPositions.front();
465 float rightStripCharge = stripAmplitudes.back();
466 float rightPos = stripPositions.back();
467 float centreCharge = (clusterCharge - leftStripCharge - rightStripCharge) / (clusterSize - 2);
468 leftStripCharge = (leftStripCharge < centreCharge) ? leftStripCharge : centreCharge;
469 rightStripCharge = (rightStripCharge < centreCharge) ? rightStripCharge : centreCharge;
470 clusterPosition = 0.5 * (leftPos + rightPos)
471 + 0.5 * (rightStripCharge - leftStripCharge) / centreCharge * pitch;
472 float sn = centreCharge / m_cutAdjacent / clusterNoise;
474 float landauHead = leftStripCharge / centreCharge;
475 double landauTail = rightStripCharge / centreCharge;
476 clusterPositionError = 0.5 * pitch * sqrt(1.0 / sn / sn
477 + 0.5 * landauHead * landauHead + 0.5 * landauTail * landauTail);
483 map<unsigned int, float> mc_relations;
484 map<unsigned int, float> truehit_relations;
485 vector<pair<unsigned int, float> > digit_weights;
486 digit_weights.reserve(clusterSize);
488 for (
size_t iDigit = clusterBounds.first; iDigit < clusterBounds.second; ++iDigit) {
490 fillRelationMap(m_mcRelation, mc_relations, iDigit);
491 fillRelationMap(m_trueRelation, truehit_relations, iDigit);
493 digit_weights.emplace_back(iDigit, stripAmplitudes[iDigit - clusterBounds.first]);
498 VxdID clusterSensorID = sensorID;
501 SVDCluster(sensorID, isU, clusterPosition, clusterPositionError, clusterTime,
502 clusterTimeErr, clusterCharge, clusterSeedCharge, clusterSize, clusterSN, clusterChi2)
506 if (!mc_relations.empty()) {
507 relClusterMCParticle.
add(clsIndex, mc_relations.begin(), mc_relations.end());
509 if (!truehit_relations.empty()) {
510 relClusterTrueHit.
add(clsIndex, truehit_relations.begin(), truehit_relations.end());
512 relClusterRecoDigit.
add(clsIndex, digit_weights.begin(), digit_weights.end());
517 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.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
OutputProbArray getProbabilities() const
Get signal time pdf.
VxdID getSensorID() const
Get the sensor ID.
short int getCellID() const
Get strip ID.
bool isUStrip() const
Get strip direction.
The SVD ShaperDigit class.
APVFloatSamples getSamples() const
Get array of samples.
std::vector< const RelationElement * > RelationLookup
Container for a RelationArray Lookup table.
Specific implementation of SensorInfo for SVD Sensors which provides additional sensor specific infor...
const TVector3 & getLorentzShift(double uCoord, double vCoord) const
Calculate Lorentz shift along a given coordinate in a magnetic field at a given position.
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.
double getVCellPosition(int vID) const
Return the position of a specific strip/pixel in v direction.
double getUPitch(double v=0) const
Return the pitch of the sensor.
double getUCellPosition(int uID, int vID=-1) const
Return the position of a specific strip/pixel in u direction.
double getVPitch(double v=0) const
Return the pitch of the sensor.
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.
Abstract base class for different kinds of events.