FEDRA emulsion software from the OPERA Collaboration
EdbTrackP Class Reference

#include <EdbPattern.h>

Inheritance diagram for EdbTrackP:
Collaboration diagram for EdbTrackP:

Public Member Functions

void AddSegment (EdbSegP *s)
 
void AddSegmentF (EdbSegP *s)
 
void AddTrack (const EdbTrackP &tr)
 
void AddVTA (EdbVTA *vta)
 
int CheckAliasSegments ()
 
int CheckMaxGap ()
 
float CHI2 ()
 
float CHI2F ()
 
void Clear ()
 
void ClearF ()
 
void ClearVTA ()
 
void ClearVTA (EdbVTA *vta)
 
void Copy (const EdbTrackP &tr)
 
Float_t DE () const
 
Int_t Dir () const
 
 EdbTrackP (EdbSegP &seg)
 
 EdbTrackP (EdbSegP *seg, float m=0.12)
 
 EdbTrackP (EdbTrackP &track)
 
 EdbTrackP (int nseg=0)
 
void FitTrack ()
 
int FitTrackKFS (bool zmax=false, float X0=5810., int design=0)
 
Double_t GetBTEfficiency ()
 
EdbSegPGetSegment (int i) const
 
EdbSegPGetSegmentF (int i) const
 
EdbSegPGetSegmentFFirst () const
 
EdbSegPGetSegmentFirst () const
 
EdbSegPGetSegmentFLast () const
 
EdbSegPGetSegmentLast () const
 
Int_t GetSegmentsAid (int &nseg) const
 
Int_t GetSegmentsFlag (int &nseg) const
 
Int_t GetSegmentsMCTrack (int &nseg) const
 
EdbSegPGetSegmentWithClosestZ (float z, float dz)
 
Float_t M () const
 
float MakePredictionTo (Float_t z, EdbSegP &ss)
 
int MakeSelector (EdbSegP &ss, bool followZ=true)
 
Int_t N () const
 
Int_t N0 () const
 
Int_t NF () const
 
Int_t Npl () const
 
Int_t PDG () const
 
Float_t PerrDown () const
 
Float_t PerrUp () const
 
void Print ()
 
void PrintNice ()
 
int RemoveAliasSegments ()
 
void RemoveSegment (EdbSegP *s)
 
void Set0 ()
 
void SetCounters ()
 
void SetDE (float de)
 
void SetM (float m)
 
void SetN0 ()
 
void SetN0 (int n0)
 
void SetNpl ()
 
void SetNpl (int npl)
 
void SetOwner ()
 
void SetPDG (int pdg)
 
void SetPerr (Float_t perrDown, Float_t perrUp)
 
void SetPerrDown (Float_t perrDown)
 
void SetPerrUp (Float_t perrUp)
 
int SetSegmentsTrack ()
 
int SetSegmentsTrack (int id)
 
void SubstituteSegment (EdbSegP *sold, EdbSegP *snew)
 
const EdbSegPTrackEnd () const
 
EdbSegPTrackExtremity (bool zpos, bool usesegpar=false) const
 
const EdbSegPTrackStart () const
 
EdbSegPTrackZmax (bool usesegpar=false) const
 
EdbSegPTrackZmin (bool usesegpar=false) const
 
void Transform (const EdbAffine2D &tr)
 
EdbVertexVertex (int zpos)
 
EdbVertexVertexE ()
 
EdbVertexVertexS ()
 
EdbVTAVTAE () const
 
EdbVTAVTAS () const
 
Float_t Wgrains () const
 
Float_t Wmean () const
 
Float_t Zend () const
 
Float_t Zmax () const
 
Float_t Zmin () const
 
Float_t Zstart () const
 
virtual ~EdbTrackP ()
 
- Public Member Functions inherited from EdbSegP
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 ()
 

Private Attributes

Float_t eDE
 total energy loss of the particle between first and last segments More...
 
Float_t eM
 invariant mass of the particle More...
 
Int_t eN0
 number of holes (if any) More...
 
Int_t eNpl
 number of plates passed through More...
 
Int_t ePDG
 particle ID from PDG More...
 
Float_t ePerrDown
 error of P() in lower direction, obtained by MCS,or shower-algorithm More...
 
Float_t ePerrUp
 error of P() in upper direction, obtained by MCS,or shower-algorithm More...
 
TSortedList * eS
 array of segments More...
 
TSortedList * eSF
 array of fitted segments More...
 
EdbVTAeVTAE
 vertex track end is attached to More...
 
EdbVTAeVTAS
 vertex track start is attached to More...
 

Additional Inherited Members

- Static Public Member Functions inherited from EdbSegP
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 inherited from EdbSegP
Float_t eTX
 
Float_t eTY
 direction tangents More...
 
Float_t eX
 
Float_t eY
 
Float_t eZ
 coordinates More...
 
- Protected Attributes inherited from EdbSegP
TMatrixD * eCOV
 covariance matrix of the parameters (x,y,tx,ty,p) More...
 

Constructor & Destructor Documentation

◆ EdbTrackP() [1/4]

EdbTrackP::EdbTrackP ( int  nseg = 0)
387{
388 Set0();
389 if(nseg>0) eS = new TSortedList();
390 if(nseg>0) { eSF = new TSortedList(); eSF->SetOwner(); }
391}
TSortedList * eS
array of segments
Definition: EdbPattern.h:117
TSortedList * eSF
array of fitted segments
Definition: EdbPattern.h:118
void Set0()
Definition: EdbPattern.cxx:432

◆ EdbTrackP() [2/4]

EdbTrackP::EdbTrackP ( EdbSegP seg)

create empty track with segment parameters

394 : EdbSegP(seg)
395{
397 eS=0;
398 eSF=0;
399 eM=0;
400 eDE=0;
401 ePDG=-999;
402 eVTAS = 0;
403 eVTAE = 0;
404 eNpl=0;
405 eN0=0;
406}
EdbSegP()
Definition: EdbSegP.h:53
Int_t ePDG
particle ID from PDG
Definition: EdbPattern.h:123
Int_t eNpl
number of plates passed through
Definition: EdbPattern.h:119
EdbVTA * eVTAE
vertex track end is attached to
Definition: EdbPattern.h:127
Float_t eM
invariant mass of the particle
Definition: EdbPattern.h:121
Int_t eN0
number of holes (if any)
Definition: EdbPattern.h:120
Float_t eDE
total energy loss of the particle between first and last segments
Definition: EdbPattern.h:122
EdbVTA * eVTAS
vertex track start is attached to
Definition: EdbPattern.h:126

