Belle II Software  release-05-02-19
NeuroTrigger.cc
1 #include <framework/logging/Logger.h>
2 #include <trg/cdc/NeuroTrigger.h>
3 
4 #include <cdc/geometry/CDCGeometryPar.h>
5 #include <framework/gearbox/Const.h>
6 #include <framework/gearbox/Unit.h>
7 #include <framework/datastore/StoreObjPtr.h>
8 #include <framework/datastore/StoreArray.h>
9 #include <trg/cdc/dbobjects/CDCTriggerNeuroConfig.h>
10 #include <trg/cdc/dataobjects/CDCTriggerTrack.h>
11 #include <string>
12 #include <cmath>
13 #include <TFile.h>
14 
15 using namespace Belle2;
16 using namespace CDC;
17 using namespace std;
18 
19 void
21 {
22  // check parameters
23  bool okay = true;
24  // ensure that length of lists matches number of sectors
25  if (p.nHidden.size() != 1 && p.nHidden.size() != p.nMLP) {
26  B2ERROR("Number of nHidden lists should be 1 or " << p.nMLP);
27  okay = false;
28  }
29  if (p.outputScale.size() != 1 && p.outputScale.size() != p.nMLP) {
30  B2ERROR("Number of outputScale lists should be 1 or " << p.nMLP);
31  okay = false;
32  }
33  bool rangeProduct = (p.phiRange.size() * p.invptRange.size() * p.thetaRange.size() * p.SLpattern.size() == p.nMLP);
34  if (!rangeProduct) {
35  if (p.phiRange.size() != 1 && p.phiRange.size() != p.nMLP) {
36  B2ERROR("Number of phiRange lists should be 1 or " << p.nMLP);
37  okay = false;
38  }
39  if (p.invptRange.size() != 1 && p.invptRange.size() != p.nMLP) {
40  B2ERROR("Number of invptRange lists should be 1 or " << p.nMLP);
41  okay = false;
42  }
43  if (p.thetaRange.size() != 1 && p.thetaRange.size() != p.nMLP) {
44  B2ERROR("Number of thetaRange lists should be 1 or " << p.nMLP);
45  okay = false;
46  }
47  if (p.SLpattern.size() != 1 && p.SLpattern.size() != p.nMLP) {
48  B2ERROR("Number of SLpattern lists should be 1 or " << p.nMLP);
49  okay = false;
50  }
51  }
52  // ensure that length of maxHitsPerSL and SLpatternMask lists matches SLpattern
53  if (p.maxHitsPerSL.size() != 1 && p.maxHitsPerSL.size() != p.SLpattern.size()) {
54  B2ERROR("Number of maxHitsPerSL lists should be 1 or " << p.SLpattern.size());
55  okay = false;
56  }
57  if (p.SLpatternMask.size() != 1 && p.SLpatternMask.size() != p.SLpattern.size()) {
58  B2ERROR("Number of SLpatternMask lists should be 1 or " << p.SLpattern.size());
59  okay = false;
60  }
61  // ensure that number of target nodes is valid
62  unsigned short nTarget = int(p.targetZ) + int(p.targetTheta);
63  if (nTarget < 1) {
64  B2ERROR("No outputs! Turn on either targetZ or targetTheta.");
65  okay = false;
66  }
67  // ensure that sector ranges are valid
68  for (unsigned iPhi = 0; iPhi < p.phiRange.size(); ++iPhi) {
69  if (p.phiRange[iPhi].size() != 2) {
70  B2ERROR("phiRange should be exactly 2 values");
71  okay = false;
72  continue;
73  }
74  if (p.phiRange[iPhi][0] >= p.phiRange[iPhi][1]) {
75  B2ERROR("phiRange[0] should be smaller than phiRange[1]");
76  okay = false;
77  }
78  if (p.phiRange[iPhi][0] < -360. || p.phiRange[iPhi][1] > 360. ||
79  (p.phiRange[iPhi][1] - p.phiRange[iPhi][0]) > 360.) {
80  B2ERROR("phiRange should be in [-360, 360], with maximal width of 360");
81  okay = false;
82  }
83  }
84  for (unsigned iPt = 0; iPt < p.invptRange.size(); ++iPt) {
85  if (p.invptRange[iPt].size() != 2) {
86  B2ERROR("invptRange should be exactly 2 values");
87  okay = false;
88  }
89  if (p.invptRange[iPt][0] >= p.invptRange[iPt][1]) {
90  B2ERROR("invptRange[0] should be smaller than invptRange[1]");
91  okay = false;
92  }
93  }
94  for (unsigned iTheta = 0; iTheta < p.thetaRange.size(); ++iTheta) {
95  if (p.thetaRange[iTheta].size() != 2) {
96  B2ERROR("thetaRange should be exactly 2 values");
97  okay = false;
98  continue;
99  }
100  if (p.thetaRange[iTheta][0] >= p.thetaRange[iTheta][1]) {
101  B2ERROR("thetaRange[0] should be smaller than thetaRange[1]");
102  okay = false;
103  }
104  if (p.thetaRange[iTheta][0] < 0. || p.thetaRange[iTheta][1] > 180.) {
105  B2ERROR("thetaRange should be in [0, 180]");
106  okay = false;
107  }
108  }
109  for (unsigned iScale = 0; iScale < p.outputScale.size(); ++iScale) {
110  if (p.outputScale[iScale].size() != 2 * nTarget) {
111  B2ERROR("outputScale should be exactly " << 2 * nTarget << " values");
112  okay = false;
113  }
114  }
115  // ensure that train sectors are valid
116  if (p.phiRange.size() != p.phiRangeTrain.size()) {
117  B2ERROR("Number of phiRange lists and phiRangeTrain lists should be equal.");
118  okay = false;
119  } else {
120  for (unsigned iPhi = 0; iPhi < p.phiRange.size(); ++iPhi) {
121  if (p.phiRangeTrain[iPhi].size() != 2) {
122  B2ERROR("phiRangeTrain should be exactly 2 values.");
123  okay = false;
124  } else if (p.phiRangeTrain[iPhi][0] > p.phiRange[iPhi][0] ||
125  p.phiRangeTrain[iPhi][1] < p.phiRange[iPhi][1]) {
126  B2ERROR("phiRangeTrain should be wider than phiRange or equal.");
127  okay = false;
128  }
129  }
130  }
131  if (p.invptRange.size() != p.invptRangeTrain.size()) {
132  B2ERROR("Number of invptRange lists and invptRangeTrain lists should be equal.");
133  okay = false;
134  } else {
135  for (unsigned iPt = 0; iPt < p.invptRange.size(); ++iPt) {
136  if (p.invptRangeTrain[iPt].size() != 2) {
137  B2ERROR("invptRangeTrain should be exactly 2 values.");
138  okay = false;
139  } else if (p.invptRangeTrain[iPt][0] > p.invptRange[iPt][0] ||
140  p.invptRangeTrain[iPt][1] < p.invptRange[iPt][1]) {
141  B2ERROR("invptRangeTrain should be wider than invptRange or equal.");
142  okay = false;
143  }
144  }
145  }
146  if (p.thetaRange.size() != p.thetaRangeTrain.size()) {
147  B2ERROR("Number of thetaRange lists and thetaRangeTrain lists should be equal.");
148  okay = false;
149  } else {
150  for (unsigned iTheta = 0; iTheta < p.thetaRange.size(); ++iTheta) {
151  if (p.thetaRangeTrain[iTheta].size() != 2) {
152  B2ERROR("thetaRangeTrain should be exactly 2 values.");
153  okay = false;
154  } else if (p.thetaRangeTrain[iTheta][0] > p.thetaRange[iTheta][0] ||
155  p.thetaRangeTrain[iTheta][1] < p.thetaRange[iTheta][1]) {
156  B2ERROR("thetaRangeTrain should be wider than thetaRange or equal.");
157  okay = false;
158  }
159  }
160  }
161 
162  if (!okay) return;
163 
164  // initialize MLPs
165  m_MLPs.clear();
166  for (unsigned iMLP = 0; iMLP < p.nMLP; ++iMLP) {
167  //get indices for sector parameters
168  vector<unsigned> indices = getRangeIndices(p, iMLP);
169  //get number of nodes for each layer
170  unsigned short maxHits = (p.maxHitsPerSL.size() == 1) ? p.maxHitsPerSL[0] : p.maxHitsPerSL[indices[3]];
171  unsigned long SLpattern = p.SLpattern[indices[3]];
172  unsigned long SLpatternMask = (p.SLpatternMask.size() == 1) ? p.SLpatternMask[0] : p.SLpatternMask[indices[3]];
173  unsigned short nInput = 27 * maxHits;
174  vector<float> nHidden = (p.nHidden.size() == 1) ? p.nHidden[0] : p.nHidden[iMLP];
175  vector<unsigned short> nNodes = {nInput};
176  for (unsigned iHid = 0; iHid < nHidden.size(); ++iHid) {
177  if (p.multiplyHidden) {
178  nNodes.push_back(nHidden[iHid] * nNodes[0]);
179  } else {
180  nNodes.push_back(nHidden[iHid]);
181  }
182  }
183  nNodes.push_back(nTarget);
184  unsigned short targetVars = int(p.targetZ) + (int(p.targetTheta) << 1);
185  //get sector ranges (initially train ranges)
186  vector<float> phiRange = p.phiRangeTrain[indices[0]];
187  vector<float> invptRange = p.invptRangeTrain[indices[1]];
188  vector<float> thetaRange = p.thetaRangeTrain[indices[2]];
189  B2DEBUG(50, "Ranges for sector " << iMLP
190  << ": phiRange [" << phiRange[0] << ", " << phiRange[1]
191  << "], invptRange [" << invptRange[0] << ", " << invptRange[1]
192  << "], thetaRange [" << thetaRange[0] << ", " << thetaRange[1]
193  << "], SLpattern " << SLpattern);
194  //get scaling values
195  vector<float> outputScale = (p.outputScale.size() == 1) ? p.outputScale[0] : p.outputScale[iMLP];
196  //convert phi and theta from degree to radian
197  phiRange[0] *= Unit::deg;
198  phiRange[1] *= Unit::deg;
199  thetaRange[0] *= Unit::deg;
200  thetaRange[1] *= Unit::deg;
201  if (p.targetTheta) {
202  outputScale[2 * int(p.targetZ)] *= Unit::deg;
203  outputScale[2 * int(p.targetZ) + 1] *= Unit::deg;
204  }
205  //create new MLP
206  m_MLPs.push_back(CDCTriggerMLP(nNodes, targetVars, outputScale,
207  phiRange, invptRange, thetaRange,
208  maxHits, SLpattern, SLpatternMask, p.tMax, p.T0fromHits,
209  p.et_option));
210  }
211  // load some values from the geometry that will be needed for the input
212  setConstants();
213 }
214 
215 vector<unsigned>
216 NeuroTrigger::getRangeIndices(const Parameters& p, unsigned isector)
217 {
218  std::vector<unsigned> indices = {0, 0, 0, 0};
219  if (p.phiRange.size() * p.invptRange.size() * p.thetaRange.size() * p.SLpattern.size() == p.nMLP) {
220  indices[0] = isector % p.phiRange.size();
221  indices[1] = (isector / p.phiRange.size()) % p.invptRange.size();
222  indices[2] = (isector / (p.phiRange.size() * p.invptRange.size())) % p.thetaRange.size();
223  indices[3] = isector / (p.phiRange.size() * p.invptRange.size() * p.thetaRange.size());
224  } else {
225  indices[0] = (p.phiRange.size() == 1) ? 0 : isector;
226  indices[1] = (p.invptRange.size() == 1) ? 0 : isector;
227  indices[2] = (p.thetaRange.size() == 1) ? 0 : isector;
228  indices[3] = (p.SLpattern.size() == 1) ? 0 : isector;
229  }
230  return indices;
231 }
232 
233 void
235 {
237  int layerId = 3;
238  int nTS = 0;
239  for (int iSL = 0; iSL < 9; ++iSL) {
240  m_TSoffset[iSL] = nTS;
241  nTS += cdc.nWiresInLayer(layerId);
242  m_TSoffset[iSL + 1] = nTS;
243  for (int priority = 0; priority < 2; ++ priority) {
244  m_radius[iSL][priority] = cdc.senseWireR(layerId + priority);
245  }
246  layerId += (iSL > 0 ? 6 : 7);
247  }
248 }
249 
250 void
251 NeuroTrigger::initializeCollections(string hitCollectionName, string eventTimeName, std::string et_option)
252 {
253  m_segmentHits.isRequired(hitCollectionName);
254  if (!((et_option == "fastestpriority") || (et_option == "zero") || (et_option == "fastest2d"))) {
255  m_eventTime.isRequired(eventTimeName);
256  }
257  m_hitCollectionName = hitCollectionName;
258 }
259 
260 vector<int>
261 NeuroTrigger::selectMLPs(float phi0, float invpt, float theta)
262 {
263  vector<int> indices = {};
264 
265  if (m_MLPs.size() == 0) {
266  B2WARNING("Trying to select MLP before initializing MLPs.");
267  return indices;
268  }
269 
270  // find all matching sectors
271  for (unsigned isector = 0; isector < m_MLPs.size(); ++isector) {
272  if (m_MLPs[isector].inPhiRange(phi0) && m_MLPs[isector].inInvptRange(invpt)
273  && m_MLPs[isector].inThetaRange(theta)) {
274  indices.push_back(isector);
275  }
276  }
277 
278  if (indices.size() == 0) {
279  B2DEBUG(150, "Track does not match any sector.");
280  B2DEBUG(150, "invpt=" << invpt << ", phi=" << phi0 * 180. / M_PI << ", theta=" << theta * 180. / M_PI);
281  }
282 
283  return indices;
284 }
285 
286 
287 int
288 NeuroTrigger::selectMLPbyPattern(std::vector<int>& MLPs, unsigned long pattern, const bool neurotrackinputmode)
289 {
290  if (MLPs.size() == 0) {
291  return -1;
292  }
293 
294  int bestIndex = -1;
295  for (unsigned i = 0; i < MLPs.size(); ++i) {
296  int isector = MLPs[i];
297  unsigned long sectorPattern = m_MLPs[isector].getSLpattern();
298  B2DEBUG(250, "hitPattern " << pattern << " sectorPattern " << sectorPattern);
299  // no hit pattern restriction -> keep looking for exact match
300  if (sectorPattern == 0) {
301  B2DEBUG(250, "found match for general sector");
302  bestIndex = isector;
303  }
304  // exact match -> keep this sector and stop searching
305  if (pattern == sectorPattern) {
306  B2DEBUG(250, "found match for hit pattern " << pattern);
307  bestIndex = isector;
308  break;
309  }
310  }
311 
312  if (bestIndex < 0) {
313  if (neurotrackinputmode) {
314  B2DEBUG(150, "No sector found to match pattern, using sector 0" << pattern << ".");
315  bestIndex = 0;
316  } else {
317  B2DEBUG(150, "No sector found to match pattern " << pattern << ".");
318  }
319  }
320  return bestIndex;
321 }
322 
323 
324 void
326 {
327  double omega = track.getOmega(); // signed track curvature
328  for (int iSL = 0; iSL < 9; ++iSL) {
329  for (int priority = 0; priority < 2; ++priority) {
330  m_alpha[iSL][priority] = (m_radius[iSL][priority] * abs(omega) < 2.) ?
331  asin(m_radius[iSL][priority] * omega / 2.) :
332  M_PI_2;
333  m_idRef[iSL][priority] = remainder(((track.getPhi0() - m_alpha[iSL][priority]) *
334  (m_TSoffset[iSL + 1] - m_TSoffset[iSL]) / 2. / M_PI),
335  (m_TSoffset[iSL + 1] - m_TSoffset[iSL]));
336  }
337  }
338 }
339 
340 void
342 {
343  unsigned precisionPhi = m_precision[0];
344  unsigned precisionAlpha = m_precision[0];
345  unsigned precisionScale = m_precision[1];
346  unsigned precisionId = m_precision[2];
347 
348  double omega = track.getOmega();
349  long phi = round(track.getPhi0() * (1 << precisionPhi));
350 
351  for (int iSL = 0; iSL < 9; ++iSL) {
352  for (int priority = 0; priority < 2; ++priority) {
353  // LUT, calculated on the fly here
354  m_alpha[iSL][priority] =
355  (m_radius[iSL][priority] * abs(omega) < 2) ?
356  round(asin(m_radius[iSL][priority] * omega / 2) * (1 << precisionAlpha)) :
357  round(M_PI_2 * (1 << precisionAlpha));
358  long dphi = (precisionAlpha >= precisionPhi) ?
359  phi - (long(m_alpha[iSL][priority]) >> (precisionAlpha - precisionPhi)) :
360  phi - (long(m_alpha[iSL][priority]) << (precisionPhi - precisionAlpha));
361  m_idRef[iSL][priority] =
362  long(dphi * round((m_TSoffset[iSL + 1] - m_TSoffset[iSL]) / 2. / M_PI * (1 << precisionScale)) /
363  (long(1) << (precisionPhi + precisionScale - precisionId)));
364  // unscale after rounding
365  m_alpha[iSL][priority] /= (1 << precisionAlpha);
366  m_idRef[iSL][priority] /= (1 << precisionId);
367  }
368  }
369 }
370 
371 double
373 {
374  int iSL = hit.getISuperLayer();
375  int priority = hit.getPriorityPosition();
376  // (((priority >> 1) & 1) - (priority & 1)) is 0, -1, 1, 0 for priority = 0, 1, 2, 3
377  double relId = hit.getSegmentID() + 0.5 * (((priority >> 1) & 1) - (priority & 1))
378  - m_TSoffset[iSL] - m_idRef[iSL][int(priority < 3)];
379  relId = remainder(relId, (m_TSoffset[iSL + 1] - m_TSoffset[iSL]));
380  return relId;
381 }
382 
383 void
384 NeuroTrigger::getEventTime(unsigned isector, const CDCTriggerTrack& track, std::string et_option, const bool neuroinputmode = false)
385 {
386  if (et_option != m_MLPs[isector].get_et_option()) {
387  B2WARNING("Used event time option is different to the one set in the MLP"
388  << LogVar("et_option", et_option) << LogVar("isector", isector)
389  << LogVar("et_option_mlp", m_MLPs[isector].get_et_option()));
390  }
391  if (et_option == "etf_or_fastestpriority") {
392  bool hasT0 = m_eventTime->hasBinnedEventT0(Const::CDC);
393  if (hasT0) {
394  m_T0 = m_eventTime->getBinnedEventT0(Const::CDC);
395  m_hasT0 = true;
396  } else {
397  m_T0 = 9999;
398  // find shortest time of related and relevant axial hits
400  track.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
401  for (unsigned ihit = 0; ihit < axialHits.size(); ++ihit) {
402  // skip hits with negative relation weight (not selected in finder)
403  if (axialHits.weight(ihit) < 0) continue;
404  unsigned short iSL = axialHits[ihit]->getISuperLayer();
405  // skip stereo hits (should not be related to track, but check anyway)
406  if ((iSL % 2 == 1) && !neuroinputmode) continue;
407  // get shortest time of relevant hits
408  double relId = getRelId(*axialHits[ihit]);
409  if (m_MLPs[isector].isRelevant(relId, iSL) &&
410  axialHits[ihit]->priorityTime() < m_T0) {
411  m_T0 = axialHits[ihit]->priorityTime();
412  }
413  }
414  if (!neuroinputmode) {
415  // find shortest time of relevant stereo hits
416  StoreArray<CDCTriggerSegmentHit> hits(m_hitCollectionName);
417  for (int ihit = 0; ihit < hits.getEntries(); ++ihit) {
418  unsigned short iSL = hits[ihit]->getISuperLayer();
419  // skip axial hits
420  if (iSL % 2 == 0) continue;
421  // get shortest time of relevant hits
422  double relId = getRelId(*hits[ihit]);
423  if (m_MLPs[isector].isRelevant(relId, iSL) && hits[ihit]->priorityTime() < m_T0) {
424  m_T0 = hits[ihit]->priorityTime();
425  }
426  }
427  }
428  if (m_T0 < 9999) {
429  m_hasT0 = true;
430  } else {
431  m_T0 = 0;
432  m_hasT0 = false;
433  }
434  }
435  } else if (et_option == "etf_or_fastest2d") {
436  bool hasT0 = m_eventTime->hasBinnedEventT0(Const::CDC);
437  if (hasT0) {
438  m_T0 = m_eventTime->getBinnedEventT0(Const::CDC);
439  m_hasT0 = true;
440  } else {
441  m_T0 = 9999;
442  // find shortest time of related and relevant axial hits
444  track.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
445  for (unsigned ihit = 0; ihit < axialHits.size(); ++ihit) {
446  // skip hits with negative relation weight (not selected in finder)
447  if (axialHits.weight(ihit) < 0) continue;
448  unsigned short iSL = axialHits[ihit]->getISuperLayer();
449  // skip stereo hits (should not be related to track, but check anyway)
450  if (iSL % 2 == 1) continue;
451  // get shortest time of relevant hits
452  double relId = getRelId(*axialHits[ihit]);
453  if (m_MLPs[isector].isRelevant(relId, iSL) &&
454  axialHits[ihit]->priorityTime() < m_T0) {
455  m_T0 = axialHits[ihit]->priorityTime();
456  }
457  }
458  if (m_T0 < 9999) {
459  m_hasT0 = true;
460  } else {
461  m_T0 = 0;
462  m_hasT0 = false;
463  }
464  }
465  } else if (et_option == "fastestpriority") {
466  B2DEBUG(200, "et_option is 'fastestpriority'");
467  m_T0 = 9999;
468  // find shortest time of related and relevant axial hits
470  track.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
471  B2DEBUG(200, "looping over axials:");
472  for (unsigned ihit = 0; ihit < axialHits.size(); ++ihit) {
473  // skip hits with negative relation weight (not selected in finder)
474  if (axialHits.weight(ihit) < 0) continue;
475  unsigned short iSL = axialHits[ihit]->getISuperLayer();
476  // skip stereo hits (should not be related to track, but check anyway)
477  if ((iSL % 2 == 1) && !neuroinputmode) continue;
478  // get shortest time of relevant hits
479  B2DEBUG(200, " check drifttime: SL" + std::to_string(iSL) + ",ID = " + std::to_string(axialHits[ihit]->getSegmentID()) + ", t = " +
480  std::to_string(axialHits[ihit]->priorityTime()));
481  double relId = getRelId(*axialHits[ihit]);
482  if (m_MLPs[isector].isRelevant(relId, iSL) &&
483  axialHits[ihit]->priorityTime() < m_T0) {
484  m_T0 = axialHits[ihit]->priorityTime();
485  B2DEBUG(200, " new t0: " << std::to_string(m_T0));
486  }
487  }
488  if (!neuroinputmode) {
489  // find shortest time of relevant stereo hits
490  StoreArray<CDCTriggerSegmentHit> hits(m_hitCollectionName);
491  B2DEBUG(200, "looping over stereos:");
492  for (int ihit = 0; ihit < hits.getEntries(); ++ihit) {
493  unsigned short iSL = hits[ihit]->getISuperLayer();
494  // skip axial hits
495  if (iSL % 2 == 0) continue;
496  // get shortest time of relevant hits
497  B2DEBUG(200, " check drifttime: SL" + std::to_string(iSL) + ",ID = " + std::to_string(hits[ihit]->getSegmentID()) + ", t = " +
498  std::to_string(hits[ihit]->priorityTime()));
499  double relId = getRelId(*hits[ihit]);
500  if (m_MLPs[isector].isRelevant(relId, iSL) && hits[ihit]->priorityTime() < m_T0) {
501  m_T0 = hits[ihit]->priorityTime();
502  B2DEBUG(200, " new t0: " << std::to_string(m_T0));
503  }
504  }
505  }
506  if (m_T0 < 9999) {
507  m_hasT0 = true;
508  } else {
509  m_T0 = 0;
510  m_hasT0 = false;
511  }
512 
513  } else if (et_option == "fastest2d") {
514  m_T0 = 9999;
515  // find shortest time of related and relevant axial hits
517  track.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
518  for (unsigned ihit = 0; ihit < axialHits.size(); ++ihit) {
519  // skip hits with negative relation weight (not selected in finder)
520  if (axialHits.weight(ihit) < 0) continue;
521  unsigned short iSL = axialHits[ihit]->getISuperLayer();
522  // skip stereo hits (should not be related to track, but check anyway)
523  if (iSL % 2 == 1) continue;
524  // get shortest time of relevant hits
525  double relId = getRelId(*axialHits[ihit]);
526  if (m_MLPs[isector].isRelevant(relId, iSL) &&
527  axialHits[ihit]->priorityTime() < m_T0) {
528  m_T0 = axialHits[ihit]->priorityTime();
529  }
530  }
531  if (m_T0 < 9999) {
532  m_hasT0 = true;
533  } else {
534  m_T0 = 0;
535  m_hasT0 = false;
536  }
537  } else if (et_option == "zero") {
538  m_hasT0 = true;
539  m_T0 = 0;
540  } else if (et_option == "etf_only") {
541  bool hasT0 = m_eventTime->hasBinnedEventT0(Const::CDC);
542  if (hasT0) {
543  m_T0 = m_eventTime->getBinnedEventT0(Const::CDC);
544  m_hasT0 = true;
545  } else {
546  B2ERROR("No ETF output, but forced to use ETF!");
547  m_hasT0 = false;
548  m_T0 = 0;
549  }
550  } else if (et_option == "etf_or_zero") {
551  bool hasT0 = m_eventTime->hasBinnedEventT0(Const::CDC);
552  if (hasT0) {
553  m_T0 = m_eventTime->getBinnedEventT0(Const::CDC);
554  m_hasT0 = true;
555  } else {
556  m_hasT0 = true;
557  m_T0 = 0;
558  }
559  } else {
560  B2ERROR("No valid parameter for et_option (" << et_option << " )!");
561  }
562 
563 }
564 void
565 NeuroTrigger::getEventTime(unsigned isector, const CDCTriggerTrack& track)
566 {
567  B2WARNING("This function is deprecated, it used the flag 'T0fromHits'. "
568  << "Give additionally et_option as argument to use the new function. "
569  << "et_option is now set for T0fromHits=True to 'etf_fastestpriority' "
570  << "and for T0fromHits=false to 'etf_or_zero'.");
571  if (m_MLPs[isector].getT0fromHits()) {
572  getEventTime(isector, track, "etf_or_fastestpriority");
573  } else {
574  getEventTime(isector, track, "etf_or_zero");
575  }
576 }
577 
578 
579 
580 unsigned long
581 NeuroTrigger::getInputPattern(unsigned isector, const CDCTriggerTrack& track, const bool neurotrackinputmode)
582 {
583  CDCTriggerMLP& expert = m_MLPs[isector];
584  unsigned long hitPattern = 0;
585  vector<unsigned> nHits;
586  nHits.assign(9, 0);
587  // loop over axial hits related to input track
589  track.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
590  for (unsigned ihit = 0; ihit < axialHits.size(); ++ ihit) {
591  // skip hits with negative relation weight (not selected in finder)
592  if (axialHits.weight(ihit) < 0) {
593  continue;
594  }
595  unsigned short iSL = axialHits[ihit]->getISuperLayer();
596  // // skip stereo hits (should not be related to track, but check anyway)
597  if ((!neurotrackinputmode) && (iSL % 2 == 1)) continue;
598  // get priority time
599  int t = (m_hasT0) ? axialHits[ihit]->priorityTime() - m_T0 : 0;
600  if (t < 0) {
601  t = 0;
602  } else if (t > expert.getTMax()) {
603  t = expert.getTMax();
604  }
605  double relId = getRelId(*axialHits[ihit]);
606 
607  if (expert.isRelevant(relId, iSL)) {
608  if (nHits[iSL] < expert.getMaxHitsPerSL()) {
609  hitPattern |= 1 << (iSL + 9 * nHits[iSL]);
610  ++nHits[iSL];
611  }
612  B2DEBUG(250, "hit in SL " << iSL);
613  } else {
614  B2DEBUG(250, "hit in SL " << iSL << " not relevant (relId = " << relId << ")");
615  }
616  }
617  if (!neurotrackinputmode) {
618  // loop over stereo hits
619  for (int ihit = 0; ihit < m_segmentHits.getEntries(); ++ ihit) {
620  unsigned short iSL = m_segmentHits[ihit]->getISuperLayer();
621  // skip axial hits
622  if (iSL % 2 == 0) continue;
623  // get priority time
624  int t = (m_hasT0) ? m_segmentHits[ihit]->priorityTime() - m_T0 : 0;
625  if (t < 0) {
626  t = 0;
627  } else if (t > expert.getTMax()) {
628  t = expert.getTMax();
629  }
630  double relId = getRelId(*m_segmentHits[ihit]);
631  if (expert.isRelevant(relId, iSL)) {
632  if (nHits[iSL] < expert.getMaxHitsPerSL()) {
633  hitPattern |= 1 << (iSL + 9 * nHits[iSL]);
634  ++nHits[iSL];
635  }
636  B2DEBUG(250, "hit in SL " << iSL);
637  } else {
638  B2DEBUG(250, "hit in SL " << iSL << " not relevant (relId = " << relId << ")");
639  }
640  }
641  }
642  B2DEBUG(250, "hitPattern " << hitPattern);
643  return hitPattern & expert.getSLpatternMask();
644 }
645 
646 vector<unsigned>
647 NeuroTrigger::selectHitsHWSim(unsigned isector, const CDCTriggerTrack& track)
648 {
649  CDCTriggerMLP& expert = m_MLPs[isector];
650  vector<unsigned> selectedHitIds = {};
651  // prepare vectors to keep best drift times, left/right and selected hit IDs
652  vector<int> tMin;
653  tMin.assign(expert.nNodesLayer(0), expert.getTMax());
654  vector<bool> LRknown;
655  LRknown.assign(expert.nNodesLayer(0), false);
656  vector<int> hitIds;
657  hitIds.assign(expert.nNodesLayer(0), -1);
658  vector<unsigned> nHits;
659  nHits.assign(9, 0);
660 
661  // loop over all hits related to input track
663  track.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
664  B2DEBUG(250, "start hit loop over all related hits");
665  for (unsigned ihit = 0; ihit < allHits.size(); ++ihit) {
666  // skip hits with negative relation weight (not selected in finder)
667  if (allHits.weight(ihit) < 0) continue;
668  unsigned short iSL = allHits[ihit]->getISuperLayer();
669  if (expert.getSLpatternUnmasked() != 0 &&
670  !((expert.getSLpatternUnmasked() >> iSL) & 1)) {
671  B2DEBUG(250, "skipping hit in SL " << iSL);
672  continue;
673  }
674  // get priority time and apply time window cut
675  int t = (m_hasT0) ? allHits[ihit]->priorityTime() - m_T0 : 0;
676  if (t < 0) {
677  t = 0;
678  } else if (t > expert.getTMax()) {
679  t = expert.getTMax();
680  }
681  double relId = getRelId(*allHits[ihit]);
682  if (expert.isRelevant(relId, iSL)) {
683  // get reference hit (worst of existing hits)
684  unsigned short iRef = iSL;
685  if (expert.getMaxHitsPerSL() > 1) {
686  if (nHits[iSL] < expert.getMaxHitsPerSL() &&
687  (expert.getSLpatternUnmasked() >> (iSL + 9 * nHits[iSL])) & 1) {
688  iRef += 9 * nHits[iSL];
689  ++nHits[iSL];
690  } else {
691  for (unsigned compare = iSL; compare < iSL + 9 * nHits[iSL]; compare += 9) {
692  if ((LRknown[iRef] && !LRknown[compare]) ||
693  (LRknown[iRef] == LRknown[compare] && tMin[iRef] < tMin[compare]))
694  iRef = compare;
695  }
696  }
697  }
698  // choose best hit (LR known before LR unknown, then shortest drift time)
699  bool useHit = false;
700  if (LRknown[iRef]) {
701  useHit = (allHits[ihit]->LRknown() && t <= tMin[iRef]);
702  } else {
703  useHit = (allHits[ihit]->LRknown() || t <= tMin[iRef]);
704  }
705  B2DEBUG(250, "relevant wire SL " << iSL << " LR " << allHits[ihit]->getLeftRight()
706  << " t " << t << " iRef " << iRef << " useHit " << useHit);
707  if (useHit) {
708  // keep drift time and LR
709  LRknown[iRef] = allHits[ihit]->LRknown();
710  tMin[iRef] = t;
711  hitIds[iRef] = allHits[ihit]->getArrayIndex();
712  }
713  }
714  }
715  // save selected hit Ids
716  for (unsigned iHit = 0; iHit < hitIds.size(); ++iHit) {
717  if (hitIds[iHit] >= 0) selectedHitIds.push_back(hitIds[iHit]);
718  }
719  return selectedHitIds;
720 }
721 vector<unsigned>
722 NeuroTrigger::selectHits(unsigned isector, const CDCTriggerTrack& track,
723  bool returnAllRelevant)
724 {
725  CDCTriggerMLP& expert = m_MLPs[isector];
726  vector<unsigned> selectedHitIds = {};
727  // prepare vectors to keep best drift times, left/right and selected hit IDs
728  vector<int> tMin;
729  tMin.assign(expert.nNodesLayer(0), expert.getTMax());
730  vector<bool> LRknown;
731  LRknown.assign(expert.nNodesLayer(0), false);
732  vector<int> hitIds;
733  hitIds.assign(expert.nNodesLayer(0), -1);
734  vector<unsigned> nHits;
735  nHits.assign(9, 0);
736 
737  // loop over axial hits related to input track
739  track.getRelationsTo<CDCTriggerSegmentHit>(m_hitCollectionName);
740  B2DEBUG(250, "start axial hit loop");
741  for (unsigned ihit = 0; ihit < axialHits.size(); ++ihit) {
742  // skip hits with negative relation weight (not selected in finder)
743  if (axialHits.weight(ihit) < 0) continue;
744  unsigned short iSL = axialHits[ihit]->getISuperLayer();
745  // skip stereo hits (should not be related to track, but check anyway)
746  if (iSL % 2 == 1) continue;
747  if (expert.getSLpatternUnmasked() != 0 &&
748  !((expert.getSLpatternUnmasked() >> iSL) & 1)) {
749  B2DEBUG(250, "skipping hit in SL " << iSL);
750  continue;
751  }
752  // get priority time and apply time window cut
753 
754  int t = (m_hasT0) ? axialHits[ihit]->priorityTime() - m_T0 : 0;
755  if (t < 0) {
756  t = 0;
757  } else if (t > expert.getTMax()) {
758  t = expert.getTMax();
759  }
760  double relId = getRelId(*axialHits[ihit]);
761  if (expert.isRelevant(relId, iSL)) {
762  // get reference hit (worst of existing hits)
763  unsigned short iRef = iSL;
764  if (expert.getMaxHitsPerSL() > 1) {
765  if (nHits[iSL] < expert.getMaxHitsPerSL() &&
766  (expert.getSLpatternUnmasked() >> (iSL + 9 * nHits[iSL])) & 1) {
767  iRef += 9 * nHits[iSL];
768  ++nHits[iSL];
769  } else {
770  for (unsigned compare = iSL; compare < iSL + 9 * nHits[iSL]; compare += 9) {
771  if ((LRknown[iRef] && !LRknown[compare]) ||
772  (LRknown[iRef] == LRknown[compare] && tMin[iRef] < tMin[compare]))
773  iRef = compare;
774  }
775  }
776  }
777  // choose best hit (LR known before LR unknown, then shortest drift time)
778  bool useHit = false;
779  if (LRknown[iRef]) {
780  useHit = (axialHits[ihit]->LRknown() && t <= tMin[iRef]);
781  } else {
782  useHit = (axialHits[ihit]->LRknown() || t <= tMin[iRef]);
783  }
784  B2DEBUG(250, "relevant wire SL " << iSL << " LR " << axialHits[ihit]->getLeftRight()
785  << " t " << t << " iRef " << iRef << " useHit " << useHit);
786  if (useHit) {
787  // keep drift time and LR
788  LRknown[iRef] = axialHits[ihit]->LRknown();
789  tMin[iRef] = t;
790  hitIds[iRef] = axialHits[ihit]->getArrayIndex();
791  }
792  if (returnAllRelevant) selectedHitIds.push_back(axialHits[ihit]->getArrayIndex());
793  }
794  }
795 
796  // loop over stereo hits, choosing up to MaxHitsPerSL per superlayer
797  B2DEBUG(250, "start stereo hit loop");
798  for (int ihit = 0; ihit < m_segmentHits.getEntries(); ++ ihit) {
799  unsigned short iSL = m_segmentHits[ihit]->getISuperLayer();
800  // skip axial hits
801  if (iSL % 2 == 0) continue;
802  if (expert.getSLpatternUnmasked() != 0 &&
803  !((expert.getSLpatternUnmasked() >> iSL) & 1)) {
804  B2DEBUG(250, "skipping hit in SL " << iSL);
805  continue;
806  }
807  // get priority time and apply time window cut
808  int t = (m_hasT0) ? m_segmentHits[ihit]->priorityTime() - m_T0 : 0;
809  if (t < 0) {
810  t = 0;
811  } else if (t > expert.getTMax()) {
812  t = expert.getTMax();
813  }
814  double relId = getRelId(*m_segmentHits[ihit]);
815  if (expert.isRelevant(relId, iSL)) {
816  // get reference hit (worst of existing hits)
817  unsigned short iRef = iSL;
818  if (expert.getMaxHitsPerSL() > 1) {
819  if (nHits[iSL] < expert.getMaxHitsPerSL() &&
820  (expert.getSLpatternUnmasked() >> (iSL + 9 * nHits[iSL])) & 1) {
821  iRef += 9 * nHits[iSL];
822  ++nHits[iSL];
823  } else {
824  for (unsigned compare = iSL; compare < iSL + 9 * nHits[iSL]; compare += 9) {
825  if ((LRknown[iRef] && !LRknown[compare]) ||
826  (LRknown[iRef] == LRknown[compare] && tMin[iRef] < tMin[compare]))
827  iRef = compare;
828  }
829  }
830  }
831  // choose best hit (LR known before LR unknown, then shortest drift time)
832  bool useHit = false;
833  if (LRknown[iRef]) {
834  useHit = (m_segmentHits[ihit]->LRknown() && t <= tMin[iRef]);
835  } else {
836  useHit = (m_segmentHits[ihit]->LRknown() || t <= tMin[iRef]);
837  }
838  B2DEBUG(250, "relevant wire SL " << iSL << " LR " << m_segmentHits[ihit]->getLeftRight()
839  << " t " << t << " iRef " << iRef << " useHit " << useHit);
840  if (useHit) {
841  // keep drift time and LR
842  LRknown[iRef] = m_segmentHits[ihit]->LRknown();
843  tMin[iRef] = t;
844  hitIds[iRef] = ihit;
845  }
846  if (returnAllRelevant) selectedHitIds.push_back(ihit);
847  }
848  }
849 
850  // save selected hit Ids
851  if (!returnAllRelevant) {
852  for (unsigned iHit = 0; iHit < hitIds.size(); ++iHit) {
853  if (hitIds[iHit] >= 0) selectedHitIds.push_back(hitIds[iHit]);
854  }
855  }
856  return selectedHitIds;
857 }
858 
859 vector<float>
860 NeuroTrigger::getInputVector(unsigned isector, const vector<unsigned>& hitIds)
861 {
862  CDCTriggerMLP& expert = m_MLPs[isector];
863  // prepare empty input vector and vectors to keep best drift times
864  vector<float> inputVector;
865  inputVector.assign(expert.nNodesLayer(0), 0.);
866  // convert hits to input values
867  vector<unsigned> nHits;
868  nHits.assign(9, 0);
869  for (unsigned ii = 0; ii < hitIds.size(); ++ii) {
870  int ihit = hitIds[ii];
871  unsigned short iSL = m_segmentHits[ihit]->getISuperLayer();
872  unsigned short iRef = iSL + 9 * nHits[iSL];
873  ++nHits[iSL];
874  int t = (m_hasT0) ? m_segmentHits[ihit]->priorityTime() - m_T0 : 0;
875  if (t < 0) {
876  t = 0;
877  } else if (t > expert.getTMax()) {
878  t = expert.getTMax();
879  }
880  int LR = m_segmentHits[ihit]->getLeftRight();
881  double relId = getRelId(*m_segmentHits[ihit]);
882  int priority = m_segmentHits[ihit]->getPriorityPosition();
883  // get scaled input values (scaling such that the absolute value of all inputs is < 1)
884  inputVector[3 * iRef] = expert.scaleId(relId, iSL);
885  float scaleT = pow(2, floor(log2(1. / expert.getTMax())));
886  inputVector[3 * iRef + 1] = (((LR >> 1) & 1) - (LR & 1)) * t * scaleT;
887  inputVector[3 * iRef + 2] = m_alpha[iSL][int(priority < 3)] * 0.5;
888  }
889  return inputVector;
890 }
891 
892 vector<float>
893 NeuroTrigger::runMLP(unsigned isector, vector<float> input)
894 {
895  const CDCTriggerMLP& expert = m_MLPs[isector];
896  vector<float> weights = expert.getWeights();
897  vector<float> layerinput = input;
898  vector<float> layeroutput = {};
899  unsigned iw = 0;
900  for (unsigned il = 1; il < expert.nLayers(); ++il) {
901  //add bias input
902  layerinput.push_back(1.);
903  //prepare output
904  layeroutput.clear();
905  layeroutput.assign(expert.nNodesLayer(il), 0.);
906  //loop over outputs
907  for (unsigned io = 0; io < layeroutput.size(); ++io) {
908  //loop over inputs
909  for (unsigned ii = 0; ii < layerinput.size(); ++ii) {
910  layeroutput[io] += layerinput[ii] * weights[iw++];
911  }
912  //apply activation function
913  layeroutput[io] = tanh(layeroutput[io] / 2.);
914  }
915  //output is new input
916  layerinput = layeroutput;
917  }
918  return expert.unscaleTarget(layeroutput);
919 }
920 
921 vector<float>
922 NeuroTrigger::runMLPFix(unsigned isector, vector<float> input)
923 {
924  unsigned precisionInput = m_precision[3];
925  unsigned precisionWeights = m_precision[4];
926  unsigned precisionLUT = m_precision[5];
927  unsigned precisionTanh = m_precision[3];
928  unsigned dp = precisionInput + precisionWeights - precisionLUT;
929 
930  const CDCTriggerMLP& expert = m_MLPs[isector];
931  // transform inputs to fixed point (cut off to input precision)
932  vector<long> inputFix(input.size(), 0);
933  for (unsigned ii = 0; ii < input.size(); ++ii) {
934  inputFix[ii] = long(input[ii] * (1 << precisionInput));
935  }
936  // transform weights to fixed point (round to weight precision)
937  vector<float> weights = expert.getWeights();
938  vector<long> weightsFix(weights.size(), 0);
939  for (unsigned iw = 0; iw < weights.size(); ++iw) {
940  weightsFix[iw] = long(round(weights[iw] * (1 << precisionWeights)));
941  }
942  // maximum input value for the tanh LUT
943  unsigned xMax = unsigned(ceil(atanh(1. - 1. / (1 << (precisionTanh + 1))) *
944  (1 << (precisionLUT + 1))));
945 
946  // run MLP
947  vector<long> layerinput = inputFix;
948  vector<long> layeroutput = {};
949  unsigned iw = 0;
950  for (unsigned il = 1; il < expert.nLayers(); ++il) {
951  // add bias input
952  layerinput.push_back(1 << precisionInput);
953  // prepare output
954  layeroutput.clear();
955  layeroutput.assign(expert.nNodesLayer(il), 0);
956  // loop over outputs
957  for (unsigned io = 0; io < layeroutput.size(); ++io) {
958  // loop over inputs
959  for (unsigned ii = 0; ii < layerinput.size(); ++ii) {
960  layeroutput[io] += layerinput[ii] * weightsFix[iw++];
961  }
962  // apply activation function -> LUT, calculated on the fly here
963  unsigned long bin = abs(layeroutput[io]) >> dp;
964  // correction to get symmetrical rounding errors
965  float x = (bin + 0.5 - 1. / (1 << (dp + 1))) / (1 << precisionLUT);
966  long tanhLUT = (bin < xMax) ? long(round(tanh(x / 2.) * (1 << precisionTanh))) : (1 << precisionTanh);
967  layeroutput[io] = (layeroutput[io] < 0) ? -tanhLUT : tanhLUT;
968  }
969  // output is new input
970  layerinput = layeroutput;
971  }
972 
973  // transform output back to float before unscaling
974  vector<float> output(layeroutput.size(), 0.);
975  for (unsigned io = 0; io < output.size(); ++io) {
976  output[io] = layeroutput[io] / float(1 << precisionTanh);
977  }
978  return expert.unscaleTarget(output);
979 }
980 
981 void
982 NeuroTrigger::save(const string& filename, const string& arrayname)
983 {
984  B2INFO("Saving networks to file " << filename << ", array " << arrayname);
985  TFile datafile(filename.c_str(), "UPDATE");
986  TObjArray* MLPs = new TObjArray(m_MLPs.size());
987  for (unsigned isector = 0; isector < m_MLPs.size(); ++isector) {
988  MLPs->Add(&m_MLPs[isector]);
989  }
990  MLPs->Write(arrayname.c_str(), TObject::kSingleKey | TObject::kOverwrite);
991  datafile.Close();
992  MLPs->Clear();
993  delete MLPs;
994 }
995 
996 bool
997 NeuroTrigger::load(const string& filename, const string& arrayname)
998 {
999  if (filename.size() < 1) {
1000  m_MLPs.clear();
1001  m_MLPs = m_cdctriggerneuroconfig->getMLPs();
1002  if (m_MLPs.size() == 0) {
1003  B2ERROR("Could not load Neurotrigger weights from database!");
1004  return false;
1005  }
1006  B2DEBUG(2, "Loaded Neurotrigger MLP weights from database: " + m_cdctriggerneuroconfig->getNNName());
1007  B2DEBUG(100, "loaded " << m_MLPs.size() << " networks from database");
1008  // load some values from the geometry that will be needed for the input
1009  setConstants();
1010  return true;
1011  } else {
1012  TFile datafile(filename.c_str(), "READ");
1013  if (!datafile.IsOpen()) {
1014  B2WARNING("Could not open file " << filename);
1015  return false;
1016  }
1017 
1018  TObjArray* MLPs = (TObjArray*)datafile.Get(arrayname.c_str());
1019  if (!MLPs) {
1020  datafile.Close();
1021  B2WARNING("File " << filename << " does not contain key " << arrayname);
1022  return false;
1023  }
1024  m_MLPs.clear();
1025  for (int isector = 0; isector < MLPs->GetEntriesFast(); ++isector) {
1026  CDCTriggerMLP* expert = dynamic_cast<CDCTriggerMLP*>(MLPs->At(isector));
1027  if (expert) m_MLPs.push_back(*expert);
1028  else B2WARNING("Wrong type " << MLPs->At(isector)->ClassName() << ", ignoring this entry.");
1029  }
1030  MLPs->Clear();
1031  delete MLPs;
1032  datafile.Close();
1033  B2DEBUG(100, "loaded " << m_MLPs.size() << " networks");
1034 
1035  // load some values from the geometry that will be needed for the input
1036  setConstants();
1037 
1038  return true;
1039  }
1040 }
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
Belle2::NeuroTrigger::updateTrack
void updateTrack(const CDCTriggerTrack &track)
Calculate 2D phi position and arclength for the given track and store them.
Definition: NeuroTrigger.cc:325
Belle2::NeuroTrigger::getEventTime
void getEventTime(unsigned isector, const CDCTriggerTrack &track, std::string et_option, const bool)
Read out the event time and store it.
Definition: NeuroTrigger.cc:384
Belle2::NeuroTrigger::setConstants
void setConstants()
Loads parameters from the geometry and precalculates some constants that will be needed.
Definition: NeuroTrigger.cc:234
Belle2::CDCTriggerTrack
Track created by the CDC trigger.
Definition: CDCTriggerTrack.h:15
Belle2::NeuroTrigger::runMLP
std::vector< float > runMLP(unsigned isector, std::vector< float > input)
Run an expert MLP.
Definition: NeuroTrigger.cc:893
Belle2::NeuroTrigger::getInputPattern
unsigned long getInputPattern(unsigned isector, const CDCTriggerTrack &track, const bool neurotrackinputmode)
Calculate input pattern for MLP.
Definition: NeuroTrigger.cc:581
Belle2::NeuroTrigger::getRangeIndices
std::vector< unsigned > getRangeIndices(const Parameters &p, unsigned isector)
Get indices for sector ranges in parameter lists.
Definition: NeuroTrigger.cc:216
Belle2::NeuroTrigger::getRelId
double getRelId(const CDCTriggerSegmentHit &hit)
Calculate phi position of a hit relative to 2D track (scaled to number of wires).
Definition: NeuroTrigger.cc:372
Belle2::NeuroTrigger::selectHits
std::vector< unsigned > selectHits(unsigned isector, const CDCTriggerTrack &track, bool returnAllRelevant=false)
Select best hits for each super layer.
Definition: NeuroTrigger.cc:722
Belle2::CDC::CDCGeometryPar
The Class for CDC Geometry Parameters.
Definition: CDCGeometryPar.h:75
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2::CDC::CDCGeometryPar::Instance
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Definition: CDCGeometryPar.cc:41
Belle2::NeuroTrigger::selectMLPbyPattern
int selectMLPbyPattern(std::vector< int > &MLPs, unsigned long pattern, const bool neurotrackinputmode)
Select one MLP from a list of sector indices.
Definition: NeuroTrigger.cc:288
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::NeuroTrigger::runMLPFix
std::vector< float > runMLPFix(unsigned isector, std::vector< float > input)
Run an expert MLP with fixed point arithmetic.
Definition: NeuroTrigger.cc:922
Belle2::Unit::deg
static const double deg
degree to radians
Definition: Unit.h:119
Belle2::NeuroTrigger::updateTrackFix
void updateTrackFix(const CDCTriggerTrack &track)
Calculate 2D phi position and arclength for the given track and store them.
Definition: NeuroTrigger.cc:341
Belle2::NeuroTrigger::selectMLPs
std::vector< int > selectMLPs(float phi0, float invpt, float theta)
Select all matching expert MLPs based on the given track parameters.
Definition: NeuroTrigger.cc:261
Belle2::RelationVector::weight
float weight(int index) const
Get weight with index.
Definition: RelationVector.h:120
Belle2::NeuroTrigger::load
bool load(const std::string &filename, const std::string &arrayname="MLPs")
Load MLPs from file.
Definition: NeuroTrigger.cc:997
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::NeuroTrigger::Parameters
Struct to keep neurotrigger parameters.
Definition: NeuroTrigger.h:37
Belle2::NeuroTrigger::initializeCollections
void initializeCollections(std::string hitCollectionName, std::string eventTimeName, std::string et_option)
set the hit collection and event time to required and store the hit collection name
Definition: NeuroTrigger.cc:251
Belle2::NeuroTrigger::save
void save(const std::string &filename, const std::string &arrayname="MLPs")
Save MLPs to file.
Definition: NeuroTrigger.cc:982
Belle2::CDCTriggerSegmentHit
Combination of several CDCHits to a track segment hit for the trigger.
Definition: CDCTriggerSegmentHit.h:16
Belle2::NeuroTrigger::getInputVector
std::vector< float > getInputVector(unsigned isector, const std::vector< unsigned > &hitIds)
Calculate input values for MLP.
Definition: NeuroTrigger.cc:860
Belle2::CDCTriggerMLP
Class to keep all parameters of an expert MLP for the neuro trigger.
Definition: CDCTriggerMLP.h:13
Belle2::NeuroTrigger::initialize
void initialize(const Parameters &p)
Set parameters and get some network independent parameters.
Definition: NeuroTrigger.cc:20
Belle2::NeuroTrigger::selectHitsHWSim
std::vector< unsigned > selectHitsHWSim(unsigned isector, const CDCTriggerTrack &track)
Select hits for each super layer from the ones related to input track.
Definition: NeuroTrigger.cc:647