BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
HoughHit Class Reference

#include <HoughHit.h>

Inheritance diagram for HoughHit:

Public Member Functions

 HoughHit ()
 HoughHit (const MdcDigi *const digi)
 HoughHit (const HoughHit &other)
HoughHitoperator= (const HoughHit &other)
void setTruthInfo (const MdcMcHit *&mcHit)
void conformalTrans (double x, double y, double r)
double getConformal_u (double, double, double)
double getConformal_v (double, double, double)
double getConformal_r (double, double, double)
void makeCir (int n, double phi_begin, double phi_last, double r)
CFCir getCir (int i) const
HoughHitType type () const
detectorType getDetectorType () const
const MdcDigidigi () const
double getBunchTime () const
HepPoint3D getMidPoint () const
HepPoint3D getEastPoint () const
HepPoint3D getWestPoint () const
double getMidX () const
double getMidY () const
int getLayerId () const
int getWireId () const
int getSlayerType () const
const MdcLayerlayer () const
const MdcSWirewire () const
double getCharge () const
double getDriftTime () const
double getDriftDist () const
double getU () const
double getV () const
double getR () const
double getDeltaD () const
double getFltLen () const
double driftTime () const
double driftTime (double tof, double z) const
double calDriftDist (double, int, double, double, double) const
double calDriftDist (double bunchTime, int ambig) const
int slayerType (int layer)
double getXTruth () const
double getYTruth () const
double getZTruth () const
double getDriftDistTruth () const
int getIdTruth () const
int getLrTruth () const
double getUTruth () const
double getVTruth () const
double getRTruth () const
HepPoint3D getPointTruth () const
void setDeltaD (double d)
void setFltLen (double flt)
int getCirList () const
int getStyle () const
void setCirList (int cir)
void setStyle (int sty)
void print () const
void printAll () const
void printTruth () const

Static Public Member Functions

static void setMdcCalibFunSvc (const IMdcCalibFunSvc *calibSvc)
static void setMdcGeomSvc (IMdcGeomSvc *geomSvc)
static void setBunchTime (double t0)

Static Public Attributes

static int _npart = 360

Detailed Description

Definition at line 24 of file HoughHit.h.

Constructor & Destructor Documentation

◆ HoughHit() [1/3]

HoughHit::HoughHit ( )

Definition at line 18 of file HoughHit.cxx.

18 {
19 _det = Global::m_gm;
20 _digiPtr = NULL;
21 _layerPtr = NULL;
22 _wirePtr = NULL;
23}
const MdcDetector * m_gm

Referenced by HoughHit(), HoughRecHit::HoughRecHit(), HoughRecHit::HoughRecHit(), HoughRecHit::HoughRecHit(), HoughRecHit::HoughRecHit(), and operator=().

◆ HoughHit() [2/3]

HoughHit::HoughHit ( const MdcDigi *const digi)

Definition at line 24 of file HoughHit.cxx.

24 {
25 _det = Global::m_gm;
26 _digiPtr = aDigi;
27 _id = aDigi->identify();
28 _layer = MdcID::layer( _id );
29 _wire = MdcID::wire( _id );
30 _layerPtr = _det->Layer( _id );
31 _wirePtr = _det->Wire( _id );
32 _rawTime = RawDataUtil::MdcTime( _digiPtr->getTimeChannel() );
33 _charge = _digiPtr->getChargeChannel();
34 _driftTime = driftTime();
35 _driftDist = calDriftDist( _bunchTime, 0 );
36 _slayerType = slayerType( _layer );
37 _lr = 0;
38 _rmid = _wirePtr->rMid();
39
40 assert( _mdcGeomSvc != NULL );
41 const MdcGeoWire* wire = _mdcGeomSvc->Wire( _layer, _wire );
42 assert( wire != NULL );
43 HepPoint3D eastP = wire->Backward() / 10.;
44 HepPoint3D westP = wire->Forward() / 10.;
45 _eastPoint = eastP;
46 _westPoint = westP;
47 _midPoint = ( eastP + westP ) / 2.;
48
49 _type = MIDPOINT;
50 conformalTrans( _midPoint.x(), _midPoint.y(), _driftDist );
51}
@ MIDPOINT
Definition HoughGlobal.h:18
HepGeom::Point3D< double > HepPoint3D
Definition HoughHit.h:18
void conformalTrans(double x, double y, double r)
Definition HoughHit.cxx:154
double driftTime() const
Definition HoughHit.cxx:170
int slayerType(int layer)
Definition HoughHit.cxx:231
double calDriftDist(double, int, double, double, double) const
Definition HoughHit.cxx:194
const MdcSWire * wire() const
Definition HoughHit.h:64
static int layer(const Identifier &id)
Values of different levels (failure returns 0).
Definition MdcID.cxx:47
static int wire(const Identifier &id)
Definition MdcID.cxx:52
static double MdcTime(int timeChannel)