◆ EdbTrackP() [3/4]

EdbTrackP::EdbTrackP ( EdbSegP seg,
float  m = 0.12 
)
409 : EdbSegP( *seg )
410{
411 eS=0;
412 eSF=0;
413 eM=0;
414 eDE=0;
415 ePDG=-999;
416 eVTAS = 0;
417 eVTAE = 0;
418 eNpl=0;
419 eN0=0;
420 AddSegment(seg);
421 SetM(m);
422}
void AddSegment(EdbSegP *s)
Definition: EdbPattern.h:214
void SetM(float m)
Definition: EdbPattern.h:154

◆ EdbTrackP() [4/4]

EdbTrackP::EdbTrackP ( EdbTrackP track)
inline
134{ Set0(); Copy(track); }
void Copy(const EdbTrackP &tr)
Definition: EdbPattern.cxx:473
Definition: bitview.h:14

◆ ~EdbTrackP()

EdbTrackP::~EdbTrackP ( )
virtual
426{
427 if(eS) { eS->Clear(); delete eS; eS=0; }
428 if(eSF) { eSF->Clear(); delete eSF; eSF=0; }
429}

Member Function Documentation

◆ AddSegment()

void EdbTrackP::AddSegment ( EdbSegP s)
inline
215 {
216 if(!eS) eS = new TSortedList();
217 eS->Add(s);
218 }
s
Definition: check_shower.C:55

◆ AddSegmentF()

void EdbTrackP::AddSegmentF ( EdbSegP s)
inline
234 {
235 if(!eSF)
236 {
237 eSF = new TSortedList();
238 eSF->SetOwner();
239 }
240 eSF->Add(s);
241 }

◆ AddTrack()

void EdbTrackP::AddTrack ( const EdbTrackP tr)
669{
670 int nseg=tr.N();
671 int nsegf=tr.NF();
672 for(int i=0; i<nseg; i++)
673 AddSegment(tr.GetSegment(i));
674 for(int i=0; i<nsegf; i++)
675 AddSegmentF(new EdbSegP(*(tr.GetSegmentF(i))));
676}
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
void AddSegmentF(EdbSegP *s)
Definition: EdbPattern.h:233

◆ AddVTA()

void EdbTrackP::AddVTA ( EdbVTA vta)
451{
452 if(vta) {
453 if(vta->Zpos()==0) eVTAE=vta;
454 else if(vta->Zpos()==1) eVTAS=vta;
455 }
456}
Int_t Zpos() const
Definition: EdbVertex.h:47

◆ CheckAliasSegments()

int EdbTrackP::CheckAliasSegments ( )
507{
508 int nalias=0;
509 for(int i=0; i<N(); i++)
510 if( GetSegment(i)->Track() != ID()) nalias++;
511 return nalias;
512}
Int_t Track() const
Definition: EdbSegP.h:150
Int_t ID() const
Definition: EdbSegP.h:147
Int_t N() const
Definition: EdbPattern.h:177
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:195

◆ CheckMaxGap()

int EdbTrackP::CheckMaxGap ( )
531{
532 int ngap=0;
533 if(N()<2) return 0;
534 EdbSegP *s1=0, *s2=0;
535 int gap = 0;
536 for(int i=0; i<N()-1; i++) {
537 s1 = GetSegment(i);
538 s2 = GetSegment(i+1);
539 gap = TMath::Abs(s2->PID()-s1->PID());
540 if( gap > ngap) ngap = gap;
541 }
542 return ngap;
543}
Definition: EdbSegP.h:21
Int_t PID() const
Definition: EdbSegP.h:148
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ CHI2()

float EdbTrackP::CHI2 ( )
1100{
1101 double dtx=0,dty=0,chi2=0;
1102 EdbSegP *seg=0;
1103 int nseg=N();
1104 for(int i=0; i<nseg; i++) {
1105 seg = GetSegment(i);
1106 dtx = seg->TX()-TX();
1107 dty = seg->TY()-TY();
1108 chi2 += TMath::Sqrt( dtx*dtx/seg->STX() +
1109 dty*dty/seg->STY() );
1110 }
1111 chi2 /= nseg;
1112 return chi2;
1113}
Float_t STX() const
Definition: EdbSegP.h:164
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
Float_t STY() const
Definition: EdbSegP.h:165
Float_t chi2
Definition: testBGReduction_By_ANN.C:14

◆ CHI2F()

float EdbTrackP::CHI2F ( )
1117{
1118 double dtx=0,dty=0,chi2=0;
1119 EdbSegP *s=0, *sf=0;
1120 int nseg=N();
1121 for(int i=0; i<nseg; i++) {
1122 s = GetSegment(i);
1123 sf = GetSegmentF(i);
1124 dtx = s->TX() - sf->TX();
1125 dty = s->TY() - sf->TY();
1126 chi2 += TMath::Sqrt( dtx*dtx/s->STX() +
1127 dty*dty/s->STY() );
1128 }
1129 chi2 /= nseg;
1130 return chi2;
1131}
EdbSegP * GetSegmentF(int i) const
Definition: EdbPattern.h:196

◆ Clear()

void EdbTrackP::Clear ( )
inline
264{ if(eS) eS->Clear(); if(eSF) eSF->Clear(); }

◆ ClearF()

void EdbTrackP::ClearF ( )
inline
265{ if(eSF) eSF->Clear(); }

◆ ClearVTA() [1/2]

void EdbTrackP::ClearVTA ( )
460{
461 eVTAE=0;
462 eVTAS=0;
463}

◆ ClearVTA() [2/2]

void EdbTrackP::ClearVTA ( EdbVTA vta)
467{
468 if(vta->Zpos()==0) eVTAE=0;
469 else if(vta->Zpos()==1) eVTAS=0;
470}

◆ Copy()

void EdbTrackP::Copy ( const EdbTrackP tr)

do the physical copy of segments

