CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Protected Types | Static Protected Member Functions | Private Attributes | Friends
reco::TrackResiduals Class Reference

#include <TrackResiduals.h>

Public Types

enum  ResidualType { X_Y_RESIDUALS, X_Y_PULLS }
 

Public Member Functions

void print (std::ostream &stream=std::cout) const
 
void print (const HitPattern &, std::ostream &stream=std::cout) const
 
double residualX (int i, const HitPattern &) const
 
double residualX (int i) const
 
double residualY (int i, const HitPattern &) const
 
double residualY (int i) const
 
void setPullXY (int idx, double pullX, double pullY)
 
void setResidualType (enum ResidualType)
 
void setResidualXY (int idx, double residualX, double residualY)
 
 TrackResiduals ()
 
 TrackResiduals (enum ResidualType)
 

Protected Types

enum  { numResiduals = 0x40 }
 number of residuals stored More...
 

Static Protected Member Functions

static unsigned char pack_pull (double)
 
static unsigned char pack_residual (double)
 
static double unpack_pull (unsigned char)
 
static double unpack_residual (unsigned char)
 

Private Attributes

unsigned char residuals_ [numResiduals]
 residuals, bitpacked two hits to a char More...
 
char residualType
 

Friends

class Trajectory
 

Detailed Description

Definition at line 15 of file TrackResiduals.h.

Member Enumeration Documentation

anonymous enum
protected

number of residuals stored

Enumerator
numResiduals 

Definition at line 45 of file TrackResiduals.h.

Enumerator
X_Y_RESIDUALS 
X_Y_PULLS 

Definition at line 22 of file TrackResiduals.h.

Constructor & Destructor Documentation

TrackResiduals::TrackResiduals ( )

Definition at line 11 of file TrackResiduals.cc.

References residuals_.

11  :
13 {
14  memset(residuals_, 0, sizeof(residuals_));
15 }
unsigned char residuals_[numResiduals]
residuals, bitpacked two hits to a char
TrackResiduals::TrackResiduals ( enum ResidualType  type)

Definition at line 17 of file TrackResiduals.cc.

References residuals_.

17  :
19 {
20  memset(residuals_, 0, sizeof(residuals_));
21 }
type
Definition: HCALResponse.h:21
unsigned char residuals_[numResiduals]
residuals, bitpacked two hits to a char

Member Function Documentation

unsigned char TrackResiduals::pack_pull ( double  pull)
staticprotected

Definition at line 132 of file TrackResiduals.cc.

References mag(), pull_char_to_double, and FWPFMaths::sgn().

Referenced by setPullXY().

133 {
134  unsigned char sgn = (pull < 0) * 0x08; // 1xxx is -abs(0xxx)
135  int mag = -1;
136  while (++mag < 8 && pull_char_to_double[mag][1] < fabs(pull));
137  return sgn + mag;
138 }
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float sgn(float val)
Definition: FWPFMaths.cc:9
static const double pull_char_to_double[8][2]
unsigned char TrackResiduals::pack_residual ( double  pull)
staticprotected

Definition at line 148 of file TrackResiduals.cc.

References mag(), pull_char_to_double, and FWPFMaths::sgn().

Referenced by setResidualXY().

149 {
150  unsigned char sgn = (pull < 0) * 0x08; // 1xxx is -abs(0xxx)
151  int mag = -1;
152  while (++mag < 8 && pull_char_to_double[mag][1] < fabs(pull)) {
153  // just counting
154  ;
155  }
156  return sgn + mag;
157 }
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float sgn(float val)
Definition: FWPFMaths.cc:9
static const double pull_char_to_double[8][2]
void TrackResiduals::print ( std::ostream &  stream = std::cout) const

Definition at line 159 of file TrackResiduals.cc.

References flags, i, numResiduals, and residuals_.

