FEDRA emulsion software from the OPERA Collaboration
Track Struct Reference

#include <Track.h>

Inheritance diagram for Track:
Collaboration diagram for Track:

Public Member Functions

void addGrain (Grain gr)
 
bool evaluateEstimators ()
 
unsigned int getAreaSum () const
 
double getBottomZ () const
 
unsigned short getCount () const
 
double getErrX () const
 
double getErrY () const
 
double getErrZ () const
 
std::vector< GraingetGrains () const
 
double getSigma () const
 
double getSX () const
 
double getSY () const
 
double getSZ () const
 
double getTopZ () const
 
double getX () const
 
double getY () const
 
double getZ () const
 
void setAreaSum (unsigned int areaSum)
 
void setBottomZ (double z)
 
void setCount (unsigned short count)
 
void setErrX (double errx)
 
void setErrY (double erry)
 
void setErrZ (double errz)
 
void setSigma (double sigma)
 
void setSX (double sx)
 
void setSY (double sy)
 
void setSZ (double sz)
 
void setTopZ (double z)
 
void setX (double x)
 
void setY (double y)
 
void setZ (double z)
 
 Track ()
 
 ~Track ()
 

Public Attributes

double dz
 
int Field
 
float FirstZ
 
TVector Intercept
 
TVector InterceptErrors
 
float LastZ
 
double meanDeltaTheta3D
 
double meanDeltaThetaXZ
 
double meanDeltaThetaYZ
 
double meanDistanceClustersTo3DLine
 
double meanGapClusterToCluster
 
double mx
 
double mxSigma
 
double my
 
double mySigma
 
BYTE * pCorrection
 
char * pCorrection
 
unsigned PointsN
 
TVectorpPoints
 
double qx
 
double qy
 
double rx
 
double ry
 
float Sigma
 
double sigmaDeltaTheta3D
 
double sigmaDeltaThetaXZ
 
double sigmaDeltaThetaYZ
 
TVector Slope
 
TVector SlopeErrors
 
boolean Valid
 
bool Valid
 
double xChi2
 
double yChi2
 

Private Member Functions

double pointTo3DlineDistance (double x0, double y0, double z0)
 

Private Attributes

unsigned int _areaSum
 
double _bottomZ
 
unsigned short _count
 
double _errx
 
double _erry
 
double _errz
 
std::vector< Grain_grains
 
double _sigma
 
double _sx
 
double _sy
 
double _sz
 
double _topZ
 
double _x
 
double _y
 
double _z
 

Constructor & Destructor Documentation

◆ Track()

Track::Track ( )
6{
9 mx = qx = rx = xChi2 = 0;
10 my = qy = ry = yChi2 = 0;
11 dz = 0;
12 // to be set as input...
13 _errx = _erry = 0.3;
14 _errz = 3.;
15}
double qy
Definition: Track.h:45
double sigmaDeltaThetaYZ
Definition: Track.h:42
double meanDeltaThetaXZ
Definition: Track.h:39
double meanDeltaTheta3D
Definition: Track.h:34
double meanGapClusterToCluster
Definition: Track.h:43
double xChi2
Definition: Track.h:44
double _errz
Definition: Track.h:19
double _erry
Definition: Track.h:18
double mx
Definition: Track.h:44
double sigmaDeltaThetaXZ
Definition: Track.h:40
double my
Definition: Track.h:45
double _errx
Definition: Track.h:17
double sigmaDeltaTheta3D
Definition: Track.h:38
double yChi2
Definition: Track.h:45
double ry
Definition: Track.h:45
double rx
Definition: Track.h:44
double meanDistanceClustersTo3DLine
Definition: Track.h:43
double dz
Definition: Track.h:46
double qx
Definition: Track.h:44
double meanDeltaThetaYZ
Definition: Track.h:41

◆ ~Track()

Track::~Track ( )
inline
34{};

Member Function Documentation

◆ addGrain()

void Track::addGrain ( Grain  gr)
inline
65{_grains.push_back(gr);}
std::vector< Grain > _grains
Definition: Track.h:26