474{
476 Clear();
477 ((EdbSegP*)(this))->Copy( *((EdbSegP*)(&tr)) );
478 SetM(tr.M());
479 SetPDG(tr.PDG());
480 SetNpl(tr.Npl());
481 SetN0(tr.N0());
482 SetDE(tr.DE());
483 AddVTA(tr.VTAS());
484 AddVTA(tr.VTAE());
485
486 int nseg=tr.N();
487 for(int i=0; i<nseg; i++)
488 AddSegment(new EdbSegP(*tr.GetSegment(i)));
489 nseg=tr.NF();
490 for(int i=0; i<nseg; i++)
491 AddSegmentF(new EdbSegP(*tr.GetSegmentF(i)));
492 if (eS) eS->SetOwner();
493 if (eSF) eSF->SetOwner();
494}
void AddVTA(EdbVTA *vta)
Definition: EdbPattern.cxx:450
void Clear()
Definition: EdbPattern.h:264
void SetNpl()
Definition: EdbPattern.h:169
void SetN0()
Definition: EdbPattern.h:162
void SetDE(float de)
Definition: EdbPattern.h:165
void SetPDG(int pdg)
Definition: EdbPattern.h:149

◆ DE()

Float_t EdbTrackP::DE ( ) const
inline
166{return eDE;}

◆ Dir()

Int_t EdbTrackP::Dir ( ) const
inline
207{return (DZ()<0) ? -1 : 1;}
Float_t DZ() const
Definition: EdbSegP.h:154

◆ FitTrack()

void EdbTrackP::FitTrack ( )

track fit by averaging of segments parameters

697{
699
700 int nseg=N();
701 float x=0,y=0,z=0,tx=0,ty=0,w=0;
702 EdbSegP *seg=0;
703 for(int i=0; i<nseg; i++) {
704 seg = GetSegment(i);
705 x += seg->X();
706 y += seg->Y();
707 z += seg->Z();
708 tx += seg->TX();
709 ty += seg->TY();
710 w += seg->W();
711 }
712 x /= nseg;
713 y /= nseg;
714 z /= nseg;
715 tx /= nseg;
716 ty /= nseg;
717 Set(ID(),x,y,tx,ty,w,Flag());
718 SetZ(z);
719}
Float_t X() const
Definition: EdbSegP.h:173
Float_t Z() const
Definition: EdbSegP.h:153
void SetZ(float z)
Definition: EdbSegP.h:125
Float_t Y() const
Definition: EdbSegP.h:174
Float_t W() const
Definition: EdbSegP.h:151
Int_t Flag() const
Definition: EdbSegP.h:149
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

◆ FitTrackKFS()

int EdbTrackP::FitTrackKFS ( bool  zmax = false,
float  X0 = 5810.,
int  design = 0 
)

if (zmax==true) track parameters are calculated at segment with max Z if (zmax==false) track parameters are calculated at segment with min Z tTODO: track parameters??

Note: momentum for track must be defined with SetP() prior to call this routine! Note: momentum dispersion for track should be defined with SetErrorP() prior to call this routine - it is necessary for vertex fit

on the output: Chi2: the full chi-square (not divided by NDF); NDF = 4 Prob: is Chi2 probability (area of the tail of Chi2-distribution) If we accept events with Prob >= ProbMin then ProbMin is the probability to reject the good event