◆ HoughHit() [3/3]

HoughHit::HoughHit ( const HoughHit & other)

Definition at line 53 of file HoughHit.cxx.

54 : _digiPtr( other._digiPtr )
55 , _det( other._det )
56 , _layerPtr( other._layerPtr )
57 , _wirePtr( other._wirePtr )
58 , _id( other._id )
59 , _layer( other._layer )
60 , _wire( other._wire )
61 , _rawTime( other._rawTime )
62 , _charge( other._charge )
63 , _eastPoint( other._eastPoint )
64 , _westPoint( other._westPoint )
65 , _midPoint( other._midPoint )
66 , _u( other._u )
67 , _v( other._v )
68 , _type( other._type )
69 , _detectorType( other._detectorType )
70 , _rmid( other._rmid )
71 , _slayerType( other._slayerType )
72 , _cirlist( other._cirlist )
73 , _style( other._style )
74 , _driftTime( other._driftTime )
75 , _driftDist( other._driftDist )
76 , CF_drift( other.CF_drift )
77 , vec_cfcir( other.vec_cfcir )
78 ,
79 // truth
80 _truthU( other._truthU )
81 , _truthV( other._truthV )
82 , _truthR( other._truthR )
83 , _truthlr( other._truthlr )
84 , _truthPoint( other._truthPoint )
85 , _truthDrift( other._truthDrift )
86 , _deltad( other._deltad )
87 , _flightLength( other._flightLength ) {}
Index other(Index i, Index j)

Member Function Documentation

◆ calDriftDist() [1/2]

double HoughHit::calDriftDist ( double bunchTime,
int ambig ) const

Definition at line 188 of file HoughHit.cxx.

188 {
189
190 // double crudeTof = 0; //FIXME
191 return calDriftDist( bunchTime + crudeTof(), ambig, 0., 0., 0. );
192}

◆ calDriftDist() [2/2]

double HoughHit::calDriftDist ( double tof,
int ambig,
double entranceAngle,
double ,
double z ) const

Definition at line 194 of file HoughHit.cxx.

195 {
196
197 double driftD;
198 // drift time ns, layer id begin with 0, entrance angle rads,
199 // lr ambig: wire ambig 1,-1,0 -> Calib 0,1,2
200 int lrCalib = 2;
201 if ( ambig == 1 ) lrCalib = 0;
202 else if ( ambig == -1 ) lrCalib = 1;
203
204 // tof in s, driftDist in cm, dirftTime in ns
205 if ( fabs( z ) > 150. || fabs( driftTime( tof, z ) ) > 1500. ) { return 9999.; }
206 driftD = 0.1 * _calibPtr->driftTimeToDist( driftTime( tof, z ), _layer, _wire, lrCalib,
207 entranceAngle ); // to cm
208 // std::cout<<"driftDist "<<"("<<_layer <<","<<_wire <<") dd "<<driftD<<" dt
209 // "<<driftTime(tof,z) <<" lr "<<lrCalib <<" eAng "<<entranceAngle <<" tof "<<tof*1.e9<<" z
210 // "<<z <<" t0walk "<<_T0Walk<<" rawT "<<_rawTime <<" tprop "<< _rawTime - driftTime(tof,z)-
211 // _T0Walk-1.e9*tof<<std::endl;
212
213 if ( abs( driftD ) < 0.00001 ) driftD = 0.00001;
214 return driftD;
215}