◆ evaluateEstimators()

bool Track::evaluateEstimators ( )
34{
35 //std::ofstream out("test.txt",std::ios::app); // used now for debugging --> later to be removed
36
37 if (_grains.size() == 0)
38 return false;
39
40 unsigned short nGrains = _count;
41 //double s = 0, sx = 0, sy = 0, sz = 0, sxx = 0, syy = 0, szz = 0, szx = 0, szy = 0;
42 //double d = _errz/_errx;
43
44 Fitter fitXZ; fitXZ.setErrXY(3.,0.33); //error on z coordinate is ~3 um, in x or y is 0.33 (pixel dimension)
45 Fitter fitYZ; fitYZ.setErrXY(3.,0.33);
46 for (unsigned short i = 0; i < nGrains; i++)
47 {
48 double x = _grains.at(i).x;
49 double y = _grains.at(i).y;
50 double z = _grains.at(i).z;
51
52 fitXZ.addPoint(z,x);
53 fitYZ.addPoint(z,y);
54
56
57 if (i < nGrains -1)
58 {
59 double dx = _grains.at(i+1).x-x;
60 double dy = _grains.at(i+1).y-y;
61 double dz = _grains.at(i+1).z-z;
62 meanGapClusterToCluster += sqrt(dx*dx+dy*dy+dz*dz);
63 }
64
65 if (i < nGrains-2)
66 {
67 double dx1 = _grains.at(i+1).x-x;
68 double dy1 = _grains.at(i+1).y-y;
69 double dz1 = _grains.at(i+1).z-z;
70 double dx2 = _grains.at(i+2).x-_grains.at(i+1).x;
71 double dy2 = _grains.at(i+2).y-_grains.at(i+1).y;
72 double dz2 = _grains.at(i+2).z-_grains.at(i+1).z;
73
74 double sx1 = dx1/dz1; // tan_x
75 double sy1 = dy1/dz1; // tan_y
76
77 double sx2 = dx2/dz2; // tan_x
78 double sy2 = dy2/dz2; // tan_y
79
80 double valueXZ = acos((dx1*dx2+dz1*dz2)/(sqrt(dx1*dx1+dz1*dz1)*sqrt(dx2*dx2+dz2*dz2)));
81 double valueYZ = acos((dy1*dy2+dz1*dz2)/(sqrt(dy1*dy1+dz1*dz1)*sqrt(dy2*dy2+dz2*dz2)));
82 double value3D = acos((dx1*dx2+dy1*dy2+dz1*dz2)/(sqrt(dx1*dx1+dy1*dy1+dz1*dz1)*sqrt(dx2*dx2+dy2*dy2+dz2*dz2)));
83 //double value3D = acos((sx1*sx2+sy1*sy2+1)/(sqrt(sx1*sx1+sy1*sy1+1)*sqrt(sx2*sx2+sy2*sy2+1)));
84
85 double signXZ = 1.;
86 double signYZ = 1.;
87
88 if (sx2<sx1)
89 signXZ*=-1;
90 if (sy2<sy1)
91 signYZ*=-1;
92
93 meanDeltaThetaXZ+= valueXZ * signXZ;
94 meanDeltaThetaYZ+= valueYZ * signYZ;
95 meanDeltaTheta3D+= value3D;
96
97 sigmaDeltaThetaXZ+=valueXZ*valueXZ;
98 sigmaDeltaThetaYZ+=valueYZ*valueYZ;
99 sigmaDeltaTheta3D+=value3D*value3D;
100
101 }
102 }
103 meanDeltaThetaXZ/=(nGrains-2.);
104 meanDeltaThetaYZ/=(nGrains-2.);
105 meanDeltaTheta3D/=(nGrains-2.);
106
107 meanGapClusterToCluster/=(nGrains-1.);
108
109 sigmaDeltaThetaXZ/=(nGrains-2.);
110 sigmaDeltaThetaYZ/=(nGrains-2.);
111 sigmaDeltaTheta3D/=(nGrains-2.);
112
116
120
121
122
123 // fit line grains on plane x-z
124 fitXZ.evaluateLinearFit();
125 mx = fitXZ.getSlope();
126 qx = fitXZ.getIntercept();
127 rx = fitXZ.getPearsonCoefficent();
128
129
130 // fit line grains on plane y-z
131 fitYZ.evaluateLinearFit();
132 my = fitYZ.getSlope();
133 qy = fitYZ.getIntercept();
134 ry = fitYZ.getPearsonCoefficent();
135
136 xChi2 = fitXZ.getChi2();
137 yChi2 = fitYZ.getChi2();
138
139 // lenght of the segment (in z).
140 dz = _grains.at(0).z-_grains.at(nGrains-1).z;
142
143 return true;
144}
Definition: fitter.h:9
void addPoint(double x, double y)
Definition: fitter.cpp:15
double getIntercept()
Definition: fitter.h:24
bool evaluateLinearFit()
Definition: fitter.cpp:28
double getSlope()
Definition: fitter.h:23
double getPearsonCoefficent()
Definition: fitter.h:25
void setErrXY(double errx, double erry)
Definition: fitter.cpp:22
double getChi2()
Definition: fitter.h:26
double pointTo3DlineDistance(double x0, double y0, double z0)
Definition: Track.cpp:17
unsigned short _count
Definition: Track.h:13