723{
736
737
738 int nseg=N();
739 if(nseg<=0) return 0;
740 if(NF()) {
741 eSF->Delete();
742 }
743
744 if(SP()<0.000001) SetErrorP(1.); // TODO: razobratsia s etimi impul'sami!
745 //printf("%d segments to fit\n",N());
746
747 //TODO - eliminate constants!!!
748
749 float dPb = 0., p = 0., m = 0.13957, e = 0.13957, de = 0., pa = 0., pn = 0.;
750 float e0 = e;
751 float eTPb = 1000./1300.;
752 float pcut = 0.050;
753 double teta0sq;
754 double dz, ptx, pty;
755 int step;
756 int istart, iend;
757
758 VtVector *par[260], *parpred[260], *pars[260], *meas[260];
759 VtSqMatrix *pred[260];
760 VtSymMatrix *cov[260], *covpred[260], *covpredinv[260], *covs[260], *dmeas[260];
761
762 int i=0;
763 if(nseg == 1) {
764 EdbSegP *s = GetSegment(0);
765 EdbSegP segf(*s);
766
767 segf.Set(s->ID(),s->X(),s->Y(),s->TX(),s->TY(),1.,s->Flag());
768 segf.SetZ(s->Z());
769 segf.SetErrorP ( SP() );
770 segf.SetChi2(0.);
771 segf.SetProb(1.);
772 segf.SetP( P() );
773 segf.SetPID( s->PID() );
774 segf.SetDZ(s->DZ());
775
776 AddSegmentF(new EdbSegP(segf));
777 return 0;
778 }
779
780 if(nseg>259) return -1;
781
782 EdbSegP *seg0=0;
783 EdbSegP *seg=0;
784
785 if( GetSegment(N()-1)->Z() < GetSegment(0)->Z() )
786 {
787 if (zmax)
788 {
789 step=-1;
790 istart=nseg-1;
791 iend=0;
792 }
793 else
794 {
795 step=1;
796 istart=0;
797 iend=nseg-1;
798 }
799 }
800 else
801 {
802 if (!zmax)
803 {
804 step=-1;
805 istart=nseg-1;
806 iend=0;
807 }
808 else
809 {
810 step=1;
811 istart=0;
812 iend=nseg-1;
813 }
814
815 }
816
817 seg0 = GetSegment(istart);
818
819 par[istart] = new VtVector( (double)(seg0->X()),
820 (double)(seg0->Y()),
821 (double)(seg0->TX()),
822 (double)(seg0->TY()) );
823 meas[istart] = new VtVector(*par[istart]);
824 pred[iend] = new VtSqMatrix(4);
825// (*pred[istart]).clear();
826// (*pred[istart])(0,0) = 1.;
827// (*pred[istart])(1,1) = 1.;
828// (*pred[istart])(2,2) = 1.;
829// (*pred[istart])(3,3) = 1.;
830 cov[istart] = new VtSymMatrix(4); // covariance matrix for seg0
831 for(int k=0; k<4; k++)
832 for(int l=0; l<4; l++) (*cov[istart])(k,l) = (seg0->COV())(k,l);
833 dmeas[istart] = new VtSymMatrix(*cov[istart]); // covariance matrix for seg0
834
835 Double_t chi2=0.;
836
837 i = istart;
838 p = P();
839 m = M();
840 if (p < pcut) p = pcut;
841 e = TMath::Sqrt((double)p*(double)p + (double)m*(double)m);
842 e0 = e;
843 while( (i+=step) != iend+step ) {
844
845 seg = GetSegment(i);
846
847 VtSymMatrix dms(4); // multiple scattering matrix
848 dms.clear();
849
850 dz = seg->Z()-seg0->Z();
851 ptx = (*par[i-step])(2); // previous tx
852 pty = (*par[i-step])(3); // previous ty
853 dPb = dz*TMath::Sqrt(1.+ptx*ptx+pty*pty); // thickness of the Pb+emulsion cell in microns
854 if ((design != 0) && (p > pcut))
855 {
856 de = EdbPhysics::DeAveragePb(p, m, TMath::Abs(eTPb*dPb));
857 if (de < 0.) de = 0.;
858 if (design < 0) de = -de;
859 if (dz >= 0.)
860 e = e - de;
861 else
862 e = e + de;
863 if (e < m) e = m;
864 pn = TMath::Sqrt((double)e*(double)e - (double)m*(double)m);
865 if (pn <= pcut) pn = pcut;
866 pa = 0.5*(p + pn);
867 p = pn;
868 }
869 else
870 {
871 pa = p;
872 }
873 teta0sq = EdbPhysics::ThetaMS2( pa, m, dPb, X0 );
874 dms(0,0) = teta0sq*dz*dz/3.;
875 dms(1,1) = dms(0,0);
876 dms(2,2) = teta0sq;
877 dms(3,3) = dms(2,2);
878// dms(2,0) = teta0sq*TMath::Abs(dz)/2.;
879 dms(2,0) = teta0sq*dz/2.;
880 dms(3,1) = dms(2,0);
881 dms(0,2) = dms(2,0);
882 dms(1,3) = dms(2,0);
883
884 pred[i-step] = new VtSqMatrix(4); //propagation matrix for track parameters (x,y,tx,ty)
885 pred[i-step]->clear();
886
887 (*pred[i-step])(0,0) = 1.;
888 (*pred[i-step])(1,1) = 1.;
889 (*pred[i-step])(2,2) = 1.;
890 (*pred[i-step])(3,3) = 1.;
891 (*pred[i-step])(0,2) = dz;
892 (*pred[i-step])(1,3) = dz;
893
894 parpred[i] = new VtVector(4); // prediction from seg0 to seg
895 *parpred[i] = (*pred[i-step])*(*par[i-step]);
896
897 covpred[i] = new VtSymMatrix(4); // covariation matrix for prediction
898 *covpred[i] = (*pred[i-step])*((*cov[i-step])*((*pred[i-step]).T()))+dms;
899
900 dmeas[i] = new VtSymMatrix(4); // original covariation matrix for seg
901 for(int k=0; k<4; k++)
902 for(int l=0; l<4; l++) (*dmeas[i])(k,l) = (seg->COV())(k,l);
903
904 covpredinv[i] = new VtSymMatrix(4);
905 (*covpredinv[i]) = (*covpred[i]).dsinv();
906 VtSymMatrix dmeasinv(4);
907 dmeasinv = (*dmeas[i]).dsinv();
908 cov[i] = new VtSymMatrix(4);
909 (*cov[i]) = (*covpredinv[i]) + dmeasinv;
910 (*cov[i]) = (*cov[i]).dsinv();
911
912 meas[i] = new VtVector( (double)(seg->X()),
913 (double)(seg->Y()),
914 (double)(seg->TX()),
915 (double)(seg->TY()) );
916
917 par[i] = new VtVector(4);
918 (*par[i]) = (*cov[i])*((*covpredinv[i])*(*parpred[i]) + dmeasinv*(*meas[i])); // new parameters for seg
919
920// chi2 += ((*par[i])-(*parpred[i]))*((*covpredinv[i])*((*par[i])-(*parpred[i]))) +
921// ((*par[i])-(*meas[i]))*(dmeasinv*((*par[i])-(*meas[i])));
922
923 VtSymMatrix dresid(4);
924 dresid = (*dmeas[i]) - (*cov[i]);
925 dresid = dresid.dsinv();
926
927 chi2 += ((*par[i])-(*meas[i]))*(dresid*((*par[i])-(*meas[i])));
928
929 seg0 = seg;
930 }
931
932 Set(ID(),(float)(*par[iend])(0),(float)(*par[iend])(1),
933 (float)(*par[iend])(2),(float)(*par[iend])(3),1.,Flag());
934 SetZ(GetSegment(iend)->Z());
935 SetPID(GetSegment(iend)->PID());
936 SetCOV( (*cov[iend]).array(), 4 );
937
938 //SetChi2((float)chi2);
939 //SetProb( (float)TMath::Prob(chi2,nseg*4));
940
941// Smoothing
942
943 double chi2p=0;
944
945 pars[iend] = new VtVector(*par[iend]);
946 covs[iend] = new VtSymMatrix(*cov[iend]);
947 VtSymMatrix dresid(4);
948 dresid = (*dmeas[iend]) - (*covs[iend]);
949 dresid = dresid.dsinv();
950 chi2p = ((*pars[iend])-(*meas[iend]))*(dresid*((*pars[iend])-(*meas[iend])));
951
952 EdbSegP segf;
953 segf.Set(ID(),(float)(*pars[iend])(0),(float)(*pars[iend])(1),
954 (float)(*pars[iend])(2),(float)(*pars[iend])(3),1.,Flag());
955 segf.SetZ(GetSegment(iend)->Z());
956 segf.SetCOV( (*covs[iend]).array(), 4 );
957 segf.SetErrorP ( SP() );
958 segf.SetChi2((float)chi2p);
959 segf.SetProb( (float)TMath::Prob(chi2p,4));
960 segf.SetW( (float)nseg );
961 segf.SetP( P() );
962 segf.SetPID( GetSegment(iend)->PID() );
963 segf.SetDZ(seg->DZ());
964
965 AddSegmentF(new EdbSegP(segf));
966
967 i=iend;
968 double DE=0.;
970
971 while( (i-=step) != istart-step ) {
972 VtSqMatrix BackTr(4);
973 BackTr = (*cov[i])*(((*pred[i]).T())*(*covpredinv[i+step]));
974 pars[i] = new VtVector(4);
975 covs[i] = new VtSymMatrix(4);
976 (*pars[i]) = (*par[i]) + BackTr*((*pars[i+step])-(*parpred[i+step]));
977 (*covs[i]) = (*cov[i]) + BackTr*(((*covs[i+step])-(*covpred[i+step]))*BackTr.T());
978 dresid = (*dmeas[i]) - (*covs[i]);
979 dresid = dresid.dsinv();
980 chi2p = ((*pars[i])-(*meas[i]))*(dresid*((*pars[i])-(*meas[i])));
981// chi2 += chi2p;
982
983 segf.Set(ID(),(float)(*pars[i])(0),(float)(*pars[i])(1),
984 (float)(*pars[i])(2),(float)(*pars[i])(3),1.,Flag());
985 segf.SetZ(GetSegment(i)->Z());
986 segf.SetCOV( (*covs[i]).array(), 4 );
987 segf.SetErrorP ( SP() );
988 segf.SetChi2((float)chi2p);
989 segf.SetProb( (float)TMath::Prob(chi2p,4));
990 segf.SetW( (float)nseg );
991 segf.SetP( P() );
992 segf.SetPID( GetSegment(i)->PID() );
993 segf.SetDZ(seg->DZ());
994 AddSegmentF(new EdbSegP(segf));
995 dz = eTPb*TMath::Abs(GetSegment(i)->Z() - GetSegment(i+step)->Z());
996 dPb = dz*TMath::Sqrt(1.+(*pars[i])(2)*(*pars[i])(2)+(*pars[i])(3)*(*pars[i])(3)); // thickness of the Pb+emulsion cell in microns
997 DE += EdbPhysics::DeAveragePbFast(P(),M(),TMath::Abs(dPb));
998 }
999 SetChi2((float)chi2);
1000 SetProb((float)TMath::Prob(chi2,4*(nseg-1)));
1001 SetW( (float)nseg );
1002 SetDE( (float)DE );
1003
1004// DEBUG
1005// printf(" e0 - e = %f, de = %f\n", e0-e, DE);
1006
1007// Delete matrixes and vectors
1008
1009 delete par[istart];
1010 par[istart] = 0;
1011 delete cov[istart];
1012 cov[istart] = 0;
1013 delete meas[istart];
1014 meas[istart] = 0;
1015 delete dmeas[istart];
1016 dmeas[istart] = 0;
1017 delete pred[istart];
1018 pred[istart] = 0;
1019 delete pars[istart];
1020 pars[istart] = 0;
1021 delete covs[istart];
1022 covs[istart] = 0;
1023 i=istart;
1024 while( (i+=step) != iend+step ) {
1025 delete pred[i];
1026 pred[i] = 0;
1027 delete parpred[i];
1028 parpred[i] = 0;
1029 delete covpred[i];
1030 covpred[i] = 0;
1031 delete covpredinv[i];
1032 covpredinv[i] = 0;
1033 delete par[i];
1034 par[i] = 0;
1035 delete cov[i];
1036 cov[i] = 0;
1037 delete meas[i];
1038 meas[i] = 0;
1039 delete dmeas[i];
1040 dmeas[i] = 0;
1041 delete pars[i];
1042 pars[i] = 0;
1043 delete covs[i];
1044 covs[i] = 0;
1045 }
1046 return 0;
1047}
T Prob(const T &rhs, int n)
Definition: Prob.hh:37
brick dz
Definition: RecDispMC.C:107
brick X0
Definition: RecDispMC.C:112
int design
Definition: RecDispMC.C:90
static double DeAveragePb(float p, float mass, float dx)
Definition: EdbPhys.cxx:132
static double DeAveragePbFast(float p, float mass, float dx)
Definition: EdbPhys.cxx:218
static double ThetaMS2(float p, float mass, float dx, float X0)
Definition: EdbPhys.cxx:51
static void DeAveragePbFastSet(float p, float mass)
Definition: EdbPhys.cxx:186
void SetPID(int pid)
Definition: EdbSegP.h:129
void SetProb(float prob)
Definition: EdbSegP.h:134
TMatrixD & COV() const
Definition: EdbSegP.h:123
Float_t P() const
Definition: EdbSegP.h:152
Float_t SP() const
Definition: EdbSegP.h:166
void SetW(float w)
Definition: EdbSegP.h:132
void SetCOV(TMatrixD &cov)
Definition: EdbSegP.h:99
void SetChi2(float chi2)
Definition: EdbSegP.h:135
void SetP(float p)
Definition: EdbSegP.h:133
void SetDZ(float dz)
Definition: EdbSegP.h:126
void SetErrorP(float sp2)
Definition: EdbSegP.h:94
Int_t NF() const
Definition: EdbPattern.h:178
Float_t M() const
Definition: EdbPattern.h:155
Float_t DE() const
Definition: EdbPattern.h:166
p
Definition: testBGReduction_AllMethods.C:8

