Belle II Software  release-06-02-00
CurlingTrackCandSplitterModule.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/CurlingTrackCandSplitterModule.h>
10 
11 // DataStore Stuff
12 #include <framework/datastore/StoreArray.h>
13 #include <framework/dataobjects/EventMetaData.h>
14 #include <framework/datastore/StoreObjPtr.h>
15 
16 // Clusters and TrueHits
17 #include <pxd/dataobjects/PXDCluster.h>
18 #include <pxd/dataobjects/PXDTrueHit.h>
19 #include <svd/dataobjects/SVDCluster.h>
20 #include <svd/dataobjects/SVDTrueHit.h>
21 
22 // SpacePoint related stuff
23 #include <tracking/spacePointCreation/SpacePoint.h>
24 #include <tracking/spacePointCreation/SpacePointTrackCand.h>
25 
26 #include <boost/format.hpp>
27 
28 #include <algorithm> // sort & unique
29 
30 using namespace std;
31 using namespace Belle2;
32 
33 REG_MODULE(CurlingTrackCandSplitter)
34 
35 // COULDDO: somehow retrieve the names under which Clusters and TrueHits are stored in the StoreArray (possible?) Or make them a parameter of the module. Otherwise it can happen, that, if more than one StoreArray of XYZCluster is present, all of them get searched, which might be unintended (or even undefined behaviour)
36 
38 {
39  setDescription("Module for checking SpacePointTrackCands for curling behaviour and (if wanted) splitting them into SpacePointTrackCands that no longer show curling behaviour. WARNING: MODULE IS DEPRECATED use SPTCReferee instead!");
40 
41  addParam("splitCurlers", m_PARAMsplitCurlers,
42  "Split curling SpacePointTrackCands into non-curling SpacePointTrackCands and store them", true);
43  addParam("nTrackStubs", m_PARAMnTrackStubs,
44  "Maximum number of SpacePointTrackCand Stubs to be created from a curling SpacePointTrackCand. Set to 0 if you want all possible TrackCand Stubs",
45  0);
46 
47  addParam("SpacePointTCName", m_PARAMsptcName, "Collection name of the SpacePointTrackCands to be analyzed for curling behaviour",
48  std::string(""));
49  addParam("curlingFirstOutName", m_PARAMcurlingOutFirstName,
50  "Collection name under which the first outgoing part of a curling TrackCand will be stored in the StoreArray. The first part of a curling Track has its origin at the interaction point.",
51  std::string(""));
52  addParam("curlingAllInName", m_PARAMcurlingAllInName,
53  "Collection name under which all ingoing parts of a curling TrackCand will be stored in the StoreArray", std::string(""));
54  addParam("curlingRestOutName", m_PARAMcurlingOutRestName,
55  "Collection name under which all but the first outgoing parts of a curling TrackCand will be stored in the StoreArray",
56  std::string(""));
57  addParam("completeCurlerName", m_PARAMcompleteCurlerName,
58  "Collection name under which all parts of a curling TrackCand will be stored in the StoreArray together. NOTE: only if this parameter is set to a non-empty string a complete (but splitted) curling TrackCand will be stored!",
59  std::string(""));
60 
61  // WARNING TODO: find out the units that are used internally!!!
62  std::vector<double> defaultOrigin = { 0., 0., 0. };
63  addParam("setOrigin", m_PARAMsetOrigin,
64  "WARNING: still need to find out the units that are used internally! Reset origin to given point. Used for determining the direction of flight of a particle for a given hit. Needs to be reset for e.g. testbeam, where origin is not at (0,0,0)",
65  defaultOrigin);
66 
67  addParam("positionAnalysis", m_PARAMpositionAnalysis,
68  "Set to true to investigate the positions of SpacePoints and TrueHits and write them to a ROOT file", false);
69 
70  std::vector<std::string> defaultRootFName;
71  defaultRootFName.push_back("PositionResiduals");
72  defaultRootFName.push_back("RECREATE");
73 
74  addParam("rootFileName", m_PARAMrootFileName,
75  "Filename and write-mode ('RECREATE' or 'UPDATE'). If given more than 2 strings this module will cause termination",
76  defaultRootFName);
77 
78 
79  addParam("useNonSingleTHinPA", m_PARAMuseNonSingleTHinPA,
80  "Switch for using SpacePoints in position Analysis that are related to more than one TrueHit", false);
81 
82  initializeCounters(); // NOTE: they get initialized in initialize again!!
83 
84  // initialize other variables to some default values to avoid unintended behaviour
85  m_saveCompleteCurler = false;
86  m_treePtr = nullptr;
87  m_rootFilePtr = nullptr;
88 }
89 
90 // ================================================= INITIALIZE =========================================================
91 void CurlingTrackCandSplitterModule::initialize()
92 {
93  initializeCounters();
94  B2INFO("CurlingTrackCandSplitter ----------------------------- initialize() -------------------------------------");
95  B2WARNING("CurlingTrackCandSplitter is deprecated and will be removed from framework in the near future! use SPTCReferee instead!");
96  // check if all necessary StoreArrays are present
97  StoreArray<SpacePointTrackCand> spacePointTCs(m_PARAMsptcName);
98  spacePointTCs.isRequired(m_PARAMsptcName);
99 
100  // count all empty input parameter strings, and issue a warning if more than one is empty (COULDDO: B2FATAL instead of warning)
101  int emptyCtr = 0;
102  if (m_PARAMcurlingOutFirstName.empty()) { emptyCtr++; }
103  if (m_PARAMcurlingAllInName.empty()) { emptyCtr++; }
104  if (m_PARAMcurlingOutRestName.empty()) { emptyCtr++; }
105 
106  if (emptyCtr > 1) {
107  B2WARNING("CurlingTrackCandSplitter::initialize: More than one of your input strings for the collection names is empty. This can lead to undeterministic behaviour since two or more collections will be stored under the same name!");
108  }
109 
110  // register new StoreArrays, and relation to original TrackCand
111  StoreArray<SpacePointTrackCand> curlingFirstOuts(m_PARAMcurlingOutFirstName);
112  curlingFirstOuts.registerInDataStore(m_PARAMcurlingOutFirstName, DataStore::c_ErrorIfAlreadyRegistered);
113  curlingFirstOuts.registerRelationTo(spacePointTCs);
114 
115  StoreArray<SpacePointTrackCand> curlingAllIns(m_PARAMcurlingAllInName);
116  curlingAllIns.registerInDataStore(m_PARAMcurlingAllInName, DataStore::c_ErrorIfAlreadyRegistered);
117  curlingAllIns.registerRelationTo(spacePointTCs);
118 
119  StoreArray<SpacePointTrackCand> curlingRestOuts(m_PARAMcurlingOutRestName);
120  curlingRestOuts.registerInDataStore(m_PARAMcurlingOutRestName, DataStore::c_ErrorIfAlreadyRegistered);
121  curlingRestOuts.registerRelationTo(spacePointTCs);
122 
123  // have to do this here, because in event() I do not want to check every time if this string is empty or not and act accordingly. If I register this with an empty string here, I can use it with an empty string in event() and only store stuff into it, when it is actually named with a non-empty string
124  StoreArray<SpacePointTrackCand> curlingCompletes(m_PARAMcompleteCurlerName);
125  curlingCompletes.registerInDataStore(m_PARAMcompleteCurlerName, DataStore::c_ErrorIfAlreadyRegistered);
126  curlingCompletes.registerRelationTo(spacePointTCs);
127 
128  if (!m_PARAMcompleteCurlerName.empty()) {
129  m_saveCompleteCurler = true;
130  B2DEBUG(1, "You put in " << m_PARAMcompleteCurlerName <<
131  " as collection name for complete curling TrackCands. Complete curling TrackCands will hence be stored.");
132  } else {
133  B2DEBUG(1,
134  "You did not put in any under which complete curling TrackCands should be stored, hence curling TrackCands will only be stored in parts.");
135  m_saveCompleteCurler = false;
136  }
137 
138  // check value for nTrackStubs and reset if necessary
139  if (m_PARAMnTrackStubs < 0) {
140  B2WARNING("CurlingTrackCandSplitter::initialize> Value of nTrackStubs is below 0: nTrackStubs = " << m_PARAMnTrackStubs <<
141  ". Resetting this value to 0 now! This means that all parts of curling TrackCands will be stored.");
142  m_PARAMnTrackStubs = 0;
143  } else { B2DEBUG(1, "Entered value for nTrackStubs = " << m_PARAMnTrackStubs); }
144 
145  B2DEBUG(1, "Entered Value for splitCurlers: " << m_PARAMsplitCurlers);
146 
147  if (m_PARAMsetOrigin.size() != 3) {
148  B2WARNING("CurlingTrackCandSplitter::initialize: Provided origin is not a 3D point! Please provide 3 values (x,y,z). Rejecting user input and setting origin to (0,0,0) for now!");
149  m_PARAMsetOrigin.clear();
150  m_PARAMsetOrigin.assign(3, 0);
151  }
152  m_origin.SetXYZ(m_PARAMsetOrigin.at(0), m_PARAMsetOrigin.at(1), m_PARAMsetOrigin.at(2));
153  B2DEBUG(10, "Set origin to (x,y,z): (" << m_origin.X() << "," << m_origin.Y() << "," << m_origin.Z() << ")");
154 
155  if (m_PARAMpositionAnalysis) {
156  // check if there are two entries and if the second value is either UPDATE or RECREATE
157  if (m_PARAMrootFileName.size() != 2 || (m_PARAMrootFileName[1] != "UPDATE" && m_PARAMrootFileName[1] != "RECREATE")) {
158  string output;
159  for (string entry : m_PARAMrootFileName) { output += "'" + entry + "' "; }
160  B2FATAL("CurlingTrackCandSplitter::initialize() : rootFileName is set wrong: entries are: " << output);
161  }
162  // create ROOT file
163  m_PARAMrootFileName[0] += ".root";
164  m_rootFilePtr = new TFile(m_PARAMrootFileName[0].c_str(), m_PARAMrootFileName[1].c_str());
165  m_treePtr = new TTree("m_treePtr", "aTree");
166 
167  // link everything to the according variables
168  for (int layer = 0; layer < c_nPlanes; ++layer) {
169  string layerString = (boost::format("%1%") % (layer +
170  1)).str(); // layer numbering starts at 1 this way (plus cppcheck complains about division by zero otherwise)
171 
172  string name = "SpacePointXGlobal_" + layerString;
173  m_treePtr->Branch(name.c_str(), &m_rootSpacePointXGlobals.at(layer));
174  name = "SpacePointYGlobal_" + layerString;
175  m_treePtr->Branch(name.c_str(), &m_rootSpacePointYGlobals.at(layer));
176  name = "SpacePointZGlobal_" + layerString;
177  m_treePtr->Branch(name.c_str(), &m_rootSpacePointZGlobals.at(layer));
178 
179  name = "SpacePointULocal_" + layerString;
180  m_treePtr->Branch(name.c_str(), &m_rootSpacePointULocals.at(layer));
181  name = "SpacePointVlocal_" + layerString;
182  m_treePtr->Branch(name.c_str(), &m_rootSpacePointVLocals.at(layer));
183 
184  name = "TrueHitXGlobal_" + layerString;
185  m_treePtr->Branch(name.c_str(), &m_rootTrueHitXGlobals.at(layer));
186  name = "TrueHitYGlobal_" + layerString;
187  m_treePtr->Branch(name.c_str(), &m_rootTrueHitYGlobals.at(layer));
188  name = "TrueHitZGlobal_" + layerString;
189  m_treePtr->Branch(name.c_str(), &m_rootTrueHitZGlobals.at(layer));
190 
191  name = "TrueHitULocal_" + layerString;
192  m_treePtr->Branch(name.c_str(), &m_rootTrueHitULocals.at(layer));
193  name = "TrueHitVLocal_" + layerString;
194  m_treePtr->Branch(name.c_str(), &m_rootTrueHitVLocals.at(layer));
195 
196  name = "PosResidualsXGlobal_" + layerString;
197  m_treePtr->Branch(name.c_str(), &m_rootPosResidueXGlobal.at(layer));
198  name = "PosResidualsYGlobal_" + layerString;
199  m_treePtr->Branch(name.c_str(), &m_rootPosResidueYGlobal.at(layer));
200  name = "PosResidualsZGlobal_" + layerString;
201  m_treePtr->Branch(name.c_str(), &m_rootPosResidueZGlobal.at(layer));
202 
203  name = "PosResidualsULocal_" + layerString;
204  m_treePtr->Branch(name.c_str(), &m_rootPosResidueULocal.at(layer));
205  name = "PosResidualsVLocal_" + layerString;
206  m_treePtr->Branch(name.c_str(), &m_rootPosResidueVLocal.at(layer));
207 
208 
209  name = "LocalPositionResiduals_" + layerString;
210  m_treePtr->Branch(name.c_str(), &m_rootLocalPosResiduals.at(layer));
211  name = "GlobalPositionResiduals_" + layerString;
212  m_treePtr->Branch(name.c_str(), &m_rootGlobalPosResiduals.at(layer));
213 
214  name = "MisMatchPosDistance_" + layerString;
215  m_treePtr->Branch(name.c_str(), &m_rootMisMatchPosDistance.at(layer));
216  name = "MisMatchPosX_" + layerString;
217  m_treePtr->Branch(name.c_str(), &m_rootMisMatchPosX.at(layer));
218  name = "MisMatchPosY_" + layerString;
219  m_treePtr->Branch(name.c_str(), &m_rootMisMatchPosY.at(layer));
220  name = "MisMatchPosZ_" + layerString;
221  m_treePtr->Branch(name.c_str(), &m_rootMisMatchPosZ.at(layer));
222 
223  name = "MisMatchPosU_" + layerString;
224  m_treePtr->Branch(name.c_str(), &m_rootMisMatchPosU.at(layer));
225  name = "MisMatchPosV_" + layerString;
226  m_treePtr->Branch(name.c_str(), &m_rootMisMatchPosV.at(layer));
227 
228  name = "MisMatchMomX_" + layerString;
229  m_treePtr->Branch(name.c_str(), &m_rootMisMatchMomX.at(layer));
230  name = "MisMatchMomY_" + layerString;
231  m_treePtr->Branch(name.c_str(), &m_rootMisMatchMomY.at(layer));
232  name = "MisMatchMomZ_" + layerString;
233  m_treePtr->Branch(name.c_str(), &m_rootMisMatchMomZ.at(layer));
234  }
235  } else {
236  m_rootFilePtr = nullptr;
237  m_treePtr = nullptr;
238  }
239 }
240 
241 // =================================================== EVENT ============================================================
242 void CurlingTrackCandSplitterModule::event()
243 {
244  StoreObjPtr<EventMetaData> eventMetaDataPtr("EventMetaData", DataStore::c_Event);
245  const int eventCounter = eventMetaDataPtr->getEvent();
246  B2DEBUG(10, "CurlingTrackCandSplitter::event(). -------------- Processing event " << eventCounter << " ----------------");
247 
248  // StoreArrays that will be used for storing
249  StoreArray<SpacePointTrackCand> outgoingFirstTCs(m_PARAMcurlingOutFirstName);
250  StoreArray<SpacePointTrackCand> ingoingAllTCs(m_PARAMcurlingAllInName);
251  StoreArray<SpacePointTrackCand> outgoingRestTCs(m_PARAMcurlingOutRestName);
252  StoreArray<SpacePointTrackCand> completeCurlingTCs(m_PARAMcompleteCurlerName);
253 
254 
255  StoreArray<SpacePointTrackCand> spacePointTCs(m_PARAMsptcName);
256  int nTCs = spacePointTCs.getEntries();
257 
258  B2DEBUG(15, "Found " << nTCs << " SpacePointTrackCands in StoreArray " << spacePointTCs.getName() << " for this event");
259 
260  RootVariables rootVariables;
261 
262  for (int iTC = 0; iTC < nTCs; ++iTC) {
263  SpacePointTrackCand* spacePointTC = spacePointTCs[iTC];
264  m_spacePointTCCtr++;
265 
266  B2DEBUG(15, "=========================== Processing SpacePointTrackCand " << iTC << " ===============================");
267  try {
268  const std::vector<int> splittingIndices = checkTrackCandForCurling(*spacePointTC, rootVariables);
269 
270  if (splittingIndices.empty()) {
271  B2DEBUG(15, "This SpacePointTrackCand shows no curling behaviour and will be added to collection: " << m_PARAMcurlingOutFirstName);
272  spacePointTC->setTrackStubIndex(0); // set TrackStubIndex to 0 (indicates, that this TrackCandidate shows no curling behaviour)
273  // add this spacePoint to the StoreArray with the first outgoing parts since the whole TC is outgoing
274  SpacePointTrackCand* newSPTC = outgoingFirstTCs.appendNew(*spacePointTC);
275  newSPTC->addRelationTo(spacePointTC);
276  m_NoCurlingTCsCtr++;
277  } else {
278  B2DEBUG(15, "This SpacePointTrackCand shows curling behaviour");
279  if (!m_PARAMsplitCurlers) {
280  B2DEBUG(15, "This SpacePointTrackCand could be split into " << splittingIndices.size() + 1 <<
281  " but will not, because splitCurlers is set to false");
282  continue; // should jump to enclosing for-loop!!! (process the next SpacePointTrackCand)
283  }
284  m_curlingTCCtr++;
285  // get the TrackCand Stubs
286  std::vector<SpacePointTrackCand> trackStubs = splitCurlingTrackCand(*spacePointTC, m_PARAMnTrackStubs, splittingIndices);
287 
288  // add the stubs to the appropriate StoreArray
289  for (SpacePointTrackCand trackStub : trackStubs) {
290  m_createdTrackStubsCtr++;
291  if (m_saveCompleteCurler) {
292  SpacePointTrackCand* newSPTC = completeCurlingTCs.appendNew(trackStub);
293  newSPTC->addRelationTo(spacePointTC);
294  B2DEBUG(25, "Added SpacePointTrackCand " << newSPTC->getArrayIndex() << " to StoreArray " << newSPTC->getArrayName());
295  }
296  if (!trackStub.isOutgoing()) {
297  SpacePointTrackCand* newSPTC = ingoingAllTCs.appendNew(trackStub);
298  newSPTC->addRelationTo(spacePointTC);
299  B2DEBUG(25, "Added SpacePointTrackCand " << newSPTC->getArrayIndex() << " to StoreArray " << newSPTC->getArrayName());
300  } else { // if not ingoing differentiate between first part and all of the rest
301  if (trackStub.getTrackStubIndex() > 1) {
302  SpacePointTrackCand* newSPTC = outgoingRestTCs.appendNew(trackStub);
303  newSPTC->addRelationTo(spacePointTC);
304  B2DEBUG(25, "Added SpacePointTrackCand " << newSPTC->getArrayIndex() << " to StoreArray " << newSPTC->getArrayName());
305  } else {
306  SpacePointTrackCand* newSPTC = outgoingFirstTCs.appendNew(trackStub);
307  newSPTC->addRelationTo(spacePointTC);
308  B2DEBUG(25, "Added SpacePointTrackCand " << newSPTC->getArrayIndex() << " to StoreArray " << newSPTC->getArrayName());
309  }
310  }
311  }
312  }
313  } catch (FoundNoTrueHit& anE) {
314  B2WARNING("Caught an exception during checking for curling behaviour: " << anE.what() <<
315  " This TrackCandidate cannot be checked for curling behaviour");
316  m_noDecisionPossibleCtr++;
317  } catch (FoundNoCluster& anE) {
318  B2WARNING("Caught an exception during checking for curling behaviour: " << anE.what() <<
319  " This TrackCandidate cannot be checked for curling behaviour");
320  m_noDecisionPossibleCtr++;
321  } catch (TrueHitsNotMatching& anE) {
322  B2WARNING("Caught an exception during checking for curling behaviour: " << anE.what() <<
323  " This TrackCandidate cannot be checked for curling behaviour");
324  m_noDecisionPossibleCtr++;
325  } catch (SpacePointTrackCand::UnsupportedDetType& anE) {
326  B2WARNING("Caught an exception during checking for curling behaviour: " << anE.what() <<
327  " This TrackCandidate cannot be checked for curling behaviour");
328  m_noDecisionPossibleCtr++;
329  }
330  }
331  // only write to root once per event
332  if (m_PARAMpositionAnalysis) { writeToRoot(rootVariables); }
333 }
334 
335 // =================================================== TERMINATE ========================================================
336 void CurlingTrackCandSplitterModule::terminate()
337 {
338  B2INFO("CurlingTrackCandSplitter::terminate(): checked " << m_spacePointTCCtr << " SpacePointTrackCands for curling behaviour. " <<
339  m_curlingTCCtr << " of them were curling and " << m_createdTrackStubsCtr << " TrackStubs were created. " << m_NoCurlingTCsCtr <<
340  " SPTCs were not curling and were merely copied into StoreArray " << m_PARAMcurlingOutFirstName << ". In " <<
341  m_noDecisionPossibleCtr << " cases no decision could be made. There were " << m_NoSingleTrueHitCtr <<
342  " SpacePoints that were related to more than one TrueHit");
343  // do ROOT file stuff
344  if (m_treePtr != nullptr) {
345  m_rootFilePtr->cd(); //important! without this the famework root I/O (SimpleOutput etc) could mix with the root I/O of this module
346  m_treePtr->Write();
347  m_rootFilePtr->Close();
348  }
349 }
350 
351 // ============================================== CHECK FOR CURLING ======================================================
352 const std::vector<int> CurlingTrackCandSplitterModule::checkTrackCandForCurling(const Belle2::SpacePointTrackCand& SPTrackCand,
353  RootVariables& rootVariables)
354 {
355  const std::vector<const Belle2::SpacePoint*>& tcSpacePoints = SPTrackCand.getHits();
356  unsigned int nHits = SPTrackCand.getNHits();
357 
358  B2DEBUG(70, "SpacePointTrackCand contains " << nHits << " SpacePoints");
359 
360  std::vector<int> returnVector; // fill this vector with indices, if no indices can be found, leave it empty
361 
362  std::pair<bool, bool>
363  directions; // only store the last two directions to decide if it has changed or not. .first is always last hit, .second is present hit.
364  directions.first =
365  true; // assume that the track points outwards from the interaction point at first. NOTE: this assumption is not dangerous here (it is in other places because it does not have to be the case). The information on the direction of flight for the first hit is stored in the returnVector itself. If the first entry is 0 (this means that the trackCand first 'pointed' towards the interaction point)
366 
367  for (unsigned int iHit = 0; iHit < nHits; ++iHit) {
368  const SpacePoint* spacePoint = tcSpacePoints[iHit];
369  auto detType = spacePoint->getType();
370 
371  B2DEBUG(100, "Now checking SpacePoint " << iHit << " in SPTC. This SpacePoint has Index " << spacePoint->getArrayIndex() <<
372  " in StoreArray " << spacePoint->getArrayName());
373 
374  // get global position and momentum for every spacePoint in the SpacePointTrackCand
375  std::pair<B2Vector3<double>, B2Vector3<double> > hitGlobalPosMom;
376 
377  if (detType == VXD::SensorInfoBase::PXD) {
378  // first get PXDCluster, from that get TrueHit
379  PXDCluster* pxdCluster =
380  spacePoint->getRelatedTo<PXDCluster>("ALL"); // COULDDO: search only certain Cluster Arrays -> get name somehow
381  // CAUTION: only looking for one TrueHit here, but there could actually be more of them 'molded' into one Cluster
382  PXDTrueHit* pxdTrueHit =
383  pxdCluster->getRelatedTo<PXDTrueHit>("ALL"); // COULDDO: search only certain PXDTrueHit arrays -> new parameter for module
384 
385  if (pxdTrueHit == nullptr) {
386  B2DEBUG(1, "Found no PXDTrueHit for PXDCluster " << pxdCluster->getArrayIndex() << " from Array " << pxdCluster->getArrayName() <<
387  ". This PXDCluster is related with SpacePoint " << spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName());
388  throw FoundNoTrueHit();
389  }
390 
391  B2DEBUG(100, "Found PXDCluster " << pxdCluster->getArrayIndex() << " and " << " PXDTrueHit " << pxdTrueHit->getArrayIndex() <<
392  " from StoreArray " << pxdTrueHit->getArrayName() << " related to this SpacePoint");
393 
394  B2DEBUG(100, "Now getting global position and momentum for PXDCluster " << pxdCluster->getArrayIndex() << " from Array " <<
395  pxdCluster->getArrayName());
396  hitGlobalPosMom = getGlobalPositionAndMomentum(pxdTrueHit);
397 
398  // if position analysis is set to true, print to root file
399  if (m_PARAMpositionAnalysis) { getValuesForRoot(spacePoint, pxdTrueHit, rootVariables); }
400 
401  } else if (detType == VXD::SensorInfoBase::SVD) {
402  // get all related SVDClusters and do some sanity checks, before getting the SVDTrueHits and then using them to get global position and momentum
403  RelationVector<SVDCluster> svdClusters =
404  spacePoint->getRelationsTo<SVDCluster>("ALL"); // COULDDO: search only certain Cluster Arrays -> get name somehow (addidional parameter?)
405  if (svdClusters.size() == 0) {
406  B2WARNING("Found no related clusters for SpacePoint " << spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName()
407  << ". With no Cluster no information if a track is curling or not can be obtained");
408  throw FoundNoCluster(); // this should also never happen, as the vice versa way is used to get to the SpacePoints in the first place (in the GFTC2SPTCConverterModule e.g.)
409  } else {
410  // collect the TrueHits, if there is more than one compare them, to see if both Clusters point to the same TrueHit
411  // WARNING there can be more! more than one TrueHit can be 'hidden' in one Cluster!!!
412  // TODO: look at this again, this seems not to work properly at the moment!!!
413  std::vector<const SVDTrueHit*> svdTrueHits;
414  for (const SVDCluster& aCluster : svdClusters) {
415  // CAUTION: there can be more than one TrueHit for a given Cluster!!!
416  RelationVector<SVDTrueHit> relTrueHits =
417  aCluster.getRelationsTo<SVDTrueHit>("ALL"); // COULDDO: search only certain SVDTrueHit arrays -> new parameter for module
418  if (relTrueHits.size() == 0) {
419  B2DEBUG(1, "Found no SVDTrueHit for SVDCluster " << aCluster.getArrayIndex() << " from Array " << aCluster.getArrayName() <<
420  ". This SVDCluster is related with SpacePoint " << spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName());
421  throw FoundNoTrueHit();
422  }
423 
424  B2DEBUG(100, "Found " << relTrueHits.size() << " TrueHits for SVDCluster " << aCluster.getArrayIndex() << " from Array " <<
425  aCluster.getArrayName());
426  for (unsigned int i = 0; i < relTrueHits.size(); ++i) { svdTrueHits.push_back(relTrueHits[i]); }
427  }
428 
429  // if there is only one cluster related to the SpacePoint simply check if one (or more TrueHits are present). Additionally checking the size for svdTrueHits again is not necessary here, because if there was only one Cluster and no TrueHits were found this part is never reached!
430  // WARNING: It is not guaranteed that this actually leads to a valid relation between SpacePoint and TrueHit!!
431  // TODO: continuing with next SpacePoint skips check for curling behavior!!! FIX this!!!
432 // if (svdClusters.size() == 1) {
433 // stringstream inds;
434 // for (const SVDTrueHit * trueHit : svdTrueHits) { inds << trueHit->getArrayIndex() << ", "; }
435 // B2DEBUG(150, "Found only one Cluster related to SpacePoint " << spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName() << ". To this Cluster " << svdTrueHits.size() << " related TrueHits were found. Indices: " << inds.str());
436 // m_NoSingleTrueHitCtr++;
437 // continue; // start over with next SpacePoint
438 // }
439 
440  // if there is at least one TrueHit in the vector check how many unique TrueHits there are
441  if (svdTrueHits.size() >= 1) {
442  B2DEBUG(150, "Found " << svdTrueHits.size() << " SVDTrueHits related to Clusters related to SpacePoint " <<
443  spacePoint->getArrayIndex() << " from Array " << spacePoint->getArrayName() << ". Now checking if they are compatible");
444 
445  // sort & unique to find the unique entries of the relation vector
446  std::sort(svdTrueHits.begin(), svdTrueHits.end());
447  unsigned int oldSize = svdTrueHits.size();
448  auto newEnd = std::unique(svdTrueHits.begin(), svdTrueHits.end());
449  svdTrueHits.resize(std::distance(svdTrueHits.begin(), newEnd));
450 
451  // If there is no overlapping TrueHit get all indices of the TrueHits and print a warning. If position Analysis is enabled calculate the position differences of the TrueHits. In the end throw an Exception (only if the TrueHits belong to a non-single Cluster SVD SpacePoint)
452  if (svdTrueHits.size() == oldSize) {
453  stringstream trueHitInds;
454  for (const SVDTrueHit* trueHit : svdTrueHits) { trueHitInds << trueHit->getArrayIndex() << ", "; }
455  B2DEBUG(1, "There is no overlapping TrueHit for SpacePoint " << spacePoint->getArrayIndex() << " from Array " <<
456  spacePoint->getArrayName() << ". The Indices of the TrueHits are: " << trueHitInds.str());
457 
458  // Only do these calculations if output to root is is enabled
459  if (m_PARAMpositionAnalysis) {
460  std::vector<B2Vector3<double> > globalPositions;
461  std::vector<B2Vector3<double> > globalMomenta;
462  // collect all values
463  for (unsigned int i = 0; i < svdTrueHits.size(); ++i) {
464  auto posMom = getGlobalPositionAndMomentum(svdTrueHits[i]);
465  globalPositions.push_back(posMom.first);
466  globalMomenta.push_back(posMom.second);
467  }
468  // WARNING: getting only layer number of first TrueHit in vector here. Although this should not change, this is never actually checked!!
469  int layer = svdTrueHits[0]->getSensorID().getLayerNumber() - 1; // layer numbering starts at 1, indexing of array at 0
470  // do the calculations (starting from one because of comparison of two elements in each run through loop)
471  for (unsigned int i = 1; i < globalPositions.size(); ++i) {
472  rootVariables.MisMatchPosResiduals.at(layer).push_back((globalPositions[i] - globalPositions[i - 1]).Mag());
473 
474  rootVariables.MisMatchPosX.at(layer).push_back((globalPositions[i] - globalPositions[i - 1]).X());
475  rootVariables.MisMatchPosY.at(layer).push_back((globalPositions[i] - globalPositions[i - 1]).Y());
476  rootVariables.MisMatchPosZ.at(layer).push_back((globalPositions[i] - globalPositions[i - 1]).Z());
477 
478  rootVariables.MisMatchPosU.at(layer).push_back((svdTrueHits[i]->getU() - svdTrueHits[i - 1]->getU()));
479  rootVariables.MisMatchPosV.at(layer).push_back((svdTrueHits[i]->getV() - svdTrueHits[i - 1]->getV()));
480 
481  B2Vector3<double> momDiff = globalMomenta[i] - globalMomenta[i - 1];
482  rootVariables.MisMatchMomX.at(layer).push_back(momDiff.X());
483  rootVariables.MisMatchMomY.at(layer).push_back(momDiff.Y());
484  rootVariables.MisMatchMomZ.at(layer).push_back(momDiff.Z());
485  }
486  }
487  // if the TrueHits are related from a singleCluster SVD SpacePoint (i.e. more than one TrueHits are molded into one Cluster) do not throw this exception but continue with the curling checking
488  if (svdClusters.size() > 1) { TrueHitsNotMatching(); }
489  }
490  }
491 
492  // if there is more than one TrueHit remaining for one SpacePoint increase the counter
493  if (svdTrueHits.size() > 1) {
494  m_NoSingleTrueHitCtr++;
495  }
496 
497  // WARNING if there are more than one matching TrueHits only the first TrueHit is used for comparison and for position analysis
498  B2DEBUG(100, "Now getting global position and momentum for SVDCluster " << svdClusters[0]->getArrayIndex() << " from Array " <<
499  svdClusters[0]->getArrayName() << " via SVDTrueHit " << svdTrueHits[0]->getArrayIndex() << " from StoreArray " <<
500  svdTrueHits[0]->getArrayName());
501  hitGlobalPosMom = getGlobalPositionAndMomentum(svdTrueHits[0]);
502 
503  // if position analysis is set to true, print to root file
504  // only do so if there is only one TrueHit (this is only for the moment!!!) OR if the switch is set to do so even if there is more than one TrueHit
505  // TODO: Decide how to handle such cases where more than one TrueHit is left and implement accordingly
506  if (m_PARAMpositionAnalysis && (svdTrueHits.size() == 1 || m_PARAMuseNonSingleTHinPA)) { getValuesForRoot(spacePoint, svdTrueHits[0], rootVariables); }
507  }
508  } else { // this should never be reached, because it should be caught in the creation of the SpacePointTrackCand which is passed to this function!
509  throw SpacePointTrackCand::UnsupportedDetType();
510  }
511 
512  // get the direction of flight for the present SpacePoint
513  directions.second = getDirectionOfFlight(hitGlobalPosMom, m_origin);
514 
515  // check if the directions have changed since the last hit, if so, add the number of the SpacePoint (inside the SpacePointTrackCand) to the returnVector
516  if (directions.first != directions.second) {
517  B2DEBUG(75, "The direction of flight has changed for SpacePoint " << iHit <<
518  " in SpacePointTrackCand. The StoreArray index of this SpacePoint is " << spacePoint->getArrayIndex() << " in " <<
519  spacePoint->getArrayName());
520  returnVector.push_back(iHit);
521  }
522  // assign old value to .first, for next comparison
523  directions.first = directions.second;
524  }
525  return returnVector;
526 }
527 
528 // ======================================= GET GLOBAL POSITION AND MOMENTUM ============================================================
529 
531 template<class TrueHit>
532 std::pair<const Belle2::B2Vector3<double>, const Belle2::B2Vector3<double> >
533 CurlingTrackCandSplitterModule::getGlobalPositionAndMomentum(TrueHit* aTrueHit)
534 {
535  // get sensor stuff (needed for pointToGlobal)
536  VxdID aVxdId = aTrueHit->getSensorID();
537 
538  B2DEBUG(100, "Getting global position and momentum vectors for TrueHit " << aTrueHit->getArrayIndex() << " from Array " <<
539  aTrueHit->getArrayName() << ". This hit has VxdID " << aVxdId);
540 
541  const VXD::GeoCache& geometry = VXD::GeoCache::getInstance();
542  const VXD::SensorInfoBase& sensorInfoBase = geometry.getSensorInfo(aVxdId);
543 
544  // get position
545  B2Vector3<double> hitLocal = B2Vector3<double>(aTrueHit->getU(), aTrueHit->getV(), 0);
546  B2Vector3<double> hitGlobal = sensorInfoBase.pointToGlobal(
547  hitLocal, true); // should work like this, since local coordinates are only 2D
548  B2DEBUG(100, "Local position of hit is (" << hitLocal.X() << "," << hitLocal.Y() << "," << hitLocal.Z() <<
549  "), Global position of hit is (" << hitGlobal.X() << "," << hitGlobal.Y() << "," << hitGlobal.Z() << ")");
550 
551  // get momentum
552  B2Vector3<double> pGlobal = sensorInfoBase.vectorToGlobal(aTrueHit->getMomentum(), true);
553  B2DEBUG(100, "Global momentum of hit is (" << pGlobal.X() << "," << pGlobal.Y() << "," << pGlobal.Z() << ")");
554 
555  return std::make_pair(hitGlobal, pGlobal);
556 }
557 
558 // ======================================= GET DIRECTION OF FLIGHT ======================================================================
559 bool CurlingTrackCandSplitterModule::getDirectionOfFlight(const
560  std::pair<const B2Vector3<double>, const B2Vector3<double>>& hitPosAndMom,
561  const B2Vector3<double>& origin)
562 {
563  B2Vector3<double> originToHit = hitPosAndMom.first - origin;
564  B2Vector3<double> momentumAtHit = hitPosAndMom.second + originToHit;
565 
566 
567  B2DEBUG(100, "Position of hit relative to origin is (" << originToHit.X() << "," << originToHit.Y() << "," << originToHit.Z() <<
568  "). Momentum relative to hit (relative to origin) (" << momentumAtHit.X() << "," << momentumAtHit.Y() << "," << momentumAtHit.Z() <<
569  ")");
570 
571  // cylindrical coordinates (?) -> use B2Vector3.Perp() to get the radial component of a given vector
572  // for spherical coordinates -> use B2Vector3.Mag() for the same purposes
573  double hitRadComp = originToHit.Perp(); // radial component of hit coordinates
574  double hitMomRadComp =
575  momentumAtHit.Perp(); // radial component of the tip of the momentum vector, when its origin would be the hit position (as it is only the direction of the momentum that matters here, units are completely ignored -> COULDDO: use the unit() method from B2Vector3
576 
577  B2DEBUG(250, " radial component of hit coordinates: " << hitRadComp <<
578  ", radial component of tip of momentum vector with its origin set to hit position: " << hitMomRadComp);
579 
580  if (hitMomRadComp < hitRadComp) {
581  B2DEBUG(100, "Direction of flight is inwards for this hit");
582  return false;
583  } else {
584  B2DEBUG(100, "Direction of flight is outwards for this hit");
585  return true;
586  }
587 }
588 
589 // ================================================ SPLIT CURLING TRACK CAND =========================================================
590 const std::vector<Belle2::SpacePointTrackCand>
591 CurlingTrackCandSplitterModule::splitCurlingTrackCand(const Belle2::SpacePointTrackCand& SPTrackCand, int NTracklets,
592  const std::vector<int>& splitIndices)
593 {
594  std::vector<SpacePointTrackCand> spacePointTCs;
595 
596  std::vector<std::pair<int, int> >
597  rangeIndices; // store pairs of Indices indicating the first and the last index of a TrackStub inside a SpacePointTrackCand
598 
599  int firstIndex = 0; // first 'first index' is 0
600  for (int index : splitIndices) {
601  rangeIndices.push_back({firstIndex, index});
602  firstIndex = index + 1; // next first index is last index + 1
603  }
604  rangeIndices.push_back({firstIndex, SPTrackCand.getNHits() - 1}); // the last TrackStub contains all hits from the last splitIndex to the end of the SpacePointTrackCand
605 
606  // NOTE: if the first index of splitIndices is zero, this means that the direction of flight for the first SpacePoint (i.e. TrueHit) of the TrackCand was not outgoing.
607  // WARNING: This is only true as long as this behaviour is not changed in checkTrackCandForCurling(...), as there the assumption is made that the direction of flight of the first hit is outgoing, and only if it is not 0 is the first entry of the returnVector of the latter.
608  bool outgoing = splitIndices[0] != 0;
609 
610  // return NTracklets (user defined) at most, or if NTracklets is 0, return all
611  // if NTracklets is 0 set it to the appropriate size for the following for-loop
612  if (NTracklets < 1) { NTracklets = rangeIndices.size(); }
613  for (unsigned iTr = 0; iTr < rangeIndices.size() && iTr < uint(NTracklets); ++iTr) {
614 
615  int lastInd = rangeIndices[iTr].second;
616  int firstInd = rangeIndices[iTr].first;
617 
618  B2DEBUG(75, "Creating Track Stub " << iTr << " of " << splitIndices.size() <<
619  " possible Track Stub for this SpacePointTrackCand. The indices for this Tracklet are (first,last): (" << firstInd << "," << lastInd
620  << "). This SpacePointTrackCand contains " << SPTrackCand.getNHits() << " SpacePoints in total.");
621 
622  // encapsulate these functions into a try-clause since both throw
623  try {
624  const std::vector<const SpacePoint*> trackletSpacePoints = SPTrackCand.getHitsInRange(firstInd, lastInd);
625  const std::vector<double> trackletSortingParams = SPTrackCand.getSortingParametersInRange(firstInd, lastInd);
626 
627  SpacePointTrackCand newSPTrackCand = SpacePointTrackCand(trackletSpacePoints, SPTrackCand.getPdgCode(), SPTrackCand.getChargeSeed(),
628  SPTrackCand.getMcTrackID());
629  newSPTrackCand.setSortingParameters(trackletSortingParams);
630 
631  // TODO: set state seed and cov seed for all but the first tracklets (first is just the seed of the original TrackCand)
632  if (iTr < 1) {
633  newSPTrackCand.set6DSeed(SPTrackCand.getStateSeed());
634  newSPTrackCand.setCovSeed(SPTrackCand.getCovSeed());
635  }
636 
637  // set direction of flight and flip it for the next track stub (track is split where the direction of flight changes so this SHOULD not introduce any errors)
638  newSPTrackCand.setFlightDirection(outgoing);
639  outgoing = !outgoing;
640 
641  // if the TrackCandidate curls this index starts at 1, if it is a curling TrackCand
642  newSPTrackCand.setTrackStubIndex(iTr + 1);
643 
644  spacePointTCs.push_back(newSPTrackCand);
645  } catch (SpacePointTrackCand::SPTCIndexOutOfBounds& anE) {
646  B2WARNING("Caught an exception while trying to split SpacePointTrackCands: " << anE.what() <<
647  " This SPTC will be skipped from splitting!");
648  }
649  }
650 
651  return spacePointTCs;
652 }
653 
654 // ======================================================= GET ROOT VALUES =========================================================
655 template <class TrueHit>
656 void CurlingTrackCandSplitterModule::getValuesForRoot(const Belle2::SpacePoint* spacePoint, const TrueHit* trueHit,
657  RootVariables& rootVariables)
658 {
659  B2DEBUG(100, "Getting positions (for ROOT output) of SpacePoint " << spacePoint->getArrayIndex() << " from Array " <<
660  spacePoint->getArrayName() << " and TrueHit " << trueHit->getArrayIndex() << " from Array " << trueHit->getArrayName());
661 
662  // get VxdIDs of spacePoint and trueHit (and their according layer numbers for storing the information in the appropriate arrays)
663  VxdID spacePointVxdId = spacePoint->getVxdID();
664  VxdID trueHitVxdId = trueHit->getSensorID();
665 
666  // get positions from SpacePoint
667  const B2Vector3<double>& spacePointGlobal =
668  spacePoint->getPosition(); // COULDDO: uneccesary, spacePoint->X(), etc. returns the same information!
669 // std::pair<double, double> spacePointUV = getUV(spacePoint);
670  TaggedUVPos spacePointUV = getUV(spacePoint);
671  const B2Vector3<double> spacePointLocal = B2Vector3<double>(spacePointUV.m_U, spacePointUV.m_V, 0);
672 
673  // get local position from TrueHit
674  const B2Vector3<double> trueHitLocal = B2Vector3<double>(trueHit->getU(), trueHit->getV(), 0);
675 
676  // get sensor Info for global position of TrueHit
677  const VXD::GeoCache& geometry = VXD::GeoCache::getInstance();
678  const VXD::SensorInfoBase& sensorInfoBase = geometry.getSensorInfo(trueHitVxdId);
679 
680  // get global position from TrueHit
681  const B2Vector3<double> trueHitGlobal = sensorInfoBase.pointToGlobal(trueHitLocal, true);
682 
683  // Layer numbering starts at 1 not at 0 so deduce layer by one to access array
684  int spLayer = spacePointVxdId.getLayerNumber() - 1;
685  int thLayer = trueHitVxdId.getLayerNumber() - 1;
686 
687  if (spLayer != thLayer) {
688  B2FATAL("Layer numbers of TrueHit and SpacePoint do not match!"); // this should never happen -> FATAL if it does, because something has gone amiss then
689  }
690 
691  // all positions collected, but only write the values to ROOT that have been set appropriately!
692  bool singleCluster = true; // for debug output
693  if (spacePointUV.m_setU) {
694  rootVariables.SpacePointULocal.at(spLayer).push_back(spacePointUV.m_U);
695  rootVariables.TrueHitULocal.at(thLayer).push_back(trueHit->getU());
696  rootVariables.PosResidueULocal.at(spLayer).push_back((spacePointUV.m_U - trueHit->getU()));
697  }
698  if (spacePointUV.m_setV) {
699  rootVariables.SpacePointVLocal.at(spLayer).push_back(spacePointUV.m_V);
700  rootVariables.TrueHitVLocal.at(thLayer).push_back(trueHit->getV());
701  rootVariables.PosResidueVLocal.at(spLayer).push_back((spacePointUV.m_V - trueHit->getV()));
702  }
703  if (spacePointUV.m_setU && spacePointUV.m_setV) {
704  rootVariables.SpacePointXGlobal.at(spLayer).push_back(spacePoint->X());
705  rootVariables.SpacePointYGlobal.at(spLayer).push_back(spacePoint->Y());
706  rootVariables.SpacePointZGlobal.at(spLayer).push_back(spacePoint->Z());
707 
708  rootVariables.TrueHitXGlobal.at(thLayer).push_back(trueHitGlobal.X());
709  rootVariables.TrueHitYGlobal.at(thLayer).push_back(trueHitGlobal.Y());
710  rootVariables.TrueHitZGlobal.at(thLayer).push_back(trueHitGlobal.Z());
711 
712  rootVariables.PosResidueXGlobal.at(spLayer).push_back((spacePointGlobal - trueHitGlobal).X());
713  rootVariables.PosResidueYGlobal.at(spLayer).push_back((spacePointGlobal - trueHitGlobal).Y());
714  rootVariables.PosResidueZGlobal.at(spLayer).push_back((spacePointGlobal - trueHitGlobal).Z());
715 
716  rootVariables.PosResiduesGlobal.at(spLayer).push_back((spacePointGlobal - trueHitGlobal).Mag());
717  rootVariables.PosResiduesLocal.at(spLayer).push_back((spacePointLocal - trueHitLocal).Mag());
718 
719  singleCluster = false;
720  }
721 
722  B2DEBUG(200, "Global (x,y,z)/Local (U,V) positions of SpacePoint: (" << spacePointGlobal.X() << "," << spacePointGlobal.Y() << ","
723  << spacePointGlobal.Z() << ")/(" << spacePointLocal.X() << "," << spacePointLocal.Y() << "). This was a singleCluster SpacePoint: "
724  << singleCluster);
725 
726  B2DEBUG(200, "Global (x,y,z)/Local (U,V) positions of TrueHit: (" << trueHitGlobal.X() << "," << trueHitGlobal.Y() << "," <<
727  trueHitGlobal.Z() << ")/(" << trueHitLocal.X() << "," << trueHitLocal.Y() << ")");
728 
729  B2DEBUG(200, "This leads to position differences global/local: " << (spacePointGlobal - trueHitGlobal).Mag() << "/" <<
730  (spacePointLocal - trueHitLocal).Mag());
731 
732 }
733 
734 // =================================== WRITE TO ROOT ===============================================================================
735 void CurlingTrackCandSplitterModule::writeToRoot(const RootVariables& rootVariables)
736 {
737  m_rootGlobalPosResiduals = rootVariables.PosResiduesGlobal;
738  m_rootLocalPosResiduals = rootVariables.PosResiduesLocal;
739 
740  m_rootSpacePointULocals = rootVariables.SpacePointULocal;
741  m_rootSpacePointVLocals = rootVariables.SpacePointVLocal;
742 
743  m_rootSpacePointXGlobals = rootVariables.SpacePointXGlobal;
744  m_rootSpacePointYGlobals = rootVariables.SpacePointYGlobal;
745  m_rootSpacePointZGlobals = rootVariables.SpacePointZGlobal;
746 
747  m_rootTrueHitULocals = rootVariables.TrueHitULocal;
748  m_rootTrueHitVLocals = rootVariables.TrueHitVLocal;
749 
750  m_rootTrueHitXGlobals = rootVariables.TrueHitXGlobal;
751  m_rootTrueHitYGlobals = rootVariables.TrueHitYGlobal;
752  m_rootTrueHitZGlobals = rootVariables.TrueHitZGlobal;
753 
754  m_rootPosResidueXGlobal = rootVariables.PosResidueXGlobal;
755  m_rootPosResidueYGlobal = rootVariables.PosResidueYGlobal;
756  m_rootPosResidueZGlobal = rootVariables.PosResidueZGlobal;
757  m_rootPosResidueULocal = rootVariables.PosResidueULocal;
758  m_rootPosResidueVLocal = rootVariables.PosResidueVLocal;
759 
760  m_rootMisMatchPosDistance = rootVariables.MisMatchPosResiduals;
761  m_rootMisMatchPosX = rootVariables.MisMatchPosX;
762  m_rootMisMatchPosY = rootVariables.MisMatchPosY;
763  m_rootMisMatchPosZ = rootVariables.MisMatchPosZ;
764 
765  m_rootMisMatchPosU = rootVariables.MisMatchPosU;
766  m_rootMisMatchPosV = rootVariables.MisMatchPosV;
767 
768  m_rootMisMatchMomX = rootVariables.MisMatchMomX;
769  m_rootMisMatchMomY = rootVariables.MisMatchMomY;
770  m_rootMisMatchMomZ = rootVariables.MisMatchMomZ;
771 
772  m_treePtr->Fill();
773 }
774 
775 CurlingTrackCandSplitterModule::TaggedUVPos CurlingTrackCandSplitterModule::getUV(const Belle2::SpacePoint* spacePoint)
776 {
777  TaggedUVPos returnVals; // initialized to: both bools false, both doubles to 0.
778 
779  // 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)
780  double normU = spacePoint->getNormalizedLocalU();
781  double normV = spacePoint->getNormalizedLocalV();
782 
783  // convert normalized coordinates to local coordinates (those are already 'unwedged')
784  std::pair<double, double> localUV = SpacePoint::convertNormalizedToLocalCoordinates(std::make_pair(normU, normV),
785  spacePoint->getVxdID());
786  //double unwedgedU = SpacePoint::getUUnwedged(localUV, spacePoint->getVxdID());
787 
788  // set values for return: set both m_U and m_V regardless if they have actually been set in the SpacePoint and simply assign the m_setX tags, as they decide if a value is actually used
789  std::pair<bool, bool> setCoords = spacePoint->getIfClustersAssigned();
790 
791  returnVals.m_setU = setCoords.first;
792  returnVals.m_setV = setCoords.second;
793  returnVals.m_U = localUV.first;
794  returnVals.m_V = localUV.second;
795 
796  return returnVals;
797 }
798 
799 // ====================================================== INITIALIZE COUNTERS =======================================================
800 void CurlingTrackCandSplitterModule::initializeCounters()
801 {
802  m_createdTrackStubsCtr = 0;
803  m_curlingTCCtr = 0;
804  m_noDecisionPossibleCtr = 0;
805  m_spacePointTCCtr = 0;
806  m_NoSingleTrueHitCtr = 0;
807  m_NoCurlingTCsCtr = 0;
808 }
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
DataType Perp() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:185
Module for checking SpacePointTrackCandidates for curling behaviour and splitting them into Track Can...
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
Class PXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: PXDTrueHit.h:31
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.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:28
Class SVDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: SVDTrueHit.h:33
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
unsigned int getNHits() const
get the number of hits (space points) in the track candidate
int getPdgCode() const
get pdg code
void setCovSeed(const TMatrixDSym &cov)
set the covariance matrix seed
int getMcTrackID() const
get the MC Track ID
const std::vector< const Belle2::SpacePoint * > & getHits() const
get hits (space points) of track candidate
void setTrackStubIndex(int trackStubInd)
set TrackStub index
const TMatrixDSym & getCovSeed() const
get the covariance matrix seed (6D).
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...
const std::vector< double > getSortingParametersInRange(int firstIndex, int lastIndex) const
get the sorting parameters in range (indices of SpacePoints inside SpacePointTrackCand) including fir...
const TVectorD & getStateSeed() const
get state seed as 6D vector
double getChargeSeed() const
get charge
void setFlightDirection(bool direction)
set the direction of flight (true is outgoing, false is ingoing).
SpacePoint typically is build from 1 PXDCluster or 1-2 SVDClusters.
Definition: SpacePoint.h:42
const B2Vector3< double > & getPosition() const
return the position vector in global coordinates
Definition: SpacePoint.h:138
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
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
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.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
T * appendNew()
Construct a new T object at the end of the array.
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:95
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
Base class to provide Sensor Information for PXD and SVD.
TVector3 pointToGlobal(const TVector3 &local, bool reco=false) const
Convert a point from local to global coordinates.
TVector3 vectorToGlobal(const TVector3 &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
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
std::array< std::vector< double >, c_nPlanes > PosResidueVLocal
V-position (local) difference between TrueHit and SpacePoint (layerwise)
std::array< std::vector< double >, c_nPlanes > MisMatchMomZ
Difference of Momentum in Z-Direction for TrueHits that do not match but are related from one SpacePo...
std::array< std::vector< double >, c_nPlanes > MisMatchPosZ
Difference of Z-positions (global) for mismatched TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > MisMatchPosU
Difference of U-positions (local) for mismatched TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > PosResidueZGlobal
Z-position (global) difference between TrueHit and SpacePoint (layerwise)
std::array< std::vector< double >, c_nPlanes > MisMatchMomX
Difference of Momentum in X-Direction for TrueHits that do not match but are related from one SpacePo...
std::array< std::vector< double >, c_nPlanes > MisMatchPosV
Difference of V-positions (local) for mismatched TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > MisMatchPosResiduals
Distance between TrueHits that do not match but are related from one SpacePoint (layerwise)
std::array< std::vector< double >, c_nPlanes > MisMatchPosY
Difference of Y-positions (global) for mismatched TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > SpacePointYGlobal
global y-positions of SpacePoints (layerwise)
std::array< std::vector< double >, c_nPlanes > MisMatchPosX
Difference of X-positions (global) for mismatched TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > SpacePointVLocal
local v-positions of SpacePoints (layerwise)
std::array< std::vector< double >, c_nPlanes > TrueHitULocal
local u-positions of TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > SpacePointZGlobal
global z-positions of SpacePoints (layerwise)
std::array< std::vector< double >, c_nPlanes > PosResidueULocal
U-position (local) differnece between TrueHit and SpacePoint (layerwise)
std::array< std::vector< double >, c_nPlanes > TrueHitXGlobal
global x-positions of TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > PosResiduesGlobal
position differences in global coordinates (layerwise)
std::array< std::vector< double >, c_nPlanes > SpacePointULocal
local u-positions of SpacePoints (layerwise)
std::array< std::vector< double >, c_nPlanes > TrueHitVLocal
local v-positions of TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > MisMatchMomY
Difference of Momentum in Y-Direction for TrueHits that do not match but are related from one SpacePo...
std::array< std::vector< double >, c_nPlanes > PosResidueYGlobal
Y-position (global) difference between TrueHit and SpacePoint (layerwise)
std::array< std::vector< double >, c_nPlanes > SpacePointXGlobal
global x-positions of SpacePoints (layerwise)
std::array< std::vector< double >, c_nPlanes > PosResiduesLocal
position differences in local coordinates (layerwise)
std::array< std::vector< double >, c_nPlanes > PosResidueXGlobal
X-position (global) difference between TrueHit and SpacePoint (layerwise)
std::array< std::vector< double >, c_nPlanes > TrueHitZGlobal
global z-positions of TrueHits (layerwise)
std::array< std::vector< double >, c_nPlanes > TrueHitYGlobal
global y-positions of TrueHits (layerwise)
struct for easier handling of getting U- & V-position of SpacePoints and some difficulties that arise...