FEDRA emulsion software from the OPERA Collaboration
EdbSegP Class Reference

#include <EdbSegP.h>

Inheritance diagram for EdbSegP:
Collaboration diagram for EdbSegP:

Public Member Functions

void addEMULDigit (TObject *a)
 
Int_t Aid (int i) const
 
bool CheckCOV () const
 
Float_t Chi2 () const
 
void Clear ()
 
Int_t Compare (const TObject *obj) const
 
void Copy (const EdbSegP &s)
 
TMatrixD & COV () const
 
Float_t DeltaR (EdbSegP *seg1) const
 
Float_t DeltaTheta (EdbSegP *seg1) const
 
Float_t DZ () const
 
Float_t DZem () const
 
 EdbSegP ()
 
 EdbSegP (const EdbSegP &s)
 
 EdbSegP (int id, float x, float y, float tx, float ty, float w=0, int flag=0)
 
TRefArray * EMULDigitArray () const
 
Int_t Flag () const
 
void ForceCOV (TMatrixD &cov)
 
Int_t ID () const
 
bool IsCompatible (EdbSegP &s, float nsigx, float nsigt) const
 
Bool_t IsEqual (const TObject *obj) const
 
Bool_t IsSortable () const
 
Int_t MCEvt () const
 
Int_t MCTrack () const
 
void MergeTo (EdbSegP &s)
 
Float_t P () const
 
Float_t Phi () const
 
Int_t PID () const
 
Int_t Plate () const
 
void Print (Option_t *opt="") const
 
void PrintNice () const
 
Float_t Prob () const
 
Float_t ProbLink (EdbSegP &s1, EdbSegP &s2)
 
void PropagateTo (float z)
 
void PropagateToCOV (float z)
 
void PropagateToDZ (float dz)
 
EdbID ScanID () const
 
void Set (int id, float x, float y, float tx, float ty, float w, int flag)
 
void Set0 ()
 
void SetAid (int a, int v, int side=0)
 
void SetChi2 (float chi2)
 
void SetCOV (double *array, int dim=5)
 
void SetCOV (TMatrixD &cov)
 
void SetDZ (float dz)
 
void SetDZem (float dz)
 
void SetErrorP (float sp2)
 
void SetErrors ()
 
void SetErrors (float sx2, float sy2, float sz2, float stx2, float sty2, float sp2=1.)
 
void SetErrors0 ()
 
void SetErrorsCOV (float sx2, float sy2, float sz2, float stx2, float sty2, float sp2=1.)
 
void SetFlag (int flag)
 
void SetID (int id)
 
void SetMC (int mEvt, int mTrack)
 
void SetP (float p)
 
void SetPID (int pid)
 
void SetPlate (int plateid)
 
void SetProb (float prob)
 
void SetProbability (float p)
 
void SetScanID (EdbID id)
 
void SetSide (int side=0)
 
void SetSZ (float sz)
 
void SetTrack (int trid)
 
void SetTX (Float_t tx)
 
void SetTY (Float_t ty)
 other functions More...
 
void SetVid (int vid, int sid)
 
void SetVolume (float w)
 
void SetW (float w)
 
void SetX (Float_t x)
 
void SetY (Float_t y)
 
void SetZ (float z)
 
Int_t Side () const
 mandatory virtual functions: More...
 
Float_t SP () const
 
Float_t STX () const
 
Float_t STY () const
 
Float_t SX () const
 
Float_t SY () const
 
Float_t SZ () const
 
Float_t Theta () const
 
Int_t Track () const
 
Float_t TX () const
 tangens = deltaX/deltaZ More...
 
Float_t TY () const
 tangens = deltaY/deltaZ More...
 
Int_t Vid (int i) const
 
Float_t Volume () const
 
Float_t W () const
 
Float_t X () const
 
Float_t Y () const
 
Float_t Z () const
 
virtual ~EdbSegP ()
 void Transform(EdbAffine2D &aff) { ((EdbTrack2D*)this)->Transform(&aff); } More...
 
- Public Member Functions inherited from EdbTrack2D
virtual void Print (Option_t *opt="") const
 
virtual void Substruct (EdbTrack2D *t)
 
virtual void Test () const
 
virtual void Transform (const EdbAffine2D *a)
 
virtual ~EdbTrack2D ()
 
- Public Member Functions inherited from EdbPoint2D
virtual void Print (Option_t *opt="") const
 
virtual void SetX (float x)=0
 
virtual void SetY (float y)=0
 
virtual void SetZ (float z)
 
virtual void Substruct (EdbPoint *p)
 
virtual void Test () const
 
virtual void TestPoint2D () const
 
virtual void Transform (const EdbAffine2D *a)
 
virtual Float_t X () const =0
 
virtual Float_t Y () const =0
 
virtual Float_t Z () const
 
virtual ~EdbPoint2D ()
 
- Public Member Functions inherited from EdbPoint
virtual void SetX (float x)=0
 
virtual void SetY (float y)=0
 
virtual void SetZ (float z)=0
 
virtual void Substruct (EdbPoint *p)=0
 
virtual void Test () const
 
virtual void Transform (const EdbAffine2D *a)
 
virtual void Transform (const EdbAffine3D *a)
 
virtual Float_t X () const =0
 
virtual Float_t Y () const =0
 
virtual Float_t Z () const =0
 
virtual ~EdbPoint ()
 
- Public Member Functions inherited from EdbAngle2D
virtual void Print (Option_t *opt="") const
 
virtual void SetTX (float x)=0
 
virtual void SetTY (float y)=0
 
virtual void Substruct (const EdbAngle2D *a)
 
virtual void Test () const
 
