Belle II Software release-09-00-00
RawSecMapMergerModule.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/VXDTFHelperTools/RawSecMapMergerModule.h>
10#include <tracking/spacePointCreation/SpacePoint.h>
11#include <tracking/trackFindingVXD/environment/VXDTFFiltersHelperFunctions.h>
12#include <tracking/trackFindingVXD/filterMap/filterFramework/SelectionVariableNamesToFunctions.h>
13#include <vxd/geometry/GeoCache.h>
14
15using namespace std;
16using namespace Belle2;
17
18//-----------------------------------------------------------------
19// Register the Module
20//-----------------------------------------------------------------
21REG_MODULE(RawSecMapMerger);
22
23//-----------------------------------------------------------------
24// Implementation
25//-----------------------------------------------------------------
26
28{
29 //Set module properties
30 setDescription("this module takes a root file containing a raw sectorMap created by the SecMapTrainerBaseModule and converts it to a sectormap which can be read by the VXDTF. Please check the parameters to be set...");
31// setPropertyFlags(c_ParallelProcessingCertified); /// WARNING this module should _not_ be used for parallel processing! Its task is to create the sector maps only once...
32
33
34 addParam("rootFileNames", m_PARAMrootFileNames,
35 "List of files (wildcards not allowed - use python glob.glob() to expand to list of files)", {"lowTestRedesign_454970355.root"});
36
37 addParam("mapNames", m_PARAMmapNames, "names of sectorMaps to be loaded.", {""});
38
39 addParam("printFullGraphs", m_PARAMprintFullGraphs,
40 "If true, the full trained graphs will be printed to screen. WARNING: produces a lot of output for full detector-cases!",
41 bool(false));
42
43 addParam("threshold", m_RelThreshold,
44 "Relative threshold (in %) used to prune the sector maps. Will remove X % of the least used subgraphs.", {0});
45}
46
47
48
49
50
51std::vector<std::string> RawSecMapMergerModule::getRootFiles(std::string mapName)
52{
53 B2INFO("RawSecMapMerger::getRootFiles(): loading mapName: " << mapName);
54
55 std::vector<std::string> files4ThisMap;
56 for (std::string& fileName : m_PARAMrootFileNames) {
57 if (fileName.find(mapName) == string::npos) {
58 B2DEBUG(20, "getRootFiles: fileName " << fileName << " was _not_ accepted for map " << mapName);
59 continue;
60 }
61 B2DEBUG(20, "getRootFiles: fileName " << fileName << " accepted for map " << mapName);
62 files4ThisMap.push_back(fileName);
63 }
64 return files4ThisMap;
65}
66
67
68
69
70
71std::unique_ptr<TChain> RawSecMapMergerModule::createTreeChain(const SectorMapConfig& configuration, const std::string& nHitString)
72{
73 B2INFO("RawSecMapMerger::createTreeChain(): loading mapName: " << configuration.secMapName << " with extension " << nHitString);
74 unique_ptr<TChain> treeChain = unique_ptr<TChain>(new TChain((configuration.secMapName + nHitString).c_str()));
75
76 std::vector<std::string> fileList = getRootFiles(configuration.secMapName);
77 for (std::string& file : fileList) { treeChain->Add(file.c_str()); }
78
79 return treeChain;
80}
81
82
83
84
85
86template<class ValueType> std::vector<BranchInterface<ValueType>> RawSecMapMergerModule::getBranches(
87 std::unique_ptr<TChain>& chain,
88 const std::vector<std::string>& branchNames)
89{
90 std::vector<BranchInterface<ValueType>> branches;
91 B2INFO("RawSecMapMerger::getBranches(): loading branches: " << branchNames.size());
92 unsigned nBranches = branchNames.size();
93
94 branches.resize(nBranches, BranchInterface<ValueType>());
95 for (unsigned fPos = 0; fPos < nBranches; fPos++) {
96 branches[fPos].name = branchNames[fPos];
97 chain->SetBranchAddress(
98 branches[fPos].name.c_str(),
99 &(branches[fPos].value),
100 &(branches[fPos].branch));
101 }
102 B2INFO("RawSecMapMerger::getBranches(): done");
103 return branches;
104}
105
106std::string
108 unsigned nHits,
109 const SectorMapConfig&,
110 std::vector<std::string>& secBranchNames,
111 std::vector<std::string>& filterBranchNames)
112{
113 if (nHits == 2) {
114 secBranchNames = { "outerSecID", "innerSecID"};
116
117 for (const auto& filterNameToFunction : twoHitsFilterNameToFunction) {
118 std::string filterName = filterNameToFunction.first ;
119 filterBranchNames.push_back(filterName);
120 }
121 return "2Hit";
122 }
123
124 if (nHits == 3) {
125 secBranchNames = { "outerSecID", "centerSecID", "innerSecID"};
127
128 for (const auto& filterNameToFunction : threeHitsFilterNameToFunction) {
129 std::string filterName = filterNameToFunction.first ;
130 filterBranchNames.push_back(filterName);
131 }
132 return "3Hit";
133 }
134
135 B2ERROR("prepareNHitSpecificStuff: wrong chainLength!");
136 return "";
137}
138
139template <class FilterType> void RawSecMapMergerModule::trainGraph(
140 SectorGraph<FilterType>& mainGraph,
141 std::unique_ptr<TChain>& chain,
142 std::vector<BranchInterface<unsigned>>& sectorBranches,
143 std::vector<BranchInterface<double>>& filterBranches)
144{
145 auto nEntries = chain->GetEntries();
146 B2DEBUG(29, "RawSecMapMerger::trainGraph(): start of " << nEntries << " entries in tree and " << sectorBranches.size() <<
147 " branches");
148 if (nEntries == 0) { B2WARNING("trainGraph: valid file but no data stored!"); return; }
149
150 auto percentMark = nEntries / 10; auto progressCounter = 0;
151 // log all sector-combinations and determine their absolute number of appearances:
152 for (auto i = 0 ; i <= nEntries; i++) {
153 if (percentMark < 1 or (i % percentMark) == 0) {
154 B2INFO("RawSecMapMerger::trainGraph(): with mark: " << percentMark << " and i=" << i << ", " << progressCounter <<
155 "% related, mainGraph has got " << mainGraph.size() << " sectors...");
156 progressCounter += 10;
157 }
158 auto thisEntry = chain->LoadTree(i);
159
160 auto ids = getSecIDs(sectorBranches, thisEntry);
161 auto currentID = SubGraphID(ids);
162
163 auto pos = mainGraph.find(currentID);
164
165 if (pos == mainGraph.end()) { B2WARNING("trainGraph: could not find subgraph " << currentID.print() << " - skipping entry..."); continue; }
166
167 for (auto& filter : filterBranches) {
168 filter.update(thisEntry);
169 pos->second.addValue(FilterType(filter.name), filter.value);
170 }
171 } // entry-loop-end
172}
173
174
175
176
177
179 std::unique_ptr<TChain>& chain,
180 std::vector<BranchInterface<unsigned>>& sectorBranches,
181 std::vector<BranchInterface<double>>& filterBranches)
182{
183 auto nEntries = chain->GetEntries();
184 B2INFO("RawSecMapMerger::buildGraph(): start of " << nEntries << " entries in tree and " << sectorBranches.size() <<
185 " branches");
186
187 // creating main graph containing all subgraphs:
188 std::vector<std::string> filterNames;
189 for (auto& entry : filterBranches) {
190 filterNames.push_back(entry.name);
191 }
192 SectorGraph<FilterType> mainGraph(filterNames);
193
194 if (nEntries == 0) { B2WARNING("buildGraph: valid file but no data stored!"); return mainGraph; }
195 auto percentMark = nEntries / 10;
196 auto progressCounter = 0;
197
198 // log all sector-combinations and determine their absolute number of appearances:
199 for (auto i = 0 ; i <= nEntries; i++) {
200 if (percentMark < 1 or (i % percentMark) == 0) {
201 B2INFO("RawSecMapMerger::buildGraph(): with mark: " << percentMark << " and i=" << i << ", " << progressCounter <<
202 "% related, mainGraph has got " << mainGraph.size() << " sectors...");
203 progressCounter += 10;
204 }
205 auto thisEntry = chain->LoadTree(i);
206
207 std::vector<unsigned> ids = getSecIDs(sectorBranches, thisEntry);
208
209 if (! good(ids))
210 continue;
211
212 auto currentID = SubGraphID(ids);
213 B2DEBUG(29, "buildGraph-SubGraphID-print: id: " << currentID.print());
214
215 auto pos = mainGraph.find(currentID);
216
217 if (pos == mainGraph.end()) { pos = mainGraph.add(currentID); }
218
219 if (pos == mainGraph.end()) { B2WARNING("could not find nor add subgraph - skipping entry..."); continue; }
220
221 pos->second.wasFound();
222
223 for (auto& filter : filterBranches) {
224 filter.update(thisEntry);
225 pos->second.checkAndReplaceIfMinMax(FilterType(filter.name), filter.value);
226 }
227 } // entry-loop-end
228
229 B2INFO("RawSecMapMerger::buildGraph(): mainGraph finished - has now size: " << mainGraph.size());
230 B2DEBUG(20, "fullGraph-Print: " << mainGraph.print());
231
232 return mainGraph;
233}
234
235
236bool RawSecMapMergerModule::good(const std::vector<unsigned>& ids)
237{
238 switch (ids.size()) {
239 case 2:
240 if (FullSecID(ids[0]).getLayerID() == FullSecID(ids[1]).getLayerID() &&
241 FullSecID(ids[0]).getLadderID() == FullSecID(ids[1]).getLadderID()
242 )
243 return false; // the ids are bad: for us a track cannot cross twice the same ladder
244 return true;
245 case 3:
246 if (FullSecID(ids[0]).getLayerID() == FullSecID(ids[1]).getLayerID() &&
247 FullSecID(ids[0]).getLadderID() == FullSecID(ids[1]).getLadderID()
248 )
249 return false; // the ids are bad: for us a track cannot cross twice the same ladder
250 if (FullSecID(ids[1]).getLayerID() == FullSecID(ids[2]).getLayerID() &&
251 FullSecID(ids[1]).getLadderID() == FullSecID(ids[2]).getLadderID()
252 )
253 return false; // the ids are bad: for us a track cannot cross twice the same ladder
254 if (FullSecID(ids[0]).getLayerID() == FullSecID(ids[2]).getLayerID() &&
255 FullSecID(ids[0]).getLadderID() == FullSecID(ids[2]).getLadderID()
256 )
257 return false; // the ids are bad: for us a track cannot cross twice the same ladder
258 return true;
259 default:
260 return true;
261 }
262}
263
264
266 std::unique_ptr<TChain>& chain,
267 std::vector<BranchInterface<unsigned>>& sectorBranches,
268 std::vector<BranchInterface<double>>& filterBranches)
269{
270 // prepare everything:
271 unsigned nEntries = chain->GetEntries();
272 unsigned percentMark = 1;
273 if (nEntries > 100) { percentMark = nEntries / 50; }
274 unsigned progressCounter = 0;
275
276 B2INFO("RawSecMapMerger::printData(): start of " << nEntries <<
277 " entries in tree and " << sectorBranches.size() <<
278 "/" << filterBranches.size() <<
279 " sector-/filter-branches");
280
281 for (unsigned i = 0 ; i < nEntries; i++) {
282 if (percentMark > 1 and (i % percentMark) != 0) { continue; }
283 progressCounter += 2;
284 B2INFO("RawSecMapMerger::printData(): entry " << i << " of " << nEntries << ": " << progressCounter << " \% of data printed.");
285
286 auto thisEntry = chain->LoadTree(i);
287
288 string out;
289 for (unsigned k = 0 ; k < sectorBranches.size(); k++) {
290 sectorBranches[k].branch->GetEntry(thisEntry);
291 out += sectorBranches[k].name + ": " + FullSecID(sectorBranches[k].value).getFullSecString() + ". ";
292 }
293 out += "\n";
294
295 for (unsigned k = 0 ; k < filterBranches.size(); k++) {
296 filterBranches[k].branch->GetEntry(thisEntry);
297 out += filterBranches[k].name + ": " + to_string(filterBranches[k].value) + ". ";
298 }
299 B2INFO(out << "\n");
300 }
301}
302
303
304
305
306
308 std::string configName, unsigned int nHitCombinations, bool print2File)
309{
310 SecMapHelper::printStaticSectorRelations<SpacePoint>(filters, configName, nHitCombinations, print2File);
311}
312
313
316//std::vector<VxdID> RawSecMapMergerModule::getCompatibleVxdIDs(const SectorMapConfig& config)
317//{
318//
319// // TODO: remove that part and use the version in the bootstrap module
320//
321// // TODO WARNING hardcoded values!
322// std::vector<unsigned> layers = { 0, 1, 2, 3, 4, 5, 6};
323// std::vector<unsigned> ladders = { 0, 8, 12, 7, 10, 12, 16};
324// std::vector<unsigned> sensors = { 0, 2, 2, 2, 3, 4, 5};
325//
326// std::vector<VxdID> vxdIDs;
327//
328// for (unsigned layerID : config.allowedLayers) {
329// for (unsigned ladderID = 0; ladderID <= ladders.at(layerID); ladderID++) {
330// if (ladderID == 0 and layerID != 0) continue; // only virtual IP (layer 0) has ladder 0
331// for (unsigned sensorID = 0; sensorID <= sensors.at(layerID); sensorID++) {
332// if (sensorID == 0 and layerID != 0) continue; // only virtual IP (layer 0) has sensor 0
333// vxdIDs.push_back(VxdID(layerID, ladderID, sensorID));
334// }
335// }
336// }
337// return vxdIDs;
338//}
339
340
341
342template <class FilterType> unsigned RawSecMapMergerModule::updateFilterSubLayerIDs(SectorGraph<FilterType>& mainGraph,
343 VXDTFFilters<SpacePoint>& segFilters)
344{
345 // get all VXD sensors in the geometry
346 // WARNING: if a different geometry in the first training step was used this may lead to difficulties
347 std::vector<VxdID> vxdIDs = VXD::GeoCache::getInstance().getListOfSensors();
348
349
350 // collect all secIDs occured in training and use them to update the sectors in the SectorID in the VXDTFFilter
351 // in particular the sublayerID which is determined from the graph
352 for (VxdID sensor : vxdIDs) {
353
354 std::vector< FullSecID> allTrainedSecIDsOfSensor = mainGraph.getAllFullSecIDsOfSensor(sensor);
355
356 // this removes all FullSecIDs which occured more than once
357 std::sort(allTrainedSecIDsOfSensor.begin(), allTrainedSecIDsOfSensor.end());
358 allTrainedSecIDsOfSensor.erase(
359 std::unique(
360 allTrainedSecIDsOfSensor.begin(),
361 allTrainedSecIDsOfSensor.end()),
362 allTrainedSecIDsOfSensor.end());
363
364 for (FullSecID sector : allTrainedSecIDsOfSensor) {
365 // the search within that function will ignore the sublayerid, the sublayer id will be set to the one in "sector"
366 bool success = segFilters.setSubLayerIDs(sector, sector.getSubLayerID());
367 // if success is false the sector was not found in the segFilters. This should not happen!
368 if (!success) B2FATAL("There is a mismatch between the FullSecIDs in the Trainings Graph and the SectorMap!");
369 }
370
371 B2DEBUG(20, "Sensor: " << sensor << " had " << allTrainedSecIDsOfSensor.size() << " trained IDs and ");
372 } // end loop sensor of vxdIDs.
373
374 return vxdIDs.size() + 1;
375}
376
377
378
379
380
381// TODO this is not yet capable of dealing with other than twoHitFilters. -> generalize!
382template <class FilterType> void RawSecMapMergerModule::getSegmentFilters(
383 const SectorMapConfig& config,
384 SectorGraph<FilterType>& mainGraph,
385 VXDTFFilters<SpacePoint>* xHitFilters,
386 int nSecChainLength)
387{
388
389 // Thomas : possible bug, the sublayer id s have been updated only for the nSecChainLength==2 case
390 /*
391 if (xHitFilters->size() == 0) {
392 unsigned nSectors = updateFilterSubLayerIDs( mainGraph, *xHitFilters);
393 B2DEBUG(20, "RawSecMapMerger::getSegmentFilters: in updateSubLayerIDs " << nSectors << " were added to secMap " <<
394 config.secMapName);
395 } else {
396 B2DEBUG(20, "RawSecMapMerger::getSegmentFilters: in given xHitFilters-container has size of " << xHitFilters->size() <<
397 " and therefore no further sectors have to be added.");
398 }
399 */
400 // after rewriting this function only updates the sublayer ids of the already existing sectors
401 // so it should only be executed once!!
402 // TODO: remove the if by a better construction!! Also what happens if for the nSecChainLength>2 case the sublayerids need updates???
403 if (nSecChainLength == 2) updateFilterSubLayerIDs(mainGraph, *xHitFilters);
404
405 B2DEBUG(20, "RawSecMapMerger::getSegmentFilters: secMap " << config.secMapName << " got the following sectors:\n" <<
406 mainGraph.print());
407
408
409 for (auto& subGraph : mainGraph) {
410
411 if (nSecChainLength == 2) {
412 add2HitFilters(*xHitFilters, subGraph.second, config);
413 } else if (nSecChainLength == 3) {
414 add3HitFilters(*xHitFilters, subGraph.second, config);
415 } else if (nSecChainLength == 4) {
416 add4HitFilters(*xHitFilters, subGraph.second, config);
417 } else { B2FATAL("nSecChainLength " << nSecChainLength << " is not within allowed range [2;4]!"); }
418 }
419}
420
421
422
423
424
426 filterContainer, SubGraph<FilterType>& subGraph, const SectorMapConfig& config)
427{
428// // // WARNING evil hack -> SelectionVariables themselves should be able to tell their own names!
429// // std::string named3D = "Distance3DSquared", namedXY = "Distance2DXYSquared", nameddZ = "Distance1DZ", namesRZ = "SlopeRZ",
430// // named3Dn = "Distance3DNormed";
431 const auto& filterCutsMap = subGraph.getFinalQuantileValues();
433
434 auto filterNameToFunctions(SelectionVariableNamesToFunctions(
436 std::string filterVals;
437 for (const auto& filterNameToFunction : filterNameToFunctions) {
438 string filterName = filterNameToFunction.first ;
439 filterVals += filterName + ": "
440 + std::to_string(filterCutsMap.at(filterName).getMin())
441 + "/"
442 + std::to_string(filterCutsMap.at(filterName).getMax()) + ", ";
443 }
444 B2DEBUG(20, "SubGraph " << subGraph.getID().print() << " - filter:min/max: " << filterVals);
445
446 VXDTFFilters<SpacePoint>::twoHitFilter_t friendSectorsSegmentFilter =
447 (
448 (
449 (filterCutsMap.at("DistanceInTimeUside").getMin() <= DistanceInTimeUside<SpacePoint>() <=
450 filterCutsMap.at("DistanceInTimeUside").getMax()) &&
451 (filterCutsMap.at("DistanceInTimeVside").getMin() <= DistanceInTimeVside<SpacePoint>() <=
452 filterCutsMap.at("DistanceInTimeVside").getMax()) &&
453 (filterCutsMap.at("Distance3DSquared").getMin() <= Distance3DSquared<SpacePoint>() <=
454 filterCutsMap.at("Distance3DSquared").getMax()) &&
455 (filterCutsMap.at("Distance2DXYSquared").getMin() <= Distance2DXYSquared<SpacePoint>() <=
456 filterCutsMap.at("Distance2DXYSquared").getMax()) &&
457 (filterCutsMap.at("Distance1DZ").getMin() <= Distance1DZ<SpacePoint>() <= filterCutsMap.at("Distance1DZ").getMax()) &&
458 (filterCutsMap.at("SlopeRZ").getMin() <= SlopeRZ<SpacePoint>() <= filterCutsMap.at("SlopeRZ").getMax()) &&
459 (filterCutsMap.at("Distance3DNormed").getMin() <= Distance3DNormed<SpacePoint>() <=
460 filterCutsMap.at("Distance3DNormed").getMax())
461 )
462 );
463
464 auto secIDs = subGraph.getID().getFullSecIDs();
465
466 // secIDs are sorted from outer to inner:
467 B2DEBUG(20, "RawSecMapMerger::add2HitFilters: now adding FriendSectorFilter for secIDs (outer/inner): " << secIDs.at(
468 0) << "/" << secIDs.at(1));
469 if (filterContainer.addTwoHitFilter(secIDs.at(0), secIDs.at(1),
470 friendSectorsSegmentFilter) == 0)
471 B2WARNING("secMap: " << config.secMapName << "Problem adding the friendship relation from the inner sector:" <<
472 secIDs.at(1) << " -> " << secIDs.at(0) << " outer sector");
473}
474
475
476
478 filterContainer, SubGraph<FilterType>& subGraph, const SectorMapConfig& config)
479{
480 const auto& filterCutsMap = subGraph.getFinalQuantileValues();
482
483 auto filterNameToFunctions(SelectionVariableNamesToFunctions(
485 std::string filterVals;
486
487 for (auto& filterNameToFunction : filterNameToFunctions) {
488 string filterName = filterNameToFunction.first ;
489 filterVals += filterName + ": "
490 + std::to_string(filterCutsMap.at(filterName).getMin())
491 + "/"
492 + std::to_string(filterCutsMap.at(filterName).getMax()) + ", ";
493 }
494 B2DEBUG(20, "SubGraph " << subGraph.getID().print() << " - filter:min/max: " << filterVals);
495
496
498 ((filterCutsMap.at("DistanceInTime").getMin() <= DistanceInTime<SpacePoint>() <= filterCutsMap.at("DistanceInTime").getMax()) &&
499 (filterCutsMap.at("Angle3DSimple").getMin() <= Angle3DSimple<SpacePoint>() <= filterCutsMap.at("Angle3DSimple").getMax()) &&
500 (filterCutsMap.at("CosAngleXY").getMin() <= CosAngleXY<SpacePoint>() <= filterCutsMap.at("CosAngleXY").getMax()) &&
501 (filterCutsMap.at("AngleRZSimple").getMin() <= AngleRZSimple<SpacePoint>() <= filterCutsMap.at("AngleRZSimple").getMax()) &&
502 (CircleDist2IP<SpacePoint>() <= filterCutsMap.at("CircleDist2IP").getMax()) &&
503 (filterCutsMap.at("DeltaSlopeRZ").getMin() <= DeltaSlopeRZ<SpacePoint>()) <= filterCutsMap.at("DeltaSlopeRZ").getMax() &&
504 (filterCutsMap.at("DeltaSlopeZoverS").getMin() <= DeltaSlopeZoverS<SpacePoint>() <=
505 filterCutsMap.at("DeltaSlopeZoverS").getMax()) &&
506 (filterCutsMap.at("DeltaSoverZ").getMin() <= DeltaSoverZ<SpacePoint>() <= filterCutsMap.at("DeltaSoverZ").getMax()) &&
507 (filterCutsMap.at("HelixParameterFit").getMin() <= HelixParameterFit<SpacePoint>() <=
508 filterCutsMap.at("HelixParameterFit").getMax()) &&
509 (filterCutsMap.at("Pt").getMin() <= Pt<SpacePoint>() <= filterCutsMap.at("Pt").getMax()) &&
510 (filterCutsMap.at("CircleRadius").getMin() <= CircleRadius<SpacePoint>() <= filterCutsMap.at("CircleRadius").getMax())
511
512 ).observe(VoidObserver());
513
514
515 auto secIDs = subGraph.getID().getFullSecIDs();
516
517 // secIDs are sorted from outer to inner:
518 B2DEBUG(20, "RawSecMapMerger::add3HitFilters: now adding FriendSectorFilter for secIDs (outer/center/inner): "
519 << secIDs.at(0) << "/"
520 << secIDs.at(1) << "/"
521 << secIDs.at(2));
522 if (filterContainer.addThreeHitFilter(secIDs.at(0), secIDs.at(1), secIDs.at(2),
523 threeHitFilter) == 0)
524 B2WARNING("secMap: " << config.secMapName << "Problem adding the friendship relation for the secIDs (outer/center/inner): "
525 << secIDs.at(0) << "/"
526 << secIDs.at(1) << "/"
527 << secIDs.at(2));
528}
529
530
531
532template <class FilterType> void RawSecMapMergerModule::add4HitFilters(
533 VXDTFFilters<SpacePoint>& /*filterContainer*/, SubGraph<FilterType>& /*subGraph*/,
534 const SectorMapConfig& /*config*/)
535{
536
537}
538
539
540
541
542
543
544
545
Class to identify a sector inside of the VXD.
Definition: FullSecID.h:33
std::string getFullSecString() const
returns the FullSecID coded as string compatible to secIDs stored in the xml-sectormaps
Definition: FullSecID.cc:116
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
unsigned updateFilterSubLayerIDs(SectorGraph< FilterType > &mainGraph, VXDTFFilters< SpacePoint > &segFilters)
updates the sublayer ID of the FullSecIDs used in the VXDTFFilters with the one used during the train...
void add3HitFilters(VXDTFFilters< SpacePoint > &filterContainer, SubGraph< FilterType > &subGraph, const SectorMapConfig &config)
WARNING TODO clean up and documentation!
SectorGraph< FilterType > buildGraph(std::unique_ptr< TChain > &chain, std::vector< BranchInterface< unsigned > > &sectorBranches, std::vector< BranchInterface< double > > &filterBranches)
build graph with secChains found in TChain.
bool m_PARAMprintFullGraphs
If true, the full trained graphs will be printed to screen.
std::vector< BranchInterface< ValueType > > getBranches(std::unique_ptr< TChain > &chain, const std::vector< std::string > &branchNames)
for given chain and names of branches: this function returns their pointers to the branch and the con...
void printVXDTFFilters(const VXDTFFilters< SpacePoint > &filters, std::string configName, unsigned int nHitCombinations, bool print2File)
for debugging purposes: print VXDTFFilters into a file of name of the sectorMapConfig.
std::unique_ptr< TChain > createTreeChain(const SectorMapConfig &configuration, const std::string &nHitString)
bundle all relevant files to a TChain
std::vector< std::string > m_PARAMrootFileNames
///////////////////////////////////////////////////////////////////////////////// member variables of...
std::vector< std::string > m_PARAMmapNames
contains names of sectorMaps to be loaded.
void printData(std::unique_ptr< TChain > &chain, std::vector< BranchInterface< unsigned > > &sectorBranches, std::vector< BranchInterface< double > > &filterBranches)
for debugging: print data for crosschecks.
std::vector< std::string > getRootFiles(std::string mapName)
returns all names of root-files fitting given parameter mapName
int m_RelThreshold
Relative threshold for pruning the sector maps (in %).
std::vector< unsigned > getSecIDs(std::vector< BranchInterface< unsigned > > &secBranches, Long64_t entry)
returns secIDs of current entry in the secBranches.
void add4HitFilters(VXDTFFilters< SpacePoint > &filterContainer, SubGraph< FilterType > &subGraph, const SectorMapConfig &config)
WARNING TODO clean up and documentation!
RawSecMapMergerModule()
Constructor of the module.
void add2HitFilters(VXDTFFilters< SpacePoint > &filterContainer, SubGraph< FilterType > &subGraph, const SectorMapConfig &config)
WARNING TODO clean up and documentation!
void getSegmentFilters(const SectorMapConfig &config, SectorGraph< FilterType > &mainGraph, VXDTFFilters< SpacePoint > *xHitFilters, int nSecChainLength)
returns all VxdIDs (sensors) compatible with given configData.
void trainGraph(SectorGraph< FilterType > &mainGraph, std::unique_ptr< TChain > &chain, std::vector< BranchInterface< unsigned > > &sectorBranches, std::vector< BranchInterface< double > > &filterBranches)
fill the graphs with raw data fitting to their filters respectively.
std::string prepareNHitSpecificStuff(unsigned nHits, const SectorMapConfig &config, std::vector< std::string > &secBranchNames, std::vector< std::string > &filterBranchNames)
sets everything which is hit-dependent.
bool good(const std::vector< unsigned > &ids)
check that the vector of FullSecIDs
contains all subgraphs.
Definition: SectorGraph.h:31
unsigned size() const
returns number of collected subgraphs so far.
Definition: SectorGraph.h:55
std::vector< FullSecID > getAllFullSecIDsOfSensor(VxdID sensor)
returns a Vector containing all FullSecIDs found for given sensor.
Definition: SectorGraph.h:300
Iterator find(SubGraphID idChain)
find entry.
Definition: SectorGraph.h:46
Iterator add(SubGraphID &newID)
add new subgraph if not added already.
Definition: SectorGraph.h:66
std::string print(bool fullPrint=true) const
returns a string giving an overview of the graph.
Definition: SectorGraph.h:76
Iterator end()
returns end of subGraphs.
Definition: SectorGraph.h:52
stores the ID of a subgraph, which is basically a chain of FullSecID coded as unsigned ints.
Definition: SubGraphID.h:24
std::vector< FullSecID > getFullSecIDs() const
returns SecIDs coded as FullSecIDs:
Definition: SubGraphID.h:185
std::string print() const
returns string of entries.
Definition: SubGraphID.h:101
contains all relevant stuff needed for dealing with a subGraph.
Definition: SubGraph.h:30
const std::unordered_map< FilterType, MinMax > & getFinalQuantileValues()
this deletes the old min and max values stored and replaces them with the quantiles to be found.
Definition: SubGraph.h:120
const SubGraphID & getID() const
returns iD of this graph
Definition: SubGraph.h:93
Class that contains all the static sectors to which the filters are attached.
Definition: VXDTFFilters.h:63
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.)) twoHitFilter_t
minimal working 2-hits-example used for redesign of VXDTF.
Definition: VXDTFFilters.h:82
bool setSubLayerIDs(FullSecID sector, int sublayer)
during the trainings phase the sublayer ids have to be updated
Definition: VXDTFFilters.h:329
int addThreeHitFilter(FullSecID outer, FullSecID center, FullSecID inner, const threeHitFilter_t &filter)
adds a three hit filter
Definition: VXDTFFilters.h:189
int addTwoHitFilter(FullSecID outer, FullSecID inner, const twoHitFilter_t &filter)
adds a two hit filter
Definition: VXDTFFilters.h:175
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.)) threeHitFilter_t
minimal working example for 3-hits:
Definition: VXDTFFilters.h:112
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:59
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
The most CPU efficient Observer for the VXDTF filter tools (even if useless).
Definition: VoidObserver.h:30
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
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
std::unordered_map< std::string, typename Variable::functionType > SelectionVariableNamesToFunctions(Belle2::Filter< Variable, Range, Options... >)
Return a map from the SelectionVariable name to the SelectionVariable function of the Variable used i...
Abstract base class for different kinds of events.
STL namespace.
simple struct for interfacing the Branch.
simple struct containing all the configuration data needed for the SecMapTrainer.