Belle II Software development
StaticSector.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#pragma once
9#include <tracking/dataobjects/FullSecID.h>
10#include <tracking/trackFindingVXD/filterMap/map/CompactSecIDs.h>
11#include <tracking/trackFindingVXD/filterMap/filterFramework/TBranchLeafType.h>
12
13#include <TF1.h>
14
15#include <unordered_map>
16#include <tuple>
17#include <utility>
18
19namespace Belle2 {
27 template < class HitType,
28 class Filter2sp, class Filter3sp, class Filter4sp >
30
31 // allows VXDTFFilters to modify private members and use private functions, which is kind of dangerous
32 // maybe better to make the one function needed to modify the filters public
33 template<class T>
34 friend class VXDTFFilters;
35
36 public:
37
39
43 {}
44
46 explicit StaticSector(FullSecID secID) :
47 m_secID(secID), m_compactSecIDsMap(nullptr)
48 {}
49
50
51 // ACCESSOR FUNCTIONS
52
56 const Filter2sp* getFilter2sp(FullSecID innerSector) const
57 {
58 auto filter = m_2spFilters.find(m_compactSecIDsMap->getCompactID(innerSector));
59 if (filter == m_2spFilters.end()) {
60 return nullptr;
61 }
62 return &(filter->second);
63 }
64
69 const Filter3sp* getFilter3sp(const FullSecID& centerID, const FullSecID& innerID) const
70 {
71 auto filter = m_3spFilters.find(m_compactSecIDsMap->getCompactID(centerID, innerID));
72 if (filter == m_3spFilters.end()) {
73 return nullptr;
74 }
75 return &(filter->second);
76 }
77
78
82 const Filter4sp* getFilter4sp(const FullSecID& /*outerCenterID*/, const FullSecID& /*innerCenterID*/,
83 const FullSecID& /*innerID*/) const
84 {
85 B2WARNING("StaticSector:getFilter4sp 4 hit, all 4 hits are yet ignored in here! TODO: implement!");
86 return nullptr;
87 }
88
89
90
93 const std::unordered_map<CompactSecIDs::sectorID_t, Filter2sp >&
95 {return m_2spFilters;}
96
99 const std::unordered_map<CompactSecIDs::secPairID_t, Filter3sp >&
101 {return m_3spFilters;}
102
105 const std::unordered_map<CompactSecIDs::secTripletID_t, Filter4sp >&
107 {return m_4spFilters;}
108
109
111 void assignCompactSecIDsMap(const CompactSecIDs& compactSecIDsMap)
112 {
113 m_compactSecIDsMap = & compactSecIDsMap;
114 }
115
116
121 void assign2spFilter(FullSecID inner, const Filter2sp filter)
122 {
124 m_inner2spSecIDs.push_back(inner);
125 }
126
129 const Filter3sp filter)
130 {
132 m_inner3spSecIDs.push_back({center, inner});
133 }
134
136 void assign4spFilter(FullSecID outerCenter, FullSecID innerCenter, FullSecID inner,
137 const Filter3sp filter)
138 {
139 m_4spFilters[ m_compactSecIDsMap->getCompactID(outerCenter, innerCenter, inner) ] =
140 filter;
141 m_inner4spSecIDs.push_back({outerCenter, innerCenter, inner});
142 }
143
144
146 // @param sublayer : the new SubLayerID will be 0 if sublayer==0, and will be 1 else
147 void setSubLayerID(int sublayer)
148 {
149 m_secID = FullSecID(m_secID.getVxdID(), (bool)sublayer, m_secID.getSecID());
150 }
151
152
154 const std::vector< FullSecID >& getInner2spSecIDs() const
155 {
156 return m_inner2spSecIDs;
157 }
158
159
161 const std::vector< std::pair<FullSecID, FullSecID> >& getInner3spSecIDs() const
162 {
163 return m_inner3spSecIDs;
164 }
165
167 const std::vector< std::tuple<FullSecID, FullSecID, FullSecID> >& getInner4spSecIDs() const
168 {
169 return m_inner4spSecIDs;
170 }
171
173 FullSecID getFullSecID() const { return m_secID; }
174
175
177
180 inline bool operator == (const StaticSector& b) const
181 {
182 return (getFullSecID() == b.getFullSecID());
183 }
184
186 inline bool operator == (const FullSecID& b) const
187 {
188 return (getFullSecID() == b);
189 }
190
192
193 private:
194
200 void modify2SPFilters(const std::vector< std::tuple<int, std::string> >& adjustFunctions)
201 {
202 // loop over all 2SP-filters
203 for (auto& filter : m_2spFilters) modifySingleFilter<Filter2sp>(filter.second, adjustFunctions);
204 }
205
211 void modify3SPFilters(const std::vector< std::tuple<int, std::string > >& adjustFunctions)
212 {
213 // loop over all 3SP-filters
214 for (auto& filter : m_3spFilters) modifySingleFilter<Filter3sp>(filter.second, adjustFunctions);
215 }
216
226 template<class FilterType>
227 void modifySingleFilter(FilterType& filter, const std::vector< std::tuple<int, std::string > >& adjustFunctions)
228 {
229 // get the "map" to the cut values, the char in the pair codes the type, and the pointer points to the value
230 std::vector< std::pair<char, void*> > accessor = {};
231 filter.getNameAndReference(&accessor);
232
233 // this will produce lots of output
234 B2DEBUG(20, std::endl << "BEFORE: " << filter.getNameAndReference() << std::endl);
235
236 // loop over all adjustfunctions
237 for (const std::tuple<int, std::string >& entry : adjustFunctions) {
238 // first is the index
239 int index = std::get<0>(entry);
240
241 if (index < 0 || index >= (int)accessor.size()) {
242 B2FATAL("Provided index is out of range! index = " << index << " number of entries = " << accessor.size());
243 }
244
245 // now do some casting magic
246 double x = 0;
247 // the secID will be treated as 0th parameter of the TF1 ([0]) if specified.
248 double y = m_secID;
249
250 char typeID = accessor[index].first;
251 void* valuePtr = accessor[index].second;
252 if (typeID == TBranchLeafType(double())) x = *((double*)valuePtr);
253 else if (typeID == TBranchLeafType(int())) x = (double)(*((int*)valuePtr));
254 else if (typeID == TBranchLeafType(float())) x = (double)(*((float*)valuePtr));
255 else if (typeID == TBranchLeafType(bool())) x = (double)(*((bool*)valuePtr));
256 else {
257 B2FATAL("Unrecognized type : " << typeID);
258 } // end else if
259
260 // create function
261 TF1 f("function", std::get<1>(entry).c_str());
262 if (!f.IsValid() || f.GetNpar() > 1) {
263 B2FATAL("No valid function provided! The provided string has to be able to be converted by TF1. Also max. 1 parameter is allowed!"
264 << " The provided string is: \"" << std::get<1>(entry).c_str() << "\"");
265 }
266
267 double result = f.EvalPar(&x, &y);
268
269 // now cast back to original type and set the value
270 if (typeID == TBranchLeafType(double())) *((double*)valuePtr) = result;
271 else if (typeID == TBranchLeafType(int())) *((int*)valuePtr) = result;
272 else if (typeID == TBranchLeafType(float())) *((float*)valuePtr) = result;
273 else if (typeID == TBranchLeafType(bool())) *((bool*)valuePtr) = result;
274 else {
275 B2FATAL("Unrecognized type : " << typeID);
276 } // end else if
277
278 }// end loop over the functions
279
280 B2DEBUG(20, "AFTER: " << filter.getNameAndReference() << std::endl);
281 }
282
285
288
290 std::vector< FullSecID > m_inner2spSecIDs;
292 std::vector< std::pair< FullSecID, FullSecID> > m_inner3spSecIDs;
294 std::vector< std::tuple< FullSecID, FullSecID, FullSecID > > m_inner4spSecIDs;
295
297 std::unordered_map<CompactSecIDs::sectorID_t, Filter2sp > m_2spFilters;
299 std::unordered_map<CompactSecIDs::secPairID_t, Filter3sp > m_3spFilters;
301 std::unordered_map<CompactSecIDs::secTripletID_t, Filter4sp > m_4spFilters;
302
303 };
304
305
307} //Belle2 namespace
This class provides a computer convenient numbering scheme for the sectors in the sector map and for ...
Definition: CompactSecIDs.h:28
sectorID_t getCompactID(const FullSecID &fullID) const
Returns the compact id of the FullSecID It does not throw exceptions (at least it should not).
Class to identify a sector inside of the VXD.
Definition: FullSecID.h:33
VxdID getVxdID() const
returns VxdID of sensor.
Definition: FullSecID.h:138
short int getSecID() const
returns SecID of current FullSecID (only unique for each sensor).
Definition: FullSecID.h:146
class to describe a static sector of the sector map.
Definition: StaticSector.h:29
const std::vector< std::tuple< FullSecID, FullSecID, FullSecID > > & getInner4spSecIDs() const
returns all IDs for inner sectors of four-sector-combinations
Definition: StaticSector.h:167
void assign2spFilter(FullSecID inner, const Filter2sp filter)
Assign the 2 space point.
Definition: StaticSector.h:121
std::vector< std::pair< FullSecID, FullSecID > > m_inner3spSecIDs
stores innerSecIDs for the attached 3-hit filters
Definition: StaticSector.h:292
const Filter3sp * getFilter3sp(const FullSecID &centerID, const FullSecID &innerID) const
Get the pointer to the 3 Space Point filter assigned to the friendship relation among this sector; wi...
Definition: StaticSector.h:69
FullSecID getFullSecID() const
returns FullSecID of this sector
Definition: StaticSector.h:173
std::vector< FullSecID > m_inner2spSecIDs
stores innerSecIDs for the attached 2-hit filters
Definition: StaticSector.h:290
const std::unordered_map< CompactSecIDs::secPairID_t, Filter3sp > & getAllFilters3sp() const
Get constant access to the whole set of 3 Space Point filters.
Definition: StaticSector.h:100
const CompactSecIDs * m_compactSecIDsMap
map from FullSecID to CompactSecID
Definition: StaticSector.h:287
const std::unordered_map< CompactSecIDs::sectorID_t, Filter2sp > & getAllFilters2sp() const
Get constant access to the whole set of 2 Space Point filters.
Definition: StaticSector.h:94
void assign4spFilter(FullSecID outerCenter, FullSecID innerCenter, FullSecID inner, const Filter3sp filter)
Parameters: pass the ID of the inner sectors (sorted from outer(left) to inner(right) and the filters...
Definition: StaticSector.h:136
FullSecID m_secID
stores its own secID
Definition: StaticSector.h:284
void assignCompactSecIDsMap(const CompactSecIDs &compactSecIDsMap)
Assign the compact sector ID to this sector.
Definition: StaticSector.h:111
StaticSector()
CONSTRUCTORS.
Definition: StaticSector.h:41
const Filter4sp * getFilter4sp(const FullSecID &, const FullSecID &, const FullSecID &) const
Get the pointer to the 4 Space Point filter assigned to the WARNING: not implemented yet.
Definition: StaticSector.h:82
std::unordered_map< CompactSecIDs::secTripletID_t, Filter4sp > m_4spFilters
stores the attached 4-hit filters
Definition: StaticSector.h:301
void modify2SPFilters(const std::vector< std::tuple< int, std::string > > &adjustFunctions)
PRIVATE MEMBERS.
Definition: StaticSector.h:200
std::unordered_map< CompactSecIDs::sectorID_t, Filter2sp > m_2spFilters
stores the attached 2-hit filters
Definition: StaticSector.h:297
const std::vector< std::pair< FullSecID, FullSecID > > & getInner3spSecIDs() const
returns all IDs for inner sectors of three-sector-combinations
Definition: StaticSector.h:161
const std::unordered_map< CompactSecIDs::secTripletID_t, Filter4sp > & getAllFilters4sp() const
Get constant access to the whole set of 4 Space Point filters.
Definition: StaticSector.h:106
void assign3spFilter(FullSecID center, FullSecID inner, const Filter3sp filter)
Parameters: pass the ID of the inner sectors (sorted from outer(left) to inner(right) and the filters...
Definition: StaticSector.h:128
void modifySingleFilter(FilterType &filter, const std::vector< std::tuple< int, std::string > > &adjustFunctions)
Function that modifies the upper and lower bounds of the Ranges contained in the filter.
Definition: StaticSector.h:227
std::vector< std::tuple< FullSecID, FullSecID, FullSecID > > m_inner4spSecIDs
stores innerSecIDs for the attached 4-hit filters
Definition: StaticSector.h:294
void modify3SPFilters(const std::vector< std::tuple< int, std::string > > &adjustFunctions)
function that modifies all 3SP-filters connected to this static sector
Definition: StaticSector.h:211
void setSubLayerID(int sublayer)
set sublayer ID, needed as it is updated in the trainings phase
Definition: StaticSector.h:147
bool operator==(const StaticSector &b) const
COMPARISON OPERATORS.
Definition: StaticSector.h:180
const std::vector< FullSecID > & getInner2spSecIDs() const
returns all IDs for inner sectors of two-sector-combinations
Definition: StaticSector.h:154
std::unordered_map< CompactSecIDs::secPairID_t, Filter3sp > m_3spFilters
stores the attached 3-hit filters
Definition: StaticSector.h:299
const Filter2sp * getFilter2sp(FullSecID innerSector) const
Get the pointer to the 2 Space Point filter assigned to the friendship relation among this sector; wi...
Definition: StaticSector.h:56
StaticSector(FullSecID secID)
constructor
Definition: StaticSector.h:46
Class that contains all the static sectors to which the filters are attached.
Definition: VXDTFFilters.h:63
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
char TBranchLeafType(const char *)
Overloading TBranchLeafType to be able to get identifier 'C' for type char*.
Abstract base class for different kinds of events.