virtual void Transform (const EdbAffine2D *a)
 
virtual Float_t TX () const =0
 tangens = deltaX/deltaZ More...
 
virtual Float_t TY () const =0
 tangens = deltaY/deltaZ More...
 
virtual ~EdbAngle2D ()
 

Static Public Member Functions

static Float_t Angle (const EdbSegP &s1, const EdbSegP &s2)
 
static Float_t Distance (const EdbSegP &s1, const EdbSegP &s2)
 
static void LinkMT (const EdbSegP *s1, const EdbSegP *s2, EdbSegP *s)
 

Public Attributes

Float_t eTX
 
Float_t eTY
 direction tangents More...
 
Float_t eX
 
Float_t eY
 
Float_t eZ
 coordinates More...
 

Protected Attributes

TMatrixD * eCOV
 covariance matrix of the parameters (x,y,tx,ty,p) More...
 

Private Attributes

Int_t eAid [2]
 [0]-AreaID, [1]-ViewID More...
 
Float_t eChi2
 chi-square More...
 
Float_t eDZ
 the length of segment along z-axis More...
 
Float_t eDZem
 the length of segment along z-axis in the emulsion More...
 
TRefArray * eEMULDigitArray
 AM+AC 27/07/07. More...
 
Int_t eFlag
 
Int_t eID
 segment id (unique in plate) More...
 
Int_t eMCEvt
 MC event number. More...
 
Int_t eMCTrack
 MC track number. More...
 
Float_t eP
 momentum of the particle More...
 
Int_t ePID
 mother pattern ID More...
 
Float_t eProb
 probability More...
 
EdbID eScanID
 Int_t ePlate; //!< plate id. More...
 
Float_t eSZ
 square of the Z-error More...
 
Int_t eTrack
 id of the track (-1) if no track, for a track object it represents the index of trk.root file (tracks tree) More...
 
Int_t eVid [2]
 [0]-view entry in the input tree, [1]-segment entry in the view More...
 
Float_t eVolume
 segment volume More...
 
Float_t eW
 weight More...
 

Detailed Description

//////////////////////////////////////////////////////////////////////// // EdbSegP // // segment with the attributes useful for processing // // ////////////////////////////////////////////////////////////////////////

Constructor & Destructor Documentation

◆ EdbSegP() [1/3]

EdbSegP::EdbSegP ( )
inline
53 {
55 Set0();
56 }
TRefArray * eEMULDigitArray
AM+AC 27/07/07.
Definition: EdbSegP.h:44
void Set0()
Definition: EdbSegP.cxx:30

◆ EdbSegP() [2/3]

EdbSegP::EdbSegP ( int  id,
float  x,
float  y,
float  tx,
float  ty,
float  w = 0,
int  flag = 0 
)
23{
25 Set0();
26 Set(id,x,y,tx,ty,w,flag);
27}
void Set(int id, float x, float y, float tx, float ty, float w, int flag)
Definition: EdbSegP.h:87
void w(int rid=2, int nviews=2)
Definition: test.C:27

◆ EdbSegP() [3/3]

EdbSegP::EdbSegP ( const EdbSegP s)
inline
58 {
60 Set0(); Copy(s);
61 }
void Copy(const EdbSegP &s)
Definition: EdbSegP.cxx:105
s
Definition: check_shower.C:55

◆ ~EdbSegP()

virtual EdbSegP::~EdbSegP ( )
inlinevirtual

void Transform(EdbAffine2D &aff) { ((EdbTrack2D*)this)->Transform(&aff); }

Member Function Documentation

◆ addEMULDigit()

void EdbSegP::addEMULDigit ( TObject *  a)
inline
77 {
78 if(!eEMULDigitArray) eEMULDigitArray = new TRefArray();
79 eEMULDigitArray->Add(a);
80 }
void a()
Definition: check_aligned.C:59

◆ Aid()

Int_t EdbSegP::Aid ( int  i) const
inline
169{if(i==0) return eAid[i]; else if(i==1) return eAid[i]%100000; else return -1;}
Int_t eAid[2]
[0]-AreaID, [1]-ViewID
Definition: EdbSegP.h:27

◆ Angle()

Float_t EdbSegP::Angle ( const EdbSegP s1,
const EdbSegP s2 
)
static
445{
446 return EdbMath::Angle3(s1.TX(),s1.TY(),s2.TX(),s2.TY());
447}
static double Angle3(float tx1, float ty1, float tx2, float ty2)
Definition: EdbMath.cxx:42
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ CheckCOV()

bool EdbSegP::CheckCOV ( ) const
inline

< Check Covariance Matrix, is there?

105 {
106 if (eCOV) return true;
107 else return false;
108 }
TMatrixD * eCOV
covariance matrix of the parameters (x,y,tx,ty,p)
Definition: EdbSegP.h:49

◆ Chi2()

Float_t EdbSegP::Chi2 ( ) const
inline
157{return eChi2;}
Float_t eChi2
chi-square
Definition: EdbSegP.h:35

◆ Clear()

void EdbSegP::Clear ( )
inline
86{ eCOV->Clear(); }

◆ Compare()

int EdbSegP::Compare ( const TObject *  obj) const

