FEDRA emulsion software from the OPERA Collaboration
TBinTracking Class Reference

#include <bitview.h>

Inheritance diagram for TBinTracking:
Collaboration diagram for TBinTracking:

Public Member Functions

EdbViewAdoptSegment (EdbView *view, float sx, float sy, float Xmin, float Xmax, float Ymin, float Ymax, float bs, int thr, float accept, int inside=1, int excludeCommonClusters=0, int rep=1)
 
void FillByteArray (TBitView *Tv)
 
TObjArray * FindBinTracks (int thr)
 
TByteMatrixGetElement (double tgx, double tgy)
 
TByteMatrixGetElement (int n)
 
int GetNelements (void)
 
TH2F * Histo (double tgx, double tgy)
 
void ImproveTracksArray (TObjArray *tracks, int thr, float accept)
 
TByteMatrixoperator[] (int n)
 
 TBinTracking ()
 
 TBinTracking (double tgxmin, double tgxmax, double tgymin, double tgymax, double tgstep)
 
 TBinTracking (double tgxmin, double tgxmax, double tgymin, double tgymax, double tgstep, double dx, double dy, double dz)
 
virtual ~TBinTracking ()
 

Public Attributes

double dX
 
double dY
 
double dZ
 
double Tgstep
 
double Tgxmax
 
double Tgxmin
 
double Tgymax
 
double Tgymin
 
float ZBase
 

Protected Attributes

TObjArray * TByteArray
 

Constructor & Destructor Documentation

◆ TBinTracking() [1/3]

TBinTracking::TBinTracking ( )
531{
532 TByteArray = new TObjArray();
533 dX = 1.0;
534 dY = 1.0;
535 dZ = 1.0;
536// TTrackArray = NULL;
537}
double dZ
Definition: bitview.h:71
TObjArray * TByteArray
Definition: bitview.h:64
double dY
Definition: bitview.h:71
double dX
Definition: bitview.h:71

◆ TBinTracking() [2/3]

TBinTracking::TBinTracking ( double  tgxmin,
double  tgxmax,
double  tgymin,
double  tgymax,
double  tgstep 
)
541{
542 Tgxmin = tgxmin;
543 Tgxmax = tgxmax;
544 Tgymin = tgymin;
545 Tgymax = tgymax;
546 Tgstep = tgstep;
547 TByteArray = new TObjArray();
548 dX = 1.0;
549 dY = 1.0;
550 dZ = 1.0;
551// TTrackArray = NULL;
552}
double Tgxmin
Definition: bitview.h:68
double Tgymax
Definition: bitview.h:68
double Tgstep
Definition: bitview.h:68
double Tgymin
Definition: bitview.h:68
double Tgxmax
Definition: bitview.h:68

◆ TBinTracking() [3/3]

TBinTracking::TBinTracking ( double  tgxmin,
double  tgxmax,
double  tgymin,
double  tgymax,
double  tgstep,
double  dx,
double  dy,
double  dz 
)
556{
557 Tgxmin = tgxmin;
558 Tgxmax = tgxmax;
559 Tgymin = tgymin;
560 Tgymax = tgymax;
561 Tgstep = tgstep;
562 TByteArray = new TObjArray();
563 dX = dx;
564 dY = dy;
565 dZ = dz;
566}
brick dz
Definition: RecDispMC.C:107

◆ ~TBinTracking()

TBinTracking::~TBinTracking ( )
virtual
570{
571 if (TByteArray) {
572 TByteArray->Delete();
573 delete TByteArray;
574 }
575/* if (TTrackArray) {
576 TTrackArray->Delete();
577 delete TTrackArray;
578 }*/
579}

Member Function Documentation

◆ AdoptSegment()

