Belle II Software  release-06-02-00
SpacePoint2TrueHitConnectorModule.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/SpacePoint2TrueHitConnectorModule.h>
10 #include <framework/datastore/StoreObjPtr.h>
11 #include <framework/dataobjects/EventMetaData.h>
12 #include <framework/core/Environment.h> // getNumberProcesses
13 
14 #include <pxd/dataobjects/PXDCluster.h>
15 #include <svd/dataobjects/SVDCluster.h>
16 #include <mdst/dataobjects/MCParticle.h>
17 
18 #include <string>
19 #include <vector>
20 #include <cmath>
21 
22 #include <algorithm>
23 #include <unordered_map>
24 #include <tracking/spacePointCreation/MapHelperFunctions.h> // map helper stuff
25 
26 #include <TFile.h>
27 #include <TTree.h>
28 #include <TVector3.h>
29 
30 using namespace Belle2;
31 using namespace std;
32 
33 
34 REG_MODULE(SpacePoint2TrueHitConnector) // register the modules
35 
37  Module()
38 {
39  setDescription("Module that tries to find the appropriate TrueHit to each SpacePoint and to register a relation between them for making MC information for SpacePoints more easily accesible for Modules that need it. Module can also be used to filter out 'fishy' SpacePoints.");
40  setPropertyFlags(c_ParallelProcessingCertified);
41 
42  addParam("storeSeperate", m_PARAMstoreSeparate,
43  "Set to false if you do not want to create seperate StoreArrays for processed SpacePoints. (i.e. a relation from SpacePoint to TrueHit will be set in the passed StoreArray. NOTE: this StoreArray will contain SpacePoints with a relation to TrueHits and such without after this module). The Names of the output StoreArrays will be the names of the input StoreArrays with 'outputSuffix' (module parameter) appended to them",
44  true);
45 
46  addParam("registerAll", m_PARAMregisterAll,
47  "If set to true, the module simply registers a relation for all TrueHits that are related to a SpacePoint (resp. its Clusters). In this way the module can be used to find all related TrueHits and then the user can decide what to do with these TrueHits (otherwise this module does some decision making). Setting this to true means that all checks (e.g. 'minWeight', 'maxPosSigma', ...) are ommitted! NOTE that some of the information is lost in this way (e.g. how many of the Clusters of a SpacePoint have been related to a TrueHit)!",
48  false);
49 
50  addParam("positionAnalysis", m_PARAMpositionAnalysis,
51  "Analyze the positions of SpacePoints and corresponding TrueHits. NOTE: if enabled a root file gets created!", false);
52 
53  addParam("requirePrimary", m_PARAMrequirePrimary,
54  "Set to true if only relations to TrueHits that are related to a primary particle should get registered.", false);
55 
56  addParam("requireProximity", m_PARAMrequireProximity,
57  "Require that the TrueHit is close to the SpacePoint (in local coordinates). The meaning of 'close' can be defined with the parameters 'maxPosSigma' and 'maxGlobalPosDiff'.",
58  true);
59 
60  std::vector<std::string> defaultInList; // default list for input StoreArrays
61  defaultInList.push_back(std::string(""));
62  addParam("TrueHitNames", m_PARAMtrueHitNames,
63  "Container names of TrueHits. NOTE: need one name per 'DetectorType' (i.e. unique entries in 'DetectorType)!", defaultInList);
64  addParam("SpacePointNames", m_PARAMspacePointNames, "Container names of SpacePoints.", {"SVDSpacePoints", "PXDSpacePoints"});
65  addParam("DetectorTypes", m_PARAMdetectorTypes,
66  "detector types to determine which entries in 'TrueHitNames' and 'SpacePointNames' belong to which detector type. Entries have to be 'SVD' or 'PXD'. NOTE: if more 'SpacePointNames' than 'DetectorTypes' get passed, the last entry in 'DetectorTypes' is assumed to be valid for all remaining 'SpacePointNames'!");
67  addParam("ClusterNames", m_PARAMclusterNames,
68  "Container names of Clusters. NOTE: need one name per 'DetectorType' (i.e. unique entries in 'DetectorType')!", defaultInList);
69 
70  std::vector<std::string> defaultRootFName = { "PositionAnalysis", "RECREATE" };
71  addParam("rootFileName", m_PARAMrootFileName,
72  "Filename and write-mode ('RECREATE' or 'UPDATE'). If given more than 2 strings this module will cause termination",
73  defaultRootFName);
74 
75  addParam("outputSuffix", m_PARAMoutputSuffix,
76  "Suffix that will be appended to the container names if 'storeSeperate' is set to true", std::string("_relTH"));
77 
78  addParam("maxGlobalPosDiff", m_PARAMmaxGlobalDiff,
79  "max difference of global position coordinates between TrueHit and SpacePoint (in each direction) in cm.", 0.05);
80 // addParam("maxLocalPosDiff", m_PARAMmaxLocalDiff, "max difference of local position coordinates between TrueHit and SpacePoint (in U & V direction) in cm. NOTE: the default value is still subject to tuning and finding the appropriate value!", 0.01);
81 
82  addParam("maxPosSigma", m_PARAMmaxPosSigma, "Define the maximum local position difference in units of PitchSize / sqrt(12).", 4.);
83  addParam("minWeight", m_PARAMminWeight,
84  "Define a minimal weight a relation between a Cluster and a TrueHit has to have for the TrueHit to be considered as possible candidate.",
85  0.);
86 
87  // initialize all counters
88  initializeCounters();
89  m_rootFilePtr = nullptr;
90  m_treePtr = nullptr;
91 
92  if (m_PARAMpositionAnalysis == true and Environment::Instance().getNumberProcesses() > 0) {
93  B2WARNING(
94  "SpacePoint2TrueHitConnector::initialize: parameter positionAnalysis (and therefore root-output) is enabled and basf2 is running in multi-threaded mode - this can cause nondeterministic behavior! "
95  << "\n"
96  << " you can suppress multi-threading for this module by writing:"
97  << "\n"
98  << "main.add_module('SpacePoint2TrueHitConnector').set_property_flags(0) "
99  << "\n"
100  << "into the steering file!");
101  }
102 }
103 
104 // ================================================================ INITIALIZE ====================================================
105 void SpacePoint2TrueHitConnectorModule::initialize()
106 {
107  B2INFO("SpacePoint2TrueHitConnector -------------------------- initialize --------------------------------");
108 
109  m_nContainers = m_PARAMspacePointNames.size(); // get number of passed arrays
110  unsigned int nTHNames = m_PARAMtrueHitNames.size();
111  unsigned int nDetTypes = m_PARAMdetectorTypes.size();
112  unsigned int nClNames = m_PARAMclusterNames.size();
113  if (m_nContainers < nTHNames || m_nContainers < nDetTypes || m_nContainers < nClNames) {
114  B2FATAL("Passed " << nTHNames << " TrueHitNames and " << nDetTypes << " DetectorTypes but number of passed SpacePointArrays is " <<
115  m_nContainers);
116  }
117  if ((nTHNames != nDetTypes) || (nClNames != nTHNames)) {
118  B2FATAL("Passed " << nTHNames << " TrueHitNames and " << nClNames << "ClusterNames but " << nDetTypes << " DetectorTypes!");
119  }
120 
121  for (unsigned int i = 0; i < m_nContainers; ++i) {
122  if (i < nDetTypes) { // handle the detector types
123  std::string detType = m_PARAMdetectorTypes.at(i);
124  if (detType.compare(std::string("SVD")) != 0 && detType.compare(std::string("PXD")) != 0) {
125  B2FATAL("Found entry " << detType << " in DetectorTypes, but only 'PXD' and 'SVD' are allowed!");
126  }
127  if (detType.compare(std::string("SVD")) == 0) {
128  m_detectorTypes.push_back(c_SVD);
133  } else {
134  m_detectorTypes.push_back(c_PXD);
139  }
140  } else { m_detectorTypes.push_back(m_detectorTypes.at(i - 1)); } // add the last entry again
141 
143  m_inputSpacePoints.at(i).first.isRequired();
144 
145  // add counters
146  m_SpacePointsCtr.push_back(0);
147  m_nRelTrueHitsCtr.push_back(array<unsigned int, 5>());
148  m_noClusterCtr.push_back(0);
149  m_noTrueHitCtr.push_back(0);
150  m_regRelationsCtr.push_back(0);
151  m_ghostHitCtr.push_back(0);
152  m_rejectedRelsCtr.push_back(0);
153  }
154 
155  if (m_PARAMstoreSeparate) {
156  if (m_PARAMoutputSuffix.empty()) {
157  B2WARNING("'outputSuffix' is empty and 'storeSeperate' is set to true. This would lead to StoreArrays with the same name. Resetting to 'outputSuffix' to '_relTH'!");
158  m_PARAMoutputSuffix = "_relTH";
159  }
160 
161  for (unsigned int i = 0; i < m_nContainers; ++i) {
162  std::string name = m_inputSpacePoints.at(i).first.getName() + m_PARAMoutputSuffix;
165 
166  // register relations (detector dependent)
167  if (m_inputSpacePoints.at(i).second == c_SVD) {
170  } else {
173  }
174  }
175  } else {
176  for (unsigned int i = 0; i < m_nContainers; ++i) {
177  if (m_inputSpacePoints.at(i).second == c_SVD) {
179  } else {
181  }
182  }
183  }
184 
185  m_maxGlobalDiff = m_PARAMmaxGlobalDiff * m_PARAMmaxGlobalDiff; // only comparing squared values in module!
186 // m_PARAMmaxLocalDiff *= m_PARAMmaxLocalDiff;
187 
188  if (m_PARAMmaxPosSigma < 0) {
189  B2WARNING("'maxPosSigma' is set to a value below 0: " << m_PARAMmaxPosSigma << "! Resetting to default (4)!");
190  m_PARAMmaxPosSigma = 4.;
191  }
192 
194 }
195 
196 // ========================================================== EVENT ===============================================================
198 {
199  StoreObjPtr<EventMetaData> eventMetaDataPtr("EventMetaData", DataStore::c_Event);
200  const int eventCtr = eventMetaDataPtr->getEvent();
201  B2DEBUG(10, "SpacePoint2TrueHitConnector::event(). Processing event " << eventCtr << " -----------------------");
202 
203  if (m_PARAMpositionAnalysis) m_rootVariables = RootVariables(); // clear the rootVariables for every event
204 
205  // loop over all containers
206  for (m_iCont = 0; m_iCont < m_nContainers; ++m_iCont) {
207  StoreArray<SpacePoint> spacePoints = m_inputSpacePoints.at(m_iCont).first;
208  const int nSpacePoints = spacePoints.getEntries();
209  e_detTypes detType = m_inputSpacePoints.at(m_iCont).second;
210  std::string detTypeStr = detType == c_SVD ? "SVD" : "PXD";
211  B2DEBUG(10, "Found " << nSpacePoints << " SpacePoints in Array " << m_inputSpacePoints.at(m_iCont).first.getName() <<
212  " for this event. detType: " << detTypeStr);
213 
214  m_SpacePointsCtr.at(m_iCont) += nSpacePoints;
215 
216  for (int iSP = 0; iSP < nSpacePoints; ++iSP) {
217  SpacePoint* spacePoint = spacePoints[iSP];
218  B2DEBUG(49, "Processing SpacePoint " << iSP << " from " << nSpacePoints);
219 
220  baseMapT trueHitMap = processSpacePoint<baseMapT>(spacePoint, detType);
221  if (trueHitMap.empty()) continue; // next SpacePoint if something went wrong
222 
223  unsigned int nUniqueTHs = getUniqueSize(trueHitMap);
224  B2DEBUG(50, "Found " << nUniqueTHs << " TrueHits (unique) related to SpacePoint " << spacePoint->getArrayIndex() << " from Array "
225  << spacePoint->getArrayName());
226  unsigned int iRels = nUniqueTHs > 4 ? 4 : nUniqueTHs - 1;
227  m_nRelTrueHitsCtr.at(m_iCont).at(iRels)++;
228 
229  // print the complete map if the debug level is set high enough
230  if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 350, PACKAGENAME())) {
231  std::string mapCont = printMap(trueHitMap);
232  B2DEBUG(250, "The TrueHits and their weights for spacePoint " << spacePoint->getArrayIndex() << ": " + mapCont);
233  }
234 
235  int thIndex = -1; // declare negative and only assign a value if ONE relation gets registered! (see doc of positionAnalysis)
236  if (m_PARAMregisterAll) {
237  registerAllRelations(spacePoint, trueHitMap, detType);
238  } else { // find THE ONE TrueHit (to rule them all, one TrueHit to find them all ...)
239  // COULDDO: wrap this up in a function
240  pair<VXDTrueHit*, double> trueHitwWeight;
241 
242  if (detType == c_PXD) trueHitwWeight = getTHwithWeight<baseMapT, PXDTrueHit>(trueHitMap, m_PXDTrueHits, spacePoint, c_PXD);
243  else trueHitwWeight = getTHwithWeight<baseMapT, SVDTrueHit>(trueHitMap, m_SVDTrueHits, spacePoint, c_SVD);
244 
245  if (trueHitwWeight.first != nullptr) {
246  registerOneRelation(spacePoint, trueHitwWeight, detType);
247  thIndex = trueHitwWeight.first->getArrayIndex();
248  } else {
249  B2DEBUG(10, "Could not relate one TrueHit to SpacePoint " << spacePoint->getArrayIndex() << ".");
251  }
252  }
253 
254  if (m_PARAMpositionAnalysis) { positionAnalysis(spacePoint, trueHitMap, thIndex, detType); }
255  } // end loop SpacePoints
256  } // end loop containers
257  if (m_PARAMpositionAnalysis) { m_treePtr->Fill(); }
258 }
259 
260 // ================================================================ TERMINATE ======================================================
262 {
263  unsigned int sumSpacePoints = accumulate(m_SpacePointsCtr.begin(), m_SpacePointsCtr.end(), 0);
264  unsigned int sumRelations = accumulate(m_regRelationsCtr.begin(), m_regRelationsCtr.end(), 0);
265 
266  B2RESULT("SpacePoint2TrueHitConnector: Got " << sumSpacePoints << " SpacePoints in " << m_nContainers <<
267  " containers and registered " << sumRelations << " relations to TrueHits");
268  if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 1, PACKAGENAME())) {
269  stringstream contSumm;
270  contSumm << "Container-wise summary: \n";
271 
272  for (unsigned int iCont = 0; iCont < m_nContainers; ++iCont) {
273  contSumm << "In Container " << iCont << " (container name: " << m_inputSpacePoints[iCont].first.getName() << ") " <<
274  m_SpacePointsCtr[iCont] << " SpacePoints were contained. " << m_regRelationsCtr[iCont] <<
275  " relations were registered.\nNumber of related TrueHits to a SpacePoint are:\n";
276  for (unsigned int i = 0; i < m_nRelTrueHitsCtr[iCont].size() - 1; ++i) { contSumm << i + 1 << " related TrueHits to a SpacePoint : " << m_nRelTrueHitsCtr[iCont].at(i) << "\n"; }
277  contSumm << " more than 4 related TrueHits to a SpacePoint: " << m_nRelTrueHitsCtr[iCont].at(4) << "\n"; // WARNING: hardcoded
278  contSumm << m_rejectedRelsCtr.at(iCont) << " SpacePoints did not get a relation, " << m_ghostHitCtr.at(
279  iCont) << " were probably ghost hits in this container!\n";
280  contSumm << m_noTrueHitCtr[iCont] << " SpacePoints had no relation to a TrueHit at all.\n";
281  }
282  B2DEBUG(1, contSumm.str());
283 
284  if (!m_PARAMregisterAll) {
285  // TODO: do this containerwise
286  stringstream furtherSummary;
287  furtherSummary << "Ommited Relations because of weight < " << m_PARAMminWeight << ": " << m_weightTooSmallCtr << "\n";
288  furtherSummary << "Rejected Relations because of non primary particle: " << m_rejectedNoPrimaryCtr;
289 // furtherSummary << "Summary for all containers:\n";
290 // furtherSummary << "possible/accepted relations for cases:\n";
291 // furtherSummary << m_all2WTHCtr << "/" << m_accAll2WTHCtr << " SP with THs (more than one) with all THs having two weights\n";
292 // furtherSummary << m_single2WTHCtr << "/" << m_accSingle2WTHCtr << " SP with THs (more than one) with only one TH having two weights\n";
293 // furtherSummary << m_nonSingle2WTHCtr << "/" << m_accNonSingle2WTHCtr << " SP with THs (more than one) with more than one but not all THs having two weights\n";
294 // furtherSummary << "In " << m_oneCluster2THCtr << " cases there was a SP with only one Cluster but more than one related TrueHits";
295  B2DEBUG(2, furtherSummary.str());
296  }
297  }
298 
300 
301 // B2INFO("total number of weights: " << m_totWeightsCtr << " of which " << m_negWeightCtr << " were negative")
302 // B2INFO("m_moreThan2Weights = " << m_moreThan2Weights);
303 }
304 
305 // ====================================================== PROCESS SPACEPOINT ======================================================
306 template<typename MapType>
308 {
309  B2DEBUG(50, "Processing SpacePoint " << spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName());
310  MapType trueHitMap;
311  try {
312  if (detType == c_PXD) {
313  trueHitMap = getRelatedTrueHits<MapType, PXDCluster, PXDTrueHit>(spacePoint, m_PXDClusters.getName(), m_PXDTrueHits.getName());
314  } else {
315  trueHitMap = getRelatedTrueHits<MapType, SVDCluster, SVDTrueHit>(spacePoint, m_SVDClusters.getName(), m_SVDTrueHits.getName());
316  }
317  } catch (NoClusterToSpacePoint& anE) {
318  B2WARNING("Caught an exception while trying to relate SpacePoints and TrueHits: " << anE.what());
319  m_noClusterCtr.at(m_iCont)++;
320  } catch (...) { // catch the rest
321  B2ERROR("Caught undefined exception while trying to relate SpacePoints and TrueHits");
322  throw; // throw further (maybe it is caught somewhere else)
323  }
324 
325  B2DEBUG(499, "trueHitMap.size() before return in processSpacePoint: " << trueHitMap.size());
326  return trueHitMap;
327 }
328 
330 template<typename MapType, typename ClusterType, typename TrueHitType>
332  std::string trueHitName)
333 {
334  MapType trueHitsMap; // map to be filled with indices (keys) and weights (values)
335 
336  RelationVector<ClusterType> spacePointClusters = spacePoint->getRelationsTo<ClusterType>(clusterName);
337  if (spacePointClusters.size() == 0) {
338  B2DEBUG(1, "Found no related Cluster for SpacePoint " << spacePoint->getArrayIndex() << " from Array " <<
339  spacePoint->getArrayIndex());
340  throw NoClusterToSpacePoint();
341  }
342  B2DEBUG(75, "Found " << spacePointClusters.size() << " related Clusters to SpacePoint " << spacePoint->getArrayIndex() <<
343  " from Array " << spacePoint->getArrayName());
344 
345  // loop over all Clusters, get all TrueHits from them and add the information to the map
346  short noTrueHits = 0;
347  for (size_t iCl = 0; iCl < spacePointClusters.size(); ++iCl) {
348 
349  const ClusterType* cluster = spacePointClusters[iCl];
350  RelationVector<TrueHitType> clusterTrueHits = cluster->template getRelationsTo<TrueHitType>(trueHitName);
351  if (clusterTrueHits.size() == 0) {
352  B2DEBUG(3, "Found no related TrueHit for Cluster " << cluster->getArrayIndex() << " contained by SpacePoint " <<
353  spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName());
354  noTrueHits++;
355  } else {
356  B2DEBUG(80, "Found " << clusterTrueHits.size() << " related TrueHits to Cluster " << cluster->getArrayIndex() << " from Array " <<
357  cluster->getArrayName());
358 
359  for (unsigned int i = 0; i < clusterTrueHits.size(); ++i) { // 'TrueHit' loop
360  if ((clusterTrueHits.weight(i) < m_PARAMminWeight) && !m_PARAMregisterAll) {
362  continue; // not take into account relations with too little weight
363  }
364  int index = clusterTrueHits[i]->getArrayIndex();
365  if (trueHitsMap.find(index) == trueHitsMap.end()) { // create entry only if it is not already in the map
366  trueHitsMap[index] = TrueHitInfo(index);
367  B2DEBUG(499, "Added new TrueHitInfo to map. Index " << index);
368  }
369  if (spacePointClusters.weight(iCl) > 0) {
370  trueHitsMap[index].setUWeight(clusterTrueHits.weight(i));
371  B2DEBUG(499, "Added UCluster to TrueHitInfo with index " << index);
372  } else {
373  trueHitsMap[index].setVWeight(clusterTrueHits.weight(i));
374  B2DEBUG(499, "Added VCluster to TrueHitInfo with index " << index);
375  }
376  }
377  }
378  }
379  // check if there was no Cluster with a related TrueHit -> then throw
380  if (noTrueHits == spacePoint->getNClustersAssigned()) {
381  m_noTrueHitCtr[m_iCont]++; // increase counter for current container
382  }
383  return trueHitsMap;
384 }
385 
386 // ================================================================= INITIALIZE COUNTERS ==========================================
388 {
389  // NEW
390  for (unsigned int i = 0; i < m_SpacePointsCtr.size(); ++i) {
391  m_SpacePointsCtr.at(i) = 0;
392  m_regRelationsCtr.at(i) = 0;
393  m_noClusterCtr.at(i) = 0;
394  m_ghostHitCtr.at(i) = 0;
395  m_noTrueHitCtr.at(i) = 0;
396  m_rejectedRelsCtr.at(i) = 0;
397 
398  for (unsigned int j = 0; j < m_nRelTrueHitsCtr.at(i).size(); ++j) { m_nRelTrueHitsCtr.at(i).at(j) = 0; }
399  }
400 
403 
404  // initialize the following here because of cppcheck complaining (not initialized in constructor)
405  m_nContainers = 0;
406  m_maxGlobalDiff = 0.;
407  m_iCont = 0;
408 
409 // m_negWeightCtr = 0;
410 // m_totWeightsCtr = 0;
411 // m_moreThan2Weights = 0;
412 //
413 // m_single2WTHCtr = 0;
414 // m_nonSingle2WTHCtr = 0;
415 // m_all2WTHCtr = 0;
416 // m_accSingle2WTHCtr = 0;
417 // m_accNonSingle2WTHCtr = 0;
418 // m_accAll2WTHCtr = 0;
419 //
420 // m_oneCluster2THCtr = 0;
421 //
422 
423 }
424 
425 // ===================================================== REGISTER ALL RELATIONS ===================================================
426 template<typename MapType>
428 {
429  B2DEBUG(50, "Registering all possible relations for SpacePoint " << spacePoint->getArrayIndex() << " from Array " <<
430  spacePoint->getArrayName() << ". storeSeparate is set to " << m_PARAMstoreSeparate);
431  SpacePoint* newSP = spacePoint; // declaring pointer here, getting new pointer if storeSeparate is true
432 
433  if (m_PARAMstoreSeparate) { // if storing in separate Array, re-register the relations to the Clusters first
434  newSP = m_outputSpacePoints.at(m_iCont).appendNew(*spacePoint);
435  B2DEBUG(50, "Added new SpacePoint to Array " << m_outputSpacePoints[m_iCont].getName() << ".");
436  if (detType == c_PXD) reRegisterClusterRelations<PXDCluster>(spacePoint, newSP, m_PXDClusters.getName());
437  else reRegisterClusterRelations<SVDCluster>(spacePoint, newSP, m_SVDClusters.getName());
438  }
439 
440  std::vector<TrueHitInfo> trueHitInfos = getAllValues(trueHitMap);
441  // sort by the number of related Clusters first and than by weight (descending order)
442  std::sort(trueHitInfos.begin(), trueHitInfos.end());
443 
444  for (const TrueHitInfo& info : trueHitInfos) {
445  if (detType == c_PXD) registerTrueHitRelation<PXDTrueHit>(newSP, info.m_Id, 1, m_PXDTrueHits); // only one Cluster for PXDs!
446  else { // for SVD relation weight is depending on which Cluster is related to TrueHit!
447  double weight = calculateRelationWeight(info, spacePoint);
448  registerTrueHitRelation<SVDTrueHit>(newSP, info.m_Id, weight, m_SVDTrueHits);
449  }
450  }
451 }
452 
453 // =============================================== REGISTER ONE RELATION ==========================================================
454 template<typename TrueHitType>
456  std::pair<TrueHitType*, double> trueHitwWeight, e_detTypes detType)
457 {
458  TrueHitType* trueHit = trueHitwWeight.first;
459  B2DEBUG(50, "Registering relation to TrueHit " << trueHit->getArrayIndex() << " from Array " << trueHit->getArrayName());
460  SpacePoint* newSP = spacePoint; // declaring pointer here, getting new pointer if storeSeparate ist true
461 
462  if (m_PARAMstoreSeparate) { // if storing in separate Array, re-register the relations to the Clusters first
463  newSP = m_outputSpacePoints.at(m_iCont).appendNew(*spacePoint);
464  B2DEBUG(50, "Added new SpacePoint to Array " << m_outputSpacePoints[m_iCont].getName() << ".");
465  if (detType == c_PXD) reRegisterClusterRelations<PXDCluster>(spacePoint, newSP, m_PXDClusters.getName());
466  else reRegisterClusterRelations<SVDCluster>(spacePoint, newSP, m_SVDClusters.getName());
467  }
468 
469  newSP->addRelationTo(trueHit, trueHitwWeight.second);
471  B2DEBUG(50, "Added Relation to TrueHit " << trueHit->getArrayIndex() << " from Array " << trueHit->getArrayName() <<
472  " for SpacePoint " << spacePoint->getArrayIndex() << " (weight = " << trueHitwWeight.second << ")");
473 }
474 
475 // ========================================================= REREGISTER CLUSTER RELATIONS =========================================
476 template<typename ClusterType>
478  Belle2::SpacePoint* newSpacePoint, std::string clusterName)
479 {
480  B2DEBUG(100, "Registering the Relations to Clusters of SpacePoint " << origSpacePoint->getArrayIndex() << " in Array " <<
481  origSpacePoint->getArrayName() << " for SpacePoint " << newSpacePoint->getArrayIndex() << " in Array " <<
482  newSpacePoint->getArrayName());
483 
484  vector<pair<ClusterType*, double> > clustersAndWeights = getRelatedClusters<ClusterType>(origSpacePoint, clusterName);
485  for (auto aCluster : clustersAndWeights) {
486  newSpacePoint->addRelationTo(aCluster.first, aCluster.second);
487  B2DEBUG(100, "Registered Relation to Cluster " << aCluster.first->getArrayIndex() << " with weight " << aCluster.second);
488  }
489 }
490 
491 // =========================================================== GET RELATED CLUSTERS ===============================================
492 template<typename ClusterType>
493 std::vector<std::pair<ClusterType*, double> > SpacePoint2TrueHitConnectorModule::getRelatedClusters(Belle2::SpacePoint* spacePoint,
494  std::string clusterName)
495 {
496  vector<pair<ClusterType*, double> > indsAndWeights;
497  RelationVector<ClusterType> relClusters = spacePoint->getRelationsTo<ClusterType>(clusterName);
498 
499  for (unsigned int iCl = 0; iCl < relClusters.size(); ++iCl) {
500  indsAndWeights.push_back(make_pair(relClusters[iCl], relClusters.weight(iCl)));
501  }
502 
503  // safety measure, should not / cannot happen (checked before)
504  if (indsAndWeights.empty()) { B2ERROR("No Clusters related to SpacePoint " << spacePoint->getArrayIndex() << "!"); }
505 
506  return indsAndWeights;
507 }
508 
509 // ====================================================== REGISTER TRUEHIT RELATIONS ==============================================
510 template<typename TrueHitType>
513 {
514  TrueHitType* trueHit = trueHits[index];
515  spacePoint->addRelationTo(trueHit, weight);
516  m_regRelationsCtr.at(m_iCont)++; // increase counter of registered relations for this container
517  B2DEBUG(50, "Added Relation to TrueHit " << index << " from Array " << trueHits.getName() << " for SpacePoint " <<
518  spacePoint->getArrayIndex() << " (weight = " << weight << ")");
519 }
520 
521 // ====================================================== POSITION ANALYSIS =======================================================
522 // TODO: debug output
523 template<typename MapType>
525  const int& index, e_detTypes detType)
526 {
527  B2DEBUG(250, "Doing position analysis for SpacePoint " << spacePoint->getArrayIndex() << " from Array " <<
528  spacePoint->getArrayName());
529  simpleBitfield<unsigned short int> relationStatus;
530 
531  // TODO TODO TODO TODO TODO TODO TODO: remove if not needed, only for tessting at the moment (i.e. do not commit)
532  pair<unsigned short int, unsigned short int> clusterSizes = getClusterSizes(spacePoint, detType);
533  pair<double, double> positionError = getLocalError(spacePoint);
534  // TODO TODO TODO TODO TODO TODO TODO: remove if not needed, only for tessting at the moment (i.e. do not commit)
535 
536  // do some checks and set the relationStatus
537  pair<bool, bool> setUV = spacePoint->getIfClustersAssigned();
538  if (setUV.first) relationStatus.addStatus(c_SpacePointU);
539  if (setUV.second) relationStatus.addStatus(c_SpacePointV);
540 
541  const vector<TrueHitInfo> trueHitInfos = getAllValues(trueHitMap);
542  unsigned int nRelations = trueHitInfos.size();
543  if (nRelations != 1) relationStatus.addStatus(c_nonUniqueRelation);
544  else relationStatus.addStatus(c_clearHit); // TODO: check if this is according to the definition of cleanHit
545 
546  B2DEBUG(999, "SpacePoint has assigned U: " << setUV.first << ", V: " << setUV.second << ". Possible TrueHits: " << nRelations);
547 
548  pair<double, double> spLocalPos = getLocalPos(spacePoint);
549  unsigned short int vxdId = spacePoint->getVxdID();
550  short nClusters = spacePoint->getNClustersAssigned();
551 
552  // loop over all TrueHitInfos and check if it is a ghostHit, and if a noise hit is contained in the SpacePoint
553  bool onceU = false, onceV = false; // check if both Clusters are at least used once
554  bool twoClusters = false;
555  for (const TrueHitInfo& info : trueHitInfos) {
556  onceU = onceU || info.m_U; // if U is set for one of the possible TrueHits, onceU is true after this loop
557  onceV = onceV || info.m_V;
558  twoClusters = twoClusters || (info.getNClusters() == 2); // if one of the TrueHits has two Clusters in the SpacePoint -> true
559  }
560  bool allSet = onceU && onceV;
561 
562  if (nClusters > 1) { // at this stage noise Clusters and ghostHits can only appear if there are two Clusters in the SpacePoint
563  if (!allSet) relationStatus.addStatus(c_noiseCluster);
564  if (allSet && !twoClusters) relationStatus.addStatus(c_ghostHit);
565  }
566 
567  // get the status that has been set until now as from now on it can differ for every relation
568  unsigned short int overAllStatus = relationStatus.getStatus();
569 
570  for (const TrueHitInfo& info : trueHitInfos) {
571  relationStatus.setStatus(overAllStatus); // reset status
572  VXDTrueHit* trueHit;
573  if (detType == c_SVD) { trueHit = m_SVDTrueHits[info.m_Id]; }
574  else { trueHit = m_PXDTrueHits[info.m_Id]; }
575  pair<TVector3, TVector3> trueHitPos = getTrueHitPositions(trueHit);
576 
577  double weightU = info.m_wU;
578  double weightV = info.m_wV;
579 
580  MCParticle* mcParticle = trueHit->getRelatedFrom<MCParticle>("ALL");
581  if (mcParticle != nullptr) {
582  if (mcParticle->hasStatus(MCParticle::c_PrimaryParticle)) relationStatus.addStatus(c_primaryParticle);
583  }
584 
585  if (info.m_Id == index) relationStatus.addStatus(c_registeredRelation);
586 
587  m_rootVariables.SpacePointULocal.push_back(spLocalPos.first);
588  m_rootVariables.SpacePointVLocal.push_back(spLocalPos.second);
589  m_rootVariables.SpacePointXGlobal.push_back(spacePoint->X());
590  m_rootVariables.SpacePointYGlobal.push_back(spacePoint->Y());
591  m_rootVariables.SpacePointZGlobal.push_back(spacePoint->Z());
592 
593  m_rootVariables.TrueHitULocal.push_back(trueHitPos.first.X());
594  m_rootVariables.TrueHitVLocal.push_back(trueHitPos.first.Y());
595  m_rootVariables.TrueHitXGlobal.push_back(trueHitPos.second.X());
596  m_rootVariables.TrueHitYGlobal.push_back(trueHitPos.second.Y());
597  m_rootVariables.TrueHitZGlobal.push_back(trueHitPos.second.Z());
598 
599  m_rootVariables.NRelations.push_back(nRelations);
600  m_rootVariables.RelationStatus.push_back(relationStatus.getStatus());
601  m_rootVariables.WeightU.push_back(weightU);
602  m_rootVariables.WeightV.push_back(weightV);
603  m_rootVariables.HitVxdID.push_back(vxdId);
604 
605  // TODO TODO TODO TODO TODO TODO TODO: remove if not needed, only for tessting at the moment (i.e. do not commit)
606  m_rootVariables.ClusterSizeU.push_back(clusterSizes.first);
607  m_rootVariables.ClusterSizeV.push_back(clusterSizes.second);
608  m_rootVariables.SpacePointErrorU.push_back(positionError.first);
609  m_rootVariables.SpacePointErrorV.push_back(positionError.second);
610  m_rootVariables.SpacePointErrorX.push_back(spacePoint->getPositionError().X());
611  m_rootVariables.SpacePointErrorY.push_back(spacePoint->getPositionError().Y());
612  m_rootVariables.SpacePointErrorZ.push_back(spacePoint->getPositionError().Z());
613  // TODO TODO TODO TODO TODO TODO TODO: remove if not needed, only for tessting at the moment (i.e. do not commit)
614 
615  B2DEBUG(999, "Branch contents of this entry:\nSPLocalU: " << spLocalPos.first << ". SPLocalV: " << spLocalPos.second << "\n" << \
616  "SPGlobalX: " << spacePoint->X() << ", SPGlobalY: " << spacePoint->Y() << ", SPGlobalZ " << spacePoint->Z() << "\n" << \
617  "THLocalU: " << trueHitPos.first.X() << ", THLocalV: " << trueHitPos.first.Y() << "\n" << \
618  "THGlobalX: " << trueHitPos.second.X() << ", THGlobalY: " << trueHitPos.second.Y() << ", THGlobalZ: " << trueHitPos.second.Z() <<
619  "\n" << \
620  "weight1: " << weightU << ", weight2: " << weightV << ", VxdID: " << vxdId << ", nRelations: " << nRelations << ", relStatus: " <<
621  relationStatus.getStatus());
622  }
623 }
624 
625 // ====================================================== GET TH WITH WEIGHT ======================================================
626 template<typename MapType, typename TrueHitType>
627 std::pair<TrueHitType*, double>
629  Belle2::SpacePoint* spacePoint, e_detTypes detType)
630 {
631  vector<TrueHitInfo> trueHitInfos = getAllValues(aMap);
632  std::pair<TrueHitType*, double> THwithWeight(nullptr, 0.0); // default return value
633  // return nullptr pointer and zero weight if there is no TrueHit (safety measure that should not actually be needed!)
634  if (trueHitInfos.empty()) return THwithWeight;
635  std::sort(trueHitInfos.begin(), trueHitInfos.end()); // sort to have best candidates at beginning
636 
637  short nClusters = spacePoint->getNClustersAssigned();
638  size_t nRelations = trueHitInfos.size(); // get the number of possible relations
639 
640  B2DEBUG(50, "Trying to select one TrueHit for SpacePoint " << spacePoint->getArrayIndex() << " from Array " <<
641  spacePoint->getArrayName() << ". SpacePoint has " << nClusters << " Clusters.");
642  B2DEBUG(150, "There are " << nRelations << " possible candidates.");
643 
644  // very verbose output only to have a look on why these TrueHits could be in the same SpacePoint, COULDDO: wrap up in function
645  if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 999, PACKAGENAME())) {
646  stringstream output;
647  output << "The candidates are: ";
648  for (auto& info : trueHitInfos) { output << info; }
649  B2DEBUG(999, output.str());
650 
651  std::pair<double, double> spacePointLocal = getLocalPos(spacePoint);
652  std::pair<bool, bool> assignedLocal = spacePoint->getIfClustersAssigned();
653 
654  B2DEBUG(999, "SpacePoint " << spacePoint->getArrayIndex() << " U: " << spacePointLocal.first << " V: " << spacePointLocal.second <<
655  " assigned: " << assignedLocal.first << ", " << assignedLocal.second);
656 
657  std::vector<int> uniqueKeys = getUniqueKeys(aMap);
658  for (int key : uniqueKeys) {
659  TrueHitType* trueHit = trueHits[key];
660  // NOTE: assuming here that there is only one MCParticle to each TrueHit
661  MCParticle* mcParticle = trueHit->template getRelatedFrom<MCParticle>("ALL");
662 
663  int mcPartId = -1, pdgCode = 0;
664  bool primary = false;
665 
666  if (mcParticle != nullptr) {
667  mcPartId = mcParticle->getArrayIndex();
668  primary = mcParticle->hasStatus(MCParticle::c_PrimaryParticle);
669  pdgCode = mcParticle->getPDG();
670  }
671 
672  B2DEBUG(999, "TrueHit " << key << " U: " << trueHit->getU() << ", V: " << trueHit->getV() << " mc Particle Id: " << mcPartId <<
673  ", primary " << primary << ", pdg: " << pdgCode);
674  }
675  }
676 
677  // loop over all TrueHitInfos and check if the number of related Clusters is equal to the number of the number of the assigned
678  // Clusters in the SpacePoint -> if so, check for compatibility and return (the first that is compatible)
679  bool ghostHit = false;
680  for (const TrueHitInfo& info : trueHitInfos) {
681  B2DEBUG(499, "Now checking trueHit: " << info);
682  if (nClusters == info.getNClusters()) {
683  TrueHitType* trueHit = trueHits[info.m_Id];
684  if (compatibleCombination(spacePoint, trueHit)) {
685  double weight = (detType == c_PXD ? 1 : calculateRelationWeight(info, spacePoint)); // if PXD weight is always 1!
686  return make_pair(trueHit, weight);
687  }
688  } else {
689  B2DEBUG(499, "The number of related Clusters ( = " << nClusters << ") and the number of associated weights ( = " <<
690  info.getNClusters() << ") do not match! This indicates a ghost hit");
691  ghostHit = true;
692  }
693  }
694  if (ghostHit) m_ghostHitCtr[m_iCont]++;
695 
696  // B2DEBUG(4999, "just before return in getTHwithWeight")
697  return THwithWeight;
698 }
699 
700 // ============================================== CHECK IF SPACEPOINT AND TRUEHIT ARE COMPATIBLE ==================================
701 template <typename TrueHitType>
703 {
704  B2DEBUG(150, "Checking TrueHit " << trueHit->getArrayIndex() << " from Array " << trueHit->getArrayName() << " and SpacePoint " <<
705  spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName() << " for compatibility");
706 
707  // check primary first
708  MCParticle* mcPart = trueHit->template getRelatedFrom<MCParticle>("ALL");
709  bool primaryPart = false; // assume secondary or background for safety
710  if (mcPart != nullptr) primaryPart = mcPart->hasStatus(MCParticle::c_PrimaryParticle);
711  if (m_PARAMrequirePrimary && !primaryPart) {
712  B2DEBUG(150, "TrueHit is not related to a primary particle but 'requirePrimary' is set to true!");
714  return false;
715  }
716 
717  // CAUTION: if further tests are added, make sure that they do not get 'over-ruled' by this one (i.e. put before this if)
719  B2DEBUG(999, "Not checking positions because 'requireProximity' is set to false");
720  return true;
721  }
722  const VxdID spacePointVxdId = spacePoint->getVxdID();
723  const VxdID trueHitVxdId = trueHit->getSensorID();
724 
725  B2DEBUG(999, "Comparing the VxdIDs, SpacePoint: " << spacePointVxdId << ", TrueHit: " << trueHitVxdId);
726  if (spacePointVxdId != trueHitVxdId) {
727  B2DEBUG(150, "SpacePoint and TrueHit do not have the same VxdID. spacePoint: " << spacePointVxdId << ", trueHit: " << trueHitVxdId);
728  return false;
729  }
730 
732  trueHitVxdId); // only have to get one, since VxdIds are already the same at this point!
733  double maxUres = SensorInfoBase.getUPitch(trueHit->getV()) / sqrt(12.) * m_PARAMmaxPosSigma;
734  double maxVres = SensorInfoBase.getVPitch(trueHit->getV()) / sqrt(12.) * m_PARAMmaxPosSigma;
735 
736  B2DEBUG(999, "maximum residual in U: " << maxUres << ", in V: " << maxVres);
737 
738  const TVector3 trueHitLocalPos = TVector3(trueHit->getU(), trueHit->getV(), 0);
739  const TVector3 trueHitGlobalPos = SensorInfoBase.pointToGlobal(trueHitLocalPos, true); // uses alignment
740 
741  std::pair<double, double> spacePointLocal = getLocalPos(spacePoint);
742  // compare only those values of the local coordinates that have been set
743  std::pair<bool, bool> setCoordinates = spacePoint->getIfClustersAssigned();
744 
745  if (setCoordinates.first) {
746  B2DEBUG(999, "Comparing the U-coordinates, SpacePoint: " << spacePointLocal.first << ", TrueHit: " << trueHitLocalPos.X() <<
747  " -> diff: " << spacePointLocal.first - trueHitLocalPos.X());
748  if (pow(spacePointLocal.first - trueHitLocalPos.X(), 2) > maxUres * maxUres) {
749  B2DEBUG(150, "The local position difference in U direction is " << spacePointLocal.first - trueHitLocalPos.X() <<
750  " but maximum local position difference is set to: " << maxUres);
751  return false;
752  }
753  }
754  if (setCoordinates.second) {
755  B2DEBUG(999, "Comparing the V-coordinates, SpacePoint: " << spacePointLocal.second << ", TrueHit: " << trueHitLocalPos.Y() <<
756  " -> diff: " << spacePointLocal.second - trueHitLocalPos.Y());
757  if (pow(spacePointLocal.second - trueHitLocalPos.Y(), 2) > maxVres * maxVres) {
758  B2DEBUG(150, "The local position difference in V direction is " << spacePointLocal.second - trueHitLocalPos.Y() <<
759  " but maximum local position difference is set to: " << maxVres);
760  return false;
761  }
762  }
763 
764  // only if both local coordinates of a SpacePoint are set, compare also the global positions!
765  if (setCoordinates.first && setCoordinates.second) {
766  B2DEBUG(999, "Comparing the global positions, SpacePoint: (" << spacePoint->X() << "," << spacePoint->Y() << "," << spacePoint->Z()
767  << "), TrueHit: (" << trueHitGlobalPos.X() << "," << trueHitGlobalPos.Y() << "," << trueHitGlobalPos.Z() << ")");
768  if (pow(spacePoint->X() - trueHitGlobalPos.X(), 2) > m_maxGlobalDiff
769  || pow(spacePoint->Y() - trueHitGlobalPos.Y(), 2) > m_maxGlobalDiff ||
770  pow(spacePoint->Z() - trueHitGlobalPos.Z(), 2) > m_maxGlobalDiff) {
771  B2DEBUG(150, "The position differences are for X: " << spacePoint->X() - trueHitGlobalPos.X() << ", Y: " << spacePoint->Y() -
772  trueHitGlobalPos.Y() << " Z: " << spacePoint->Z() - trueHitGlobalPos.Z() << " but the maximum position difference is set to: " <<
773  sqrt(m_PARAMmaxGlobalDiff));
774  return false;
775  }
776  } else {
777  B2DEBUG(5, "For SpacePoint " << spacePoint->getArrayIndex() <<
778  " one of the local coordinates was not assigned. The global positions and the un-assigned local coordinate were not compared!");
779  }
780 
781  return true;
782 }
783 
784 // ===================================================== GET LOCAL SPACEPOINT COORDINATES =========================================
786 {
787  // get the normalized local coordinates from SpacePoint and convert them to local coordinates (have to do so because at the slanted parts the local U-position is dependant on the local V-position)
788  // NOTE: second way is to convert the global position of the SpacePoint to local position via the sensorInfoBase (yields same results)
789  double normU = spacePoint->getNormalizedLocalU();
790  double normV = spacePoint->getNormalizedLocalV();
791  return spacePoint->convertNormalizedToLocalCoordinates(std::make_pair(normU, normV), spacePoint->getVxdID());
792 }
793 
794 // =============================================================== GET TRUEHIT POSITIONS ==========================================
795 template<typename TrueHitType>
796 std::pair<TVector3, TVector3> SpacePoint2TrueHitConnectorModule::getTrueHitPositions(TrueHitType* trueHit)
797 {
798 // TrueHitType* trueHit = trueHits[index];
799  const TVector3 localPos = TVector3(trueHit->getU(), trueHit->getV(), 0);
800 
801  const VxdID trueHitVxdId = trueHit->getSensorID();
802  VXD::SensorInfoBase SensorInfoBase = VXD::GeoCache::getInstance().getSensorInfo(trueHitVxdId);
803  const TVector3 globalPos = SensorInfoBase.pointToGlobal(localPos, true); // uses alignment
804 
805  return make_pair(localPos, globalPos);
806 }
807 
808 // ================================================================= INITIALIZE ROOT FILE =========================================
810 {
811  if (m_PARAMrootFileName.size() != 2 || (m_PARAMrootFileName[1] != "UPDATE" && m_PARAMrootFileName[1] != "RECREATE")) {
812  string output;
813  for (string entry : m_PARAMrootFileName) {
814  output += "'" + entry + "' ";
815  }
816  B2FATAL("CurlingTrackCandSplitter::initialize() : rootFileName is set wrong: entries are: " << output);
817  }
818 
819  string fileName = m_PARAMrootFileName[0] + ".root";
820  m_rootFilePtr = new TFile(fileName.c_str(), m_PARAMrootFileName[1].c_str());
821  m_treePtr = new TTree("PosAnaTree", "Position Analysis");
822 
823  // link the variables
824  m_treePtr->Branch("SPLocalU", &m_rootVariables.SpacePointULocal);
825  m_treePtr->Branch("SPLocalV", &m_rootVariables.SpacePointVLocal);
826  m_treePtr->Branch("SPGlobalX", &m_rootVariables.SpacePointXGlobal);
827  m_treePtr->Branch("SPGlobalY", &m_rootVariables.SpacePointYGlobal);
828  m_treePtr->Branch("SPGlobalZ", &m_rootVariables.SpacePointZGlobal);
829 
830  m_treePtr->Branch("THLocalU", &m_rootVariables.TrueHitULocal);
831  m_treePtr->Branch("THLocalV", &m_rootVariables.TrueHitVLocal);
832  m_treePtr->Branch("THGlobalX", &m_rootVariables.TrueHitXGlobal);
833  m_treePtr->Branch("THGlobalY", &m_rootVariables.TrueHitYGlobal);
834  m_treePtr->Branch("THGlobalZ", &m_rootVariables.TrueHitZGlobal);
835 
836  m_treePtr->Branch("WeightU", &m_rootVariables.WeightU);
837  m_treePtr->Branch("WeightV", &m_rootVariables.WeightV);
838 
839  m_treePtr->Branch("relStatus", &m_rootVariables.RelationStatus);
840  m_treePtr->Branch("VxdID", &m_rootVariables.HitVxdID);
841  m_treePtr->Branch("nRelations", &m_rootVariables.NRelations);
842 
843  m_treePtr->Branch("clusterSizeU", &m_rootVariables.ClusterSizeU);
844  m_treePtr->Branch("clusterSizeV", &m_rootVariables.ClusterSizeV);
845  m_treePtr->Branch("SPErrorU", &m_rootVariables.SpacePointErrorU);
846  m_treePtr->Branch("SPErrorV", &m_rootVariables.SpacePointErrorV);
847  m_treePtr->Branch("SPErrorX", &m_rootVariables.SpacePointErrorX);
848  m_treePtr->Branch("SPErrorY", &m_rootVariables.SpacePointErrorY);
849  m_treePtr->Branch("SPErrorZ", &m_rootVariables.SpacePointErrorZ);
850 }
851 
852 // ========================================================= CLOSE ROOT FILE ======================================================
854 {
855  if (m_treePtr != nullptr && m_rootFilePtr != nullptr) {
856  m_rootFilePtr->cd(); //important! without this the famework root I/O (SimpleOutput etc) could mix with the root I/O of this module
857 // m_treePtr->Write(); // using TTree::Write() instead, which calls this any way
858  m_rootFilePtr->Write();
859  m_rootFilePtr->Close();
860  }
861 }
862 
863 // TODO TODO TODO TODO TODO TODO TODO: remove if not needed, only for tessting at the moment (i.e. do not commit)
864 std::pair<unsigned short int, unsigned short int> SpacePoint2TrueHitConnectorModule::getClusterSizes(Belle2::SpacePoint* spacePoint,
865  e_detTypes detType)
866 {
867  std::pair<unsigned short int, unsigned short int> clusterSizes = { 0, 0 };
868  if (detType == c_PXD) {
869  vector<pair<PXDCluster*, double> > relClusters = getRelatedClusters<PXDCluster>(spacePoint, m_PXDClusters.getName());
870  return make_pair(relClusters[0].first->getUSize(), relClusters[0].first->getVSize());
871  } else {
872  vector<pair<SVDCluster*, double> > relClusters = getRelatedClusters<SVDCluster>(spacePoint, m_SVDClusters.getName());
873  for (auto cluster : relClusters) {
874  if (cluster.first->isUCluster()) clusterSizes.first = cluster.first->getSize();
875  else clusterSizes.second = cluster.first->getSize();
876  }
877  }
878  return clusterSizes;
879 }
880 
882 {
883  auto detType = spacePoint->getType();
884  pair<double, double> errors = { -1., -1. };
885  if (detType == VXD::SensorInfoBase::PXD) {
886  vector<pair<PXDCluster*, double> > relClusters = getRelatedClusters<PXDCluster>(spacePoint, m_PXDClusters.getName());
887  errors.first = relClusters[0].first->getUSigma();
888  errors.second = relClusters[0].first->getVSigma();
889  } else if (detType == VXD::SensorInfoBase::SVD) {
890  vector<pair<SVDCluster*, double> > relClusters = getRelatedClusters<SVDCluster>(spacePoint, m_SVDClusters.getName());
891  for (auto cluster : relClusters) {
892  if (cluster.first->isUCluster()) { errors.first = cluster.first->getPositionSigma(); }
893  else { errors.second = cluster.first->getPositionSigma(); }
894  }
895  } else {
896  B2ERROR("Detector type not known in SpacePoint2TrueHitConnector::getLocalError() !");
897  }
898 
899  return errors;
900 }
901 // TODO TODO TODO TODO TODO TODO TODO: remove if not needed, only for tessting at the moment (i.e. do not commit)
902 
903 // ==================================================== CALCULATE RELATION WEIGHT ==================================================
905 {
906  bool isUAssigned = spacePoint->getIfClustersAssigned().first; // get if the Cluster in the SpacePoint is a U-Cluster
907  // get if two Clusters are related to the Truehit (always false for single Cluster SPs
908  bool bothClusters = trueHitInfo.m_U && trueHitInfo.m_V;
909  // calculate the additional weight that identifies which Clusters are related to this TrueHit
910  // 0 -> both Clusters have Relation to TrueHit, 10 -> only U-Cluster, 20 -> only V-Cluster has Relation to TrueHit
911  short addWeight = 20 - (bothClusters ? 20 : isUAssigned ? 10 : 0);
912 
913  return addWeight + (bothClusters ? 2 : 1);
914 }
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:420
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:416
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:418
@ 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
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
bool hasStatus(unsigned short int bitmask) const
Return if specific status bit is set.
Definition: MCParticle.h:129
int getArrayIndex() const
Get 0-based index of the particle in the corresponding MCParticle list.
Definition: MCParticle.h:244
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
const std::string & getName() const
Returns the name of the module.
Definition: Module.h:187
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
float weight(int index) const
Get weight with index.
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.
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
std::vector< std::string > m_PARAMclusterNames
names of containers of Clusters
std::pair< TrueHitType *, double > getTHwithWeight(const MapType &aMap, Belle2::StoreArray< TrueHitType > trueHits, Belle2::SpacePoint *spacePoint, e_detTypes detType)
get the TrueHit from information that is stored in the map (conditions are checked in the following o...
std::vector< std::string > m_PARAMtrueHitNames
names of containers of TrueHits
std::vector< unsigned int > m_SpacePointsCtr
Number of SpacePoints presented to the module.
std::vector< std::string > m_PARAMdetectorTypes
detector type names as strings to determine which name belongs to which detector type
unsigned int m_nContainers
number of passed containers -> storing the size of an input vector for not having to obtain it every ...
void reRegisterClusterRelations(Belle2::SpacePoint *origSpacePoint, Belle2::SpacePoint *newSpacePoint, std::string clusterName="ALL")
register all the relations to Clusters that origSpacePoint had for newSpacePoint
RootVariables m_rootVariables
Root variables used for collecting data eventwise.
std::vector< unsigned int > m_rejectedRelsCtr
Number of SpacePoints that were not related to a TrueHit (i.e.
void event() override
event: try to find the appropriate TrueHit to all SpacePoints
MapType getRelatedTrueHits(Belle2::SpacePoint *spacePoint, std::string clusterName="ALL", std::string trueHitName="ALL")
get all the related TrueHits to the SpacePoint, including their weights in a map (multimap!...
std::pair< double, double > getLocalError(Belle2::SpacePoint *spacePoint)
get the position error of SpacePoints in local coordinates @retuns .first is U position error,...
void registerAllRelations(Belle2::SpacePoint *spacePoint, MapType trueHitMap, e_detTypes detType)
register a Relation to all the TrueHits in the trueHitMap for the passed SpacePoint
std::vector< unsigned int > m_noClusterCtr
Number of SpacePoints without relation to a Cluster (i.e.
std::vector< std::pair< Belle2::StoreArray< Belle2::SpacePoint >, e_detTypes > > m_inputSpacePoints
StoreArray of all input SpacePoints.
double m_PARAMmaxGlobalDiff
maximum difference of global position coordinates for each direction between TrueHit and SpacePoint
void initializeRootFile()
initialize the root file that is used for output
std::unordered_map< int, TrueHitInfo > baseMapT
typedef for shorter notation throughout the module
std::vector< unsigned int > m_regRelationsCtr
Number of registered relations.
unsigned int m_weightTooSmallCtr
Count the omitted relations because of a too small weight.
void terminate() override
terminate: print some summary information
unsigned int m_iCont
'helper variable' needed to not have to pass one integer down to processSpacePoint only to have a han...
std::vector< std::pair< ClusterType *, double > > getRelatedClusters(Belle2::SpacePoint *spacePoint, std::string clusterName="ALL")
get the pointers to the related Clusters and the weight of the original Relations between the spacePo...
std::vector< std::array< unsigned int, 5 > > m_nRelTrueHitsCtr
counting different numbers of related TrueHits (to a SpacePoint) with one variable
void positionAnalysis(Belle2::SpacePoint *spacePoint, const MapType &trueHitMap, const int &index, e_detTypes detType)
Analyze the position of SpacePoints and corresponding TrueHits.
bool m_PARAMstoreSeparate
switch for storing the SpacePoints that can be related to a TrueHit into separate StoreArrays,...
unsigned int m_rejectedNoPrimaryCtr
Count how many times a relation was rejected because TrueHit was not related to primary.
std::vector< Belle2::StoreArray< Belle2::SpacePoint > > m_outputSpacePoints
StoreArray of all output SpacePoints.
std::pair< unsigned short int, unsigned short int > getClusterSizes(Belle2::SpacePoint *spacePoint, e_detTypes detType)
get the sizes of the related Clusters of a SpacePoint
MapType processSpacePoint(Belle2::SpacePoint *spacePoint, e_detTypes detType)
process a SpacePoint.
std::pair< TVector3, TVector3 > getTrueHitPositions(TrueHitType *trueHit)
get the local (.first) and global (.second) position of a TrueHit (passed by index)
std::vector< e_detTypes > m_detectorTypes
storing the detector types for each container in vector, needed in initialize
std::vector< std::string > m_PARAMrootFileName
name and update status of root file
double m_maxGlobalDiff
storing the squared value of m_PARAMmaxGlobalDiff here to not alter the parameter input
void registerTrueHitRelation(Belle2::SpacePoint *spacePoint, int index, double weight, Belle2::StoreArray< TrueHitType > trueHits)
register the relation between a SpacePoint and the TrueHit (passed by index in the correspoinding Tru...
std::vector< unsigned int > m_noTrueHitCtr
Number of SpacePoints that contained a Cluster to which no TrueHit could be found (i....
std::string m_PARAMoutputSuffix
suffix that will be appended to the StoreArray names of the output StoreArrays
Belle2::StoreArray< Belle2::SVDTrueHit > m_SVDTrueHits
SVDTrueHits StoreArray used throughout the module.
bool m_PARAMrequirePrimary
require the TrueHit to be related to a primary particle in order for the relation to get registered!
double m_PARAMmaxPosSigma
defining th maximum difference of local coordinates in units of PitchSize / sqrt(12)
Belle2::StoreArray< Belle2::PXDTrueHit > m_PXDTrueHits
PXDTrueHits StoreArray used throughout the module.
Belle2::StoreArray< Belle2::PXDCluster > m_PXDClusters
PXDTClusters StoreArray used throughout the module.
bool m_PARAMregisterAll
switch for registereing all relations for all TrueHits for all SpacePoints (there can be more than 1 ...
Belle2::StoreArray< Belle2::SVDCluster > m_SVDClusters
PXDClusters StoreArray used throughout the module.
void registerOneRelation(Belle2::SpacePoint *spacePoint, std::pair< TrueHitType *, double > trueHitwWeight, e_detTypes)
register Relation between SpacePoint and TrueHit (and if neccessary also between SpacePoint and Clust...
double calculateRelationWeight(const TrueHitInfo &trueHitInfo, Belle2::SpacePoint *spacePoint)
calculate the Relation weight to be used (for SVD only, although method works with PXD as well!...
bool m_PARAMrequireProximity
require the TrueHit to be close to the SpacePoint.
std::vector< std::string > m_PARAMspacePointNames
names of containers of SpacePoints
std::pair< double, double > getLocalPos(Belle2::SpacePoint *spacePoint)
get the local position of a SpacePoint
void initializeCounters()
initialize all counters to 0 WARNING: only call in constructor of module!
bool compatibleCombination(Belle2::SpacePoint *spacePoint, TrueHitType *trueHit)
compares the TrueHit and the SpacePoint positions (global) to decide whether they are compatible NOTE...
bool m_PARAMpositionAnalysis
switch for doing the analysis of positions of SpacePoints and TrueHits
double m_PARAMminWeight
define a minimal weight a relation between Cluster and TrueHit.
std::vector< unsigned int > m_ghostHitCtr
Number of SpacePoints that are considered ghost hits.
SpacePoint typically is build from 1 PXDCluster or 1-2 SVDClusters.
Definition: SpacePoint.h:42
unsigned short getNClustersAssigned() const
Returns the number of Clusters assigned to this SpacePoint.
Definition: SpacePoint.h:165
double getNormalizedLocalV() const
Return normalized local coordinates of the cluster in v (0 <= posV <= 1).
Definition: SpacePoint.h:154
VxdID getVxdID() const
Return the VxdID of the sensor on which the the cluster of the SpacePoint lives.
Definition: SpacePoint.h:148
const B2Vector3< double > & getPositionError() const
return the hitErrors in sigma of the global position
Definition: SpacePoint.h:141
std::pair< bool, bool > getIfClustersAssigned() const
Returns, if u(v)-coordinate is based on cluster information.
Definition: SpacePoint.h:162
Belle2::VXD::SensorInfoBase::SensorType getType() const
Return SensorType (PXD, SVD, ...) on which the SpacePoint lives.
Definition: SpacePoint.h:145
static std::pair< double, double > convertNormalizedToLocalCoordinates(const std::pair< double, double > &hitNormalized, Belle2::VxdID vxdID, const Belle2::VXD::SensorInfoBase *aSensorInfo=nullptr)
converts a hit in sensor-independent relative coordinates into local coordinate of given sensor.
Definition: SpacePoint.cc:179
double Z() const
return the z-value of the global position of the SpacePoint
Definition: SpacePoint.h:129
double X() const
return the x-value of the global position of the SpacePoint
Definition: SpacePoint.h:123
double getNormalizedLocalU() const
Return normalized local coordinates of the cluster in u (0 <= posU <= 1).
Definition: SpacePoint.h:151
double Y() const
return the y-value of the global position of the SpacePoint
Definition: SpacePoint.h:126
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
const std::string & getName() const
Return name under which the object is saved in the DataStore.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
Class VXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: VXDTrueHit.h:36
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:66
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:213
Base class to provide Sensor Information for PXD and SVD.
double getUPitch(double v=0) const
Return the pitch of the sensor.
double getVPitch(double v=0) const
Return the pitch of the sensor.
TVector3 pointToGlobal(const TVector3 &local, bool reco=false) const
Convert a point from local to global coordinates.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
helper class for setting up a bitfield that can be used to store several flags in one variable TODO: ...
void setStatus(T __statusBits)
set the status of the bitfield (CAUTION: overwrites any previously defined status!...
void addStatus(T __statusBits)
add a status to the bitfield (if it has not already been added)
const T getStatus() const
get the status of the bitfield
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
std::vector< typename MapType::key_type > getUniqueKeys(const MapType &aMap)
get the unique keys of a map (i.e.
std::string printMap(const MapType &aMap)
get the contents of the map as string.
unsigned int getUniqueSize(const MapType &aMap)
get the number of unique keys inside the map NOTE: for non-multimap this is the same as ....
std::vector< typename MapType::mapped_type > getAllValues(const MapType &aMap)
get all values in the map (i.e.
Abstract base class for different kinds of events.
helper struct to access root variables inside the module
std::vector< unsigned int > ClusterSizeV
size of the v-cluster (resp.
std::vector< double > WeightV
weight of the relation between the V-Cluster of the SpacePoint and the TrueHit
std::vector< double > SpacePointErrorZ
positiion error of SpacePoint in Z direction (global)
std::vector< double > TrueHitXGlobal
TrueHit global X-position.
std::vector< double > SpacePointErrorX
positiion error of SpacePoint in X direction (global)
std::vector< double > SpacePointVLocal
SpacePoint local V-position.
std::vector< unsigned int > ClusterSizeU
size of the u-cluster (resp.
std::vector< double > SpacePointErrorY
positiion error of SpacePoint in Y direction (global)
std::vector< unsigned int > NRelations
Number of related TrueHits to a SpacePoint.
std::vector< double > SpacePointXGlobal
SpacePoint global X-position.
std::vector< unsigned short int > HitVxdID
VxdID of the SpacePoint/TrueHit.
std::vector< double > SpacePointErrorU
position error of SpacePoint in U direction
std::vector< unsigned short int > RelationStatus
different flags of the relation stored in here (see c_relationStatus)
std::vector< double > SpacePointYGlobal
SpacePoint global Y-position.
std::vector< double > SpacePointULocal
SpacePoint local U-position.
std::vector< double > SpacePointZGlobal
SpacePoint global Z-position.
std::vector< double > TrueHitYGlobal
TrueHit global Y-position.
std::vector< double > SpacePointErrorV
position error of SpacePoint in V direction
std::vector< double > WeightU
weight of the relation between the U-Cluster of the SpacePoint and the TrueHit
std::vector< double > TrueHitZGlobal
TrueHit global Z-position.
helper struct that holds information that is needed for the registration of the relation between Spac...
bool m_U
if true, U-Cluster is used by SpacePoint
bool m_V
if true, V-Cluster is used by SpacePoint