421{
422 const EdbSegP *s=(EdbSegP*)obj;
423 double f1=0, f2=0;
424
425 if ( Abs(s->Z()- Z()) >0.000001 ) { f1 = Z(); f2 = s->Z(); } // usual case for tracks
426 else if( Abs(s->X()- X()) >0.000001 ) { f1 = X(); f2 = s->X(); }
427 else if( Abs(s->Y()- Y()) >0.000001 ) { f1 = Y(); f2 = s->Y(); }
428 else if( Abs(s->TX()- TX())>0.000001 ) { f1 = TX(); f2 = s->TX(); }
429 else if( Abs(s->TY()- TY())>0.000001 ) { f1 = TY(); f2 = s->TY(); }
430 else if( Abs(s->W()- W()) >0.000001 ) { f1 = W(); f2 = W(); }
431
432 if (f1>f2) return 1;
433 else if (f1<f2) return -1;
434 else return 0;
435}
Definition: EdbSegP.h:21
Float_t X() const
Definition: EdbSegP.h:173
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174
Float_t W() const
Definition: EdbSegP.h:151

◆ Copy()

void EdbSegP::Copy ( const EdbSegP s)

106{
107 SetPID(s.PID());
108 Set(s.ID(),s.X(),s.Y(),s.TX(),s.TY(),s.W(),s.Flag());
109 SetZ(s.Z());
110 if(s.CheckCOV()) SetCOV(s.COV());
111 SetSZ(s.SZ());
112 SetVid(s.Vid(0),s.Vid(1));
113 SetAid(s.Aid(0),s.Aid(1),s.Side());
114 SetProb(s.Prob());
115 SetVolume(s.Volume());
116 SetDZ(s.DZ());
117 SetDZem(s.DZem());
118 SetP(s.P());
119 SetChi2(s.Chi2());
120 SetTrack(s.Track());
121 SetMC( s.MCEvt(), s.MCTrack() );
122 SetScanID(s.ScanID());
123 if(s.eEMULDigitArray) {
124 eEMULDigitArray = new TRefArray(*s.eEMULDigitArray);
125 printf("copy: %d %d\n",s.eEMULDigitArray->GetEntries(),eEMULDigitArray->GetEntries());
126 }
127}
void SetVolume(float w)
Definition: EdbSegP.h:136
void SetPID(int pid)
Definition: EdbSegP.h:129
void SetProb(float prob)
Definition: EdbSegP.h:134
void SetScanID(EdbID id)
Definition: EdbSegP.h:143
void SetDZem(float dz)
Definition: EdbSegP.h:127
void SetTrack(int trid)
Definition: EdbSegP.h:131
void SetZ(float z)
Definition: EdbSegP.h:125
void SetCOV(TMatrixD &cov)
Definition: EdbSegP.h:99
void SetChi2(float chi2)
Definition: EdbSegP.h:135
void SetSZ(float sz)
Definition: EdbSegP.h:124
void SetP(float p)
Definition: EdbSegP.h:133
void SetMC(int mEvt, int mTrack)
Definition: EdbSegP.h:141
void SetDZ(float dz)
Definition: EdbSegP.h:126
void SetVid(int vid, int sid)
Definition: EdbSegP.h:137
void SetAid(int a, int v, int side=0)
Definition: EdbSegP.h:138

◆ COV()

TMatrixD & EdbSegP::COV ( ) const
inline
123{return *eCOV;}

◆ DeltaR()

Float_t EdbSegP::DeltaR ( EdbSegP seg1) const
inline
191 {
192 return TMath::Sqrt( TMath::Power(X()-seg1->X(),2)+TMath::Power(Y()-seg1->Y(),2) );
193 }

◆ DeltaTheta()

Float_t EdbSegP::DeltaTheta ( EdbSegP seg1) const
inline
188 {
189 return TMath::Sqrt( TMath::Power(TX()-seg1->TX(),2)+TMath::Power(TY()-seg1->TY(),2) );
190 }

◆ Distance()

Float_t EdbSegP::Distance ( const EdbSegP s1,
const EdbSegP s2 
)
static
437 {
438 double dx=s1.X()-s2.X();
439 double dy=s1.Y()-s2.Y();
440 double dz=s1.Z()-s2.Z();
441 return TMath::Sqrt(dx*dx+dy*dy+dz*dz);
442}
brick dz
Definition: RecDispMC.C:107

◆ DZ()

Float_t EdbSegP::DZ ( ) const
inline
154{return eDZ;}
Float_t eDZ
the length of segment along z-axis
Definition: EdbSegP.h:39

◆ DZem()

Float_t EdbSegP::DZem ( ) const
inline
155{return eDZem;}
Float_t eDZem
the length of segment along z-axis in the emulsion
Definition: EdbSegP.h:40

◆ EMULDigitArray()

TRefArray * EdbSegP::EMULDigitArray ( ) const
inline
82{ return eEMULDigitArray;}

◆ Flag()

Int_t EdbSegP::Flag ( ) const
inline
149{return eFlag;}
Int_t eFlag
Definition: EdbSegP.h:28

◆ ForceCOV()

void EdbSegP::ForceCOV ( TMatrixD &  cov)
inline

< to correctly copy COV matrices in MakeTracksTree

110 {
111 if(!(&cov)) return;
112 if(eCOV) *eCOV = cov;
113 else eCOV = new TMatrixD(cov);
114 }

◆ ID()

Int_t EdbSegP::ID ( ) const
inline
147{return eID;}
Int_t eID
segment id (unique in plate)
Definition: EdbSegP.h:25

◆ IsCompatible()

bool EdbSegP::IsCompatible ( EdbSegP s,
float  nsigx,
float  nsigt 
) const

return true if segments are closer then nsig sigma in all coordinates assumed that z is the same