◆ GetBTEfficiency()

Double_t EdbTrackP::GetBTEfficiency ( )

Returns Basetrack Scanning efficiency, estimated by BT in the Tracks and the number of holes. Formula: efficiency_per_BT = ( N()-2 ) / ( Npl()-2 ) This makes sense when having a track with at least 3 plates crossing.

1135{
1141 Double_t nseg=(Double_t)N();
1142 Double_t npl=(Double_t)Npl();
1143
1144 if (nseg==2 && npl==2) return 1;
1145 Double_t BTEfficiency=(nseg-2)/(npl-2);
1146
1147 return BTEfficiency;
1148}
Int_t Npl() const
Definition: EdbPattern.h:171

◆ GetSegment()

EdbSegP * EdbTrackP::GetSegment ( int  i) const
inline
195{return (eS) ? (EdbSegP*)(eS->At(i)) : 0; }

◆ GetSegmentF()

EdbSegP * EdbTrackP::GetSegmentF ( int  i) const
inline
196{return (eSF) ? (EdbSegP*)(eSF->At(i)) : 0;}

◆ GetSegmentFFirst()

EdbSegP * EdbTrackP::GetSegmentFFirst ( ) const
inline
192{return (eSF) ? (EdbSegP*)(eSF->First()) : 0;}

◆ GetSegmentFirst()

EdbSegP * EdbTrackP::GetSegmentFirst ( ) const
inline
189{return (eS) ? (EdbSegP*)(eS->First()) : 0;}

◆ GetSegmentFLast()

EdbSegP * EdbTrackP::GetSegmentFLast ( ) const
inline
193{return (eSF) ? (EdbSegP*)(eSF->Last()) : 0;}

◆ GetSegmentLast()

EdbSegP * EdbTrackP::GetSegmentLast ( ) const
inline
190{return (eS) ? (EdbSegP*)(eS->Last()) : 0;}

◆ GetSegmentsAid()