◆ getAreaSum()

unsigned int Track::getAreaSum ( ) const
inline
67{return _areaSum;}
unsigned int _areaSum
Definition: Track.h:12

◆ getBottomZ()

double Track::getBottomZ ( ) const
inline
80{return _bottomZ;}
double _bottomZ
Definition: Track.h:25

◆ getCount()

unsigned short Track::getCount ( ) const
inline
68{return _count;}

◆ getErrX()

double Track::getErrX ( ) const
inline
72{return _errx;}

◆ getErrY()

double Track::getErrY ( ) const
inline
73{return _erry;}

◆ getErrZ()

double Track::getErrZ ( ) const
inline
74{return _errz;}

◆ getGrains()

std::vector< Grain > Track::getGrains ( ) const
inline
81{return _grains;}

◆ getSigma()

double Track::getSigma ( ) const
inline
78{return _sigma;}
double _sigma
Definition: Track.h:23

◆ getSX()

double Track::getSX ( ) const
inline
75{return _sx;}
double _sx
Definition: Track.h:20

◆ getSY()

double Track::getSY ( ) const
inline
76{return _sy;}
double _sy
Definition: Track.h:21

◆ getSZ()

double Track::getSZ ( ) const
inline
77{return _sz;}
double _sz
Definition: Track.h:22

◆ getTopZ()

double Track::getTopZ ( ) const
inline
79{return _topZ;}
double _topZ
Definition: Track.h:24

◆ getX()

double Track::getX ( ) const
inline
69{return _x;}
double _x
Definition: Track.h:14

◆ getY()

double Track::getY ( ) const
inline
70{return _y;}
double _y
Definition: Track.h:15

◆ getZ()

double Track::getZ ( ) const
inline
71{return _z;}
double _z
Definition: Track.h:16

◆ pointTo3DlineDistance()

double Track::pointTo3DlineDistance ( double  x0,
double  y0,
double  z0 
)
private
18{
19 double dist = 0;
20 double qp_x = _x-x0;
21 double qp_y = _y-y0;
22 double qp_z = _z-z0;
23 double vx = _sx;
24 double vy = _sy;
25 double vz = 1;
26 double qr = (qp_x*vx+qp_y*vy+qp_z*vz)/sqrt(vx*vx+vy*vy+vz*vz);
27 double qp =sqrt(qp_x*qp_x+qp_y*qp_y+qp_z*qp_z);
28 dist = sqrt(qp*qp-qr*qr);
29 //std::cout << qr << "\t" << qp << "\t" << dist << std::endl;
30 return dist;
31}
brick z0
Definition: RecDispMC.C:106