301{
304 float dtx=TX()-s.TX();
305 if( dtx*dtx > STX()*nsigt*nsigt ) return false;
306 float dty=TY()-s.TY();
307 if( dty*dty > STY()*nsigt*nsigt ) return false;
308 float dz=s.Z()-Z();
309 float dx=X()+TX()*dz-s.X();
310 if( dx*dx > SX()*nsigx*nsigx ) return false;
311 float dy=Y()+TY()*dz-s.Y();
312 if( dy*dy > SY()*nsigx*nsigx ) return false;
313 return true;
314}
Float_t STX() const
Definition: EdbSegP.h:164
Float_t SY() const
Definition: EdbSegP.h:163
Float_t STY() const
Definition: EdbSegP.h:165
Float_t SX() const
Definition: EdbSegP.h:162

◆ IsEqual()

Bool_t EdbSegP::IsEqual ( const TObject *  obj) const

400{
401 const EdbSegP *s=(EdbSegP*)obj;
402 if(s->ID() != ID()) return false;
403 if(s->PID() != PID()) return false;
404 if(s->Vid(0) != Vid(0)) return false;
405 if(s->Vid(1) != Vid(1)) return false;
406 if(s->Aid(0) != Aid(0)) return false;
407 if(s->Aid(1) != Aid(1)) return false;
408 if(s->Side() != Side()) return false;
409 if( Abs(s->Z()- Z()) >0.000001 ) return false;
410 if( Abs(s->X()- X()) >0.000001 ) return false;
411 if( Abs(s->Y()- Y()) >0.000001 ) return false;
412 if( Abs(s->TX()- TX())>0.000001 ) return false;
413 if( Abs(s->TY()- TY())>0.000001 ) return false;
414 if( Abs(s->W()- W()) >0.000001 ) return false;
415 if(s->Plate() != Plate()) return false;
416 return true;
417}
Int_t ID() const
Definition: EdbSegP.h:147
Int_t Plate() const
Definition: EdbSegP.h:159
Int_t Side() const
mandatory virtual functions:
Definition: EdbSegP.h:170
Int_t PID() const
Definition: EdbSegP.h:148
Int_t Aid(int i) const
Definition: EdbSegP.h:169
Int_t Vid(int i) const
Definition: EdbSegP.h:168

◆ IsSortable()

Bool_t EdbSegP::IsSortable ( ) const
inline
198{ return kTRUE; }

◆ LinkMT()

void EdbSegP::LinkMT ( const EdbSegP s1,
const EdbSegP s2,
EdbSegP s 
)
static

Segments fit by Andrey Aleksandrov (Jul-2003)

131{
133
134 Double_t dz = s2->Z() - s1->Z();
135 Double_t dz2 = dz*dz;
136
137 Double_t q1,q2,w1,w2;
138 Double_t d1,d2,dxx11,dxx22;
139 Double_t dtt01,dtt02,dtx01,dtx02;
140 Double_t dxx01,dxx02,dxt01,dxt02;
141 Double_t xm1,xm2,sx0,sy0,stx0,sty0;
142
143 Double_t q;
144
145 if(dz==0.0) {
146 s->SetZ(s1->Z());
147 s->SetID(s1->ID());
148
149 q1 = s1->SX();
150 q2 = s2->SX();
151 w1 = s1->STX();
152 w2 = s2->STX();
153
154 sx0 = q1*q2/(q1+q2);
155 q = (s1->X()/q1+s2->X()/q2)*sx0;
156 s->SetX(q);
157 stx0 = w1*w2/(w1+w2);
158 q = (s1->TX()/w1+s2->TX()/w2)*stx0;
159 s->SetTX(q);
160
161 q1 = s1->SY();
162 q2 = s2->SY();
163 w1 = s1->STY();
164 w2 = s2->STY();
165
166 sy0 = q1*q2/(q1+q2);
167 q = (s1->Y()/q1+s2->Y()/q2)*sy0;
168 s->SetY(q);
169 sty0 = w1*w2/(w1+w2);
170 q = (s1->TY()/w1+s2->TY()/w2)*sty0;
171 s->SetTY(q);
172
173 s->SetErrors(sx0,sy0,0.0,stx0,sty0);
174 s->SetW( s1->W()+s2->W() );
175 return;
176 }
177
178 q = 0.5*(s1->Z()+s2->Z());
179 Double_t dzr = 1.0/dz;
180
181 s->SetZ(q);
182 s->SetID(s1->ID());
183
184 q1 = s1->SX();
185 q2 = s2->SX();
186 w1 = s1->STX();
187 w2 = s2->STX();
188
189 q = dz2*w2+q2;
190 d1 = 1.0/(q+q1);
191 xm1 = (q*s1->X()+(s2->X()-dz*s2->TX())*q1)*d1;
192
193 q = dz2*w1+q1;
194 d2 = 1.0/(q+q2);
195 xm2 = (q*s2->X()+(s1->X()+dz*s1->TX())*q2)*d2;
196
197
198 dtt01 = q2*d2;
199 dtt02 = q1*d1;
200 dxx11 = 1.0-dtt02;
201 dxx22 = 1.0-dtt01;
202 dxx01 = 0.5*(dxx11+dtt01);
203 dxx02 = 0.5*(dxx22+dtt02);
204 dxt01 = 0.5*dz*dtt01;
205 dxt02 = -0.5*dz*dtt02;
206 dtx01 = dzr*(dtt01-dxx11);
207 dtx02 = dzr*(dxx22-dtt02);
208
209 q = (xm1+xm2)*0.5;
210 s->SetX(q);
211 q = (xm2-xm1)*dzr;
212 s->SetTX(q);
213 sx0 = dxx01*dxx01*q1+dxx02*dxx02*q2+dxt01*dxt01*w1+dxt02*dxt02*w2;
214 stx0 = dtx01*dtx01*q1+dtx02*dtx02*q2+dtt01*dtt01*w1+dtt02*dtt02*w2;
215
216 q1 = s1->SY();
217 q2 = s2->SY();
218 w1 = s1->STY();
219 w2 = s2->STY();
220
221 q = dz2*w2+q2;
222 d1 = 1.0/(q+q1);
223 xm1 = (q*s1->Y()+(s2->Y()-dz*s2->TY())*q1)*d1;
224
225 q = dz2*w1+q1;
226 d2 = 1.0/(q+q2);
227 xm2 = (q*s2->Y()+(s1->Y()+dz*s1->TY())*q2)*d2;
228
229 dtt01 = q2*d2;
230 dtt02 = q1*d1;
231 dxx11 = 1.0-dtt02;
232 dxx22 = 1.0-dtt01;
233 dxx01 = 0.5*(dxx11+dtt01);
234 dxx02 = 0.5*(dxx22+dtt02);
235 dxt01 = 0.5*dz*dtt01;
236 dxt02 = -0.5*dz*dtt02;
237 dtx01 = dzr*(dtt01-dxx11);
238 dtx02 = dzr*(dxx22-dtt02);
239
240 q = (xm1+xm2)*0.5;
241 s->SetY(q);
242 q = (xm2-xm1)*dzr;
243 s->SetTY(q);
244 sy0 = dxx01*dxx01*q1+dxx02*dxx02*q2+dxt01*dxt01*w1+dxt02*dxt02*w2;
245 sty0 = dtx01*dtx01*q1+dtx02*dtx02*q2+dtt01*dtt01*w1+dtt02*dtt02*w2;
246
247 s->SetErrors(sx0,sy0,0.0,stx0,sty0);
248 s->SetW( s1->W()+s2->W() );
249 s->SetDZ(dz);
250}
q
Definition: testBGReduction_AllMethods.C:55

