ATLAS Offline Software
FPGATrackSimFitConstantBank.cxx
Go to the documentation of this file.
1 // Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
2 
3 #include <Eigen/StdVector>
7 #include "limits.h"
8 #include <TMath.h>
9 #include <Eigen/Dense>
10 #include <iostream>
11 #include <iomanip>
12 #include <fstream>
13 #include <string>
14 #include <cassert>
15 #include <cstdio>
16 #include <cmath>
17 #include <bitset>
18 #include <ieee754.h>
19 
20 #include <sstream>
21 
22 FPGATrackSimFitConstantBank::FPGATrackSimFitConstantBank(FPGATrackSimPlaneMap const * pmap, int ncoords, std::string const & fname, bool /*isFirstStage*/, int missingPlane) :
23  AthMessaging ("FPGATrackSimFitConstantBank"),
24  m_pmap(pmap),
25  m_bankID(0),
26  m_nsectors(0),
27  m_ncoords(ncoords),
28  m_nconstr(0),
29  m_npixcy(0),
30  m_missingPlane(missingPlane),
31 // m_isFirstStage(isFirstStage),
32  m_isIdealCoordFit(true)
33 {
34  std::ifstream geocfile(fname);
35  if (not geocfile.is_open()) {
36  ATH_MSG_WARNING("FitConstants file: " << fname << " cannot be opened");
37  }
38  else {
39  ATH_MSG_INFO("Reading " << fname);
40  // Read the file header
41  readHeader(geocfile);
42  ATH_MSG_INFO("Settings: m_ncoords="<<m_ncoords<<" m_npars="<<m_npars);
43  // Read the sector constants
44  readSectorInfo(geocfile);
45  // Pre-calculate the majority logic elements
46  if (m_missingPlane == -1)
48 
49  if (sizeof(float) * CHAR_BIT != 32)
50  ATH_MSG_WARNING("Floating points on this computer are not 32 bit. This may cause a problem for the hardware agreement. Be careful!");
51 
52  setIdealCoordFit(true);
54  }
55 }
56 
57 
58 // Reads the file header, sets configuration and sizes the sector variables
59 void FPGATrackSimFitConstantBank::readHeader(std::ifstream & geocfile)
60 {
61  std::string key;
62  int ival;
63  int nplanes;
64 
65  for (int i = 0; i < 5; ++i) getline(geocfile, key); // skip 5 lines
66 
67  geocfile >> key;
68  if (key != "NPLANES") ATH_MSG_ERROR("Invalid file format");
69  geocfile >> nplanes;
70 
71  geocfile >> key;
72  if (key!="NSECTORS") ATH_MSG_ERROR("Invalid file format");
73  geocfile >> m_nsectors;
74 
75  geocfile >> key;
76  if (key!="NDIM") ATH_MSG_ERROR("Invalid file format");
77  geocfile >> ival;
78 
79  // Set derived configuration variables
80  if (ival==1) m_npars = 3;
81  else if (ival==2) m_npars = 5;
82  else ATH_MSG_ERROR("Number of dimensions invalid");
83 
84  m_npixcy = 2*(m_ncoords - nplanes); // difference between number of coordinates and planes is number of 2d measurements
86 
87  // Allocate the block of pointer per sector
88  m_sector_good.resize(m_nsectors);
93 
95 
99 }
100 
101 
102 // Reads the main body of the constants file
103 void FPGATrackSimFitConstantBank::readSectorInfo(std::ifstream & geocfile)
104 {
105  std::string key;
106  int ival;
107  double dval; // temp value used to conver from double to float
108 
109  for (int isec = 0; isec < m_nsectors; ++isec)
110  {
111  // Read sector number
112  geocfile >> key >> ival;
113  if (ival != isec) ATH_MSG_WARNING("Reached sector " << ival << " but expected sector " << isec);
114 
115  // Read fit parameters
116  for (int ipar = 0; ipar < m_npars; ++ipar)
117  {
118  geocfile >> key;
119  // no check on key?
120  for (int icoord = 0; icoord < m_ncoords; ++icoord)
121  {
122  geocfile >> dval;
123  if (geocfile.fail())
124  ATH_MSG_WARNING("par loop (1) key=\"" << key << "\" ipar,icoord=" << ipar << "," << icoord);
125 
126  if (dval != 0) m_sector_good[isec] = true;
127  m_fit_pars(isec, ipar, icoord) = dval;
128  }
129  }
130 
131  // Read k averages
132  geocfile >> key;
133  assert(key=="kaverages");
134  for (int ik=0;ik<m_nconstr;++ik)
135  {
136  geocfile >> dval;
137  if (dval != 0) m_sector_good[isec] = true;
138  m_kaverage(isec,ik) = dval;
139  }
140 
141  // Read kernel
142  geocfile >> key;
143  assert(key=="kernel");
144  for (int ik=0;ik<m_nconstr;++ik)
145  {
146  for (int ik2=0;ik2<m_ncoords;++ik2)
147  {
148  geocfile >> dval;
149  if (dval != 0) m_sector_good[isec] = true;
150  m_kernel(isec, ik, ik2) = dval;
151  }
152  }
153 
154  // Read fit const
155  for (int ipar=0;ipar<m_npars;++ipar)
156  {
157  geocfile >> key;
158  geocfile >> dval;
159  if (dval != 0) m_sector_good[isec] = true;
160  m_fit_const(isec, ipar) = dval;
161  }
162 
163  }
164 }
165 
166 
167 // Majority fits are those for which we have missing hits or when we drop a hit on 8/8 to improve chi2
169 {
170 
171  for (int isec=0;isec!=m_nsectors;++isec)
172  {
173  // Fill m_maj_a and m_maj_kk
174  for (int ix=0;ix!=m_ncoords;++ix)
175  {
176  for (int row=0;row!=m_nconstr;++row)
177  m_maj_a(isec, ix) += m_kaverage(isec, row) * m_kernel(isec, row, ix);
178  for (int jx=0;jx!=m_ncoords;++jx)
179  {
180  for (int row=0;row!=m_nconstr;++row)
181  m_maj_kk(isec, ix, jx) += m_kernel(isec, row, ix) * m_kernel(isec, row, jx);
182  }
183  }
184 
185  // Fill m_maj_invkk
186 
187  // Pixel layers (2 coordinates each), assume these come first!
188  for (int ix=0;ix!=m_npixcy;ix+=2)
189  {
190  float det = m_maj_kk(isec, ix, ix) * m_maj_kk(isec, ix+1, ix+1)
191  - m_maj_kk(isec, ix+1, ix) * m_maj_kk(isec, ix, ix+1);
192  if (det == 0) continue;
193 
194  m_maj_invkk(isec, ix, ix) = m_maj_kk(isec, ix+1, ix+1)/det;
195  m_maj_invkk(isec, ix, ix+1) = -m_maj_kk(isec, ix+1, ix)/det;
196  m_maj_invkk(isec, ix+1, ix) = -m_maj_kk(isec, ix, ix+1)/det;
197  m_maj_invkk(isec, ix+1, ix+1) = m_maj_kk(isec, ix, ix)/det;
198  }
199 
200  // Strip layers (1 coordinate)
201  for (int ix=m_npixcy;ix!=m_ncoords;++ix) {
202  if (m_maj_kk(isec, ix, ix) == 0) continue;
203  m_maj_invkk(isec, ix, ix) = 1./m_maj_kk(isec, ix, ix);
204  }
205  }
206 }
207 
208 
209 // Prepares the inverted constants used in the invlinfit()
211 {
212  m_invfit_consts.resize(m_nsectors);
213 
214  for (int isec=0;isec!=m_nsectors;++isec)
215  {
216  Eigen::MatrixXf thisMatrix(m_ncoords, m_ncoords);
217 
218  // Parameters, by row
219  for (int ip=0;ip!=m_npars;++ip)
220  for (int ix=0;ix!=m_ncoords;++ix)
221  thisMatrix(ip,ix)= (m_fit_pars(isec, ip, ix));
222 
223  // Constraints, by row
224  for (int ic=0;ic!=m_nconstr;++ic)
225  for (int ix=0;ix!=m_ncoords;++ix)
226  thisMatrix(ic+m_npars,ix) = m_kernel(isec, ic, ix);
227 
228  // reset usable coordinate
229  for (int ix = 0; ix != m_ncoords; ++ix)
230  m_WCs(isec,ix) = 1; // start with all WCs
231 
232  // now let's check if we have any "null" row or column, and if so, we will get rid of them temporarily
233  // we don't have access to module/WC information at this stage should just look for non-zero values
234  for (int ix=0;ix!=m_ncoords;++ix) {
235  bool foundrow = false;
236  bool foundcolumn = false;
237  for (int iy=0;iy!=m_ncoords;++iy) {
238  if (abs(thisMatrix(ix,iy)) > MTX_TOLERANCE) foundrow = true;
239  if (abs(thisMatrix(iy,ix)) > MTX_TOLERANCE) foundcolumn = true;
240  }
241  if (foundrow && foundcolumn) m_WCs(isec,ix) = 0;
242  }
243 
244  // now count the number of usable entries
245  size_t nDimToUse = 0; for (int ix=0;ix!=m_ncoords;++ix) if(!m_WCs(isec,ix)) ++nDimToUse;
246 
247  Eigen::MatrixXf tempMatrix(nDimToUse, nDimToUse); // this is the matrix we will invert
248 
249  size_t counteri = 0, counterj = 0;
250  for (int i = 0; i < m_ncoords; i++) {
251  counterj = 0;
252  if (!m_WCs(isec,i)) {
253  for (int j = 0; j < m_ncoords; j++) {
254  if (!m_WCs(isec,j)) {
255  tempMatrix(counteri,counterj) = thisMatrix(i,j); // if we use this coordinate, copy it over
256  counterj++; // iterate counter
257  }
258  }
259  counteri++;
260  }
261  }
262 
263  // Invert the matrix
264  Eigen::MatrixXf inverted = tempMatrix.inverse();
265 
266  // this is the full matrix that includes inverted coordinates plus the padded zeros
267  Eigen::MatrixXf fullInverted(m_ncoords, m_ncoords);
268 
269  // now we put things back into the full matrix, adding back the padded zeros
270  counteri = 0; counterj = 0;
271  for (int i = 0; i < m_ncoords; i++) {
272  counterj = 0;
273  if (!m_WCs(isec,i)) {
274  for (int j = 0; j < m_ncoords; j++) {
275  if (!m_WCs(isec,j)) {
276  fullInverted(i,j) = inverted(counteri,counterj); // if we use this coordinate, copy it over
277  counterj++; // iterate counter
278  }
279  else {
280  fullInverted(i,j) = 0;
281  }
282  }
283  counteri++;
284  }
285  else {
286  for (int j = 0; j < m_ncoords; j++) {
287  fullInverted(i,j) = 0;
288  }
289  }
290  }
291  m_invfit_consts[isec]= fullInverted;
292 
293  }
294 }
295 
296 
297 
299 // Main Interface Functions
300 // This functions takes in hit coordinates from 'track' and does the linear fit
301 // to calculate the track parameters, which are populated back into 'track'.
302 // Returns true on success.
303 bool FPGATrackSimFitConstantBank::linfit(sector_t sector, FPGATrackSimTrack & track, bool isSecondStage) const
304 {
305  // Do majority guess if it's needed
306  if (m_missingPlane == -1 && !m_isIdealCoordFit)
307  {
308  int nmissing = track.getNMissing();
309  if (nmissing > 0)
310  {
311  int guess_res = missing_point_guess(sector, track, !isSecondStage, false);
312  if (nmissing != guess_res) return false; // majority failed
313  }
314  }
315 
316  // evaluate the chi2
317  linfit_chisq(sector, track);
318 
319  // Do the fit
320  linfit_pars_eval(sector, track);
321 
322  return true;
323 }
324 
325 
326 int FPGATrackSimFitConstantBank::missing_point_guess(sector_t sector, FPGATrackSimTrack &track, bool isFirstStage, bool doExtrapolation) const
327 {
328 
329  std::vector<int> coordsmask(m_ncoords);
330  std::vector<int> missid;
331 
332  // Keep track of which hits are missing
333  int nmissing = 0;
334  std::vector<int> missing_planes;
335  for (int j = 0; j < m_ncoords ; ++j)
336  {
337  coordsmask[j] = (track.getHitMap()>>j) & 1;
338  if (!coordsmask[j]) {
339  int plane = m_pmap->getCoordLayer(j);
340  if(missing_planes.size() == 0) {
341  missing_planes.push_back(plane);
342  }
343  else {
344  for (unsigned k = 0; k < missing_planes.size(); k++) {
345  if(plane == missing_planes[k]) break;
346  if(k == missing_planes.size() - 1) missing_planes.push_back(plane);
347  }
348  }
349  missid.push_back(j);
350  nmissing++;
351  }
352  }
353 
354  if(!doExtrapolation){
355  if(isFirstStage){
356  if(missing_planes.size() > 1){
357  ATH_MSG_WARNING("missing_point_guess() Can only guess 1 layer in the first stage");
358  return 0;
359  }
360  if (nmissing > 2){
361  ATH_MSG_WARNING("missing_point_guess() Cannot guess more than two coordinates in the first stage");
362  return 0;
363  }
364  }
365  else{
366  if(missing_planes.size() > 2){
367  ATH_MSG_WARNING("missing_point_guess() Can only guess 2 layers in the second stage");
368  return 0;
369  }
370  if (nmissing > 4){
371  ATH_MSG_WARNING("missing_point_guess() Cannot guess more than four coordinates in the second stage");
372  return 0;
373  }
374  }
375  }
376 
377  Eigen::MatrixXd coef(nmissing,nmissing);
378  Eigen::VectorXd a(nmissing);
379 
380  for (int i=0; i<nmissing; ++i)
381  {
382  a[i] = -m_maj_a(sector, missid[i]);
383  // complete the calculation with the terms that depend on the known hits
384  for (int col=0; col<m_ncoords; ++col)
385  {
386  if (!coordsmask[col]) continue;
387 
388  a[i] -= m_maj_kk(sector, col, missid[i])*track.getPhiCoord(m_pmap->getCoordLayer(col));
389 
390  if (m_pmap->getDim(m_pmap->getCoordLayer(col)) == 2) { // do two at a time if 2d, then skip ahead
391  a[i] -= m_maj_kk(sector, col+1, missid[i])*track.getEtaCoord(m_pmap->getCoordLayer(col));
392  ++col;
393  }
394  }
395 
396  for (int j=0; j<nmissing; ++j)
397  {
398  // here for missing pieces need to know the coordinate and sector number only
399  coef(i,j) = m_maj_kk(sector, missid[i], missid[j]);
400  }
401  }
402 
403  if (coef.determinant()==0) return -1;
404  Eigen::VectorXd missing_hits = coef.inverse()*a; // our missing hits!
405 
406  for(int m=0; m<nmissing; m++){
407 
408  int missedplane = m_pmap->getCoordLayer(missid[m]);
409 
410  if(m_pmap->getDim(m_pmap->getCoordLayer(missid[m])) == 1){
411 
412  FPGATrackSimHit newhit;
413  if(!doExtrapolation)
415  else
418  newhit.setLayer(missedplane);
419  newhit.setSection(0);
420  if (m_isIdealCoordFit) {
421  // TODO: all the missing hit logic will eventually need spacepoint updates.
422  double target_r = track.getIdealRadius(missedplane);
423  newhit.setX(target_r*TMath::Cos(missing_hits[m]));
424  newhit.setY(target_r*TMath::Sin(missing_hits[m]));
425  }
426  else{
427  newhit.setPhiCoord(missing_hits[m]);
428  }
429 
430  track.setFPGATrackSimHit(missedplane, newhit);
431  }
432  else if (m_pmap->getDim(m_pmap->getCoordLayer(missid[m])) == 2){
433 
434  FPGATrackSimHit newhit;
435  if(!doExtrapolation)
437  else
440  newhit.setLayer(missedplane);
441  newhit.setSection(0);
442  if (m_isIdealCoordFit) {
443  // TODO This will also eventually ned spacepoint updates.
444  double target_r = track.getIdealRadius(missedplane);
445  newhit.setX(target_r*TMath::Cos(missing_hits[m]));
446  newhit.setY(target_r*TMath::Sin(missing_hits[m]));
447  newhit.setZ(missing_hits[m+1]);
448  }
449  else {
450  newhit.setPhiCoord(missing_hits[m]);
451  newhit.setEtaCoord(missing_hits[m+1]);
452  }
453  m++; //skip ahead
454 
455  track.setFPGATrackSimHit(missedplane, newhit);
456  }
457  }
458 
459  return nmissing;
460 }
461 
463 {
464  float chi2 = 0;
465 
466  for (int i = 0; i < m_nconstr; i++)
467  {
468  float chi_component = m_kaverage(sector, i);
469  for (int ix = 0; ix != m_npixcy/2; ix++) // pxl plane loop (divide by two to get number of pix planes from pix coords)
470  {
471  chi_component += m_kernel(sector, i, 2*ix) * trk.getPhiCoord(m_pmap->getCoordLayer(2*ix));
472  chi_component += m_kernel(sector, i, 2*ix+1) * trk.getEtaCoord(m_pmap->getCoordLayer(2*ix));
473  }
474 
475  for (int ix = m_npixcy; ix != m_ncoords; ix++) // strip coords, easier
476  {
477  chi_component += m_kernel(sector, i, ix) * trk.getPhiCoord(m_pmap->getCoordLayer(ix));
478  }
479 
480  chi2 += chi_component * chi_component;
481  }
482 
483  trk.setChi2(chi2);
484 }
485 
486 
487 
488 /* This function takes as input the pointer to the geometrical constants
489  structure, the coords array and and array that will contain the
490  results of the fit.
491 
492  The fits return the number of useful parameter for the fits.
493  */
495 {
496  std::vector<float> pars(m_npars);
497 
498  for (int ip = 0; ip < m_npars; ip++)
499  {
500  pars[ip] = m_fit_const(sector, ip);
501 
502  for (int coord = 0; coord < m_ncoords; coord++) {
503  pars[ip] += m_fit_pars(sector, ip, coord) * trk.getPhiCoord(m_pmap->getCoordLayer(coord));
504  if (m_pmap->getDim(m_pmap->getCoordLayer(coord)) == 2) { // do two at a time if 2d, then skip ahead
505  pars[ip] += m_fit_pars(sector, ip, coord+1) * trk.getEtaCoord(m_pmap->getCoordLayer(coord));
506  ++coord;
507  }
508  }
509  }
510 
511  trk.setQOverPt(pars[0]);
512  trk.setD0(pars[1]);
513  trk.setPhi(pars[2]); // angle is moved within -pi to +pi
514  trk.setZ0(pars[3]);
515  trk.setEta(pars[4]);
516  if (trk.getDoDeltaGPhis()) {
517  trk.setQOverPt(trk.getHoughY()/1000.0 - trk.getQOverPt()); // final q/pT = q/pT from HT + delta q/pT
518  trk.setPhi(trk.getHoughX() - trk.getPhi()); // final phi_0 = phi_0 from HT + delta phi_0
519  }
520 }
521 
522 
523 
527 void FPGATrackSimFitConstantBank::invlinfit(sector_t sector, FPGATrackSimTrack &track, double const *constr) const
528 {
529  // vector of the acutal parameters, it is zeroed
530  Eigen::VectorXf pars(m_ncoords);
531 
532  for (int ip = 0; ip < m_npars; ++ip)
533  {
534  // The first elements are the track parameters. The track are shifted using the sector constants
535  pars(ip) = track.getParameter(ip) - m_fit_const(sector, ip);
536  if (ip == 2) pars(ip) = remainder(pars(ip), 2.*M_PI);
537  }
538  for (int ic = 0; ic < m_nconstr; ++ic)
539  {
540  // The rest of the paramaters are the external constraints.
541  // The external constraint is also shifted by the kAverage value
542  if (!constr)
543  pars(m_npars+ic) = -m_kaverage(sector, ic);
544  else
545  pars(m_npars+ic) = constr[ic] - m_kaverage(sector, ic);
546  }
547 
548  // The raw hits are obtained multiplying the parameters to the inverted constants
549  Eigen::VectorXf rawhits = (m_invfit_consts[sector]) * pars;
550 
551  // The hits are assigned to the original track
552  for (int j = 0; j < m_ncoords ; ++j) {
553  int plane = m_pmap->getCoordLayer(j);
554  SiliconTech Tech = m_pmap->getDetType(plane);
555  FPGATrackSimHit hit;
557  hit.setDetType(Tech);
558  hit.setLayer(plane);
559  if (m_WCs(sector,j))
560  hit.setPhiCoord(-1); // to keep track later on, set hit position negative
561  else
562  hit.setPhiCoord(rawhits(j));
563 
564  if (Tech == SiliconTech::pixel) {
565  if (m_WCs(sector,j))
566  hit.setEtaCoord(-1); // to keep track later on, set hit position negative
567  else
568  hit.setEtaCoord(rawhits(j+1));
569 
570  ++j; // skip a coordinate if doing two at once
571  }
572  track.setFPGATrackSimHit(plane, hit);
573  }
574 }
FPGATrackSimTrack::getPhiCoord
float getPhiCoord(int ilayer) const
make_hlt_rep.pars
pars
Definition: make_hlt_rep.py:90
query_example.row
row
Definition: query_example.py:24
FPGATrackSimHit::setSection
void setSection(unsigned v)
Definition: FPGATrackSimHit.h:88
FPGATrackSimTrack::setD0
void setD0(float v)
Definition: FPGATrackSimTrack.h:84
SiliconTech::strip
@ strip
python.SystemOfUnits.m
int m
Definition: SystemOfUnits.py:91
HitType::extrapolated
@ extrapolated
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
FPGATrackSimTrack::getDoDeltaGPhis
bool getDoDeltaGPhis() const
Definition: FPGATrackSimTrack.h:22
FPGATrackSimFitConstantBank::readSectorInfo
void readSectorInfo(std::ifstream &geocfile)
Definition: FPGATrackSimFitConstantBank.cxx:103
FPGATrackSimTrack
Definition: FPGATrackSimTrack.h:16
FPGATrackSimFitConstantBank::m_kernel
vector3D< float > m_kernel
Definition: FPGATrackSimFitConstantBank.h:85
FPGATrackSimHit::setPhiCoord
void setPhiCoord(float v)
Definition: FPGATrackSimHit.h:99
FPGATrackSimTrack::getPhi
float getPhi() const
Definition: FPGATrackSimTrack.h:35
FPGATrackSimFitConstantBank::invlinfit
void invlinfit(sector_t sector, FPGATrackSimTrack &track, double const *constr) const
This method uses the track parameters and additional constraints to use the constants to calculate th...
Definition: FPGATrackSimFitConstantBank.cxx:527
FPGATrackSimFitConstantBank::readHeader
void readHeader(std::ifstream &geocfile)
Definition: FPGATrackSimFitConstantBank.cxx:59
FPGATrackSimFitConstantBank.h
M_PI
#define M_PI
Definition: ActiveFraction.h:11
FPGATrackSimFitConstantBank::m_maj_a
vector2D< float > m_maj_a
Definition: FPGATrackSimFitConstantBank.h:88
FPGATrackSimFitConstantBank::m_fit_const
vector2D< float > m_fit_const
Definition: FPGATrackSimFitConstantBank.h:84
FPGATrackSimFitConstantBank::m_missingPlane
int m_missingPlane
Definition: FPGATrackSimFitConstantBank.h:71
FPGATrackSimFitConstantBank::m_nconstr
int m_nconstr
Definition: FPGATrackSimFitConstantBank.h:69
FPGATrackSimHit::setY
void setY(float v)
Definition: FPGATrackSimHit.h:123
m_nsectors
G4int m_nsectors
GU add phi min and phi span of overall volume.
Definition: PresParameterDef.h:60
FPGATrackSimFitConstantBank::m_fit_pars
vector3D< float > m_fit_pars
Definition: FPGATrackSimFitConstantBank.h:83
FPGATrackSimFitConstantBank::m_nsectors
int m_nsectors
Definition: FPGATrackSimFitConstantBank.h:66
FPGATrackSimFitConstantBank::m_ncoords
int m_ncoords
Definition: FPGATrackSimFitConstantBank.h:68
FPGATrackSimTrack::getHoughY
float getHoughY() const
Definition: FPGATrackSimTrack.h:31
HitType::guessed
@ guessed
FPGATrackSimPlaneMap::getCoordLayer
uint32_t getCoordLayer(uint32_t coord) const
Definition: FPGATrackSimPlaneMap.h:101
FPGATrackSimConstants.h
FPGATrackSimTrack::setZ0
void setZ0(float v)
Definition: FPGATrackSimTrack.h:86
FPGATrackSimTrack::getQOverPt
float getQOverPt() const
Definition: FPGATrackSimTrack.h:32
FPGATrackSimHit
Definition: FPGATrackSimHit.h:38
FPGATrackSimPlaneMap::getDim
uint32_t getDim(size_t logiLayer) const
Definition: FPGATrackSimPlaneMap.h:81
FPGATrackSimTrack::setPhi
void setPhi(float v, bool ForceRange=true)
Definition: FPGATrackSimTrack.cxx:159
FPGATrackSimFitConstantBank::FPGATrackSimFitConstantBank
FPGATrackSimFitConstantBank(FPGATrackSimPlaneMap const *pmap, int ncoords, std::string const &fname, bool isFirstStage, int missingPlane=-1)
Definition: FPGATrackSimFitConstantBank.cxx:22
FPGATrackSimFitConstantBank::m_npixcy
int m_npixcy
Definition: FPGATrackSimFitConstantBank.h:70
FPGATrackSimHit::setDetType
void setDetType(SiliconTech detType)
Definition: FPGATrackSimHit.h:52
FPGATrackSimHit::setX
void setX(float v)
Definition: FPGATrackSimHit.h:122
WriteCellNoiseToCool.ival
ival
Definition: WriteCellNoiseToCool.py:337
vector3D::resize
void resize(size_t x1, size_t x2, size_t x3, T const &t=T())
Definition: FPGATrackSimVectors.h:166
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
FPGATrackSimFitConstantBank::linfit_chisq
void linfit_chisq(sector_t sector, FPGATrackSimTrack &trk) const
Definition: FPGATrackSimFitConstantBank.cxx:462
SiliconTech
SiliconTech
Definition: FPGATrackSimTypes.h:25
lumiFormat.i
int i
Definition: lumiFormat.py:92
FPGATrackSimPlaneMap::getDetType
SiliconTech getDetType(size_t logiLayer) const
Definition: FPGATrackSimPlaneMap.h:82
WritePulseShapeToCool.det
det
Definition: WritePulseShapeToCool.py:204
FPGATrackSimHit::setEtaCoord
void setEtaCoord(float v)
Definition: FPGATrackSimHit.h:100
chi2
double chi2(TH1 *h0, TH1 *h1)
Definition: comparitor.cxx:522
AthMessaging
Class to provide easy MsgStream access and capabilities.
Definition: AthMessaging.h:55
find_tgc_unfilled_channelids.ip
ip
Definition: find_tgc_unfilled_channelids.py:3
FPGATrackSimTrack::getEtaCoord
float getEtaCoord(int ilayer) const
FPGATrackSimFitConstantBank::m_npars
int m_npars
Definition: FPGATrackSimFitConstantBank.h:67
FPGATrackSimPlaneMap
Definition: FPGATrackSimPlaneMap.h:62
grepfile.ic
int ic
Definition: grepfile.py:33
sector_t
int32_t sector_t
Definition: FPGATrackSimTypes.h:21
FPGATrackSimFitConstantBank::m_maj_kk
vector3D< float > m_maj_kk
Definition: FPGATrackSimFitConstantBank.h:89
remainder
std::vector< std::string > remainder(const std::vector< std::string > &v1, const std::vector< std::string > &v2)
list of entries in a vector that are not in another
Definition: compareFlatTrees.cxx:44
MTX_TOLERANCE
#define MTX_TOLERANCE
Definition: FPGATrackSimConstants.h:10
FPGATrackSimFitConstantBank::prepareInvFitConstants
void prepareInvFitConstants()
Definition: FPGATrackSimFitConstantBank.cxx:210
FPGATrackSimFitConstantBank::linfit_pars_eval
void linfit_pars_eval(sector_t sector, FPGATrackSimTrack &trk) const
Definition: FPGATrackSimFitConstantBank.cxx:494
FPGATrackSimFitConstantBank::missing_point_guess
int missing_point_guess(sector_t sector, FPGATrackSimTrack &track, bool isFirstStage, bool doExtrapolation) const
Definition: FPGATrackSimFitConstantBank.cxx:326
FPGATrackSimTrack::setQOverPt
void setQOverPt(float v)
Definition: FPGATrackSimTrack.h:83
query_example.col
col
Definition: query_example.py:7
JetVoronoiDiagramHelpers::coord
double coord
Definition: JetVoronoiDiagramHelpers.h:45
FPGATrackSimTrack::setEta
void setEta(float v)
Definition: FPGATrackSimTrack.h:87
FPGATrackSimFitConstantBank::m_pmap
FPGATrackSimPlaneMap const * m_pmap
Definition: FPGATrackSimFitConstantBank.h:60
FPGATrackSimHit::setZ
void setZ(float v)
Definition: FPGATrackSimHit.h:124
FPGATrackSimFitConstantBank::linfit
bool linfit(sector_t sector, FPGATrackSimTrack &track, bool isSecondStage) const
Definition: FPGATrackSimFitConstantBank.cxx:303
HitType::mapped
@ mapped
FPGATrackSimFitConstantBank::m_invfit_consts
std::vector< Eigen::MatrixXf, Eigen::aligned_allocator< Eigen::MatrixXf > > m_invfit_consts
Definition: FPGATrackSimFitConstantBank.h:92
python.AthDsoLogger.fname
string fname
Definition: AthDsoLogger.py:67
FPGATrackSimFitConstantBank::setIdealCoordFit
void setIdealCoordFit(bool v)
Definition: FPGATrackSimFitConstantBank.h:56
FPGATrackSimFitConstantBank::m_sector_good
std::vector< bool > m_sector_good
Definition: FPGATrackSimFitConstantBank.h:82
a
TList * a
Definition: liststreamerinfos.cxx:10
FPGATrackSimHit::setLayer
void setLayer(unsigned v)
Definition: FPGATrackSimHit.h:87
FPGATrackSimTrack::setChi2
void setChi2(float v)
Definition: FPGATrackSimTrack.h:88
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
FPGATrackSimFitConstantBank::calculateMajority
void calculateMajority()
Definition: FPGATrackSimFitConstantBank.cxx:168
vector2D::resize
void resize(size_t x1, size_t x2, T const &t=T())
Definition: FPGATrackSimVectors.h:51
FPGATrackSimFitConstantBank::m_WCs
vector2D< int > m_WCs
Definition: FPGATrackSimFitConstantBank.h:94
xAOD::track
@ track
Definition: TrackingPrimitives.h:512
FPGATrackSimFitConstantBank::m_maj_invkk
vector3D< float > m_maj_invkk
Definition: FPGATrackSimFitConstantBank.h:90
FPGATrackSimFitConstantBank::m_kaverage
vector2D< float > m_kaverage
Definition: FPGATrackSimFitConstantBank.h:86
FPGATrackSimTrack::getHoughX
float getHoughX() const
Definition: FPGATrackSimTrack.h:30
FPGATrackSimFitConstantBank::m_isIdealCoordFit
bool m_isIdealCoordFit
Definition: FPGATrackSimFitConstantBank.h:73
WriteCalibToCool.coef
coef
Definition: WriteCalibToCool.py:582
SiliconTech::pixel
@ pixel
FPGATrackSimTrack.h
fitman.k
k
Definition: fitman.py:528
FPGATrackSimHit::setHitType
void setHitType(HitType type)
Definition: FPGATrackSimHit.h:51
mapkey::key
key
Definition: TElectronEfficiencyCorrectionTool.cxx:37