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