Belle II Software  release-05-02-19
VXDTFFilters.h
1 /********************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2015 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Eugenio Paoloni, *
7  * Jakob Lettenbichler *
8  * *
9  * This software is provided "as is" without any warranty. *
10  *******************************************************************************/
11 #ifndef VXDTFFILTERS_HH
12 #define VXDTFFILTERS_HH
13 
14 #include <tracking/dataobjects/FullSecID.h>
15 
16 #include <tracking/spacePointCreation/SpacePoint.h>
17 
18 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/Distance1DZ.h>
19 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/Distance3DNormed.h>
20 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/DistanceInTimeUside.h>
21 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/DistanceInTimeVside.h>
22 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/SlopeRZ.h>
23 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/Distance1DZSquared.h>
24 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/Distance2DXYSquared.h>
25 #include <tracking/trackFindingVXD/filterMap/twoHitVariables/Distance3DSquared.h>
26 
27 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/DistanceInTime.h>
28 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/Angle3DSimple.h>
29 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/CosAngleXY.h>
30 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/AngleRZSimple.h>
31 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/CircleDist2IP.h>
32 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/DeltaSlopeRZ.h>
33 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/DeltaSlopeZoverS.h>
34 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/DeltaSoverZ.h>
35 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/HelixParameterFit.h>
36 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/Pt.h>
37 #include <tracking/trackFindingVXD/filterMap/threeHitVariables/CircleRadius.h>
38 
39 #include <tracking/trackFindingVXD/filterMap/filterFramework/Shortcuts.h>
40 #include <tracking/trackFindingVXD/filterTools/ObserverPrintResults.h>
41 #include <tracking/trackFindingVXD/filterMap/filterFramework/VoidObserver.h> // empty observer
42 
43 #include <tracking/dataobjects/SectorMapConfig.h>
44 
45 #include <framework/logging/Logger.h>
46 
47 #include <vxd/dataobjects/VxdID.h>
48 #include <tracking/trackFindingVXD/filterMap/map/CompactSecIDs.h>
49 #include <tracking/trackFindingVXD/segmentNetwork/StaticSector.h>
50 
51 #include <TFile.h>
52 #include <TString.h>
53 
54 namespace Belle2 {
65  template<class point_t>
66  class VXDTFFilters final {
67  public:
68 
69 
71 // typedef decltype((0. <= Distance3DSquared<Belle2::SpacePoint>() <= 0.).observe(ObserverPrintResults())) twoHitFilter_t;
72 
73 
75  typedef decltype(
76  (
77  0. <= DistanceInTimeUside<point_t>() <= 0. &&
78  0. <= DistanceInTimeVside<point_t>() <= 0. &&
79  0. <= Distance3DSquared<point_t>() <= 0.&&
80  0. <= Distance2DXYSquared<point_t>() <= 0.&&
81  0. <= Distance1DZ<point_t>() <= 0.&&
82  0. <= SlopeRZ<point_t>() <= 0.&&
83  0. <= Distance3DNormed<point_t>() <= 0.
84  )
86 
87  // March9th2016: TODO: we want to use a big observer observing everything - Working title: MegaObserver.
88 // typedef decltype(
89 // ((0. <= Distance3DSquared<Belle2::SpacePoint>() <= 0.).observe(VoidObserver())&&
90 // (0. <= Distance2DXYSquared<Belle2::SpacePoint>() <= 0.).observe(VoidObserver())&&
91 // (0. <= Distance1DZ<Belle2::SpacePoint>() <= 0.).observe(VoidObserver())&&
92 // (0. <= SlopeRZ<Belle2::SpacePoint>() <= 0.).observe(VoidObserver())&&
93 // (0. <= Distance3DNormed<Belle2::SpacePoint>() <= 0.).observe(VoidObserver())).observe(VoidObserver())
94 // ) twoHitFilter_t;
95 
97 // typedef decltype((0. <= Angle3DSimple<point_t>() <= 0.).observe(Observer3HitPrintResults())) threeHitFilter_t;
98 
99 
101  typedef decltype(
102  (
103  0. <= DistanceInTime<point_t>() <= 0. &&
104  0. <= Angle3DSimple<point_t>() <= 0.&&
105  0. <= CosAngleXY<point_t>() <= 0.&&
106  0. <= AngleRZSimple<point_t>() <= 0.&&
107  CircleDist2IP<point_t>() <= 0.&&
108  0. <= DeltaSlopeRZ<point_t>() <= 0.&&
109  0. <= DeltaSlopeZoverS<point_t>() <= 0.&&
110  0. <= DeltaSoverZ<point_t>() <= 0.&&
111  0. <= HelixParameterFit<point_t>() <= 0.&&
112  0. <= Pt<point_t>() <= 0.&&
113  0. <= CircleRadius<point_t>() <= 0.
114  )
116 
117 
119  typedef StaticSector< point_t, twoHitFilter_t, threeHitFilter_t , int >
121 
124  {
125  m_staticSectors.resize(2);
126  // The first static sector is not used and will never be since the first
127  // compact id is 1 and compact id = 0 is reserved to signal an error.
128  m_staticSectors[0] = nullptr;
129  // initialize the first slot of the Static sector vector
130  m_staticSectors[1] = nullptr;
131  }
132 
134  ~VXDTFFilters()
135  {
136  // delete the static sectors
137  for (staticSector_t* aSector : m_staticSectors) if (aSector != nullptr) delete aSector;
138  }
139 
145  int addSectorsOnSensor(const std::vector<double>& normalizedUsup,
146  const std::vector<double>& normalizedVsup,
147  const std::vector< std::vector<FullSecID> >& sectorIds)
148  {
149 
150  auto addedSectors = m_compactSecIDsMap.addSectors(normalizedUsup,
151  normalizedVsup,
152  sectorIds);
153 
154  if ((int) addedSectors != ((int) normalizedVsup.size() + 1) *
155  ((int) normalizedUsup.size() + 1))
156  return addedSectors ;
157 
158  addedSectors = 0;
159  try {
160  for (auto secIDrow : sectorIds)
161  for (auto secID : secIDrow) {
162  auto compactID = m_compactSecIDsMap.getCompactID(secID);
163 
164  if ((int) m_staticSectors.size() <= compactID)
165  m_staticSectors.resize(compactID + 1);
166 
167  m_staticSectors[ compactID ] = new staticSector_t(secID) ;
168  m_staticSectors[ compactID ]->assignCompactSecIDsMap(m_compactSecIDsMap);
169 
170  addedSectors ++;
171  }
172  } catch (...) { return addedSectors ; }
173 
174  return addedSectors;
175  }
176 
178  int addTwoHitFilter(FullSecID outer,
179  FullSecID inner,
180  const twoHitFilter_t& filter)
181  {
182  if (m_staticSectors.size() <= m_compactSecIDsMap[ outer ] ||
183  m_compactSecIDsMap[ outer ] == 0)
184  return 0;
185 
186  m_staticSectors[m_compactSecIDsMap[outer]]->assign2spFilter(inner, filter);
187  return 1;
188  }
189 
190 
192  int addThreeHitFilter(FullSecID outer,
193  FullSecID center,
194  FullSecID inner,
195  const threeHitFilter_t& filter)
196  {
197 
198  if (m_staticSectors.size() <= m_compactSecIDsMap[ outer ] ||
199  m_compactSecIDsMap[ outer ] == 0 ||
200  m_compactSecIDsMap[ center ] == 0 ||
201  m_compactSecIDsMap[ inner ] == 0)
202  return 0;
203 
204  m_staticSectors[m_compactSecIDsMap[outer]]->assign3spFilter(center,
205  inner,
206  filter);
207  return 1;
208  }
209 
210 
213  { return m_compactSecIDsMap[outer]; }
214 
215 
216  // const StaticSectorType& getSector(VxdID aSensorID,
217  // std::pair<float, float> normalizedLocalCoordinates) const
218  // { }
219 
223  const staticSector_t* getStaticSector(const FullSecID secID) const
224  {
225  auto sectorPosition = m_compactSecIDsMap[ secID ];
226  if (sectorPosition == 0) return nullptr;
227  return m_staticSectors[ sectorPosition ];
228  }
229 
234  const twoHitFilter_t getTwoHitFilters(const FullSecID& outer,
235  const FullSecID& inner) const
236  {
237  // TODO: sanity checks
238  static twoHitFilter_t just_in_case;
239  const auto staticSector = m_staticSectors[ m_compactSecIDsMap[ outer ] ];
240  // catch case when sector is not part of the sectorMap:
241  if (staticSector == nullptr)
242  return just_in_case;
243  const auto* filterPtr = staticSector->getFilter2sp(inner);
244  if (filterPtr == nullptr)
245  return just_in_case;
246  return *filterPtr;
247  }
248 
249 
251  bool areCoordinatesValid(VxdID aSensorID,
252  double normalizedU, double normalizedV) const
253  {
254  return m_compactSecIDsMap.areCoordinatesValid(aSensorID,
255  normalizedU, normalizedV);
256  }
257 
258 
260  // JKL: what happens here if no FullSecID could be found?
261  // EP: you get an exception
262  FullSecID getFullID(VxdID aSensorID, double normalizedU, double normalizedV) const
263  {
264  // TODO WARNING how to catch bad cases?
265  return m_compactSecIDsMap.getFullSecID(aSensorID, normalizedU, normalizedV);
266  }
267 
269  FullSecID getFullID(CompactSecIDs::sectorID_t compactSecID) const
270  {
271  return m_staticSectors.at(compactSecID)->getFullSecID();
272  }
273 
275  const SectorMapConfig& getConfig(void) const { return m_testConfig; }
279  void setConfig(const SectorMapConfig& config) { m_testConfig = config; }
280 
281 
283  const CompactSecIDs& getCompactIDsMap() const { return m_compactSecIDsMap; }
284 
285 
287  const std::vector< staticSector_t*>& getStaticSectors() const { return m_staticSectors; }
288 
290  unsigned size() const { return m_compactSecIDsMap.getSize(); }
291 
293  bool persistOnRootFile(void) const
294  {
295 
296  if (! m_testConfig.Write("config"))
297  return false;
298  if (! persistSectors())
299  return false;
300 
301  if (! persistFilters())
302  return false;
303 
304  return true;
305  };
306 
308  bool retrieveFromRootFile(const TString* dirName)
309  {
310  // locked filters cannot be modified
311  if (m_preventModification) {
312  B2FATAL("Trying to modify a locked filter! A locked filter is not supposed to be changed anymore!");
313  }
314 
315  if (! m_testConfig.Read("config"))
316  return false;
317 
318  if (! retrieveSectors(dirName))
319  return false;
320 
321  if (! retrieveFilters(dirName))
322  return false;
323 
324  return true;
325  };
326 
330  bool setSubLayerIDs(FullSecID sector, int sublayer)
331  {
332  // locked filters cannot be modified
333  if (m_preventModification) {
334  B2FATAL("Trying to modify a locked filter! A locked filter is not supposed to be changed anymore!");
335  }
336 
337  // first update the static sector
338  // the static sector is retrieved from the compactid which automatically ignores the sublayer id
339  auto sectorPosition = m_compactSecIDsMap[ sector ];
340  if (sectorPosition == 0) return false;
341  staticSector_t* staticsector = m_staticSectors[ sectorPosition ];
342  if (!staticsector) return false;
343  staticsector->setSubLayerID(sublayer);
344 
345  // then update the fullsectorid in the compactsectoridmap, the sublayerid of sector will be ignored when searching for sector to update
346  return m_compactSecIDsMap.setSubLayerID(sector, sublayer);
347  }
348 
352  void lockFilters() { m_preventModification = true; };
353 
357  void modify2SPFilters(const std::vector< std::tuple<int, std::string > >& adjustFunctions)
358  {
359 
360  // locked filters cannot be modified
362  B2FATAL("Trying to modify a locked filter! A locked filter is not supposed to be changed anymore!");
363  }
364 
365  for (auto staticSector : m_staticSectors) {
366  if (staticSector == nullptr) continue;
367  staticSector->modify2SPFilters(adjustFunctions);
368  }
369  };
370 
374  void modify3SPFilters(const std::vector< std::tuple<int, std::string > >& adjustFunctions)
375  {
376 
377  // locked filters cannot be modified
378  if (m_preventModification) {
379  B2FATAL("Trying to modify a locked filter! A locked filter is not supposed to be changed anymore!");
380  }
381 
382  for (auto staticSector : m_staticSectors) {
383  if (staticSector == nullptr) continue;
384  staticSector->modify3SPFilters(adjustFunctions);
385  }
386  };
387 
388  private:
389 
391  bool persistSectors(void) const
392  {
393  TTree* tree = new TTree(c_CompactSecIDstreeName, c_CompactSecIDstreeName);
394  UInt_t layer, ladder, sensor;
395  tree->Branch("layer" , & layer , "layer/i");
396  tree->Branch("ladder", & ladder, "ladder/i");
397  tree->Branch("sensor", & sensor, "sensor/i");
398 
399  std::vector< double >* normalizedUsup = new std::vector< double> ();
400  tree->Branch("normalizedUsup", & normalizedUsup);
401 
402  std::vector< double >* normalizedVsup = new std::vector< double> ({1., 2., 3., 4.});
403  tree->Branch("normalizedVsup", & normalizedVsup);
404 
405  std::vector< std::vector< unsigned int > >* fullSecIDs =
406  new std::vector< std::vector< unsigned int > > ();
407  tree->Branch("fullSecID", & fullSecIDs);
408 
409  unsigned nOfLayers = m_compactSecIDsMap.nOfLayers();
410  for (layer = 0 ; layer < nOfLayers ; layer ++) {
411  unsigned nOfLadders = m_compactSecIDsMap.nOfLadders(layer);
412  for (ladder = 0; ladder < nOfLadders ; ladder ++) {
413  unsigned nOfSensors = m_compactSecIDsMap.nOfSensors(layer, ladder);
414  for (sensor = 0; sensor < nOfSensors ; sensor ++) {
415  normalizedUsup->clear();
416  normalizedVsup->clear();
417  fullSecIDs->clear();
418  auto sectorsOnSensor =
419  m_compactSecIDsMap.getSectorsOnSensor(layer, ladder, sensor);
420  sectorsOnSensor.get(normalizedUsup, normalizedVsup, fullSecIDs);
421  tree->Fill();
422  }
423  }
424  }
425  delete normalizedVsup;
426  delete normalizedUsup;
427  delete fullSecIDs;
428  return true;
429  }
430 
432  bool retrieveSectors(const TString* dirName)
433  {
434  TString treeName = *dirName;
435  treeName.Append("/");
436  treeName.Append(c_CompactSecIDstreeName);
437  TTree* tree = (TTree*) gFile->Get(treeName);
438  UInt_t layer, ladder, sensor;
439  if (tree->SetBranchAddress("layer" , & layer) < 0) B2FATAL("VXDTFFilters: invalid branch address");
440  if (tree->SetBranchAddress("ladder", & ladder) < 0) B2FATAL("VXDTFFilters: invalid branch address");
441  if (tree->SetBranchAddress("sensor", & sensor) < 0) B2FATAL("VXDTFFilters: invalid branch address");
442 
443  std::vector< double >* normalizedUsup = new std::vector< double> ();
444  if (tree->SetBranchAddress("normalizedUsup", & normalizedUsup) < 0) B2FATAL("VXDTFFilters: invalid branch address");
445 
446  std::vector< double >* normalizedVsup = new std::vector< double> ({1., 2., 3., 4.});
447  if (tree->SetBranchAddress("normalizedVsup", & normalizedVsup) < 0) B2FATAL("VXDTFFilters: invalid branch address");
448 
449  std::vector< std::vector< unsigned int > >* fullSecIDs =
450  new std::vector< std::vector< unsigned int > > ();
451  if (tree->SetBranchAddress("fullSecID", & fullSecIDs) < 0) B2FATAL("VXDTFFilters: invalid branch address");
452 
453 
454  for (Long64_t i = 0; i < tree->GetEntries() ; i++) {
455  tree->GetEntry(i);
456  this->addSectorsOnSensor(* normalizedUsup,
457  * normalizedVsup,
458  * fullSecIDs);
459  }
460 
461  delete normalizedVsup;
462  delete normalizedUsup;
463  delete fullSecIDs;
464  return true;
465  }
466 
468  bool persistFilters(void) const
469  {
470 
471  TTree* sp2tree = new TTree("SegmentFilters", "SegmentFilters");
472  twoHitFilter_t twoHitFilter;
473  twoHitFilter.persist(sp2tree, "filter");
474 
475  unsigned int outerFullSecID2sp, innerFullSecID2sp;
476  sp2tree->Branch("outerFullSecID", & outerFullSecID2sp);
477  sp2tree->Branch("innerFullSecID", & innerFullSecID2sp);
478 
479 
480  TTree* sp3tree = new TTree("TripletsFilters", "TripletFilters");
481  threeHitFilter_t threeHitFilter;
482  threeHitFilter.persist(sp3tree, "filter");
483 
484  unsigned int outerFullSecID3sp, centerFullSecID3sp,
485  innerFullSecID3sp;
486  sp3tree->Branch("outerFullSecID", & outerFullSecID3sp);
487  sp3tree->Branch("centerFullSecID", & centerFullSecID3sp);
488  sp3tree->Branch("innerFullSecID", & innerFullSecID3sp);
489 
490  for (auto staticSector : m_staticSectors) {
491  if (staticSector == nullptr)
492  // Why there is an empty sector per layer?
493  continue;
494 
495 
496  outerFullSecID3sp = outerFullSecID2sp = staticSector->getFullSecID();
497  auto segmentFilters = staticSector->getAllFilters2sp();
498  for (auto compactIdFilterPair : segmentFilters) {
499  auto innerCompactId = compactIdFilterPair.first;
500  innerFullSecID2sp = getFullID(innerCompactId);
501  twoHitFilter = compactIdFilterPair.second;
502  sp2tree->Fill();
503  }
504 
505  auto tripletFilters = staticSector->getAllFilters3sp();
506  for (auto compactIdFilterPair : tripletFilters) {
507  CompactSecIDs::sectorID_t id_center, id_inner;
508  CompactSecIDs::extractCompactID(compactIdFilterPair.first, id_center, id_inner);
509  centerFullSecID3sp = getFullID(id_center);
510  innerFullSecID3sp = getFullID(id_inner);
511  threeHitFilter = compactIdFilterPair.second;
512  sp3tree->Fill();
513  }
514 
515 
516  }
517 
518  return true;
519  }
520 
522  bool retrieveFilters(const TString* dirName)
523  {
524  TString sp2treeName = *dirName;
525  sp2treeName.Append("/SegmentFilters");
526  TTree* sp2tree = (TTree*) gFile->Get(sp2treeName);
527  if (!sp2tree)
528  return false;
529 
530  twoHitFilter_t twoHitFilter;
531  twoHitFilter.setBranchAddress(sp2tree, "filter");
532 
533  unsigned int outerFullSecID2sp, innerFullSecID2sp;
534  if (sp2tree->SetBranchAddress("outerFullSecID", & outerFullSecID2sp) < 0) B2FATAL("VXDTFFilters: invalid branch address");
535  if (sp2tree->SetBranchAddress("innerFullSecID", & innerFullSecID2sp) < 0) B2FATAL("VXDTFFilters: invalid branch address");
536 
537  for (Long64_t i = 0 ; i < sp2tree->GetEntries() ; i++) {
538  sp2tree->GetEntry(i);
539 
540  // cross check to only put filters into the map which outer sector is also on outer layer!
541  FullSecID outer_secid_2sp(outerFullSecID2sp);
542  FullSecID inner_secid_2sp(innerFullSecID2sp);
543  // equal layer numbers are allowed!
544  if (outer_secid_2sp.getLayerNumber() < inner_secid_2sp.getLayerNumber()) {
545  B2WARNING("Outer sector is not on outer layer! Not adding this filter. \"Outer\" layer number: "
546  << outer_secid_2sp.getLayerNumber() << " \"Inner\" layer number " << inner_secid_2sp.getLayerNumber());
547  continue;
548  }
549 
550  // add filter to the map
551  if (!addTwoHitFilter(outer_secid_2sp, inner_secid_2sp,
552  twoHitFilter))
553  return false;
554 
555  }
556 
557  TString sp3treeName = *dirName;
558  sp3treeName.Append("/TripletsFilters");
559  TTree* sp3tree = (TTree*) gFile->Get(sp3treeName);
560  if (! sp3tree)
561  return false;
562  threeHitFilter_t threeHitFilter;
563  threeHitFilter.setBranchAddress(sp3tree, "filter");
564 
565  unsigned int outerFullSecID3sp, centerFullSecID3sp,
566  innerFullSecID3sp;
567  if (sp3tree->SetBranchAddress("outerFullSecID", & outerFullSecID3sp) < 0) B2FATAL("VXDTFFilters: invalid branch address");
568  if (sp3tree->SetBranchAddress("centerFullSecID", & centerFullSecID3sp) < 0) B2FATAL("VXDTFFilters: invalid branch address");
569  if (sp3tree->SetBranchAddress("innerFullSecID", & innerFullSecID3sp) < 0) B2FATAL("VXDTFFilters: invalid branch address");
570 
571  for (Long64_t i = 0 ; i < sp3tree->GetEntries() ; i++) {
572  sp3tree->GetEntry(i);
573 
574  // cross check to only put filters which layers have correct order
575  FullSecID outer_secid_3sp(outerFullSecID3sp);
576  FullSecID center_secid_3sp(centerFullSecID3sp);
577  FullSecID inner_secid_3sp(innerFullSecID3sp);
578  // equal layer numbers are allowed
579  if (outer_secid_3sp.getLayerNumber() < center_secid_3sp.getLayerNumber() or
580  center_secid_3sp.getLayerNumber() < inner_secid_3sp.getLayerNumber()) {
581  B2WARNING("Layers not in the correct order for Triplet filter! Will not add filter! Outer layer number: " <<
582  outer_secid_3sp.getLayerNumber() << " center layer number: " << center_secid_3sp.getLayerNumber() <<
583  " inner layer number: " << inner_secid_3sp.getLayerNumber());
584  continue;
585  }
586 
587  // add the filter to the map
588  if (!addThreeHitFilter(outer_secid_3sp, center_secid_3sp,
589  inner_secid_3sp,
590  threeHitFilter))
591  return false;
592 
593  }
594 
595  return true;
596  }
597 
601  int addSectorsOnSensor(const std::vector< double>& normalizedUsup,
602  const std::vector< double>& normalizedVsup,
603  const std::vector< std::vector< unsigned int >>&
604  fullSecIDsBaseType)
605  {
606  std::vector< std::vector< FullSecID >> fullSecIDs;
607 
608  for (auto col : fullSecIDsBaseType) {
609  std::vector< FullSecID > tmp_col;
610  for (auto id : col)
611  tmp_col.push_back(FullSecID(id));
612  fullSecIDs.push_back(tmp_col);
613  }
614 
615  return addSectorsOnSensor(normalizedUsup, normalizedVsup, fullSecIDs);
616  }
617 
623  CompactSecIDs m_compactSecIDsMap;
624 
628  std::vector< staticSector_t* > m_staticSectors;
629 
633 
635  const char* c_CompactSecIDstreeName = "CompactSecIDs";
636 
641  bool m_preventModification = false;
642  };
643 
645 }
646 
647 
648 #endif
Belle2::VXDTFFilters::size
unsigned size() const
returns number of compact secIDs stored for this filter-container.
Definition: VXDTFFilters.h:299
Belle2::VXDTFFilters::modify2SPFilters
void modify2SPFilters(const std::vector< std::tuple< int, std::string > > &adjustFunctions)
modifies the 2SP-filters according to the functions given
Definition: VXDTFFilters.h:366
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::CompactSecIDs::nOfSensors
int nOfSensors(int layer, int ladder) const
get the number of sensors on
Definition: CompactSecIDs.h:334
Belle2::VXDTFFilters::getConfig
const SectorMapConfig & getConfig(void) const
returns the configuration settings for this VXDTFFilters.
Definition: VXDTFFilters.h:284
Belle2::CompactSecIDs::nOfLadders
int nOfLadders(int layer) const
get the number of ladders on
Definition: CompactSecIDs.h:325
Belle2::VXDTFFilters
Class that contains all the static sectors to which the filters are attached.
Definition: VXDTFFilters.h:75
Belle2::VXDTFFilters::retrieveFromRootFile
bool retrieveFromRootFile(const TString *dirName)
Retrieves from the current TDirectory all the VXDTFFilters.
Definition: VXDTFFilters.h:317
Belle2::VXDTFFilters::staticSector_t
StaticSector< point_t, twoHitFilter_t, threeHitFilter_t, int > staticSector_t
typedef to make a static sector type more readable.
Definition: VXDTFFilters.h:129
Belle2::CompactSecIDs::sectorID_t
uint16_t sectorID_t
Typedef of the compact Id for a single sector.
Definition: CompactSecIDs.h:45
Belle2::CompactSecIDs::extractCompactID
static void extractCompactID(secPairID_t pair_id, sectorID_t &id1, sectorID_t &id2)
Sets the two compact Sector id.
Definition: CompactSecIDs.h:210
Belle2::VXDTFFilters::retrieveSectors
bool retrieveSectors(const TString *dirName)
Read the whole CompactSecIDs from the current TDirectory.
Definition: VXDTFFilters.h:441
Belle2::CompactSecIDs::nOfLayers
int nOfLayers(void) const
get the number of layers in this CompactSecIDs
Definition: CompactSecIDs.h:320
Belle2::filter
std::map< ExpRun, std::pair< double, double > > filter(const std::map< ExpRun, std::pair< double, double >> &runs, double cut, std::map< ExpRun, std::pair< double, double >> &runsRemoved)
filter events to remove runs shorter than cut, it stores removed runs in runsRemoved
Definition: Splitter.cc:43
Belle2::VXDTFFilters::addTwoHitFilter
int addTwoHitFilter(FullSecID outer, FullSecID inner, const twoHitFilter_t &filter)
adds a two hit filter
Definition: VXDTFFilters.h:187
Belle2::VXDTFFilters::twoHitFilter_t
decltype((0.<=DistanceInTimeUside< point_t >()<=0. &&0.<=DistanceInTimeVside< point_t >()<=0. &&0.<=Distance3DSquared< point_t >()<=0.&&0.<=Distance2DXYSquared< point_t >()<=0.&&0.<=Distance1DZ< point_t >()<=0.&&0.<=SlopeRZ< point_t >()<=0.&&0.<=Distance3DNormed< point_t >()<=0.)) typedef twoHitFilter_t
minimal working 2-hits-example used for redesign of VXDTF.
Definition: VXDTFFilters.h:94
Belle2::CompactSecIDs::getCompactID
sectorID_t getCompactID(const FullSecID &fullID) const
Returns the compact id of the FullSecID.
Definition: CompactSecIDs.h:124
Belle2::VXDTFFilters::m_staticSectors
std::vector< staticSector_t * > m_staticSectors
This vector contains all the static sectors on a sector map.
Definition: VXDTFFilters.h:637
Belle2::VXDTFFilters::m_preventModification
bool m_preventModification
The filters are not supposed to be altered after initialization (typically in the Module::initialize(...
Definition: VXDTFFilters.h:650
Belle2::VXDTFFilters::areCoordinatesValid
bool areCoordinatesValid(VxdID aSensorID, double normalizedU, double normalizedV) const
check if using getFullID() would be safe (true if it is safe):
Definition: VXDTFFilters.h:260
Belle2::CompactSecIDs::setSubLayerID
bool setSubLayerID(FullSecID &sector, int sublayer)
set the SublayerID of the sector
Definition: CompactSecIDs.h:345
Belle2::VXDTFFilters::getStaticSector
const staticSector_t * getStaticSector(const FullSecID secID) const
returns pointer to static sector for given fullSecID.
Definition: VXDTFFilters.h:232
Belle2::VXDTFFilters::persistOnRootFile
bool persistOnRootFile(void) const
Persists (i.e.: writes) on the current TDirectory the whole object.
Definition: VXDTFFilters.h:302
Belle2::FullSecID
Class to identify a sector inside of the VXD.
Definition: FullSecID.h:43
Belle2::VXDTFFilters::getStaticSectors
const std::vector< staticSector_t * > & getStaticSectors() const
JKL: intended for some checks only - returns CompactIDsMap storing the static sectors.
Definition: VXDTFFilters.h:296
Belle2::CompactSecIDs::getSectorsOnSensor
SectorsOnSensor< sectorID_t > getSectorsOnSensor(unsigned layer, unsigned ladder, unsigned sensor) const
Getter for IDs of all sectors on a sensor.
Definition: CompactSecIDs.h:141
Belle2::CompactSecIDs::areCoordinatesValid
bool areCoordinatesValid(VxdID aSensorID, double normalizedU, double normalizedV) const
JKL: returns true if operator() will not throw an exception.
Definition: CompactSecIDs.h:293
Belle2::VXDTFFilters::lockFilters
void lockFilters()
This function should be called only AFTER all adjustments to the filters have been performed.
Definition: VXDTFFilters.h:361
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::VXDTFFilters::persistSectors
bool persistSectors(void) const
Persists all the sectors on the current TDirectory.
Definition: VXDTFFilters.h:400
Belle2::VXDTFFilters::modify3SPFilters
void modify3SPFilters(const std::vector< std::tuple< int, std::string > > &adjustFunctions)
modifies the 3SP-filters according to the functions given
Definition: VXDTFFilters.h:383
Belle2::VXDTFFilters::getCompactID
CompactSecIDs::sectorID_t getCompactID(FullSecID outer) const
returns compactSecID for given FullSecID, == 0 if not found.
Definition: VXDTFFilters.h:221
Belle2::CompactSecIDs::getSize
int getSize() const
Returns the number of sectors defined so far.
Definition: CompactSecIDs.h:63
Belle2::FullSecID::getLayerNumber
short int getLayerNumber() const
returns LayerID compatible with basf2 standards.
Definition: FullSecID.h:132
Belle2::VXDTFFilters::addSectorsOnSensor
int addSectorsOnSensor(const std::vector< double > &normalizedUsup, const std::vector< double > &normalizedVsup, const std::vector< std::vector< FullSecID > > &sectorIds)
To add an array of sectors on a sensor.
Definition: VXDTFFilters.h:154
Belle2::StaticSector
class to describe a static sector of the sector map.
Definition: StaticSector.h:39
Belle2::VXDTFFilters::setSubLayerIDs
bool setSubLayerIDs(FullSecID sector, int sublayer)
during the trainings phase the sublayer ids have to be updated
Definition: VXDTFFilters.h:339
Belle2::VXDTFFilters::setConfig
void setConfig(const SectorMapConfig &config)
set the configuration which is used to create this filter
Definition: VXDTFFilters.h:288
Belle2::VXDTFFilters::m_compactSecIDsMap
CompactSecIDs m_compactSecIDsMap
This member takes care of converting the [layer][ladder] [sensor][sector] multi index into a linear i...
Definition: VXDTFFilters.h:632
Belle2::CompactSecIDs::addSectors
int addSectors(const std::vector< double > &normalizedUsup, const std::vector< double > &normalizedVsup, const std::vector< std::vector< FullSecID >> &fullSecIDs)
This method defines all the sectors on a given sensor.
Definition: CompactSecIDs.h:79
Belle2::VXDTFFilters::addThreeHitFilter
int addThreeHitFilter(FullSecID outer, FullSecID center, FullSecID inner, const threeHitFilter_t &filter)
adds a three hit filter
Definition: VXDTFFilters.h:201
Belle2::VXDTFFilters::getTwoHitFilters
const twoHitFilter_t getTwoHitFilters(const FullSecID &outer, const FullSecID &inner) const
retrieves a two hit filter:
Definition: VXDTFFilters.h:243
Belle2::SectorMapConfig
simple struct containing all the configuration data needed for the SecMapTrainer.
Definition: SectorMapConfig.h:36
Belle2::VXDTFFilters::c_CompactSecIDstreeName
const char * c_CompactSecIDstreeName
name of the tree the SecIDs are stored in when persisted
Definition: VXDTFFilters.h:644
Belle2::VXDTFFilters::m_testConfig
SectorMapConfig m_testConfig
Configuration: i.e.
Definition: VXDTFFilters.h:641
Belle2::VXDTFFilters::persistFilters
bool persistFilters(void) const
Persists on the current TDirectory the StaticSectors.
Definition: VXDTFFilters.h:477
Belle2::VXDTFFilters::threeHitFilter_t
decltype((0.<=DistanceInTime< point_t >()<=0. &&0.<=Angle3DSimple< point_t >()<=0.&&0.<=CosAngleXY< point_t >()<=0.&&0.<=AngleRZSimple< point_t >()<=0.&&CircleDist2IP< point_t >()<=0.&&0.<=DeltaSlopeRZ< point_t >()<=0.&&0.<=DeltaSlopeZoverS< point_t >()<=0.&&0.<=DeltaSoverZ< point_t >()<=0.&&0.<=HelixParameterFit< point_t >()<=0.&&0.<=Pt< point_t >()<=0.&&0.<=CircleRadius< point_t >()<=0.)) typedef threeHitFilter_t
minimal working example for 3-hits:
Definition: VXDTFFilters.h:124
Belle2::CompactSecIDs::getFullSecID
FullSecID getFullSecID(VxdID aSensorID, double normalizedU, double normalizedV) const
Returns a fullSecID for given sensor and pair of coordinates.
Definition: CompactSecIDs.h:161
Belle2::VXDTFFilters::retrieveFilters
bool retrieveFilters(const TString *dirName)
Retrieves from the current TDirectory the StaticSectors.
Definition: VXDTFFilters.h:531
Belle2::VXDTFFilters::~VXDTFFilters
~VXDTFFilters()
Destructor.
Definition: VXDTFFilters.h:143
Belle2::VXDTFFilters::getFullID
FullSecID getFullID(VxdID aSensorID, double normalizedU, double normalizedV) const
returns fullSecID for given sensorID and local coordinates.
Definition: VXDTFFilters.h:271
Belle2::VXDTFFilters::getCompactIDsMap
const CompactSecIDs & getCompactIDsMap() const
JKL: intended for some checks only - returns CompactIDsMap storing the static sectors.
Definition: VXDTFFilters.h:292
Belle2::StaticSector::setSubLayerID
void setSubLayerID(int sublayer)
set sublayer ID, needed as it is updated in the trainings phase
Definition: StaticSector.h:155