EdbView * TBinTracking::AdoptSegment ( EdbView view,
float  sx,
float  sy,
float  Xmin,
float  Xmax,
float  Ymin,
float  Ymax,
float  bs,
int  thr,
float  accept,
int  inside = 1,
int  excludeCommonClusters = 0,
int  rep = 1 
)
736{
737 EdbView *vmod = view;
738// EdbView *vmod = new EdbView(*view);
740// TClonesArray *SAr = new TClonesArray("EdbSegment");
741 TObjArray *SAr = new TObjArray();
742 TBitView* Tv = new TBitView();
743 Tv->bitSize = bs;
744 Tv->FillBitView(vmod,cfx,cfy,Xmin,Xmax,Ymin,Ymax,inside,rep);
745 FillByteArray(Tv);
746 TObjArray* tracks = FindBinTracks(thr);
748 ImproveTracksArray(tracks,thr,bs); // was accept instead of bs
750 float X[200], Y[200] ,Z[200];
751 int ntracks = tracks->GetEntries();
752 track* tr = 0;
753 EdbCluster* cl = 0;
754 int nclusters = vmod->Nclusters();
755 int sid = 0;
756 int cyclelimit = 10;
757// cout<<"tracks="<<ntracks<<endl;
758 for (long t=0; t<ntracks; t++) {
759 tr = (track *)(tracks->At(t));
760 tr->x = tr->x * bs + Xmin;
761 tr->y = tr->y * bs + Ymin;
762 EdbSegment* s = new EdbSegment();
763 float resinit = 0;
764 int pulsinit = 0;
765 int cyc = 0;
766CYCLE:
767 cyc++;
768 int n = 0;
770/* float teta,phi;
771 TrackAngle(*tr,&teta,&phi);
772*/
773 for (long i=0; i<nclusters; i++) {
774 int flag = false;
775 cl = vmod->GetCluster(i);
776 float xc = cl->GetX();
777 float yc = cl->GetY();
778 float zc = cl->GetZ();
779 if (In(xc,yc,Xmin,Xmax,Ymin,Ymax)&&SatCondAccept(tr,xc,yc,zc,dX,dY,dZ,accept)) {
780/* (!(inside)||
781 ((zc>(view->GetNframesTop()?Tv->GetZ2():Tv->GetZ4()))
782 &&(zc<(view->GetNframesTop()?Tv->GetZ1():Tv->GetZ3()))))) {
783*/
784 //printf("distance=%f\n",Distance(cl->GetX(), cl->GetY(),cl->GetZ(),*tr));
785 //if (rep!=1) {
786 if (inside) {
787 for (int j=0; j<Tv->GetNumberOfLayers(); j++) {
788// cout<<c->GetZ()<<" "<<GetLayer(j)->GetZ()<<endl;
789 if (zc==Tv->GetLayer(j)->GetZ()) {
790 flag = true;
791 break;
792 }
793 }
794 if (!flag) continue;
795 }
797 if (int last = s->GetNelements()) {
798 EdbCluster *clsucc = ((EdbCluster *)(s->GetElements()->Last()));
799 if (zc==clsucc->GetZ()) {
800 if (Distance(clsucc->GetX(),
801 clsucc->GetY(),
802 clsucc->GetZ(),
803 *tr) <= Distance(xc,yc,zc,*tr)) continue;
804 else {
805 (*(s->GetElements()))[last-1] = cl;
806 X[n-1] = xc;
807 Y[n-1] = yc;
808 Z[n-1] = zc;
809 continue;
810 }
811 }
812 }
814 s->AddElement(cl);
815 X[n] = xc;
816 Y[n] = yc;
817 Z[n] = zc;
818 n++;
819 }
820 }
821 float Ax,Bx,erAx,erBx,sx,Ay,By,erAy,erBy,sy;;
822 int flagx = MY_LFIT(Z,X,n,1,Ax,Bx,erAx,erBx,sx);
823 int flagy = MY_LFIT(Z,Y,n,1,Ay,By,erAy,erBy,sy);
824 if (flagx&&flagy) {
825 int puls = s->GetNelements();
826 int side = ((EdbCluster *)(s->GetElements()->At(0)))->GetSide();
827 float maxz = ((EdbCluster *)(s->GetElements()->At(0)))->GetZ();
828 float minz = maxz;
829 for (long j=0; j<puls; j++) {
830 cl = (EdbCluster *)(s->GetElements()->At(j));
831 if (minz>=cl->GetZ()) minz = cl->GetZ();
832 if (maxz<cl->GetZ()) maxz = cl->GetZ();
833 }
834 float res = sqrt((sx*sx)+(sy*sy));
836 if (((pulsinit!=puls)||(res!=resinit))&&cyc<cyclelimit) {
837 pulsinit = puls;
838 resinit = res;
839 float zl = ZBase; //side?maxz:minz;
840 tr->tgx = Ax;
841 tr->tgy = Ay;
842 tr->z = zl;
843 tr->x = Ax*zl+Bx;
844 tr->y = Ay*zl+By;
845 s->GetElements()->Clear();
846 goto CYCLE;
847 }
849// cout<<"cycle="<<cyc<<endl;
850 float dz = maxz-minz;
851 s->Set(tr->x,tr->y,tr->z,Ax,Ay,dz,side,puls,sid);
853/* for (int nclust=0; nclust<puls; nclust++) {
854 cl = (EdbCluster *)(s->GetElements()->At(nclust));
855 printf("< clnum=%d ", nclust);
856 cl->Print();
857 }
858*/
859// new ((*SAr)[sid]) EdbSegment(*s);
860 SAr->Add(s);
861 sid++;
862 }
863// delete s;
864 }
865 delete Tv;
866 delete tracks;
867
869 if (sid==0) return vmod; //SAr;
870// int newsid=sid;
871 EdbSegment* seg0 = 0;
872 EdbSegment* seg = 0;
873 for (int s0=0; s0<sid-1; s0++) {
874//REPEAT:
875 seg0 = (EdbSegment *)(SAr->At(s0));
876 int newid=seg0->GetID();
877 int numofdel = 0;
878 int rep = 0;
879
880 for (int s=s0+1; s<sid; s++) {
881 seg = (EdbSegment *)(SAr->At(s));
882/* if ((seg->GetX0()==seg0->GetX0())&&(seg->GetY0()==seg0->GetY0())&&
883 (seg->GetZ0()==seg0->GetZ0())&&
884 (seg->GetTx()==seg0->GetTx())&&(seg->GetTy()==seg0->GetTy())&&
885 (seg->GetDz()==seg0->GetDz())&&(seg->GetPuls()==seg0->GetPuls())) {
886*/
887 if (((fabs(seg->GetX0()-seg0->GetX0())<bs&&fabs(seg->GetY0()-seg0->GetY0())<bs)&&
888// (seg->GetZ0()==seg0->GetZ0())&&
889 fabs(seg->GetTx()-seg0->GetTx())<Tgstep&&fabs(seg->GetTy()-seg0->GetTy())<Tgstep)
890// &&(seg->GetDz()==seg0->GetDz())&&(seg->GetPuls()==seg0->GetPuls())
891// added at 13.11.2004 ///////////////////////////////
892 ||(excludeCommonClusters
893 &&CommonClusters(seg0,seg))
895 )
896 {
898 if (seg->GetPuls()>seg0->GetPuls()) {
899 (*(SAr))[s] = seg0;
900 (*(SAr))[s0] = seg;
901 ((EdbSegment *)(SAr->At(s0)))->SetID(seg0->GetID());
902 seg0 = seg;
903 rep=1;
904 }
906 delete (SAr->RemoveAt(s));
907 numofdel++;
908// sid--;
909 }
910 else {
911 newid++;
912 seg->SetID(newid);
913// seg0 = seg;
914 }
915 }
916 SAr->Compress();
917// printf("sid=%d, del=%d\n",sid,numofdel);
918 sid -= numofdel;
919 if (rep) s0--;
920 }
923
924 vmod->GetHeader()->SetNsegments(sid);
925// printf("sid=%d\n",sid);
926// EdbSegment* seg = 0;
927 for (int ss=0; ss<sid; ss++) {
928 seg = (EdbSegment *)(SAr->At(ss));
929 EdbSegment *segnew = new EdbSegment(
930 seg->GetX0(),
931 seg->GetY0(),
932 seg->GetZ0(),
933 seg->GetTx(),
934 seg->GetTy(),
935 seg->GetDz(),
936 seg->GetSide(),
937 seg->GetPuls(),
938 seg->GetID());
939 vmod->AddSegment(segnew);
940 delete segnew;
941 seg0 = ((EdbSegment *)(vmod->GetSegments()->At(ss)));
942// cout<<seg->GetNelements()<<" "<<seg->GetID()<<endl;
943 for (int cc=0; cc<seg->GetNelements(); cc++) {
944 cl = (EdbCluster *)(seg->GetElements()->At(cc));
945 seg0->AddElement(cl);
946 cl->SetSegment(seg->GetID());
947 }
948// vmod->AddSegment(seg);
949 }
950 SAr->Delete();
951 delete SAr;
953 return vmod;
954}
float zl[dim]
Definition: RecDispMC_Profiles.C:61
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
Expr< UnaryOp< Fabs< T >, Expr< A, T, D >, T >, T, D > fabs(const Expr< A, T, D > &rhs)
Definition: UnaryOperators.hh:96
int In(float X, float Y, float Xmin, float Xmax, float Ymin, float Ymax)
Definition: bitview.cxx:18
int CommonClusters(EdbSegment *seg0, EdbSegment *seg1)
Definition: bitview.cxx:28
float Distance(float xp, float yp, float zp, float x0, float y0, float z0, float teta, float phi)
Definition: bitview.cxx:88
int SatCondAccept(track *tr, float xc, float yc, float zc, double dx, double dy, double dz, float accept)
Definition: bitview.cxx:153
int MY_LFIT(float *X, float *Y, int L, int KEY, float &A, float &B, float &erA, float &erB, float &s)
Definition: bitview.cxx:180
Definition: EdbCluster.h:19
void SetSegment(int seg)
Definition: EdbCluster.h:49
Float_t GetX() const
Definition: EdbCluster.h:51
Float_t GetY() const
Definition: EdbCluster.h:52
Float_t GetZ() const
Definition: EdbCluster.h:53
virtual Float_t GetDz() const
Definition: EdbSegment.h:43
virtual Float_t GetX0() const
Definition: EdbSegment.h:38
virtual Float_t GetTx() const
Definition: EdbSegment.h:41
virtual Float_t GetZ0() const
Definition: EdbSegment.h:40
virtual Float_t GetY0() const
Definition: EdbSegment.h:39
virtual Float_t GetTy() const
Definition: EdbSegment.h:42
segment of the track
Definition: EdbSegment.h:63
Int_t GetPuls() const
Definition: EdbSegment.h:90
void AddElement(TObject *element)
Definition: EdbSegment.cxx:150
Int_t GetNelements() const
Definition: EdbSegment.h:114
TObjArray * GetElements() const
Definition: EdbSegment.h:117
Int_t GetID() const
Definition: EdbSegment.h:92
void SetID(int id)
Definition: EdbSegment.h:96
Int_t GetSide() const
Definition: EdbSegment.h:89
void SetNsegments(int nseg)
Definition: EdbView.h:86
Base scanning data object: entry into Run tree.
Definition: EdbView.h:134
TClonesArray * GetSegments() const
Definition: EdbView.h:165
EdbViewHeader * GetHeader() const
Definition: EdbView.h:163
EdbSegment * AddSegment(float x, float y, float z, float tx, float ty, float dz=0, int side=0, int puls=0, int id=-1)
Definition: EdbView.h:231
Int_t Nclusters() const
Definition: EdbView.h:215
EdbCluster * GetCluster(int i) const
Definition: EdbView.h:218
float ZBase
Definition: bitview.h:69
void ImproveTracksArray(TObjArray *tracks, int thr, float accept)
Definition: bitview.cxx:655
void FillByteArray(TBitView *Tv)
Definition: bitview.cxx:604
TObjArray * FindBinTracks(int thr)
Definition: bitview.cxx:633
float GetZ(void)
Definition: bitmatrix.h:37
Definition: bitview.h:27
void FillBitView(EdbView *view, float sx, float sy, float Xmin, float Xmax, float Ymin, float Ymax, int inside=1, int rep=1)
Definition: bitview.cxx:407
TBitMatrix * GetLayer(int n)
Definition: bitview.cxx:401
int GetNumberOfLayers(void)
Definition: bitview.h:47
float bitSize
Definition: bitview.h:35
Definition: bitview.h:14
TTree * t
Definition: check_shower.C:4
s
Definition: check_shower.C:55
TTree * tracks
Definition: check_tr.C:19
ss
Definition: energy.C:62
Double_t X
Definition: tlg2couples.C:76
Double_t Y
Definition: tlg2couples.C:76
Double_t Z
Definition: tlg2couples.C:104