◆ MCEvt()

Int_t EdbSegP::MCEvt ( ) const
inline
145{return eMCEvt;}
Int_t eMCEvt
MC event number.
Definition: EdbSegP.h:43

◆ MCTrack()

Int_t EdbSegP::MCTrack ( ) const
inline
146{return eMCTrack;}
Int_t eMCTrack
MC track number.
Definition: EdbSegP.h:42

◆ MergeTo()

void EdbSegP::MergeTo ( EdbSegP s)

create linked segment at Z of s2 TODO - navesti nauku covariantnuiu

340{
343
344 PropagateTo( s.Z() );
345
346 float wx1,wx2, wy1,wy2;
347 float wtx1,wtx2, wty1,wty2;
348
349 wx1 = 1/SX();
350 wx2 = 1/s.SX();
351 wy1 = 1/SY();
352 wy2 = 1/s.SY();
353 wtx1 = 1/STX();
354 wtx2 = 1/s.STX();
355 wty1 = 1/STY();
356 wty2 = 1/s.STY();
357
358 eX = (X()*wx1 + s.X()*wx2)/(wx1+wx2);
359 eY = (Y()*wy1 + s.Y()*wy2)/(wy1+wy2);
360 (*eCOV)(0,0)= 1./(wx1+wx2);
361 (*eCOV)(1,1)= 1./(wy1+wy2);
362
363 eTX = (TX()*wtx1 + s.TX()*wtx2)/(wtx1+wtx2);
364 eTY = (TY()*wty1 + s.TY()*wty2)/(wty1+wty2);
365 (*eCOV)(2,2)= 1./(wtx1+wtx2);
366 (*eCOV)(3,3)= 1./(wty1+wty2);
367
368 eZ = s.Z();
369 eSZ = TMath::Sqrt(( SZ() + s.SZ())/2);
370
371 eW = W()+s.W();
372
373 eID = s.ID();
374 ePID = s.PID();
375 eFlag = s.Flag();
376}
Float_t eX
Definition: EdbSegP.h:31
Float_t eTX
Definition: EdbSegP.h:32
Float_t eZ
coordinates
Definition: EdbSegP.h:31
Float_t eW
weight
Definition: EdbSegP.h:37
void PropagateTo(float z)
Definition: EdbSegP.cxx:293
Float_t SZ() const
Definition: EdbSegP.h:167
Float_t eSZ
square of the Z-error
Definition: EdbSegP.h:34
Float_t eTY
direction tangents
Definition: EdbSegP.h:32
Float_t eY
Definition: EdbSegP.h:31
Int_t ePID
mother pattern ID
Definition: EdbSegP.h:24

◆ P()

Float_t EdbSegP::P ( ) const
inline
152{return eP;}
Float_t eP
momentum of the particle
Definition: EdbSegP.h:41

◆ Phi()

Float_t EdbSegP::Phi ( ) const
inline
183{return TMath::ATan2(eTY,eTX);}

◆ PID()

Int_t EdbSegP::PID ( ) const
inline
148{return ePID;}

◆ Plate()

Int_t EdbSegP::Plate ( ) const
inline
159{return eScanID.ePlate;}
Int_t ePlate
Definition: EdbID.h:11
EdbID eScanID
Int_t ePlate; //!< plate id.
Definition: EdbSegP.h:45

◆ Print()

void EdbSegP::Print ( Option_t *  opt = "") const
virtual

Reimplemented from EdbTrack2D.