int EdbTrackP::GetSegmentsAid ( int &  nseg) const
625{
626 int nseg = N();
627 if (nseg < 2) return -1;
628 if (nseg > 300) nseg = 300;
629 int count[300];
630 int aid[300], ai;
631 int i, f = 0;
632 for(i=0; i<nseg; i++)
633 {
634 count[i] = 0;
635 aid[i] = -1;
636 ai = GetSegment(i)->Aid(1);
637 f = ai + (GetSegment(i)->Aid(0))*1000;
638 if (f < 0) continue;
639 aid[i] = ai;
640 for (int j=0; j<nseg; j++)
641 {
642 if ( f == (GetSegment(j)->Aid(1) + (GetSegment(i)->Aid(0))*1000)) count[i]++;
643 }
644 }
645 int cmax = 0;
646 int aidmax = -1;
647 for(i=0; i<nseg; i++)
648 {
649 if (count[i] > cmax)
650 {
651 cmax = count[i];
652 aidmax = aid[i];
653 }
654 }
655 nsegf = cmax;
656 return aidmax;
657}
FILE * f
Definition: RecDispMC.C:150
Int_t Aid(int i) const
Definition: EdbSegP.h:169

◆ GetSegmentsFlag()

int EdbTrackP::GetSegmentsFlag ( int &  nseg) const
592{
593 printf(" EdbTrackP::GetSegmentsFlag: obsolete, to be deleted! now replaced by GetSegmentsMCTrack\n");
594
595 int nseg = N();
596 if (nseg < 2) return -1;
597 if (nseg > 300) nseg = 300;
598 int count[300];
599 int i, f = 0;
600 for(i=0; i<nseg; i++)
601 {
602 count[i] = 0;
603 f = GetSegment(i)->Flag();
604 if (f < 0) continue;
605 for (int j=0; j<nseg; j++)
606 {
607 if ( f == GetSegment(j)->Flag()) count[i]++;
608 }
609 }
610 int cmax = 0;
611 int flagmax = -1;
612 for(i=0; i<nseg; i++)
613 {
614 if (count[i] > cmax)
615 {
616 cmax = count[i];
617 flagmax = GetSegment(i)->Flag();
618 }
619 }
620 nsegf = cmax;
621 return flagmax;
622}

◆ GetSegmentsMCTrack()

int EdbTrackP::GetSegmentsMCTrack ( int &  nseg) const
546{
547 // return: ID of the MC track contributed to this reconstructed track
548 // with the maximal number of segments
549 // nsegf - return number of segments originated from the selected MC track
550
551 int nseg = N();
552 if (nseg < 2) return -1;
553 if (nseg > 300) nseg = 300;
554 int count[300];
555 int i, f = 0;
556 for(i=0; i<nseg; i++)
557 {
558 count[i] = 0;
559 f = GetSegment(i)->MCTrack();
560 if (f < 0) continue;
561 for (int j=0; j<nseg; j++)
562 {
563 if ( f == GetSegment(j)->MCTrack()) count[i]++;
564 }
565 }
566 int cmax = 0;
567 int mctrackmax = -1;
568 for(i=0; i<nseg; i++)
569 {
570 if (count[i] > cmax)
571 {
572 cmax = count[i];
573 mctrackmax = GetSegment(i)->MCTrack();
574 }
575 }
576 nsegf = cmax;
577 return mctrackmax;
578}
Int_t MCTrack() const
Definition: EdbSegP.h:146

◆ GetSegmentWithClosestZ()

EdbSegP * EdbTrackP::GetSegmentWithClosestZ ( float  z,
float  dz 
)
680{
681 float dzmin=dzMax;
682 EdbSegP *sbest=0;
683 int nseg=N();
684 for(int i=0; i<nseg; i++) {
685 EdbSegP *s = GetSegment(i);
686 float dz = Abs(s->eZ-z);
687 if( dz < dzmin ) {
688 dzmin=dz;
689 sbest=s;
690 }
691 }
692 return sbest;
693}

◆ M()

Float_t EdbTrackP::M ( ) const
inline
155{return eM;}

◆ MakePredictionTo()

float EdbTrackP::MakePredictionTo ( Float_t  z,
EdbSegP ss 
)

1051{
1052 float dz = Zmax()-Zmin();
1053 const EdbSegP *tr=0;
1054 if ( z <= Zmin() ) tr = TrackZmin(); //TODO: fitted - not fitted
1055 else if( z >= Zmax() ) tr = TrackZmax();
1056 else {
1057 for(int i=0; i<N(); i++) {
1058 EdbSegP *s = GetSegment(i);
1059 if( Abs(s->Z()-z) < dz ) { dz = Abs(s->Z()-z); tr=s; } // select nearest segment (TODO: correct interpolation)
1060 }
1061 }
1062 if(!tr) // no segments in this track: use track body
1063 tr = this;
1064 dz = Abs(tr->Z()-z);
1065 ss.Copy(*tr);
1066 ss.PropagateTo(z);
1067 return dz;
1068}
EdbSegP * TrackZmax(bool usesegpar=false) const
Definition: EdbPattern.h:199
Float_t Zmin() const
Definition: EdbPattern.cxx:1205
EdbSegP * TrackZmin(bool usesegpar=false) const
Definition: EdbPattern.h:198
Float_t Zmax() const
Definition: EdbPattern.cxx:1213
ss
Definition: energy.C:62

◆ MakeSelector()

int EdbTrackP::MakeSelector ( EdbSegP ss,
bool  followZ = true 
)

1072{
1073 if(N()<2) return 0;
1074 ss.SetCOV( GetSegment(0)->COV() ); // TODO ?
1075 const EdbSegP *tr;
1076 if (NF())
1077 {
1078 if( followZ ) tr = TrackZmax();
1079 else tr = TrackZmin();
1080 }
1081 else
1082 {
1083 if( followZ ) tr = GetSegmentLast();
1084 else tr = GetSegmentFirst();
1085 }
1086 ss.SetTX(tr->TX());
1087 ss.SetTY(tr->TY());
1088 ss.SetX(tr->X());
1089 ss.SetY(tr->Y());
1090 ss.SetZ(tr->Z());
1091 ss.SetPID(tr->PID());
1092
1093 if( tr->PID() > GetSegment(0)->PID() ) return 1;
1094 if( tr->PID() > GetSegment(N()-1)->PID() ) return 1;
1095 return -1;
1096}
EdbSegP * GetSegmentLast() const
Definition: EdbPattern.h:190
EdbSegP * GetSegmentFirst() const
Definition: EdbPattern.h:189