Referenced by calDriftDist(), and HoughHit().

◆ conformalTrans()

void HoughHit::conformalTrans ( double x,
double y,
double r )

Definition at line 154 of file HoughHit.cxx.

154 {
155 _u = 2 * x / ( x * x + y * y - r * r );
156 _v = 2 * y / ( x * x + y * y - r * r );
157 CF_drift = 2 * r / ( x * x + y * y - r * r );
158}
Double_t x[10]

Referenced by HoughHit().

◆ digi()

const MdcDigi * HoughHit::digi ( ) const
inline

Definition at line 53 of file HoughHit.h.

53{ return _digiPtr; }

Referenced by HoughTrack::cald_layer(), HoughRecHit::HoughRecHit(), printAll(), and HoughHitList::remove().

◆ driftTime() [1/2]

double HoughHit::driftTime ( ) const

Definition at line 170 of file HoughHit.cxx.

170 {
171 double tprop = _calibPtr->getTprop( _layer, 0 );
172 double T0Walk =
173 _calibPtr->getT0( _layer, _wire ) + _calibPtr->getTimeWalk( _layer, _charge );
174 // tof in ns, driftTime in ns, _T0Walk in ns
175 double driftT = _rawTime - T0Walk - 1.e9 * _bunchTime - tprop;
176 return driftT;
177}

Referenced by HoughHitList::addMdcDigiList(), calDriftDist(), HoughTrack::find_pair_hit(), HoughTrack::find_stereo_hit(), and HoughHit().

◆ driftTime() [2/2]

double HoughHit::driftTime ( double tof,
double z ) const

Definition at line 178 of file HoughHit.cxx.

178 {
179 double tprop = _calibPtr->getTprop( _layer, z * 10. );
180 double T0Walk =
181 _calibPtr->getT0( _layer, _wire ) + _calibPtr->getTimeWalk( _layer, _charge );
182 // tof in ns, driftTime in ns, _T0Walk in ns
183 double driftT = _rawTime - T0Walk - 1.e9 * tof - tprop;
184
185 return driftT;
186}

◆ getBunchTime()

double HoughHit::getBunchTime ( ) const
inline

Definition at line 54 of file HoughHit.h.

54{ return _bunchTime; }

◆ getCharge()

double HoughHit::getCharge ( ) const
inline

Definition at line 65 of file HoughHit.h.

65{ return _charge; }

◆ getCir()

CFCir HoughHit::getCir ( int i) const
inline

Definition at line 47 of file HoughHit.h.

47{ return vec_cfcir[i]; }

◆ getCirList()

int HoughHit::getCirList ( ) const
inline

Definition at line 102 of file HoughHit.h.

102{ return _cirlist; }

Referenced by HoughTrack::cald_layer(), and HoughMap::select_slant().

◆ getConformal_r()

double HoughHit::getConformal_r ( double x,
double y,
double r )

Definition at line 166 of file HoughHit.cxx.

166 {
167 return 2 * r / ( x * x + y * y - r * r );
168}

Referenced by setTruthInfo().

◆ getConformal_u()

double HoughHit::getConformal_u ( double x,
double y,
double r )

Definition at line 160 of file HoughHit.cxx.

160 {
161 return 2 * x / ( x * x + y * y - r * r );
162}

Referenced by setTruthInfo().

◆ getConformal_v()

double HoughHit::getConformal_v ( double x,
double y,
double r )

Definition at line 163 of file HoughHit.cxx.

163 {
164 return 2 * x / ( x * x + y * y - r * r );
165}

Referenced by setTruthInfo().

◆ getDeltaD()

double HoughHit::getDeltaD ( ) const
inline

Definition at line 73 of file HoughHit.h.

73{ return _deltad; } // l1l2-R+-r

◆ getDetectorType()

detectorType HoughHit::getDetectorType ( ) const
inline

Definition at line 52 of file HoughHit.h.

52{ return _detectorType; }

◆ getDriftDist()

double HoughHit::getDriftDist ( ) const
inline

Definition at line 67 of file HoughHit.h.

67{ return _driftDist; }