380{
381 printf("EdbSegP: ID= %d PID = %d Vid = %d %d \t Aid = %d %d Flag = %d\n",
382 ID(),PID(),eVid[0],eVid[1],eAid[0],eAid[1], Flag() );
383 printf("x,y,z,dz,tx,ty,w,flag = %f %f %f %f %f %f\n",
384 X(),Y(),Z(),DZ(),TX(),TY() );
385 printf("W = %f P= %f \t Prob = %f \t Chi2 = %f \t Track = %d \n",
386 W(),P(),Prob(),Chi2(), Track() );
387
388 if(eCOV) eCOV->Print();
389}
Int_t Track() const
Definition: EdbSegP.h:150
Float_t Prob() const
Definition: EdbSegP.h:156
Float_t DZ() const
Definition: EdbSegP.h:154
Float_t Chi2() const
Definition: EdbSegP.h:157
Float_t P() const
Definition: EdbSegP.h:152
Int_t eVid[2]
[0]-view entry in the input tree, [1]-segment entry in the view
Definition: EdbSegP.h:26
Int_t Flag() const
Definition: EdbSegP.h:149

◆ PrintNice()

void EdbSegP::PrintNice ( ) const

393{
394 printf( "EdbSegP: ID= %8d PID = %8d, x= %14.3f, y= %14.3f, z= %14.3f, tx= %7.4f, ty= %7.4f, w= %7.4f, chi2= %7.4f Flag= %6d P= %6.1f MCEvt= %6d \n", ID(),PID(),X(),Y(),Z(),TX(),TY(),W(),Chi2(), Flag(), P(), MCEvt() );
395 return;
396}
Int_t MCEvt() const
Definition: EdbSegP.h:145

◆ Prob()

Float_t EdbSegP::Prob ( ) const
inline
156{return eProb;}
Float_t eProb
probability
Definition: EdbSegP.h:36

◆ ProbLink()

float EdbSegP::ProbLink ( EdbSegP s1,
EdbSegP s2 
)

return probability of the correct link in case Up/Down - specifics is that the position errors are neglected)

318{
321
322 double dz = s2.Z() - s1.Z();
323
324 double tx = (s2.X() - s1.X()) / dz;
325 double ty = (s2.Y() - s1.Y()) / dz;
326
327 double dtx1 = (s1.TX()-tx)*(s1.TX()-tx)/s1.STX();
328 double dty1 = (s1.TY()-ty)*(s1.TY()-ty)/s1.STY();
329 double dtx2 = (s2.TX()-tx)*(s2.TX()-tx)/s2.STX();
330 double dty2 = (s2.TY()-ty)*(s2.TY()-ty)/s2.STY();
331
332 double chi2 = TMath::Sqrt(dtx1 + dty1 + dtx2 + dty2);
333 double p3 = TMath::Prob(chi2*chi2,4);
334
335 return s1.Prob()*s2.Prob()*p3;
336}
T Prob(const T &rhs, int n)
Definition: Prob.hh:37
Float_t chi2
Definition: testBGReduction_By_ANN.C:14

◆ PropagateTo()

void EdbSegP::PropagateTo ( float  z)
294{
295 float dz = z-Z();
297}
void PropagateToDZ(float dz)
Definition: EdbSegP.cxx:281

◆ PropagateToCOV()

void EdbSegP::PropagateToCOV ( float  z)
254{
255 float dz = z-Z();
256 eX = X() + TX()*dz;
257 eY = Y() + TY()*dz;
258 eZ = z;
259
260 VtSqMatrix pred(4); //propagation matrix for track parameters (x,y,tx,ty)
261 pred.clear();
262 pred(0,0) = 1.;
263 pred(1,1) = 1.;
264 pred(2,2) = 1.;
265 pred(3,3) = 1.;
266 pred(0,2) = dz;
267 pred(1,3) = dz;
268
269 VtSymMatrix cov(4); // covariance matrix for seg0
270 for(int k=0; k<4; k++)
271 for(int l=0; l<4; l++) cov(k,l) = (COV())(k,l);
272
273 VtSymMatrix covpred(4); // covariation matrix for prediction
274 covpred = pred*(cov*(pred.T()));
275
276 for(int k=0; k<4; k++)
277 for(int l=0; l<4; l++) (COV())(k,l) = covpred(k,l);
278}
TMatrixD & COV() const
Definition: EdbSegP.h:123

◆ PropagateToDZ()

void EdbSegP::PropagateToDZ ( float  dz)
282{
283 eX = eX + eTX*dz;
284 eY = eY + eTY*dz;
285 eZ += dz;
286 if(eCOV) {
287 (*eCOV)(0,0) = SX() + STX()*dz*dz;
288 (*eCOV)(1,1) = SY() + STY()*dz*dz;
289 }
290}

◆ ScanID()

EdbID EdbSegP::ScanID ( ) const
inline
160{return eScanID;}

◆ Set()

void EdbSegP::Set ( int  id,
float  x,
float  y,
float  tx,
float  ty,
float  w,
int  flag 
)
inline
88 { eID=id; eX=x; eY=y; eTX=tx; eTY=ty; eW=w; eFlag=flag; }
UInt_t id
Definition: tlg2couples.C:117

◆ Set0()

void EdbSegP::Set0 ( )
31{
32 ePID=0;
33 eID=0;
34 eVid[0]=eVid[1]=0;
35 eAid[0]=eAid[1]=0;
36 eFlag=0;
37 eTrack=-1;
38 eX=eY=eZ=eTX=eTY=eSZ=0;
39 eProb=0;
40 eW=0;
41 eVolume=0;
42 eDZ=0;
43 eP=-999.;
44 eChi2=0;
45 eCOV=0;
46 eMCTrack=-999;
47 eMCEvt=-999;
48}
Int_t eTrack
id of the track (-1) if no track, for a track object it represents the index of trk....
Definition: EdbSegP.h:29
Float_t eVolume
segment volume
Definition: EdbSegP.h:38

◆ SetAid()

void EdbSegP::SetAid ( int  a,
int  v,
int  side = 0 
)
inline
138{ eAid[0]=a; eAid[1]=100000*side+v; }