160 {
161  stream << "TrackResiduals" << std::endl;
162  std::ios_base::fmtflags flags = stream.flags();
163  stream.setf (std::ios_base::hex, std::ios_base::basefield);
164  stream.setf (std::ios_base::showbase);
165 
166  for (int i = 0; i < numResiduals; i++) {
167  unsigned char residual = residuals_[i];
168  printf("0x%x\n", residual);
169  //stream << residual << std::endl;
170  }
171 
172  stream.flags(flags);
173 }
int i
Definition: DBlmapReader.cc:9
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
unsigned char residuals_[numResiduals]
residuals, bitpacked two hits to a char
void TrackResiduals::print ( const HitPattern h,
std::ostream &  stream = std::cout 
) const

Definition at line 175 of file TrackResiduals.cc.

References reco::HitPattern::getHitPattern(), i, reco::HitPattern::numberOfHits(), residualX(), residualY(), reco::HitPattern::TRACK_HITS, and reco::HitPattern::validHitFilter().

176 {
177  stream << "TrackResiduals" << std::endl;
178  for (int i = 0; i < h.numberOfHits(HitPattern::TRACK_HITS); i++) {
179  stream << (h.validHitFilter(h.getHitPattern(HitPattern::TRACK_HITS, i)) ? "valid hit: " : "invalid hit: ")
180  << "( " << residualX(i, h) << " , " << residualY(i, h) << " )"
181  << std::endl;
182  }
183 }
int i
Definition: DBlmapReader.cc:9
double residualY(int i, const HitPattern &) const
static bool validHitFilter(uint16_t pattern)
Definition: HitPattern.h:787
double residualX(int i, const HitPattern &) const
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:515
int numberOfHits(HitCategory category) const
Definition: HitPattern.h:807
double TrackResiduals::residualX ( int  i,
const HitPattern h 
) const

get the residual of the ith hit (needs the hit pattern to figure out which hits are valid)

Definition at line 83 of file TrackResiduals.cc.

References customizeTrackingMonitorSeedNumber::idx, and index_to_hitpattern().

Referenced by FWTrackResidualDetailView::prepareData(), and print().

84 {
85  int idx = index_to_hitpattern(i, h);
86  if (idx == -999) {
87  return -999;
88  }
89 
90  return residualX(idx);
91 }
int i
Definition: DBlmapReader.cc:9
static int index_to_hitpattern(int i_hitpattern, const HitPattern &h)
double residualX(int i, const HitPattern &) const
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double TrackResiduals::residualX ( int  i) const

get the residual of the ith valid hit, with no regard for alignment with the HitPattern

Definition at line 38 of file TrackResiduals.cc.

References assert(), residuals_, residualType, unpack_pull(), unpack_residual(), X_Y_PULLS, and X_Y_RESIDUALS.

39 {
40  switch (residualType) {
41  case X_Y_RESIDUALS:
42  return unpack_residual(residuals_[i] >> 4);
43  case X_Y_PULLS:
44  return unpack_pull(residuals_[i] >> 4);
45  default:
46  assert(0);
47  }
48  return 0;
49 }
int i
Definition: DBlmapReader.cc:9
assert(m_qm.get())
static double unpack_residual(unsigned char)
unsigned char residuals_[numResiduals]
residuals, bitpacked two hits to a char
static double unpack_pull(unsigned char)
double TrackResiduals::residualY ( int  i,
const HitPattern h 
) const

Definition at line 93 of file TrackResiduals.cc.

References customizeTrackingMonitorSeedNumber::idx, and index_to_hitpattern().

Referenced by FWTrackResidualDetailView::prepareData(), and print().

94 {
95  int idx = index_to_hitpattern(i, h);
96  if (idx == -999) {
97  return -999;
98  }
99 
100  return residualY(idx);
101 }
int i
Definition: DBlmapReader.cc:9
double residualY(int i, const HitPattern &) const
static int index_to_hitpattern(int i_hitpattern, const HitPattern &h)
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double TrackResiduals::residualY ( int  i) const

Definition at line 51 of file TrackResiduals.cc.

References assert(), residuals_, residualType, unpack_pull(), unpack_residual(), X_Y_PULLS, and X_Y_RESIDUALS.