Referenced by HoughTrack::cald_layer(), and HoughPeak::collectHits().

◆ getDriftDistTruth()

double HoughHit::getDriftDistTruth ( ) const
inline

Definition at line 87 of file HoughHit.h.

87{ return _truthDrift; }

◆ getDriftTime()

double HoughHit::getDriftTime ( ) const
inline

Definition at line 66 of file HoughHit.h.

66{ return _driftTime; }

◆ getEastPoint()

HepPoint3D HoughHit::getEastPoint ( ) const
inline

Definition at line 56 of file HoughHit.h.

56{ return _eastPoint; }

◆ getFltLen()

double HoughHit::getFltLen ( ) const
inline

Definition at line 74 of file HoughHit.h.

74{ return _flightLength; } // R*dtheta

◆ getIdTruth()

int HoughHit::getIdTruth ( ) const
inline

Definition at line 88 of file HoughHit.h.

88{ return _truthId; }

◆ getLayerId()

int HoughHit::getLayerId ( ) const
inline

◆ getLrTruth()

int HoughHit::getLrTruth ( ) const
inline

Definition at line 89 of file HoughHit.h.

89{ return _truthlr; }

◆ getMidPoint()

HepPoint3D HoughHit::getMidPoint ( ) const
inline

Definition at line 55 of file HoughHit.h.

55{ return _midPoint; }

◆ getMidX()

double HoughHit::getMidX ( ) const
inline

Definition at line 58 of file HoughHit.h.

58{ return _midPoint.x(); }

Referenced by HoughTrack::cald_layer(), and HoughPeak::collectHits().

◆ getMidY()

double HoughHit::getMidY ( ) const
inline

Definition at line 59 of file HoughHit.h.

59{ return _midPoint.y(); }

Referenced by HoughTrack::cald_layer(), and HoughPeak::collectHits().

◆ getPointTruth()

HepPoint3D HoughHit::getPointTruth ( ) const
inline

Definition at line 93 of file HoughHit.h.

93{ return _truthPoint; }

◆ getR()

double HoughHit::getR ( ) const
inline

Definition at line 72 of file HoughHit.h.

72{ return CF_drift; } // drift in CFS

◆ getRTruth()

double HoughHit::getRTruth ( ) const
inline

Definition at line 92 of file HoughHit.h.

92{ return _truthR; }

◆ getSlayerType()

int HoughHit::getSlayerType ( ) const
inline

◆ getStyle()

int HoughHit::getStyle ( ) const
inline

Definition at line 103 of file HoughHit.h.

103{ return _style; }

Referenced by HoughTrack::cald_layer(), and HoughMap::select_slant().

◆ getU()

double HoughHit::getU ( ) const
inline

Definition at line 70 of file HoughHit.h.

70{ return _u; }

Referenced by HoughMap::select_slant().

◆ getUTruth()

double HoughHit::getUTruth ( ) const
inline

Definition at line 90 of file HoughHit.h.

90{ return _truthU; }

Referenced by HoughTrack::cald_layer().

◆ getV()

double HoughHit::getV ( ) const
inline

Definition at line 71 of file HoughHit.h.

71{ return _v; }

Referenced by HoughMap::select_slant().

◆ getVTruth()

double HoughHit::getVTruth ( ) const
inline

Definition at line 91 of file HoughHit.h.

91{ return _truthV; }

Referenced by HoughTrack::cald_layer().

◆ getWestPoint()

HepPoint3D HoughHit::getWestPoint ( ) const
inline

Definition at line 57 of file HoughHit.h.

57{ return _westPoint; }

◆ getWireId()

int HoughHit::getWireId ( ) const
inline

Definition at line 61 of file HoughHit.h.

61{ return _wire; }

Referenced by HoughTrack::find_pair_hit(), HoughTrack::find_stereo_hit(), and small_layer().

◆ getXTruth()

double HoughHit::getXTruth ( ) const
inline

Definition at line 84 of file HoughHit.h.

84{ return _truthPoint.x(); }

◆ getYTruth()

double HoughHit::getYTruth ( ) const
inline

Definition at line 85 of file HoughHit.h.