◆ SetChi2()

void EdbSegP::SetChi2 ( float  chi2)
inline
135{ eChi2=chi2; }

◆ SetCOV() [1/2]

void EdbSegP::SetCOV ( double *  array,
int  dim = 5 
)
inline
116 {
117 if(!array) return;
118 if(!eCOV) eCOV = new TMatrixD(5,5);
119 for(int k=0; k<dim; k++)
120 for(int l=0; l<dim; l++) (*eCOV)(k,l) = array[k*dim + l];
121 }

◆ SetCOV() [2/2]

void EdbSegP::SetCOV ( TMatrixD &  cov)
inline
99 {
100 if(!(&cov)) return;
101 if(eCOV) eCOV->Copy(cov);
102 else eCOV = new TMatrixD(cov);
103 }

◆ SetDZ()

void EdbSegP::SetDZ ( float  dz)
inline
126{ eDZ=dz; }

◆ SetDZem()

void EdbSegP::SetDZem ( float  dz)
inline
127{ eDZem=dz; }

◆ SetErrorP()

void EdbSegP::SetErrorP ( float  sp2)
inline
94 {
95 if(!eCOV) eCOV = new TMatrixD(5,5);
96 (*eCOV)(4,4) = (double)sp2;
97 }

◆ SetErrors() [1/2]

void EdbSegP::SetErrors ( )
inline
90{SetErrors( 1.,1.,0.,.0001,.0001,1.);}
void SetErrors()
Definition: EdbSegP.h:90

◆ SetErrors() [2/2]

void EdbSegP::SetErrors ( float  sx2,
float  sy2,
float  sz2,
float  stx2,
float  sty2,
float  sp2 = 1. 
)

setting the diagonal elements of covariance matrix

60{
62 SetErrors0();
63 (*eCOV)(0,0) = (double)sx2;
64 (*eCOV)(1,1) = (double)sy2;
65 (*eCOV)(2,2) = (double)stx2;
66 (*eCOV)(3,3) = (double)sty2;
67 (*eCOV)(4,4) = (double)sp2;
68 eSZ = sz2;
69}
void SetErrors0()
Definition: EdbSegP.cxx:51

◆ SetErrors0()

void EdbSegP::SetErrors0 ( )

just init covariance matrix

52{
54 if(!eCOV) eCOV = new TMatrixD(5,5);
55 else eCOV->Zero();
56}

◆ SetErrorsCOV()

void EdbSegP::SetErrorsCOV ( float  sx2,
float  sy2,
float  sz2,
float  stx2,
float  sty2,
float  sp2 = 1. 
)

calculation of the non-diagonal covariance matrix considering the input sigma being in the track plane (Y - is transversal axis)

73{
76
77 SetErrors( sx2, sy2, sz2, stx2, sty2, sp2 );
78
79 //double Phi = -(TMath::ATan2((double)TY(),(double)TX()));
80 double Phi = (TMath::ATan2((double)TY(),(double)TX()));
81 TMatrixD t(5,5);
82 TMatrixD tt(5,5);
83 t(0,0) = TMath::Cos(Phi);
84 t(0,1) = -(TMath::Sin(Phi));
85 t(1,0) = TMath::Sin(Phi);
86 t(1,1) = TMath::Cos(Phi);
87 t(2,2) = TMath::Cos(Phi);
88 t(2,3) = -(TMath::Sin(Phi));
89 t(3,2) = TMath::Sin(Phi);
90 t(3,3) = TMath::Cos(Phi);
91 t(4,4) = 1.;
92 tt(0,0) = t(0,0);
93 tt(1,0) = t(0,1);
94 tt(0,1) = t(1,0);
95 tt(1,1) = t(1,1);
96 tt(2,2) = t(2,2);
97 tt(3,2) = t(2,3);
98 tt(2,3) = t(3,2);
99 tt(3,3) = t(3,3);
100 tt(4,4) = t(4,4);
101 (*eCOV) = t*((*eCOV)*tt);
102}
Float_t Phi() const
Definition: EdbSegP.h:183
TTree * t
Definition: check_shower.C:4

◆ SetFlag()

void EdbSegP::SetFlag ( int  flag)
inline
130{ eFlag=flag; }

◆ SetID()

void EdbSegP::SetID ( int  id)
inline
128{ eID=id; }

◆ SetMC()

void EdbSegP::SetMC ( int  mEvt,
int  mTrack 
)
inline
141{ eMCEvt=mEvt; eMCTrack=mTrack; }

◆ SetP()

void EdbSegP::SetP ( float  p)
inline
133{ eP=p; }
p
Definition: testBGReduction_AllMethods.C:8

◆ SetPID()

void EdbSegP::SetPID ( int  pid)
inline
129{ ePID=pid; }
int pid[1000]
Definition: m2track.cpp:13

◆ SetPlate()

void EdbSegP::SetPlate ( int  plateid)
inline
142{ eScanID.ePlate = plateid; }

◆ SetProb()

void EdbSegP::SetProb ( float  prob)
inline
134{ eProb=prob; }

◆ SetProbability()

void EdbSegP::SetProbability ( float  p)
inline
140{ eProb=p; }

◆ SetScanID()

void EdbSegP::SetScanID ( EdbID  id)
inline
143{ eScanID = id; }

◆ SetSide()

void EdbSegP::SetSide ( int  side = 0)
inline
139{ int v=eAid[1]%100000; eAid[1]=100000*side+v; }

◆ SetSZ()

void EdbSegP::SetSZ ( float  sz)
inline
124{ eSZ=sz; }

◆ SetTrack()

