Belle II Software development
SPTCRefereeModule.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 <tracking/modules/spacePointCreator/SPTCRefereeModule.h>
10#include <tracking/dataobjects/RecoTrack.h>
11
12#include <framework/datastore/StoreArray.h>
13#include <framework/datastore/StoreObjPtr.h>
14#include <framework/dataobjects/EventMetaData.h>
15
16#include <tracking/spacePointCreation/SpacePoint.h>
17#include <pxd/dataobjects/PXDTrueHit.h>
18#include <svd/dataobjects/SVDTrueHit.h>
19#include <vxd/dataobjects/VxdID.h>
20
21#include <framework/geometry/B2Vector3.h>
22#include <boost/range/adaptor/reversed.hpp> // for ranged based loops in reversed order
23
24using namespace Belle2;
25
26REG_MODULE(SPTCReferee); // register the module
27
29{
30 setDescription("Module that does some sanity checks on SpacePointTrackCands to prevent some problematic cases to be "
31 "forwarded to other modules that rely on 'unproblematic' cases (e.g. FilterCalculator). "
32 "Different checks can be enabled by setting the according flags. Using MC information for "
33 "the tests can also be switched on/off for tests where MC information can be helpful.");
35
36 // names
37 addParam("sptcName", m_PARAMsptcName,
38 "Container name of the SpacePointTrackCands to be checked (input)",
40 addParam("newArrayName", m_PARAMnewArrayName,
41 "Container name of SpacePointTrackCands if 'storeNewArray' is set to true",
43 addParam("curlingSuffix", m_PARAMcurlingSuffix,
44 "Suffix that will be used to get a name for the StoreArray in which the trackStubs that are obtained by "
45 "splitting a curling SPTC get stored. NOTE: If 'storeNewArray' is set to true, "
46 "this will not be used and all output SPTCs will be in the same Array!",
48
49 // flags
50 addParam("checkSameSensor", m_PARAMcheckSameSensor,
51 "Check if two subsequent SpacePoints are on the same sensor",
53 addParam("checkMinDistance", m_PARAMcheckMinDistance,
54 "Check if two subsequent SpacePoints are separated by more than 'minDistance'",
56 addParam("checkCurling", m_PARAMcheckCurling,
57 "Check the SpacePointTrackCand for curling behaviour and mark it as curling if it does",
59 addParam("splitCurlers", m_PARAMsplitCurlers,
60 "Split curling SpacePointTrackCands and save the TrackStubs in separate StoreArrays",
62 addParam("keepOnlyFirstPart", m_PARAMkeepOnlyFirstPart,
63 "Keep only the first part of a curling SpacePointTrackCand (e.g. when only this is needed)",
65 addParam("useMCInfo", m_PARAMuseMCInfo,
66 "Set to true if the use of MC information (e.g. from underlying TrueHits) for the checks is wanted, "
67 "and to false if the checks should all be done with information that can be obtained from "
68 "SpacePoints directly. NOTE: the tests without MC information have to be developed first!",
70 addParam("kickSpacePoint", m_PARAMkickSpacePoint,
71 "Set to true if only the 'problematic' SpacePoint shall be kicked and not the whole SpacePointTrackCand",
73 addParam("storeNewArray", m_PARAMstoreNewArray,
74 "Set to true if the checked SpacePointTrackCands should be stored in a new StoreArray."
75 "WARNING: all previously registered relations get lost in this way!",
77
78 // other
79 addParam("minDistance", m_PARAMminDistance,
80 "Minimal Distance [cm] that two subsequent SpacePoints have to be separated if 'checkMinDistance' is enabled",
82 addParam("setOrigin", m_PARAMsetOrigin, "WARNING: still need to find out the units that are used internally! "
83 "Reset origin to given point. Used for determining the direction of flight of a particle for a "
84 "given hit. Needs to be reset for e.g. testbeam, where origin is not at (0,0,0)",
86
87 addParam("minNumSpacePoints", m_PARAMminNumSpacePoints,
88 "minimum number of space points that a track candidate has to "
89 "contain (added later, set to 0 to reproduce old behavior",
91
92 addParam("checkIfFitted", m_PARAMcheckIfFitted,
93 "If true a flag is set in the SpacePointTrackCandidate if any related RecoTrack "
94 "with successful track fit is found",
96
97 // initialize counters (cppcheck)
99}
100
101// ======================================================================= INITIALIZE =============================================
103{
104 B2INFO("SPTCReferee::initialize(): ------------------------------------------------ ");
105 // check if StoreArray of SpacePointTrackCands is her
107
108 // register new StoreArray
110// StoreArray<SpacePointTrackCand> m_optionalOutputSpacePointTrackCands(m_PARAMnewArrayName);
114 } else {
116 B2DEBUG(20, "StoreArray name of the curling parts: " << m_curlingArrayName);
117// StoreArray<SpacePointTrackCand> m_curlingSpacePointTrackCands(m_curlingArrayName);
121 }
122
123 // sanity checks on the other parameters
125 if (m_PARAMminDistance < 0) {
126 B2WARNING("minDistance set to value below 0: " << m_PARAMminDistance <<
127 ", Taking the absolute value and resetting 'minDistance' to that!");
129 }
130 }
131
132 B2DEBUG(20, "Provided Parameters: checkSameSensor - " << m_PARAMcheckSameSensor << ", checkMinDistance - " <<
134 << ", checkCurling - " << m_PARAMcheckCurling << ", splitCurlers - " << m_PARAMsplitCurlers << ", keepOnlyFirstPart - " <<
135 m_PARAMkeepOnlyFirstPart << ", useMCInfo - " << m_PARAMuseMCInfo << ", kickSpacePoint - " << m_PARAMkickSpacePoint);
136 if (m_PARAMsetOrigin.size() != 3) {
137 B2WARNING("CurlingTrackCandSplitter::initialize: Provided origin is not a 3D point! Please provide 3 values (x,y,z). "
138 "Rejecting user input and setting origin to (0,0,0) for now!");
139 m_PARAMsetOrigin.clear();
140 m_PARAMsetOrigin.assign(3, 0);
141 }
143 B2DEBUG(20, "Set origin to (x,y,z): (" << m_origin.X() << "," << m_origin.Y() << "," << m_origin.Z() << ")");
144
146}
147
148// ======================================================================== EVENT =================================================
150{
151 StoreObjPtr<EventMetaData> eventMetaDataPtr("EventMetaData", DataStore::c_Event);
152 const int eventCtr = eventMetaDataPtr->getEvent();
153 B2DEBUG(20, "Processing event " << eventCtr << " -----------------------");
154
155 const int nTCs = m_inputSpacePointTrackCands.getEntries();
156
157 m_totalTrackCandCtr += nTCs;
158
159 B2DEBUG(20, "Found " << nTCs << " SpacePointTrackCands in Array " << m_inputSpacePointTrackCands.getName() << " for this event");
160
161 for (int iTC = 0; iTC < nTCs; ++iTC) { // loop over all TrackCands
163 B2DEBUG(20, "Processing SpacePointTrackCand " << iTC << ": It has " << trackCand->getNHits() << " SpacePoints in it");
164
165 if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 200, PACKAGENAME())) { trackCand->print(); }
166 B2DEBUG(20, "refereeStatus of TrackCand before tests: " << trackCand->getRefereeStatus() << " -> " <<
167 trackCand->getRefereeStatusString());
168
169 // if all tests will be performed -> add checkedByReferee status to the SPTC,
170 // CAUTION: if there are new tests this has to be updated!!!, WARNING: if curling check fails,
171 // checkedForCurling will return false but hasRefereeStatus(c_checkedByReferee) will return true after this module!
174 }
175 bool allChecksClean = true; // assume that all tests will be passed, change to false if one of them fails
176 CheckInfo prevChecksInfo;
177
178
179 // added check for the number of space points in the track candidate
180 if ((int)(trackCand->getNHits()) < m_PARAMminNumSpacePoints) {
181 allChecksClean = false;
182 }
183
184
185 // set a flag if a fitted recotrack is found for that trackcand
187 // take any related recotrack
188 RelationVector<RecoTrack> relatedRecoTracks = trackCand->getRelationsTo<RecoTrack>("ALL");
189 if (relatedRecoTracks.size() >= 1) {
190 // assume that there is only one!
191 if (relatedRecoTracks[0]->wasFitSuccessful()) {
193 } else {
194 allChecksClean = false;
195 B2DEBUG(20, "Found RecoTrack was not fitted! Will not use this track candidate for training.");
196 }
197 } else {
198 allChecksClean = false;
199 B2DEBUG(20, "No related RecoTrack found. Will not use that track candidate for training");
200 }
201 }
202
203
204 // check same sensors if desired
206 const std::vector<int> sameSensorInds = checkSameSensor(trackCand);
207 std::get<0>(prevChecksInfo) = sameSensorInds;
208 if (!sameSensorInds.empty()) {
210 allChecksClean = false;
211 // assign the actually removed indices to the prevChecksInfo
213 std::get<0>(prevChecksInfo) = removeSpacePoints(trackCand, sameSensorInds);
214 } else {
215 // only add status if the SpacePoints on the same sensors have not been removed!
217 }
218 } else {
219 B2DEBUG(20, "Found no two subsequent SpacePoints on the same sensor for this SpacePointTrackCand ("
220 << iTC << " in Array " << m_inputSpacePointTrackCands.getName() << ")");
221 }
223 B2DEBUG(20, "refereeStatus of TrackCand after checkSameSensor " << trackCand->getRefereeStatus() << " -> " <<
224 trackCand->getRefereeStatusString());
225 }
226
227
228 // check min distance if desired
230 const std::vector<int> lowDistanceInds = checkMinDistance(trackCand, m_PARAMminDistance);
231 std::get<1>(prevChecksInfo) = lowDistanceInds;
232 if (!lowDistanceInds.empty()) {
234 allChecksClean = false;
235 // assign the actually removed indices to the prevChecksInfo
237 std::get<1>(prevChecksInfo) = removeSpacePoints(trackCand, lowDistanceInds);
238 } else {
239 // only add status if the SpacePoints not far enough apart have not been removed!
241 }
242 } else {
243 B2DEBUG(20, "Found no two subsequent SpacePoints that were closer than " << m_PARAMminDistance <<
244 " cm together for this SpacePointTrackCand (" << iTC << " in Array " << m_inputSpacePointTrackCands.getName() << ")");
245 }
247 B2DEBUG(20, "refereeStatus of TrackCand after checkMinDistance " << trackCand->getRefereeStatus() << " -> " <<
248 trackCand->getRefereeStatusString());
249 }
250
251 // vector of TrackStubs that shall be saved to another StoreArray
252 std::vector<SpacePointTrackCand> curlingTrackStubs;
253 // check curling if desired
255 // setting the TrackStubIndex to 0 implies that this trackCand has been checked for curling.
256 // (If there is something wrong in the curling check this value is reset to -1!)
257 trackCand->setTrackStubIndex(0);
258 const std::vector<int> curlingSplitInds = checkCurling(trackCand, m_PARAMuseMCInfo);
259 if (!curlingSplitInds.empty()) {
260 if (!(curlingSplitInds.at(0) == 0 && curlingSplitInds.size() == 1)) {
261 // this means essentially that the direction of flight for this SPTC is inwards for all SpacePoints!
263 allChecksClean = false;
265 curlingTrackStubs = splitTrackCand(trackCand, curlingSplitInds, m_PARAMkeepOnlyFirstPart, prevChecksInfo, m_PARAMkickSpacePoint);
266 if (curlingTrackStubs.empty()) {
267 B2ERROR("The vector returned by splitTrackCand is empty!");
268 } // safety measure
269 }
270 // set this to the original SPTC only after splitting to avoid having this status in the trackStubs
272 } else {
273 B2DEBUG(20, "The only entry in the return vector of checkCurling is 0! The direction of flight is inwards for the whole SPTC!");
274 trackCand->setFlightDirection(false);
276 }
277 } else {
278 B2DEBUG(20, "SpacePointTrackCand " << trackCand->getArrayIndex() << " is not curling!");
279 }
280 B2DEBUG(20, "refereeStatus of TrackCand after checkCurling " << trackCand->getRefereeStatus() << " -> " <<
281 trackCand->getRefereeStatusString());
282 }
283
284 // PROCESSING AFTER CHECKS
285 if (allChecksClean) trackCand->addRefereeStatus(SpacePointTrackCand::c_checkedClean);
286
287 B2DEBUG(20, "referee Status of SPTC after referee module: " << trackCand->getRefereeStatus() << " -> " <<
288 trackCand->getRefereeStatusString());
289 if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 200, PACKAGENAME())) { trackCand->print();}
290
291 // store in appropriate StoreArray
293 if (!trackCand->isCurling()) { copyToNewStoreArray(trackCand, m_optionalOutputSpacePointTrackCands); }
294 else {
295 for (const SpacePointTrackCand& trackStub : curlingTrackStubs) { addToStoreArray(trackStub, m_optionalOutputSpacePointTrackCands, trackCand); }
296 }
297 } else {
298 if (trackCand->isCurling()) {
299 for (const SpacePointTrackCand& trackStub : curlingTrackStubs) { addToStoreArray(trackStub, m_curlingSpacePointTrackCands, trackCand); }
300 }
301 }
302 }
303}
304
305// ============================================================================= TERMINATE ========================================
307{
308 // TODO: info output more sophisticated
309 std::stringstream summary;
311 summary << "Checked for consecutive SpacePoints on same sensor and found "
312 << m_SameSensorCtr << " TrackCands showing this behavior.\n";
313 }
315 summary << "Checked for minimal distance between two consecutive SpacePoints and found "
316 << m_minDistanceCtr << " TrackCands with SpacePoints not far enough apart.\n";
317 }
319 summary << m_kickedSpacePointsCtr << " SpacePoints have been removed from SpacePointTrackCands\n";
320 }
322 summary << m_curlingTracksCtr << " SPTCs were curling. Registered "
323 << m_regTrackStubsCtr << " track stubs. 'splitCurlers' was set to "
324 << m_PARAMsplitCurlers << ", 'keepOnlyFirstPart' was set to "
325 << m_PARAMkeepOnlyFirstPart << ". There were "
326 << m_allInwardsCtr << " SPTCs that had flight direction 'inward' for all SpacePoints in them";
327 }
328
329 B2INFO("SPTCRefere::terminate(): Module got " << m_totalTrackCandCtr << " SpacePointTrackCands. \n" << summary.str());
330
332 B2WARNING("The curling checking without MC Information is at the moment at a very crude and unsophisticated state. "
333 "If you have MC information available you should use it to do this check!");
334 }
335}
336
337// ====================================================================== CHECK SAME SENSORS ======================================
339{
340 B2DEBUG(20, "Checking SpacePointTrackCand " << trackCand->getArrayIndex() << " from Array " << trackCand->getArrayName() <<
341 " for consecutive SpacePoints on the same sensor");
342 std::vector<int> sameSensorInds; // return vector
343
344 std::vector<const SpacePoint*> spacePoints = trackCand->getHits();
345
346 // catch cases where the TC has no space points! (Yes that happens!)
347 if (spacePoints.size() == 0) return sameSensorInds;
348
349 VxdID lastSensorId = spacePoints.at(0)->getVxdID();
350
351 for (unsigned int iSp = 1; iSp < spacePoints.size(); ++iSp) {
352 VxdID sensorId = spacePoints.at(iSp)->getVxdID();
353 B2DEBUG(20, "Checking SpacePoint " << iSp << ". (ArrayIndex " << spacePoints.at(iSp)->getArrayIndex() <<
354 ") SensorId of this SpacePoint: " << sensorId << ", SensorId of last SpacePoint: " << lastSensorId);
355 if (sensorId == lastSensorId) {
356 // push back the index of the first SpacePoint (50:50 chance of getting the right one without further testing) -> retrieving the other index is no big science from this index!!
357 sameSensorInds.push_back(iSp - 1);
358 B2DEBUG(20, "SpacePoint " << iSp << " and " << iSp - 1 << " are on the same sensor: " << sensorId);
359 }
360 lastSensorId = sensorId;
361 }
362
363 return sameSensorInds;
364}
365
366// ========================================================================= CHECK MIN DISTANCE ===================================
367const std::vector<int> SPTCRefereeModule::checkMinDistance(Belle2::SpacePointTrackCand* trackCand, double minDistance)
368{
369 B2DEBUG(20, "Checking the distances between consecutive SpacePoints for SpacePointTrackCand " << trackCand->getArrayIndex() <<
370 " from Array " << trackCand->getArrayIndex());
371 std::vector<int> lowDistanceInds; // return vector
372
373 std::vector<const SpacePoint*> spacePoints = trackCand->getHits();
374
375 // catch case where the track candidate has no spacepoints
376 if (spacePoints.size() == 0) return lowDistanceInds;
377
378 B2Vector3F oldPosition = spacePoints.at(0)->getPosition();
379
380 for (unsigned int iSp = 1; iSp < spacePoints.size(); ++iSp) {
381 B2Vector3F position = spacePoints.at(iSp)->getPosition();
382 B2Vector3F diffPos = oldPosition - position;
383 B2DEBUG(20, "Position of SpacePoint " << iSp << " (ArrayIndex " << spacePoints.at(iSp)->getArrayIndex() << "): (" << position.X() <<
384 "," << position.Y() << "," << position.Z() << "), Position of SpacePoint " << iSp - 1 << ": (" << oldPosition.X() << "," <<
385 oldPosition.Y() << "," << oldPosition.Z() << ") --> old - new = (" << diffPos.X() << "," << diffPos.Y() << "," << diffPos.Z() <<
386 ")");
387
388 if (diffPos.Mag() <= minDistance) {
389 B2DEBUG(20, "Position difference is " << diffPos.Mag() << " but minDistance is set to " << minDistance << ". SpacePoints: " << iSp
390 << " and " << iSp - 1);
391 // push back the index of the first SpacePoint (50:50 chance of getting the right one without further testing)
392 lowDistanceInds.push_back(iSp);
393 }
394 oldPosition = position;
395 }
396
397 return lowDistanceInds;
398}
399
400// ============================================================= REMOVE SPACEPOINTS ===============================================
401const std::vector<int>
402SPTCRefereeModule::removeSpacePoints(Belle2::SpacePointTrackCand* trackCand, const std::vector<int>& indsToRemove)
403{
404 std::vector<int> removedInds; // return vector
405 try {
406 unsigned int nInds = indsToRemove.size();
407 B2DEBUG(20, "Got " << nInds << " indices to remove from SPTC " << trackCand->getArrayIndex());
408
409 int nRemoved = 0;
410 for (int index : boost::adaptors::reverse(indsToRemove)) { // reverse iteration as trackCand gets 'resized' with every remove
411 B2DEBUG(20, "Removing " << nRemoved + 1 << " from " << nInds << ". index = " << index); // +1 only for better readability
412 trackCand->removeSpacePoint(index);
413 nRemoved++;
415 B2DEBUG(20, "Removed SpacePoint " << index << " from SPTC " << trackCand->getArrayIndex());
416 // NOTE: this way if a removed SpacePoint is "at the edge" between two trackStubs the status will be assigned to the second of those!
417 removedInds.push_back(index - (nInds - nRemoved));
418 }
420 } catch (SpacePointTrackCand::SPTCIndexOutOfBounds& anE) {
421 B2WARNING("Caught an Exception while trying to remove a SpacePoint from a SpacePointTrackCand: " << anE.what());
422 }
423
424 return removedInds;
425}
426// =========================================================== CHECK CURLING ======================================================
427const std::vector<int> SPTCRefereeModule::checkCurling(Belle2::SpacePointTrackCand* trackCand, bool useMCInfo)
428{
429 std::vector<int> splitInds; // return vector
430
431 //catch cases where there are no space points in the trackCand!
432 if (trackCand->getHits().size() == 0) return splitInds;
433
434 // Only do curling checking if useMCInfo is false, OR if useMCInfo is true if the SPTCs SpacePoints have been checked for a relation to TrueHits!
436
437 std::string mcInfoStr = useMCInfo ? std::string("with") : std::string("without");
438 B2DEBUG(20, "Checking SpacePointTrackCand " << trackCand->getArrayIndex() << " from Array " << trackCand->getArrayName() <<
439 " for curling behavior " << mcInfoStr << " MC Information");
440
441 // get the SpacePoints of the TrackCand
442 const std::vector<const SpacePoint*>& tcSpacePoints = trackCand->getHits();
443 B2DEBUG(20, "SPTC has " << tcSpacePoints.size() << " SpacePoints");
444
445 // get the directions of flight for every SpacePoint
446 const std::vector<bool> dirsOfFlight = getDirectionsOfFlight(tcSpacePoints, useMCInfo);
447
448 // if(trackCand->getNHits() != dirsOfFlight.size()) B2FATAL("did not get a direction of flight for every SpacePoint"); // should not /cannot happen
449
450 // loop over all entries of dirsOfFlight and compare them pair-wise. If they change -> add Index to splitInds.
451 if (!dirsOfFlight.at(0)) {
452 // if the direction of flight is inwards for the first hit, push_back 0 -> make information accessible from outside this function
453 splitInds.push_back(0);
454 B2DEBUG(20, "Direction of flight was inwards for first SpacePoint of this SPTC");
455 }
456 // DEBUG output
457 B2DEBUG(20, "Direction of flight is " << dirsOfFlight.at(0) << " for SpacePoint " << 0 << " of this SPTC");
458 for (unsigned int i = 1; i < dirsOfFlight.size(); ++i) {
459 B2DEBUG(20, "Direction of flight is " << dirsOfFlight.at(i) << " for SpacePoint " << i << " of this SPTC");
460 if (dirsOfFlight.at(i) ^ dirsOfFlight.at(i - 1)) {
461 splitInds.push_back(i); // NOTE: using the bitoperator for XOR here to determine if the bools differ!
462 B2DEBUG(20, "Direction of flight has changed from SpacePoint " << i - 1 << " to " << i << ".");
463 }
464 } // END DEBUG output
465 } else {
466 B2ERROR("'useMCInfo' is set to true, but SpacePoints of SPTC have not been checked for relations to TrueHits! Not Checking this SPTC for curling!");
467 trackCand->setTrackStubIndex(-1); // reset to not being checked for curling
468 }
469 return splitInds;
470}
471
472// ============================================================ SPLIT CURLING TRACK CAND ==========================================
473std::vector<Belle2::SpacePointTrackCand>
474SPTCRefereeModule::splitTrackCand(const Belle2::SpacePointTrackCand* trackCand, const std::vector<int>& splitIndices,
475 bool onlyFirstPart, const CheckInfo& prevChecksInfo, bool removedHits)
476{
477 std::vector<SpacePointTrackCand> trackStubs; // return vector
478
479 B2DEBUG(20, "Splitting SpacePointTrackCand " << trackCand->getArrayIndex() << " from Array " << trackCand->getArrayName() <<
480 ": number of entries in splitIndices " << splitIndices.size());
481 // int trackStub = 0;
482 bool dirOfFlight = splitIndices.at(0) != 0; // if first entry is zero the direction of flight is false (= ingoing)
483
484 B2DEBUG(20, "first entry of passed vector<int> is " << splitIndices.at(0) << " --> direction of flight is " << dirOfFlight);
485 // if the first entry of splitIndices is zero the first TrackStub is from 0 to second entry instead of from 0 to first entry
486 int firstLast = dirOfFlight ? splitIndices.at(0) : splitIndices.at(1);
487 std::vector<std::pair<int, int> >
488 rangeIndices; // .first is starting, .second is final index for each TrackStub. Store them in vector to be able to easily loop over them
489 rangeIndices.push_back(std::make_pair(0, firstLast));
490
491 if (!onlyFirstPart) { // if more than the first part is desired push_back the other ranges too
492 unsigned int iStart = dirOfFlight ? 1 : 2;
493 for (unsigned int i = iStart; i < splitIndices.size(); ++i) {
494 rangeIndices.push_back(std::make_pair(splitIndices.at(i - 1), splitIndices.at(i)));
495 }
496 // last TrackStub is from last split index to end of TrackCand
497 rangeIndices.push_back(std::make_pair(splitIndices.at(splitIndices.size() - 1), trackCand->getNHits()));
498 }
499 B2DEBUG(20, "There will be " << rangeIndices.size() << " TrackStubs created for this TrackCand. (size of the passed splitIndices: "
500 << splitIndices.size() << ", onlyFirstPart " << onlyFirstPart);
501
502 if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 999, PACKAGENAME())) {
503 std::stringstream dbOutput;
504 dbOutput << "The indices that will be used for splitting the SPTC: ";
505 for (auto entry : rangeIndices) { dbOutput << "[" << entry.first << "," << entry.second << ") "; }
506 B2DEBUG(20, dbOutput.str());
507 }
508
509 // loop over all entries in range indices and create a SpacePointTrackCand from it
510 for (unsigned int iTs = 0; iTs < rangeIndices.size(); ++iTs) {
511 int firstInd = rangeIndices.at(iTs).first;
512 int lastInd = rangeIndices.at(iTs).second;
513
514 unsigned short int refStatus = getCheckStatus(trackCand);
515
516 B2DEBUG(20, "Trying to create TrackStub from SPTC " << trackCand->getArrayIndex() << " with indices [" << firstInd << "," <<
517 lastInd << ")");
518 // encapsulate in try block to catch indices out of range
519 try {
520 const std::vector<const SpacePoint*> spacePoints = trackCand->getHitsInRange(firstInd, lastInd);
521 const std::vector<double> sortingParams = trackCand->getSortingParametersInRange(firstInd, lastInd);
522
523 // create new TrackCand
524 SpacePointTrackCand trackStub(spacePoints, trackCand->getPdgCode(), trackCand->getChargeSeed(), trackCand->getMcTrackID());
525 trackStub.setSortingParameters(sortingParams);
526
527 // set the state seed and the cov seed only for the first trackStub of the TrackCand
528 if (iTs < 1) {
529 trackStub.set6DSeed(trackCand->getStateSeed());
530 trackStub.setCovSeed(trackCand->getCovSeed());
531 }
532
533 // set the direction of flight and flip it afterwards, because next trackCand hs changed direction of flight
534 trackStub.setFlightDirection(dirOfFlight);
535 dirOfFlight = !dirOfFlight;
536
537 // trackStub index starts at 1 for curling SPTCs. NOTE: this might be subject to changes with the new bitfield in SpacePointTrackCand
538 trackStub.setTrackStubIndex(iTs + 1);
539
540 // determine and set the referee status of this trackStub based upon the information from the previous tests
541 const std::vector<int>& sameSensInds = std::get<0>(prevChecksInfo);
542 const std::vector<int>& lowDistInds = std::get<0>(prevChecksInfo);
543 bool hasSameSens = vectorHasValueBetween(sameSensInds, rangeIndices.at(iTs));
544 bool hasLowDist = vectorHasValueBetween(lowDistInds, rangeIndices.at(iTs));
545 if ((hasSameSens || hasLowDist) && removedHits) refStatus += SpacePointTrackCand::c_removedHits;
546 if (hasSameSens && !removedHits) refStatus += SpacePointTrackCand::c_hitsOnSameSensor;
547 if (hasLowDist && !removedHits) refStatus += SpacePointTrackCand::c_hitsLowDistance;
548
549 trackStub.setRefereeStatus(refStatus);
550 B2DEBUG(20, "Set TrackStubIndex " << iTs + 1 << " and refereeStatus " << trackStub.getRefereeStatus() <<
551 " for this trackStub (refStatus string: " << trackStub.getRefereeStatusString());
552
553 trackStubs.push_back(trackStub);
554 if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 499, PACKAGENAME())) { trackStub.print(); }
555 } catch (SpacePointTrackCand::SPTCIndexOutOfBounds& anE) {
556 B2WARNING("Caught an exception while trying to split a curling SpacePointTrackCand: " << anE.what() <<
557 " This trackStub will not be created!");
558 }
559 }
560
561 return trackStubs;
562}
563
564// ========================================================= GET DIRECTIONS OF FLIGHT =============================================
565const std::vector<bool>
566SPTCRefereeModule::getDirectionsOfFlight(const std::vector<const Belle2::SpacePoint*>& spacePoints, bool useMCInfo)
567{
568 std::vector<bool> dirsOfFlight; // return vector
569
570 if (useMCInfo) {
571 try {
572 for (const SpacePoint* spacePoint : spacePoints) { // loop over all SpacePoints
573 if (spacePoint->getType() == VXD::SensorInfoBase::PXD) {
574 dirsOfFlight.push_back(getDirOfFlightTrueHit<PXDTrueHit>(spacePoint, m_origin));
575 } else if (spacePoint->getType() == VXD::SensorInfoBase::SVD) {
576 dirsOfFlight.push_back(getDirOfFlightTrueHit<SVDTrueHit>(spacePoint, m_origin));
577 } else throw
578 SpacePointTrackCand::UnsupportedDetType(); // NOTE: should never happen, because SpacePointTrackCand can only handle PXD and SVD at the moment!
579 }
580 } catch (SpacePointTrackCand::UnsupportedDetType& anE) {
581 B2FATAL("Caught a fatal exception while checking if a SpacePointTrackCand curls: " <<
582 anE.what()); // FATAL because if this happens this needs some time to implement and it affects more than only this module!
583 }
584 } else {
585 dirsOfFlight = getDirsOfFlightSpacePoints(spacePoints, m_origin);
586 }
587
588 return dirsOfFlight;
589}
590
591// ================================================ GET DIRECTION OF FLIGHT FROM TRUEHIT ==========================================
592template <typename TrueHitType>
594{
595 TrueHitType* trueHit = spacePoint->template getRelatedTo<TrueHitType>("ALL"); // COULDDO: search only certain arrays
596
597 if (trueHit == nullptr) {B2FATAL("Found no TrueHit to SpacePoint " << spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName()); }
598
599 // get SensorId - needed for transforming local to global coordinates
600 // cppcheck-suppress nullPointerRedundantCheck
601 VxdID vxdID = trueHit->getSensorID();
602
603 const VXD::SensorInfoBase& sensorInfoBase = VXD::GeoCache::getInstance().getSensorInfo(vxdID);
604 // cppcheck-suppress nullPointerRedundantCheck
605 B2Vector3F position = sensorInfoBase.pointToGlobal(B2Vector3F(trueHit->getU(), trueHit->getV(), 0), true); // global position
606 // cppcheck-suppress nullPointerRedundantCheck
607 B2Vector3F momentum = sensorInfoBase.vectorToGlobal(trueHit->getMomentum(), true); // global momentum
608
609 B2DEBUG(20, "Getting the direction of flight for SpacePoint " << spacePoint->getArrayIndex() << ", related to TrueHit " <<
610 // cppcheck-suppress nullPointerRedundantCheck
611 trueHit->getArrayIndex() << ". Both are on Sensor " << vxdID << ". (TrueHit) Position: (" << position.x() << "," << position.y() <<
612 "," << position.z() << "), (TrueHit) Momentum: (" << momentum.x() << "," << momentum.y() << "," << momentum.z() << ")");
613
614 return getDirOfFlightPosMom(position, momentum, origin);
615}
616
617// ==================================================== GET DIRECTION OF FLIGHT FROM SPACEPOINT ===================================
618std::vector<bool>
619SPTCRefereeModule::getDirsOfFlightSpacePoints(const std::vector<const Belle2::SpacePoint*>& spacePoints, B2Vector3F origin)
620{
621 std::vector<bool> dirsOfFlight; // return vector
622
623 B2Vector3F oldPosition = origin; // assumption: first position is origin
624 for (unsigned int iSP = 0; iSP < spacePoints.size(); ++iSP) {
625 B2Vector3F position = spacePoints.at(iSP)->getPosition();
626 // estimate momentum by linearizing between old position and new position -> WARNING: not a very good estimate!!!
627 B2Vector3F momentumEst = position - oldPosition;
628 B2DEBUG(20, "Getting the direction of flight for SpacePoint " << spacePoints.at(iSP)->getArrayIndex() << ". Position: (" <<
629 position.x() << "," << position.y() << "," << position.z() << "), estimated momentum: (" << momentumEst.x() << "," <<
630 momentumEst.y() << "," << momentumEst.z() << ")");
631 dirsOfFlight.push_back(getDirOfFlightPosMom(position, momentumEst, origin));
632 oldPosition = position; // reassign for next round
633 }
634
635 return dirsOfFlight;
636}
637
638// =============================================== GET DIRECTION OF FLIGHT FROM POSITION AND MOMENTUM =============================
640{
641 // calculate the position relative to the set origin, and add the momentum to the position to get the direction of flight
642 B2Vector3F originToHit = position - origin;
643
644 B2DEBUG(20, "Position relative to origin: (" << originToHit.x() << "," << originToHit.y() << "," << originToHit.z() <<
645 "). Momentum : (" << momentum.x() << "," << momentum.y() << "," <<
646 momentum.z() << ").");
647
648 // get dot product of momentum and hit position for the perpendicular component only!
649 float dot_xy = originToHit.x() * momentum.x() + originToHit.y() * momentum.y();
650
651 B2DEBUG(20, "result dot product xy component between position and momentum: " << dot_xy);
652
653 if (dot_xy < 0) {
654 B2DEBUG(20, "Direction of flight is inwards for this hit");
655 return false;
656 } else {
657 B2DEBUG(20, "Direction of flight is outwards for this hit");
658 return true;
659 }
660}
661
662// ============================================ COPY TO NEW STORE ARRAY ===========================================================
665{
666 SpacePointTrackCand* newTC = newStoreArray.appendNew(*trackCand);
667 newTC->addRelationTo(trackCand);
668 B2DEBUG(20, "Added new SPTC to StoreArray " << newStoreArray.getName() << " and registered relation to SPTC " <<
669 trackCand->getArrayIndex() << " from Array " << trackCand->getArrayName());
670}
671
672// =================================================== ADD TO STORE ARRAY =========================================================
675 const Belle2::SpacePointTrackCand* origTrackCand)
676{
677 SpacePointTrackCand* newTC = storeArray.appendNew(trackCand);
678 newTC->addRelationTo(origTrackCand);
679 B2DEBUG(20, "Added new SPTC to StoreArray " << storeArray.getName() << " and registered relation to SPTC " <<
680 origTrackCand->getArrayIndex() << " from Array " << origTrackCand->getArrayName());
682}
683
684// ======================================================== GET CHECK STATUS ======================================================
686{
687 unsigned short int status = trackCand->getRefereeStatus();
691 return status;
692}
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
DataType y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:427
DataType z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:429
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:159
DataType x() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:425
void SetXYZ(DataType x, DataType y, DataType z)
set all coordinates using data type
Definition: B2Vector3.h:464
@ c_DontWriteOut
Object/array should be NOT saved by output modules.
Definition: DataStore.h:71
@ c_ErrorIfAlreadyRegistered
If the object/array was already registered, produce an error (aborting initialisation).
Definition: DataStore.h:72
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
@ c_Debug
Debug: for code development.
Definition: LogConfig.h:26
static LogSystem & Instance()
Static method to get a reference to the LogSystem instance.
Definition: LogSystem.cc:31
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
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
std::string getArrayName() const
Get name of array this object is stored in, or "" if not found.
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
unsigned int m_kickedSpacePointsCtr
counter of kicked SpacePoints
bool m_PARAMkeepOnlyFirstPart
parameter for keeping only the first part of a curling SpacePointTrackCand
std::vector< bool > getDirsOfFlightSpacePoints(const std::vector< const Belle2::SpacePoint * > &spacePoints, B2Vector3F origin)
get the directions of flight for a vector of SpacePoints using only information from SpacePoints (i....
std::string m_curlingArrayName
name of the StoreArray in which the trackStubs from a curling SPTC are stored
const std::vector< int > checkCurling(Belle2::SpacePointTrackCand *trackCand, bool useMCInfo)
Check if the SpacePointTrackCand shows curling behavior.
StoreArray< SpacePointTrackCand > m_inputSpacePointTrackCands
Input SpacePointTrackCands StoreArray.
unsigned int m_regTrackStubsCtr
counter for the number of track stubs that were registered by this module
StoreArray< SpacePointTrackCand > m_curlingSpacePointTrackCands
Curling SpacePointTrackCands StoreArray.
unsigned int m_SameSensorCtr
counter for TrackCands with SpacePoints on the same sensor
void copyToNewStoreArray(const Belle2::SpacePointTrackCand *trackCand, Belle2::StoreArray< Belle2::SpacePointTrackCand > newStoreArray)
copy the SpacePointTrackCand to a new StoreArray and register a relation to the original trackCand
double m_PARAMminDistance
minimal distance two subsequent SpacePoints have to be separated
void initialize() override
initialize: check StoreArrays, initialize counters, ...
unsigned int m_minDistanceCtr
counter for TrackCands with SpacePoints not far enough apart
const std::vector< int > checkMinDistance(Belle2::SpacePointTrackCand *trackCand, double minDistance)
Check if two subsequent SpacePoints are separated by at least the provided minDistance.
unsigned short int getCheckStatus(const Belle2::SpacePointTrackCand *trackCand)
get the checked referee status of a SPTC (i.e.
std::vector< double > m_PARAMsetOrigin
assumed interaction point from which the SpacePointTrackCands emerge.
void event() override
event: check SpacePointTrackCands
bool getDirOfFlightTrueHit(const Belle2::SpacePoint *spacePoint, B2Vector3F origin)
get the direction of flight for a SpacePoint by using information from the underlying TrueHit NOTE: t...
void terminate() override
terminate: print some summary information
std::tuple< std::vector< int >, std::vector< int > > CheckInfo
typedef for storing the outcome of previously done checks to have them available later.
unsigned int m_totalTrackCandCtr
counter for the total number of TrackCands
std::string m_PARAMsptcName
Name of input container of SpacePointTrackCands.
B2Vector3F m_origin
origin used internally.
std::vector< Belle2::SpacePointTrackCand > splitTrackCand(const Belle2::SpacePointTrackCand *trackCand, const std::vector< int > &splitIndices, bool onlyFirstPart, const CheckInfo &prevChecksInfo, bool removedHits)
split a curling SpacePointTrackCand into TrackStubs.
StoreArray< SpacePointTrackCand > m_optionalOutputSpacePointTrackCands
Optional output SpacePointTrackCands StoreArray.
int m_PARAMminNumSpacePoints
only keep track candidates which have at least m_PARAMminNumSpacePoints space points
bool getDirOfFlightPosMom(B2Vector3F position, B2Vector3F momentum, B2Vector3F origin)
get the direction of flight provided the global position and momentum of a SpacePoint/TrueHit for the...
unsigned int m_allInwardsCtr
counter for the number of SPTCs which have direction of flight inward for all SpacePoints in them
bool m_PARAMkickSpacePoint
parameter for indicating if only the 'problematic' SpacePoint shall be removed from the SPTC or if th...
bool m_PARAMstoreNewArray
parameter for indicating if all checked SpacePointTrackCands should be stored in a new StoreArray NOT...
unsigned int m_curlingTracksCtr
counter for tracks that curl
void addToStoreArray(const Belle2::SpacePointTrackCand &trackCand, Belle2::StoreArray< Belle2::SpacePointTrackCand > storeArray, const Belle2::SpacePointTrackCand *origTrackCand)
register the SpacePointTrackCand (i.e.
bool vectorHasValueBetween(std::vector< T > V, std::pair< T, T > P)
function to determine if any of the values in vector V are between the values of P (i....
const std::vector< int > removeSpacePoints(Belle2::SpacePointTrackCand *trackCand, const std::vector< int > &indsToRemove)
remove the SpacePoint passed to this function from the SpacePointTrackCand
bool m_PARAMcheckSameSensor
parameter for indicating if the check for subsequent SpacePoints being on the same sensor should be d...
bool m_PARAMcheckMinDistance
parameter for indicating if the check for the minimal distance between two subsequent SpacePoints sho...
std::string m_PARAMcurlingSuffix
Suffix that will be used to get a name for the StoreArray that holds the trackStubs that were obtaine...
std::string m_PARAMnewArrayName
Name of the output container of SpacePointTrackCands if 'storeNewArray' is set to true.
const std::vector< int > checkSameSensor(Belle2::SpacePointTrackCand *trackCand)
Check if two subsequent SpacePoints are on the same sensor.
void initializeCounters()
initialize all counters to 0
bool m_PARAMuseMCInfo
parameter for indicating if MC information should be used or not
const std::vector< bool > getDirectionsOfFlight(const std::vector< const Belle2::SpacePoint * > &spacePoints, bool useMCInfo)
get the directions of Flight for every SpacePoint in the passed vector.
bool m_PARAMcheckIfFitted
if true it is looked for any related RecoTrack and if that RecoTrack has a valid fit.
bool m_PARAMsplitCurlers
parameter for switching on/off the splitting of curling SpacePointTrackCands
bool m_PARAMcheckCurling
parameter for indicating if the SpacePointTrackCand should be checked for curling
Storage for (VXD) SpacePoint-based track candidates.
void setSortingParameters(const std::vector< double > &sortParams)
set the sorting parameters
void set6DSeed(const TVectorD &state6D)
set the 6D state seed
void print(int debuglevel=25, const Option_t *="") const
print the Track Candidate in its "full beauty".
unsigned int getNHits() const
get the number of hits (space points) in the track candidate
int getPdgCode() const
get pdg code
const TMatrixDSym & getCovSeed() const
get the covariance matrix seed (6D).
const std::vector< const Belle2::SpacePoint * > & getHits() const
get hits (space points) of track candidate
void removeSpacePoint(int indexInTrackCand)
remove a SpacePoint (and its sorting parameter) from the SpacePointTrackCand
unsigned short int getRefereeStatus(unsigned short int bitmask=USHRT_MAX) const
Return the referee status code of the SpacePointTrackCand.
void setCovSeed(const TMatrixDSym &cov)
set the covariance matrix seed
int getMcTrackID() const
get the MC Track ID
@ c_curlingTrack
bit 8: SPTC is curling (resp.
@ c_checkedTrueHits
bit 5: All SpacePoints of the SPTC have a relation to at least one TrueHit.
@ c_removedHits
bit 4: SpacePoints were removed from this SPTC.
@ c_hitsLowDistance
bit 3: SPTC has two (or more) SpacePoints that are not far enough apart.
@ c_checkedClean
bit 1: SPTC shows no 'problematic' behaviour.
@ c_checkedMinDistance
bit 7: It has been checked if two consecutive SpacePoints are far enough apart.
@ c_hasFittedRecoTrack
bit 13: SPTC is related to a RecoTrack which has a successful fit.
@ c_checkedSameSensors
bit 6: It has been checked, if two consecutive SpacePoints are on the same sensor for this SPTC.
@ c_checkedByReferee
bit 0: SPTC has been checked by a Referee (all possible tests).
@ c_hitsOnSameSensor
bit 2: SPTC has two (or more) SpacePoints on same sensor.
bool isCurling() const
get if the TrackCand is curling.
void setRefereeStatus(unsigned short int bitmask)
set referee status (resets the complete to the passed status!)
void setTrackStubIndex(int trackStubInd)
set TrackStub index
const std::vector< const Belle2::SpacePoint * > getHitsInRange(int firstInd, int lastInd) const
get hits (SpacePoints) in range (indices of SpacePoint inside SpacePointTrackCand) including first in...
bool hasRefereeStatus(unsigned int short bitmask) const
Check if the SpacePointTrackCand has the status characterized by the bitmask.
const std::vector< double > getSortingParametersInRange(int firstIndex, int lastIndex) const
get the sorting parameters in range (indices of SpacePoints inside SpacePointTrackCand) including fir...
std::string getRefereeStatusString(std::string delimiter=" ") const
get the refereeStatus as a string (easier to read than an unsigned short int)
double getChargeSeed() const
get charge
const TVectorD & getStateSeed() const
get state seed as 6D vector
void setFlightDirection(bool direction)
set the direction of flight (true is outgoing, false is ingoing).
void addRefereeStatus(unsigned short int bitmask)
add a referee status
SpacePoint typically is build from 1 PXDCluster or 1-2 SVDClusters.
Definition: SpacePoint.h:42
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 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
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
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
Base class to provide Sensor Information for PXD and SVD.
ROOT::Math::XYZVector pointToGlobal(const ROOT::Math::XYZVector &local, bool reco=false) const
Convert a point from local to global coordinates.
ROOT::Math::XYZVector vectorToGlobal(const ROOT::Math::XYZVector &local, bool reco=false) const
Convert a vector from local to global coordinates.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
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
B2Vector3< float > B2Vector3F
typedef for common usage with float
Definition: B2Vector3.h:519
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
DataType at(unsigned i) const
safe member access (with boundary check!)
Definition: B2Vector3.h:751
Abstract base class for different kinds of events.