85{ return _truthPoint.y(); }

◆ getZTruth()

double HoughHit::getZTruth ( ) const
inline

Definition at line 86 of file HoughHit.h.

86{ return _truthPoint.z(); }

◆ layer()

const MdcLayer * HoughHit::layer ( ) const
inline

Definition at line 63 of file HoughHit.h.

63{ return _layerPtr; }

Referenced by slayerType().

◆ makeCir()

void HoughHit::makeCir ( int n,
double phi_begin,
double phi_last,
double r )

Definition at line 243 of file HoughHit.cxx.

243 {
244 vec_cfcir.clear();
245 vector<CFCir>().swap( vec_cfcir );
246 for ( int i = 0; i < n; i++ )
247 {
248 double phi_slice = ( phi_last - phi_begin ) / n;
249 double phi = phi_begin + ( 1 / 2. + i ) * phi_slice;
250 double x = _u + r * cos( phi );
251 double y = _v + r * sin( phi );
252 vec_cfcir.push_back( CFCir( x, y, phi, n, _u, _v, r ) );
253 }
254}
const Int_t n

◆ operator=()

HoughHit & HoughHit::operator= ( const HoughHit & other)

Definition at line 89 of file HoughHit.cxx.

89 {
90 _digiPtr = ( other._digiPtr );
91 _det = ( other._det );
92 _layerPtr = ( other._layerPtr );
93 _wirePtr = ( other._wirePtr );
94 _id = ( other._id );
95 _layer = ( other._layer );
96 _wire = ( other._wire );
97 _rawTime = ( other._rawTime );
98 _charge = ( other._charge );
99 _eastPoint = ( other._eastPoint );
100 _westPoint = ( other._westPoint );
101 _midPoint = ( other._midPoint );
102 _u = ( other._u );
103 _v = ( other._v );
104 _type = ( other._type );
105 _detectorType = ( other._detectorType );
106 _rmid = ( other._rmid );
107 _slayerType = ( other._slayerType );
108 _cirlist = ( other._cirlist );
109 _style = ( other._style );
110 _driftTime = ( other._driftTime );
111 _driftDist = ( other._driftDist );
112 CF_drift = ( other.CF_drift );
113 vec_cfcir = ( other.vec_cfcir );
114 // truth
115 _truthU = ( other._truthU );
116 _truthV = ( other._truthV );
117 _truthR = ( other._truthR );
118 _truthlr = ( other._truthlr );
119 _truthDrift = ( other._truthDrift );
120 _truthPoint = ( other._truthPoint );
121 _deltad = ( other._deltad );
122 _flightLength = ( other._flightLength );
123
124 return *this;
125}

Referenced by HoughRecHit::operator=().

◆ print()

void HoughHit::print ( ) const

Definition at line 217 of file HoughHit.cxx.

217 {
218 std::cout << "(" << _layer << "," << _wire << ") " << std::endl;
219}

◆ printAll()

void HoughHit::printAll ( ) const

Definition at line 221 of file HoughHit.cxx.

221 {
222 std::cout << "(" << _layer << "," << _wire << ") id " << this->digi()->getTrackIndex()
223 << " xyz " << _midPoint << endl;
224}
const MdcDigi * digi() const
Definition HoughHit.h:53
int getTrackIndex() const
Definition RawData.cxx:38

◆ printTruth()

void HoughHit::printTruth ( ) const

Definition at line 225 of file HoughHit.cxx.

225 {
226 std::cout << "(" << _layer << "," << _wire << ") truth " << _truthId << " xyz "
227 << _truthPoint << " uv " << _truthU << "," << _truthV << ") "
228 << " cir list " << _cirlist << " ambig " << _truthlr << std::endl;
229}

◆ setBunchTime()

void HoughHit::setBunchTime ( double t0)
inlinestatic

Definition at line 37 of file HoughHit.h.

37{ _bunchTime = t0; }

Referenced by MdcHoughFinder::execute(), and MdcHoughFinder::initialize().

◆ setCirList()

void HoughHit::setCirList ( int cir)
inline

Definition at line 104 of file HoughHit.h.