void EdbSegP::SetTrack ( int  trid)
inline
131{ eTrack=trid; }

◆ SetTX()

void EdbSegP::SetTX ( Float_t  tx)
inline
179{ eTX=tx; }

◆ SetTY()

void EdbSegP::SetTY ( Float_t  ty)
inline

other functions

◆ SetVid()

void EdbSegP::SetVid ( int  vid,
int  sid 
)
inline
137{ eVid[0]=vid; eVid[1]=sid; }

◆ SetVolume()

void EdbSegP::SetVolume ( float  w)
inline
136{ eVolume=w; }

◆ SetW()

void EdbSegP::SetW ( float  w)
inline
132{ eW=w; }

◆ SetX()

void EdbSegP::SetX ( Float_t  x)
inline
177{ eX=x; }

◆ SetY()

void EdbSegP::SetY ( Float_t  y)
inline
178{ eY=y; }

◆ SetZ()

void EdbSegP::SetZ ( float  z)
inlinevirtual

Reimplemented from EdbPoint2D.

125{ eZ=z; }

◆ Side()

Int_t EdbSegP::Side ( ) const
inline

mandatory virtual functions:

◆ SP()

Float_t EdbSegP::SP ( ) const
inline
166{ if(!eCOV) return 0; return (Float_t)(*eCOV)(4,4); }

◆ STX()

Float_t EdbSegP::STX ( ) const
inline
164{ if(!eCOV) return 0; return (Float_t)(*eCOV)(2,2); }

◆ STY()

Float_t EdbSegP::STY ( ) const
inline
165{ if(!eCOV) return 0; return (Float_t)(*eCOV)(3,3); }

◆ SX()

Float_t EdbSegP::SX ( ) const
inline
162{ if(!eCOV) return 0; return (Float_t)(*eCOV)(0,0); }

◆ SY()

Float_t EdbSegP::SY ( ) const
inline
163{ if(!eCOV) return 0; return (Float_t)(*eCOV)(1,1); }

◆ SZ()

Float_t EdbSegP::SZ ( ) const
inline
167{ return eSZ; }

◆ Theta()

Float_t EdbSegP::Theta ( ) const
inline
184{return TMath::Sqrt(eTY*eTY+eTX*eTX);}

◆ Track()

Int_t EdbSegP::Track ( ) const
inline
150{return eTrack;}

◆ TX()

Float_t EdbSegP::TX ( ) const
inlinevirtual

tangens = deltaX/deltaZ

Implements EdbAngle2D.

175{return eTX;}

◆ TY()

Float_t EdbSegP::TY ( ) const
inlinevirtual

tangens = deltaY/deltaZ

Implements EdbAngle2D.

176{return eTY;}

◆ Vid()

Int_t EdbSegP::Vid ( int  i) const
inline
168{return (i < 0 || i > 1) ? -1 : eVid[i];}

◆ Volume()

Float_t EdbSegP::Volume ( ) const
inline
158{return eVolume;}

◆ W()

Float_t EdbSegP::W ( ) const
inline
151{return eW;}

◆ X()

Float_t EdbSegP::X ( ) const
inlinevirtual

Implements EdbPoint2D.

173{return eX;}

◆ Y()

Float_t EdbSegP::Y ( ) const
inlinevirtual

Implements EdbPoint2D.

174{return eY;}

◆ Z()

Float_t EdbSegP::Z ( ) const
inlinevirtual

Reimplemented from EdbPoint2D.

153{return eZ;}

Member Data Documentation

◆ eAid

Int_t EdbSegP::eAid[2]
private

[0]-AreaID, [1]-ViewID

◆ eChi2

Float_t EdbSegP::eChi2
private

chi-square

◆ eCOV

TMatrixD* EdbSegP::eCOV
protected

covariance matrix of the parameters (x,y,tx,ty,p)

◆ eDZ

Float_t EdbSegP::eDZ
private

the length of segment along z-axis

◆ eDZem

Float_t EdbSegP::eDZem
private

the length of segment along z-axis in the emulsion

◆ eEMULDigitArray

TRefArray* EdbSegP::eEMULDigitArray
private

AM+AC 27/07/07.

◆ eFlag

Int_t EdbSegP::eFlag
private

◆ eID

Int_t EdbSegP::eID
private

segment id (unique in plate)

◆ eMCEvt

Int_t EdbSegP::eMCEvt
private

MC event number.

◆ eMCTrack

Int_t EdbSegP::eMCTrack
private

MC track number.

◆ eP

Float_t EdbSegP::eP
private

momentum of the particle

◆ ePID

Int_t EdbSegP::ePID
private

mother pattern ID

◆ eProb

Float_t EdbSegP::eProb
private

probability

◆ eScanID

EdbID EdbSegP::eScanID
private

Int_t ePlate; //!< plate id.

brick:plate:major:minor

◆ eSZ

Float_t EdbSegP::eSZ
private

square of the Z-error

◆ eTrack

Int_t EdbSegP::eTrack
private

id of the track (-1) if no track, for a track object it represents the index of trk.root file (tracks tree)

◆ eTX

Float_t EdbSegP::eTX

◆ eTY

Float_t EdbSegP::eTY

direction tangents

◆ eVid

Int_t EdbSegP::eVid[2]
private

[0]-view entry in the input tree, [1]-segment entry in the view

◆ eVolume

Float_t EdbSegP::eVolume
private

segment volume

◆ eW

Float_t EdbSegP::eW
private

weight

◆ eX

Float_t EdbSegP::eX

◆ eY

Float_t EdbSegP::eY

◆ eZ

Float_t EdbSegP::eZ

coordinates


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