◆ N()

Int_t EdbTrackP::N ( ) const
inline
177{ return (eS) ? eS->GetSize() : 0; }

◆ N0()

Int_t EdbTrackP::N0 ( ) const
inline
163{return eN0;}

◆ NF()

Int_t EdbTrackP::NF ( ) const
inline
178{ return (eSF)? eSF->GetSize() : 0; }

◆ Npl()

Int_t EdbTrackP::Npl ( ) const
inline
171{return eNpl;}

◆ PDG()

Int_t EdbTrackP::PDG ( ) const
inline
150{return ePDG;}

◆ PerrDown()

Float_t EdbTrackP::PerrDown ( ) const
inline
261{return ePerrDown;}
Float_t ePerrDown
error of P() in lower direction, obtained by MCS,or shower-algorithm
Definition: EdbPattern.h:125

◆ PerrUp()

Float_t EdbTrackP::PerrUp ( ) const
inline
260{return ePerrUp;}
Float_t ePerrUp
error of P() in upper direction, obtained by MCS,or shower-algorithm
Definition: EdbPattern.h:124

◆ Print()

void EdbTrackP::Print ( )
1153{
1154 int nseg=0, nsegf=0;
1155 nseg = N();
1156 nsegf = NF();
1157
1158 printf("EdbTrackP with %d segments and %d fitted segments\n", nseg, nsegf );
1159 printf("particle mass = %f\n", M() );
1160 ((EdbSegP*)this)->Print();
1161
1162 if(nseg)
1163 for(int i=0; i<nseg; i++)
1164 ((EdbSegP*)(eS->At(i)))->Print();
1165
1166}

◆ PrintNice()

void EdbTrackP::PrintNice ( )
1170{
1171 int nseg=0, nsegf=0;
1172 nseg = N();
1173 nsegf = NF();
1174
1175 printf("EdbTrackP with %d segments and %d fitted segments. Estimated BT Scanning efficiency = %.02f \n", nseg, nsegf, GetBTEfficiency() );
1176 printf(" N PID ID X Y Z TX TY W P Flag MC Track Chi2 Prob Mass\n");
1177 printf(" %3d %8d %13.2f %13.2f %13.2f %7.4f %7.4f %5.1f %7.2f %7d %7d %7d %7.4f %7.4f %5.3f\n",
1178 PID(), ID(),X(),Y(),Z(), TX(), TY(), W(), P(), Flag(), MCEvt(), Track(), Chi2(), Prob(), M());
1179
1180 EdbSegP *s=0;
1181 if(nseg)
1182 for(int i=0; i<nseg; i++) {
1183 s = GetSegment(i);
1184 printf("%3d %3d %8d %13.2f %13.2f %13.2f %7.4f %7.4f %5.1f %7.2f %7d %7d %7d %7.4f %7.4f\n",
1185 i, s->PID(), s->ID(),s->X(),s->Y(),s->Z(), s->TX(),s->TY(), s->W(), s->P(), s->Flag(), s->MCEvt(), s->Track(),s->Chi2(),s->Prob());
1186 }
1187}
Float_t Prob() const
Definition: EdbSegP.h:156
Float_t Chi2() const
Definition: EdbSegP.h:157
Int_t MCEvt() const
Definition: EdbSegP.h:145
Double_t GetBTEfficiency()
Definition: EdbPattern.cxx:1134

◆ RemoveAliasSegments()

int EdbTrackP::RemoveAliasSegments ( )
516{
517 int nalias=0;
518 EdbSegP *s=0;
519 for(int i=0; i<N(); i++) {
520 s = GetSegment(i);
521 if( s->Track() != ID()) {
522 nalias++;
524 }
525 }
526 return nalias;
527}
void RemoveSegment(EdbSegP *s)
Definition: EdbPattern.h:219

◆ RemoveSegment()

void EdbTrackP::RemoveSegment ( EdbSegP s)
inline
220 {
221 if(!eS) return;
222 eS->Remove(s);
223 s->SetTrack(-1);
224 SetCounters();
225 }
void SetCounters()
Definition: EdbPattern.h:159

◆ Set0()

void EdbTrackP::Set0 ( )
433{
434 ((EdbSegP*)this)->Set0();
435 eS=0;
436 eSF=0;
437 eM=0;
438 eDE=0;
439 ePDG=-999;
440 eVTAS = 0;
441 eVTAE = 0;
442 eNpl=0;
443 eN0=0;
444 ePerrUp=0;
445 ePerrDown=0;
446 return;
447}

◆ SetCounters()

void EdbTrackP::SetCounters ( )
inline
159{ SetNpl(); SetN0(); }

◆ SetDE()

void EdbTrackP::SetDE ( float  de)
inline
165{ eDE=de; }

◆ SetM()

void EdbTrackP::SetM ( float  m)
inline
154{ eM=m; }

◆ SetN0() [1/2]

void EdbTrackP::SetN0 ( )
inline
162{ eN0 = eNpl-N(); }

◆ SetN0() [2/2]

void EdbTrackP::SetN0 ( int  n0)
inline
161{ eN0 = n0; }

◆ SetNpl() [1/2]

void EdbTrackP::SetNpl ( )
inline
170 { if(eS) eNpl = 1+TMath::Abs(GetSegment(0)->PID() - GetSegment(N()-1)->PID()); }

◆ SetNpl() [2/2]

void EdbTrackP::SetNpl ( int  npl)
inline
168{ eNpl=npl; }

◆ SetOwner()

void EdbTrackP::SetOwner ( )
inline
139{ if(eS) eS->SetOwner(true); }

◆ SetPDG()

void EdbTrackP::SetPDG ( int  pdg)
inline
149{ ePDG=pdg; }

◆ SetPerr()

void EdbTrackP::SetPerr ( Float_t  perrDown,
Float_t  perrUp 
)
1231{
1232 ePerrUp=perrUp;
1233 ePerrDown=perrDown;
1234}

◆ SetPerrDown()

void EdbTrackP::SetPerrDown ( Float_t  perrDown)
1226{
1227 ePerrDown=perrDown;
1228}

◆ SetPerrUp()

void EdbTrackP::SetPerrUp ( Float_t  perrUp)
1221{
1222 ePerrUp=perrUp;
1223}

◆ SetSegmentsTrack() [1/2]