104{ _cirlist = cir; }

◆ setDeltaD()

void HoughHit::setDeltaD ( double d)
inline

Definition at line 98 of file HoughHit.h.

98{ _deltad = d; } // l1l2-R+-r

Referenced by HoughTrack::cald_layer().

◆ setFltLen()

void HoughHit::setFltLen ( double flt)
inline

Definition at line 99 of file HoughHit.h.

99{ _flightLength = flt; } // R*dtheta

Referenced by HoughTrack::cald_layer().

◆ setMdcCalibFunSvc()

void HoughHit::setMdcCalibFunSvc ( const IMdcCalibFunSvc * calibSvc)
inlinestatic

Definition at line 35 of file HoughHit.h.

35{ _calibPtr = calibSvc; }

Referenced by MdcHoughFinder::initialize().

◆ setMdcGeomSvc()

void HoughHit::setMdcGeomSvc ( IMdcGeomSvc * geomSvc)
inlinestatic

Definition at line 36 of file HoughHit.h.

36{ _mdcGeomSvc = geomSvc; }

Referenced by MdcHoughFinder::initialize().

◆ setStyle()

void HoughHit::setStyle ( int sty)
inline

Definition at line 105 of file HoughHit.h.

105{ _style = sty; }

◆ setTruthInfo()

void HoughHit::setTruthInfo ( const MdcMcHit *& mcHit)

Definition at line 127 of file HoughHit.cxx.

127 {
128 if ( !aMcHit ) return;
129 _truthDrift = aMcHit->getDriftDistance() / 10.; // mm to cm
130 _truthPoint.setX( aMcHit->getPositionX() / 10. ); // mm to cm
131 _truthPoint.setY( aMcHit->getPositionY() / 10. );
132 _truthPoint.setZ( aMcHit->getPositionZ() / 10. );
133
134 int mcLR = aMcHit->getPositionFlag();
135 // if (mcLR == 0) mcLR = -1;//FIXME
136 // for truth reserve from digi
137 if ( mcLR == 1 ) mcLR = -1; //
138 if ( mcLR == 0 ) mcLR = 1; //
139 _truthId = aMcHit->getTrackIndex();
140 _truthlr = mcLR;
141
142 // fix
143 _truthU = getConformal_u( _truthPoint.x(), _truthPoint.y(), _truthDrift );
144 _truthV = getConformal_v( _truthPoint.y(), _truthPoint.x(), _truthDrift );
145 _truthR = getConformal_r( _truthPoint.x(), _truthPoint.y(), _truthDrift );
146
147 // std::cout<<__FILE__<<" "<<_layer<<","<<_wire<<" "<<_truthId<<" truth "<<_truthPoint<<"
148 // truthU
149 // "<<_truthU<<std::endl;
150
151 _type = DIGIWITHTRUTH;
152}
@ DIGIWITHTRUTH
Definition HoughGlobal.h:18
double getConformal_v(double, double, double)
Definition HoughHit.cxx:163
double getConformal_r(double, double, double)
Definition HoughHit.cxx:166
double getConformal_u(double, double, double)
Definition HoughHit.cxx:160

◆ slayerType()

int HoughHit::slayerType ( int layer)

Definition at line 231 of file HoughHit.cxx.

231 {
232 // get nominal shift cell of this layer
233 double nomShift = _mdcGeomSvc->Layer( layer )->nomShift();
234
235 if ( nomShift > 0 ) return 1;
236 else if ( nomShift < 0 ) return -1;
237 else return 0;
238
239 return -999;
240}
const MdcLayer * layer() const
Definition HoughHit.h:63

Referenced by HoughHit().

◆ type()

HoughHitType HoughHit::type ( ) const
inline

Definition at line 51 of file HoughHit.h.

51{ return _type; }

◆ wire()

const MdcSWire * HoughHit::wire ( ) const
inline

Definition at line 64 of file HoughHit.h.

64{ return _wirePtr; }

Referenced by HoughHit().

Member Data Documentation

◆ _npart

int HoughHit::_npart = 360
static

Definition at line 48 of file HoughHit.h.

Referenced by MdcHoughFinder::initialize().


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