52 {
53  switch (residualType) {
54  case X_Y_RESIDUALS:
55  return unpack_residual(residuals_[i] & 0x0f);
56  case X_Y_PULLS:
57  return unpack_pull(residuals_[i] & 0x0f);
58  default:
59  assert(0);
60  }
61  return 0;
62 }
int i
Definition: DBlmapReader.cc:9
assert(m_qm.get())
static double unpack_residual(unsigned char)
unsigned char residuals_[numResiduals]
residuals, bitpacked two hits to a char
static double unpack_pull(unsigned char)
void TrackResiduals::setPullXY ( int  idx,
double  pullX,
double  pullY 
)

Definition at line 103 of file TrackResiduals.cc.

References assert(), customizeTrackingMonitorSeedNumber::idx, numResiduals, pack_pull(), residuals_, residualType, and X_Y_PULLS.

Referenced by trajectoryToResiduals().

104 {
106  if (idx >= numResiduals) {
107  edm::LogWarning("TrackResiduals") << " setting pulls over the array size.";
108  return;
109  }
110 
111  residuals_[idx] = (pack_pull(pullX) << 4) | pack_pull(pullY);
112 }
assert(m_qm.get())
unsigned char residuals_[numResiduals]
residuals, bitpacked two hits to a char
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
static unsigned char pack_pull(double)
void TrackResiduals::setResidualType ( enum ResidualType  type)

Definition at line 23 of file TrackResiduals.cc.

References residualType.

24 {
26 }
type
Definition: HCALResponse.h:21
void TrackResiduals::setResidualXY ( int  idx,
double  residualX,
double  residualY 
)

Definition at line 28 of file TrackResiduals.cc.

References assert(), customizeTrackingMonitorSeedNumber::idx, numResiduals, pack_residual(), residuals_, residualType, and X_Y_RESIDUALS.

Referenced by trajectoryToResiduals().

29 {
31  if (idx >= numResiduals) {
32  edm::LogWarning("TrackResiduals") << " setting residual over the array size.";
33  return;
34  }
36 }
assert(m_qm.get())
double residualY(int i, const HitPattern &) const
unsigned char residuals_[numResiduals]
residuals, bitpacked two hits to a char
double residualX(int i, const HitPattern &) const
static unsigned char pack_residual(double)
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double TrackResiduals::unpack_pull ( unsigned char  pull)
staticprotected

Definition at line 125 of file TrackResiduals.cc.

References mag(), pull_char_to_double, and FWPFMaths::sgn().

Referenced by residualX(), and residualY().

126 {
127  int sgn = 1 - 2 * ((pull & 0x08) >> 3);
128  unsigned char mag = pull & 0x07;
129  return sgn * (pull_char_to_double[mag][0] + pull_char_to_double[mag][1]) / 2;
130 }
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float sgn(float val)
Definition: FWPFMaths.cc:9
static const double pull_char_to_double[8][2]
double TrackResiduals::unpack_residual ( unsigned char  pull)
staticprotected

Definition at line 141 of file TrackResiduals.cc.

References mag(), pull_char_to_double, and FWPFMaths::sgn().

Referenced by residualX(), and residualY().

142 {
143  int sgn = 1 - 2 * ((pull & 0x08) >> 3);
144  unsigned char mag = pull & 0x07;
145  return sgn * (pull_char_to_double[mag][0] + pull_char_to_double[mag][1]) / 2;
146 }
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float sgn(float val)
Definition: FWPFMaths.cc:9
static const double pull_char_to_double[8][2]

Friends And Related Function Documentation

friend class Trajectory
friend

Definition at line 18 of file TrackResiduals.h.

Member Data Documentation

unsigned char reco::TrackResiduals::residuals_[numResiduals]
private

residuals, bitpacked two hits to a char

Definition at line 53 of file TrackResiduals.h.

Referenced by print(), residualX(), residualY(), setPullXY(), setResidualXY(), and TrackResiduals().

char reco::TrackResiduals::residualType
private

Definition at line 54 of file TrackResiduals.h.

Referenced by residualX(), residualY(), setPullXY(), setResidualType(), and setResidualXY().