◆ FillByteArray()

void TBinTracking::FillByteArray ( TBitView Tv)
605{
606 ZBase = Tv->ZBase;
607 long tgxi = (long)((Tgxmax-Tgxmin)/Tgstep +1);
608 long tgyi = (long)((Tgymax-Tgymin)/Tgstep +1);
609 long xs, ys;
610 TBitMatrix *TBit = 0;
611 Tv->GetLayer(0)->GetSize(&xs,&ys);
612 for (long j=0; j<tgyi; j++) {
613 for (long i=0; i<tgxi; i++) {
614 double tgx = Tgxmin + (i+0.5)*Tgstep;
615 double tgy = Tgymin + (j+0.5)*Tgstep;
616 TByteMatrix *TByte = new TByteMatrix(xs,ys,tgx,tgy);
617 for (long k=0; k<Tv->GetNumberOfLayers(); k++) {
618 TBit = Tv->GetLayer(k);
619 long shiftX = (long)((TBit->GetZ()-ZBase)*tgx/Tv->bitSize);
620 long shiftY = (long)((TBit->GetZ()-ZBase)*tgy/Tv->bitSize);
621 TBit->MoveMatrix(-shiftY,-shiftX);
622 *TByte += *TBit;
623 }
624// TByte->Print(0);
625// _getch();
626 TByteArray->Add(TByte);
627// delete TByte;
628 }
629 }
630}
long GetSize(long *ix_size=NULL, long *iy_size=NULL)
Definition: basematrix.cxx:190
void MoveMatrix(long ud, long lr)
Definition: basematrix.cxx:313
Definition: bitmatrix.h:13
float ZBase
Definition: bitview.h:34
Definition: bytematrix.h:12
MFTYPE32 long(MFTYPE MPTYPE *MOCRHOOKFCTPTR)(long HookType
Definition: milocr.h:32

◆ FindBinTracks()

TObjArray * TBinTracking::FindBinTracks ( int  thr)
634{
635 TByteMatrix* TByte = 0;
636 TObjArray* Tob = new TObjArray();
637 for (long i=0; i<GetNelements(); i++) {
638 TByte = GetElement(i);
639 double tgx, tgy;
640 TByte->GetAngle(&tgx, &tgy);
641 TObjArray* points = TByte->Pick8(0,thr);
642 MyPoint* point = 0;
643 for (long ii=0; ii<points->GetEntries(); ii++) {
644 point = (MyPoint *)(points->At(ii));
645 Tob->Add(new track(tgx, tgy, point->x, point->y, ZBase,point->height));
647// TByte->Print(0);
649 }
650 delete points;
651 }
652 return Tob;
653}
Definition: basematrix.h:9
long y
Definition: basematrix.h:12
int height
Definition: basematrix.h:13
long x
Definition: basematrix.h:11
TByteMatrix * GetElement(int n)
Definition: bitview.cxx:582
int GetNelements(void)
Definition: bitview.h:80
void GetAngle(double *tgx, double *tgy)
Definition: bytematrix.cxx:99
TObjArray * Pick8(int type, int thr)
Definition: bytematrix.cxx:340

◆ GetElement() [1/2]

TByteMatrix * TBinTracking::GetElement ( double  tgx,
double  tgy 
)
590{
591 if ((tgx>=Tgxmin)&&(tgx<Tgxmax)&&(tgy>=Tgymin)&&(tgy<Tgymax)) {
592 TByteMatrix *TByte = 0;
593 double Tgx, Tgy;
594 for (long i=0; i<GetNelements(); i++) {
595 TByte = GetElement(i);
596 TByte->GetAngle(&Tgx,&Tgy);
597 if ((fabs(Tgx-tgx)<=Tgstep)&&(fabs(Tgy-tgy)<=Tgstep)) return TByte;
598 }
599 }
600 return NULL;
601}
#define NULL
Definition: nidaqmx.h:84

◆ GetElement() [2/2]

TByteMatrix * TBinTracking::GetElement ( int  n)
583{
584 if (n<TByteArray->GetEntries()) return (TByteMatrix *)TByteArray->At(n);
585 else return NULL;
586}
cout<< tr-> GetEntries()<< endl

◆ GetNelements()

int TBinTracking::GetNelements ( void  )
inline
80{return TByteArray->GetEntries();}

◆ Histo()

TH2F * TBinTracking::Histo ( double  tgx,
double  tgy 
)
958{
959 TByteMatrix *tbyte = 0;
960 if ((tbyte = GetElement(tgx,tgy)) == NULL) return 0;
961 long xs, ys;
962 tbyte->GetSize(&xs,&ys);
963 printf("xs=%ld, ys=%ld\n",xs,ys);
964 TH2F* h2 = new TH2F("h2","h2",xs,0.,xs,ys,0., ys);
965 for (long j=0; j<ys; j++) {
966 for (long i=0; i<xs; i++) {
967 h2->Fill(i+.5,j+.5,tbyte->GetByte(0,i,j));
968 }
969 }
970 return h2;
971}
int GetByte(int type, long i, long j)
Definition: bytematrix.cxx:106
TH1F * h2
Definition: energy.C:19

◆ ImproveTracksArray()

void TBinTracking::ImproveTracksArray ( TObjArray *  tracks,
int  thr,
float  accept 
)
656{
657 track* tr = 0;
658 track* tr1 = 0;
659
660 int xn[8] = { 1, -1, 0, 0, 1, 1, -1, -1 }; //pixel 3x3 suburbs
661 int yn[8] = { 0, 0, 1, -1, 1, -1, 1, -1 };
662
663 long nc = (long)((Tgymax-Tgymin)/Tgstep);
664 long nr = (long)((Tgxmax-Tgxmin)/Tgstep);
665
666 int pic=0,pic1=0;
667
668// printf("hist: %d %d\n",nc,nr);
669// int npeaks=0;
670 int ic,ir,in; // declare here to fix marazmatic eror messages of VC++
671 for(ic=0; ic<nc; ic++) {
672 for(ir=0; ir<nr; ir++) {
673 float tgx0 = Tgxmin + ir*Tgstep;
674 float tgx1 = Tgxmin + (ir+1)*Tgstep;
675 float tgy0 = Tgymin + ic*Tgstep;
676 float tgy1 = Tgymin + (ic+1)*Tgstep;
677 int nmax = -1;
678 for (int i=0; i<tracks->GetEntries(); i++) {
679 if ( (tr = (track *)(tracks->At(i))) )
680 if (In(tr->tgx,tr->tgy,tgx0,tgx1,tgy0,tgy1)) {
681 nmax = i;
682 pic = tr->height;
683 pic1 = pic+1;
684 break;
685 }
686 }
687 if (nmax==-1) continue;
688 if (pic<thr) {
689 delete (track *)(tracks->RemoveAt(nmax));
690 tracks->Compress();
691 continue;
692 }
693// tr->height = pic1;
694// pic = GetByte(type,ic,ir);
695// if( pix < thr) goto NEXT;
696// pix1=pix+1;
697 int n = 0;
698 int nn[8];
699 for(in=0; in<8; in++) {
700 int icc = ic+xn[in];
701 int irr = ir+yn[in];
702 tgx0 = Tgxmin + irr*Tgstep;
703 tgx1 = Tgxmin + (irr+1)*Tgstep;
704 tgy0 = Tgymin + icc*Tgstep;
705 tgy1 = Tgymin + (icc+1)*Tgstep;
706 for (int i=0; i<tracks->GetEntries(); i++) {
707 if ((tr1 = (track *)(tracks->At(i))) == tr) continue;
708 if (In(tr1->tgx,tr1->tgy,tgx0,tgx1,tgy0,tgy1)&&
709 (fabs(tr->x-tr1->x)<accept&&fabs(tr->y-tr1->y)<accept)) {
710 if (pic1<=tr1->height) goto NEXT;
711 nn[n] = i;
712 n++;
713 }
714 }
715 }
716 if (n) {
717//printf(" %d \n",n);
718 for (int i=0; i<n; i++) delete (track *)(tracks->RemoveAt(nn[i]));
719 tr->height = pic1;
720 tracks->Compress();
721 }
722// pp = GetByte(type, ic+xn[in], ir+yn[in] );
723// if( pix1 <= pp ) goto NEXT;
724// SetByte(type,ic,ir,pix1);
725// printf("peak8(%d,%d) = %d\n",ic,ir,pix);
726// Ta->Add(new MyPoint(ic,ir,pix));
727NEXT:
728 continue;
729 }
730 }
731}
double tgy
Definition: bitview.h:17
double tgx
Definition: bitview.h:16
double x
Definition: bitview.h:18
double y
Definition: bitview.h:19

◆ operator[]()

TByteMatrix * TBinTracking::operator[] ( int  n)
inline
79{return GetElement(n);}

Member Data Documentation

◆ dX

double TBinTracking::dX

◆ dY

double TBinTracking::dY

◆ dZ

double TBinTracking::dZ

◆ TByteArray

TObjArray* TBinTracking::TByteArray
protected

◆ Tgstep

double TBinTracking::Tgstep

◆ Tgxmax

double TBinTracking::Tgxmax

◆ Tgxmin

double TBinTracking::Tgxmin

◆ Tgymax

double TBinTracking::Tgymax

◆ Tgymin

double TBinTracking::Tgymin

◆ ZBase

float TBinTracking::ZBase

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