int EdbTrackP::SetSegmentsTrack ( )
inline
247{return SetSegmentsTrack(ID());}
int SetSegmentsTrack()
Definition: EdbPattern.h:247

◆ SetSegmentsTrack() [2/2]

int EdbTrackP::SetSegmentsTrack ( int  id)
inline
246{for(int i=0; i<N(); i++) GetSegment(i)->SetTrack(id); return N();}
void SetTrack(int trid)
Definition: EdbSegP.h:131

◆ SubstituteSegment()

void EdbTrackP::SubstituteSegment ( EdbSegP sold,
EdbSegP snew 
)
inline
227 {
228 if(!eS) return;
229 eS->Remove(sold);
230 eS->Add(snew);
231 sold->SetTrack(-1);
232 }

◆ TrackEnd()

const EdbSegP * EdbTrackP::TrackEnd ( ) const
1247{
1248 if(!N()) return (EdbSegP*)this;
1249 if(Dir()>0) return GetSegmentLast();
1250 if(Dir()<0) return GetSegmentFirst();
1251 return (EdbSegP*)this;
1252}
Int_t Dir() const
Definition: EdbPattern.h:207

◆ TrackExtremity()

EdbSegP * EdbTrackP::TrackExtremity ( bool  zpos,
bool  usesegpar = false 
) const
inline
202 { return zpos? TrackZmin(usesegpar) : TrackZmax(usesegpar); } // 0-end, 1-start (as in vertex class)

◆ TrackStart()

const EdbSegP * EdbTrackP::TrackStart ( ) const
1238{
1239 if(!N()) return (EdbSegP*)this;
1240 if(Dir()<0) return GetSegmentLast();
1241 if(Dir()>0) return GetSegmentFirst();
1242 return (EdbSegP*)this;
1243}

◆ TrackZmax()

EdbSegP * EdbTrackP::TrackZmax ( bool  usesegpar = false) const
inline
199{ if(usesegpar || (!eSF)) return GetSegmentLast(); else return GetSegmentFLast(); }
EdbSegP * GetSegmentFLast() const
Definition: EdbPattern.h:193

◆ TrackZmin()

EdbSegP * EdbTrackP::TrackZmin ( bool  usesegpar = false) const
inline
198{ if(usesegpar || (!eSF)) return GetSegmentFirst(); else return GetSegmentFFirst(); }
EdbSegP * GetSegmentFFirst() const
Definition: EdbPattern.h:192

◆ Transform()

void EdbTrackP::Transform ( const EdbAffine2D tr)

apply transformation to all track elements

498{
500 ((EdbSegP*)(this))->Transform(&aff);
501 for(int i=0; i<N(); i++) GetSegment(i)->Transform(&aff);
502 for(int i=0; i<NF(); i++) GetSegmentF(i)->Transform(&aff);
503}
virtual void Transform(const EdbAffine2D *a)
void Transform(const EdbAffine2D &tr)
Definition: EdbPattern.cxx:497

◆ Vertex()

EdbVertex * EdbTrackP::Vertex ( int  zpos)
inline
147{return zpos? VertexS(): VertexE();}
EdbVertex * VertexS()
Definition: EdbPattern.cxx:1191
EdbVertex * VertexE()
Definition: EdbPattern.cxx:1198

◆ VertexE()

EdbVertex * EdbTrackP::VertexE ( )
1199{
1200 if(eVTAE) return eVTAE->GetVertex();
1201 return 0;
1202}
EdbVertex * GetVertex() const
Definition: EdbVertex.h:52

◆ VertexS()

EdbVertex * EdbTrackP::VertexS ( )
1192{
1193 if(eVTAS) return eVTAS->GetVertex();
1194 return 0;
1195}

◆ VTAE()

EdbVTA * EdbTrackP::VTAE ( ) const
inline
144{return eVTAE;}

◆ VTAS()

EdbVTA * EdbTrackP::VTAS ( ) const
inline
143{return eVTAS;}

◆ Wgrains()

float EdbTrackP::Wgrains ( ) const
660{
661 float w=0.;
662 int nseg=N();
663 for(int i=0; i<nseg; i++) w+=GetSegment(i)->W();
664 return w;
665}

◆ Wmean()

Float_t EdbTrackP::Wmean ( ) const
582{
583 int n = N();
584 double wtot=0;
585 for(int i=0; i<n; i++) wtot+=GetSegment(i)->W();
586 wtot/=n;
587 return (Float_t)wtot;
588}

◆ Zend()

Float_t EdbTrackP::Zend ( ) const
inline
211{return TrackEnd()->Z();}
const EdbSegP * TrackEnd() const
Definition: EdbPattern.cxx:1246

◆ Zmax()

Float_t EdbTrackP::Zmax ( ) const
1214{
1215 Float_t zmax = Z();
1216 if (N() && GetSegmentLast()->Z() > zmax) zmax = GetSegmentLast()->Z();
1217 return zmax;
1218}

◆ Zmin()

Float_t EdbTrackP::Zmin ( ) const
1206{
1207 Float_t zmin = Z();
1208 if (N() && GetSegmentFirst()->Z() < zmin) zmin = GetSegmentFirst()->Z();
1209 return zmin;
1210}

◆ Zstart()

Float_t EdbTrackP::Zstart ( ) const
inline
210{return TrackStart()->Z();}
const EdbSegP * TrackStart() const
Definition: EdbPattern.cxx:1237

Member Data Documentation

◆ eDE

Float_t EdbTrackP::eDE
private

total energy loss of the particle between first and last segments

◆ eM

Float_t EdbTrackP::eM
private

invariant mass of the particle

◆ eN0

Int_t EdbTrackP::eN0
private

number of holes (if any)

◆ eNpl

Int_t EdbTrackP::eNpl
private

number of plates passed through

◆ ePDG

Int_t EdbTrackP::ePDG
private

particle ID from PDG

◆ ePerrDown

Float_t EdbTrackP::ePerrDown
private

error of P() in lower direction, obtained by MCS,or shower-algorithm

◆ ePerrUp

Float_t EdbTrackP::ePerrUp
private

error of P() in upper direction, obtained by MCS,or shower-algorithm

◆ eS

TSortedList* EdbTrackP::eS
private

array of segments

◆ eSF

TSortedList* EdbTrackP::eSF
private

array of fitted segments

◆ eVTAE

EdbVTA* EdbTrackP::eVTAE
private

vertex track end is attached to

◆ eVTAS

EdbVTA* EdbTrackP::eVTAS
private

vertex track start is attached to


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