◆ setAreaSum()

void Track::setAreaSum ( unsigned int  areaSum)
inline
49{_areaSum = areaSum;}

◆ setBottomZ()

void Track::setBottomZ ( double  z)
inline
64{_bottomZ = z;}

◆ setCount()

void Track::setCount ( unsigned short  count)
inline
50{_count = count;}

◆ setErrX()

void Track::setErrX ( double  errx)
inline
55{_errx = errx;}

◆ setErrY()

void Track::setErrY ( double  erry)
inline
56{_erry = erry;}

◆ setErrZ()

void Track::setErrZ ( double  errz)
inline
57{_errz = errz;}

◆ setSigma()

void Track::setSigma ( double  sigma)
inline
62{_sigma = sigma;}

◆ setSX()

void Track::setSX ( double  sx)
inline
59{_sx = sx;}

◆ setSY()

void Track::setSY ( double  sy)
inline
60{_sy = sy;}

◆ setSZ()

void Track::setSZ ( double  sz)
inline
61{_sz = sz;}

◆ setTopZ()

void Track::setTopZ ( double  z)
inline
63{_topZ = z;}

◆ setX()

void Track::setX ( double  x)
inline
51{_x = x;}

◆ setY()

void Track::setY ( double  y)
inline
52{_y = y;}

◆ setZ()

void Track::setZ ( double  z)
inline
53{_z = z;}

Member Data Documentation

◆ _areaSum

unsigned int Track::_areaSum
private

◆ _bottomZ

double Track::_bottomZ
private

◆ _count

unsigned short Track::_count
private

◆ _errx

double Track::_errx
private

◆ _erry

double Track::_erry
private

◆ _errz

double Track::_errz
private

◆ _grains

std::vector<Grain> Track::_grains
private

◆ _sigma

double Track::_sigma
private

◆ _sx

double Track::_sx
private

◆ _sy

double Track::_sy
private

◆ _sz

double Track::_sz
private

◆ _topZ

double Track::_topZ
private

◆ _x

double Track::_x
private

◆ _y

double Track::_y
private

◆ _z

double Track::_z
private

◆ dz

double Track::dz

◆ Field

int Track::Field

◆ FirstZ

float Track::FirstZ

◆ Intercept

TVector Track::Intercept

◆ InterceptErrors

TVector Track::InterceptErrors

◆ LastZ

float Track::LastZ

◆ meanDeltaTheta3D

double Track::meanDeltaTheta3D

◆ meanDeltaThetaXZ

double Track::meanDeltaThetaXZ

◆ meanDeltaThetaYZ

double Track::meanDeltaThetaYZ

◆ meanDistanceClustersTo3DLine

double Track::meanDistanceClustersTo3DLine

◆ meanGapClusterToCluster

double Track::meanGapClusterToCluster

◆ mx

double Track::mx

◆ mxSigma

double Track::mxSigma

◆ my

double Track::my

◆ mySigma

double Track::mySigma

◆ pCorrection [1/2]

BYTE* Track::pCorrection

◆ pCorrection [2/2]

char* Track::pCorrection

◆ PointsN

unsigned Track::PointsN

◆ pPoints

TVector * Track::pPoints

◆ qx

double Track::qx

◆ qy

double Track::qy

◆ rx

double Track::rx

◆ ry

double Track::ry

◆ Sigma

float Track::Sigma

◆ sigmaDeltaTheta3D

double Track::sigmaDeltaTheta3D

◆ sigmaDeltaThetaXZ

double Track::sigmaDeltaThetaXZ

◆ sigmaDeltaThetaYZ

double Track::sigmaDeltaThetaYZ

◆ Slope

TVector Track::Slope

◆ SlopeErrors

TVector Track::SlopeErrors

◆ Valid [1/2]

boolean Track::Valid

◆ Valid [2/2]

bool Track::Valid

◆ xChi2

double Track::xChi2

◆ yChi2

double Track::yChi2

The documentation for this struct was generated from the following files: