FEDRA emulsion software from the OPERA Collaboration
EdbView Class Reference

Base scanning data object: entry into Run tree. More...

#include <EdbView.h>

Inheritance diagram for EdbView:
Collaboration diagram for EdbView:

Public Member Functions

EdbClusterAddCluster (EdbCluster *c)
 
EdbClusterAddCluster (float x, float y, float z, float a, float v, int f, int s, int seg=-1)
 
void AddFrame (EdbFrame *frame)
 
void AddFrame (int id, float z, int ncl=0, int npix=0)
 
EdbSegmentAddSegment (EdbSegment &s)
 
void AddSegment (EdbSegment *s)
 
EdbSegmentAddSegment (float x, float y, float z, float tx, float ty, float dz=0, int side=0, int puls=0, int id=-1)
 
void AddTrack (EdbTrack *t)
 
int AttachClustersToSegments ()
 
int AttachClustersToSegmentsFast ()
 
int AttachClustersToSegmentsSlow ()
 
int AttachSegmentsToTracks ()
 
void Clear ()
 
void DeleteClustersFog ()
 
void Draw (Option_t *option="")
 
TPolyMarker3D * DrawClustersFog (Option_t *opt=0) const
 
TPolyMarker3D * DrawClustersSegments (Option_t *opt=0) const
 
 EdbView ()
 
void GenerateClustersFog (float density)
 
void GenerateClustersSegment (EdbSegment *segment, int n0=25, float sigma=.1)
 
void GenerateClustersTrack (EdbTrack *track, int n0=25, float sigma=.1)
 
void GenerateFrames (int n=32)
 
Int_t GetAreaID () const
 
EdbClusterGetCluster (int i) const
 
TClonesArray * GetClusters () const
 
voidGetClustersAddr ()
 
TList * GetClustersFrame (int frame) const
 
Float_t GetDZbot () const
 
Float_t GetDZtop () const
 
EdbFrameGetFrame (int i) const
 
TClonesArray * GetFrames () const
 
voidGetFramesAddr ()
 
EdbViewHeaderGetHeader () const
 
voidGetHeaderAddr ()
 
Long_t GetLastSystemTime () const
 
Int_t GetNframes () const
 
Int_t GetNframesBot () const
 
Int_t GetNframesTop () const
 
EdbSegmentGetSegment (int i) const
 
TClonesArray * GetSegments () const
 
voidGetSegmentsAddr ()
 
TObjArray * GetSegmentsClusters () const
 
Int_t GetTime () const
 
EdbTrackGetTrack (int i) const
 
TClonesArray * GetTracks () const
 
voidGetTracksAddr ()
 
Int_t GetViewID () const
 
Float_t GetXview () const
 
Float_t GetYview () const
 
Float_t GetZ0bot () const
 
Float_t GetZ0top () const
 
Float_t GetZ1 () const
 
Float_t GetZ2 () const
 
Float_t GetZ3 () const
 
Float_t GetZ4 () const
 
Int_t Nclusters () const
 
Int_t Nsegments () const
 
Int_t Ntracks () const
 
void Print (Option_t *opt=0) const
 
void PrintClusters (Option_t *opt=0) const
 
Int_t ReadView (char *fname)
 
void Scale (float zscale)
 
void SetAreaID (int id)
 
void SetCoordXY (float x, float y)
 
void SetCoordZ (float z1, float z2, float z3, float z4)
 
void SetLastSystemTime (Long_t time)
 
void SetNframes (int top, int bot)
 
void Shift (float zshift)
 
void Transform (const EdbAffine2D *aff)
 
float Xmin () const
 View coordinates are in pixels. More...
 
float Ymin () const
 starting from 0 More...
 
float ZFrameMax () const
 
float ZFrameMin () const
 
float Zmax () const
 
float Zmin () const
 
virtual ~EdbView ()
 

Private Attributes

TClonesArray * eClusters
 array of Clusters More...
 
TClonesArray * eFrames
 array of Frames (images) More...
 
EdbViewHeadereHeader
 View header. More...
 
Long_t eLastSystemTime
 system time when view was saved More...
 
TClonesArray * eSegments
 array of Segments More...
 
TClonesArray * eTracks
 array of Tracks More...
 

Detailed Description

Base scanning data object: entry into Run tree.

Constructor & Destructor Documentation

◆ EdbView()

EdbView::EdbView ( )

◆ ~EdbView()

EdbView::~EdbView ( )
virtual
41{
42 // Clear();
43
44 if(eClusters) { eClusters->Delete(); SafeDelete(eClusters); }
45 if(eSegments) { eSegments->Delete(); SafeDelete(eSegments); }
46 if(eTracks) { eTracks->Delete(); SafeDelete(eTracks); }
47 if(eHeader) SafeDelete(eHeader);
48 if(eFrames) { eFrames->Delete(); SafeDelete(eFrames); }
49}
TClonesArray * eClusters
array of Clusters
Definition: EdbView.h:140
TClonesArray * eSegments
array of Segments
Definition: EdbView.h:141
EdbViewHeader * eHeader
View header.
Definition: EdbView.h:138
TClonesArray * eFrames
array of Frames (images)
Definition: EdbView.h:143
TClonesArray * eTracks
array of Tracks
Definition: EdbView.h:142

Member Function Documentation

◆ AddCluster() [1/2]

EdbCluster * EdbView::AddCluster ( EdbCluster c)
inline
226 {return (EdbCluster*)(new((*eClusters)[eClusters->GetLast()+1]) EdbCluster( *c )); }
Definition: EdbCluster.h:19

◆ AddCluster() [2/2]

EdbCluster * EdbView::AddCluster ( float  x,
float  y,
float  z,
float  a,
float  v,
int  f,
int  s,
int  seg = -1 
)
inline
229 {return (EdbCluster*)(new((*eClusters)[eClusters->GetLast()+1]) EdbCluster(x,y,z, a,v,f,s,seg)); }
FILE * f
Definition: RecDispMC.C:150
void a()
Definition: check_aligned.C:59
s
Definition: check_shower.C:55

◆ AddFrame() [1/2]

void EdbView::AddFrame ( EdbFrame frame)
276{
277 int i = eFrames->GetLast()+1;
278 new((*eFrames)[i++]) EdbFrame( *frame );
279
280 //printf("EdbView::AddFrame: \n");
281 //frame->Print();
282}
Definition: EdbFrame.h:20

◆ AddFrame() [2/2]

void EdbView::AddFrame ( int  id,
float  z,
int  ncl = 0,
int  npix = 0 
)
269{
270 int i = eFrames->GetLast()+1;
271 new((*eFrames)[i++]) EdbFrame( id,z,ncl,npix );
272}

◆ AddSegment() [1/3]

EdbSegment * EdbView::AddSegment ( EdbSegment s)

return directly pointer to the segment

247{
249 int i = eSegments->GetLast()+1;
250 return new((*eSegments)[i++]) EdbSegment( s );
251}
segment of the track
Definition: EdbSegment.h:63

◆ AddSegment() [2/3]

void EdbView::AddSegment ( EdbSegment s)
255{
256 int i = eSegments->GetLast()+1;
257 new((*eSegments)[i++]) EdbSegment( *s );
258}

◆ AddSegment() [3/3]

EdbSegment * EdbView::AddSegment ( float  x,
float  y,
float  z,
float  tx,
float  ty,
float  dz = 0,
int  side = 0,
int  puls = 0,
int  id = -1 
)
inline
233 {return (EdbSegment*)(new((*eSegments)[eSegments->GetLast()+1]) EdbSegment(x,y,z, tx,ty, dz, side, puls, id)); }
brick dz
Definition: RecDispMC.C:107

◆ AddTrack()

void EdbView::AddTrack ( EdbTrack t)
262{
263 int i = eTracks->GetLast()+1;
264 new((*eTracks)[i++]) EdbTrack( *t );
265}
Track linked from segments.
Definition: EdbSegment.h:128
TTree * t
Definition: check_shower.C:4

◆ AttachClustersToSegments()

int EdbView::AttachClustersToSegments ( )
361{
362 int nca=-1;
364 if(nca<0) {
366 printf("AttachClustersToSegments: segments do not sorted - use slow algorithm!\n");
367 }
368 return nca;
369}
int AttachClustersToSegmentsFast()
Definition: EdbView.cxx:372
int AttachClustersToSegmentsSlow()
Definition: EdbView.cxx:409

◆ AttachClustersToSegmentsFast()

int EdbView::AttachClustersToSegmentsFast ( )

assume that cluster->GetSegment() == segment entry :

373{
375
376 EdbCluster *cl=0;
377 int ncl = eClusters->GetLast()+1;
378 EdbSegment *seg=0;
379 int nseg = eSegments->GetLast()+1;
380
381 // TODO: why elements not cleared???
382
383// for(int is=0; is<nseg; is++){
384// seg = (EdbSegment*)(eSegments->UncheckedAt(is));
385// if(seg->GetNelements()>0) seg->GetElements()->Clear();
386// }
387
388 int iseg;
389 int nca=0;
390 for(int i=0; i<ncl; i++){
391 cl = (EdbCluster*)(eClusters->UncheckedAt(i));
392 iseg=cl->GetSegment();
393 if( iseg < 0 ) continue;
394 if( iseg > nseg-1 ) return -1;
395 seg = (EdbSegment*)(eSegments->At(iseg));
396 seg->AddElement(cl);
397 nca++;
398 if(seg->GetNelements()>seg->GetPuls()) {
399 printf("AttachClustersToSegmentsFast: ncl(%d) > puls(%d) %d\n",
400 seg->GetNelements(),seg->GetPuls(), seg->GetID());
401 printf("Area: %d View: %d Segment: %d, i,iseg = %d %d\n",GetAreaID(), GetViewID(), seg->GetID(), i,iseg);
402 }
403
404 }
405 return nca;
406}
Int_t GetSegment() const
Definition: EdbCluster.h:58
Int_t GetPuls() const
Definition: EdbSegment.h:90
void AddElement(TObject *element)
Definition: EdbSegment.cxx:150
Int_t GetNelements() const
Definition: EdbSegment.h:114
Int_t GetID() const
Definition: EdbSegment.h:92
Int_t GetViewID() const
Definition: EdbView.h:190
Int_t GetAreaID() const
Definition: EdbView.h:191

◆ AttachClustersToSegmentsSlow()

int EdbView::AttachClustersToSegmentsSlow ( )

if segments do not ordered :

410{
412
413 EdbCluster *cl=0;
414 int ncl = eClusters->GetLast()+1;
415
416 EdbSegment *seg=0;
417 int nseg = eSegments->GetLast()+1;
418 int id,nca=0;
419 for(int is=0; is<nseg; is++){
420 seg = (EdbSegment*)(eSegments->UncheckedAt(is));
421 if(seg->GetNelements()>0) seg->GetElements()->Clear();
422 id = seg->GetID();
423 for(int i=0; i<ncl; i++){
424 cl = (EdbCluster*)(eClusters->UncheckedAt(i));
425 if( cl->GetSegment() != id ) continue;
426 seg->AddElement(cl);
427 nca++;
428 }
429 }
430 return nca;
431}
TObjArray * GetElements() const
Definition: EdbSegment.h:117
UInt_t id
Definition: tlg2couples.C:117

◆ AttachSegmentsToTracks()

int EdbView::AttachSegmentsToTracks ( )

Assumed that:
tracks.eID - number of segments in a track
Segments are written consecutively in the order they appear in microtracks (i.e. m0c0,m0c1,m0c2,m1c0,m1c1,m2c0,m2c1...)

333{
337
338 int nseg = eSegments->GetLast()+1; if(nseg<1) return 0;
339 int ntr = eTracks->GetLast()+1; if(ntr<1) return 0;
340
341 int icnt=0;
342 for(int itr=0; itr<ntr; itr++)
343 {
344 EdbTrack *t = (EdbTrack*)(eTracks->UncheckedAt(itr));
345 if(!t) continue;
346 int nst = t->GetID(); if(nst<1) continue;
347 if( t->GetElements() ) t->GetElements()->Clear();
348 for(int ist=0; ist<nst; ist++)
349 {
350 if(icnt >= nseg) break;
351 EdbSegment *s = (EdbSegment*)(eSegments->UncheckedAt(icnt++));
352 if(s) t->AddElement(s);
353 }
354 }
355 Log(3,"EdbView::AttachSegmentsToTracks", "ntr = %d nseg = %d icnt = %d",ntr,nseg,icnt);
356 return icnt;
357}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75

◆ Clear()

void EdbView::Clear ( )
80{
81 if(eClusters) { eClusters->Clear(); }
82 if(eSegments) {
83 TObjArray *e=0;
84 for(int i=0; i<Nsegments(); i++) {
85 e=((EdbSegment*)eSegments->At(i))->GetElements();
86 if( e ) e->Clear();
87 }
88 eSegments->Clear();
89 }
90 if(eTracks) { eTracks->Clear(); }
91 if(eFrames) { eFrames->Clear(); }
92}
Int_t Nsegments() const
Definition: EdbView.h:216

◆ DeleteClustersFog()

void EdbView::DeleteClustersFog ( )
209{
210 EdbCluster *cluster =0;
211
212 int ncl = Nclusters();
213 for( int icl=0; icl<ncl; icl++) {
214 cluster = (EdbCluster*)(eClusters->UncheckedAt(icl));
215 if( cluster->GetSegment() < 0 ) eClusters->Remove(cluster);
216 }
217
218 eClusters->Compress();
219}
Int_t Nclusters() const
Definition: EdbView.h:215

◆ Draw()

void EdbView::Draw ( Option_t *  option = "")
467{
468 float xmin=-200, xmax=200;
469 float ymin=-150, ymax=150;
470 //float xmin=GetXview()-500, xmax=xmin+1000.;
471 //float ymin=GetYview()-500, ymax=ymin+1000.;
472
473 float zmin = TMath::Min( GetZ1(), GetZ4() );
474 float zmax = TMath::Max( GetZ1(), GetZ4() );
475
476 TH3F *hist = new TH3F("hist","title",1,xmin,xmax,1,ymin,ymax,1,zmin,zmax);
477 //TH3F *hist = new TH3F("hist","title",1,0,512,1,0,512,1,zmin,zmax);
478
479 hist->Draw();
480
481 if(eClusters) {
482 TPolyMarker3D *marksFog = DrawClustersFog();
483 marksFog->SetMarkerColor(kBlue);
484 marksFog->SetMarkerStyle(21);
485 marksFog->SetMarkerSize(.3);
486 marksFog->Draw();
487
488 TPolyMarker3D *marksSeg = DrawClustersSegments();
489 marksSeg->SetMarkerColor(kRed);
490 marksSeg->SetMarkerStyle(21);
491 marksSeg->SetMarkerSize(.3);
492 marksSeg->Draw();
493 }
494
495 if(eTracks) {
496 // Draw tracks ....
497
498 EdbTrack *tr = 0;
499 int ntr = Ntracks();
500
501 for(int itr=0; itr<ntr; itr++) {
502 tr = (EdbTrack*)eTracks->At(itr);
503
504 TPolyLine3D *line = new TPolyLine3D(2);
505 line->SetPoint(0, tr->GetX0(), tr->GetY0(), tr->GetZ0() );
506 line->SetPoint(1,
507 tr->GetX0() + tr->GetTx()*tr->GetDz(),
508 tr->GetY0() + tr->GetTy()*tr->GetDz(),
509 tr->GetZ0()+ tr->GetDz() );
510
511 line->SetLineColor(kGreen);
512 line->SetLineWidth(10);
513 line->Draw();
514 }
515 }
516
517 if(eSegments) {
518 // Draw segments ....
519
520 EdbSegment *seg = 0;
521 int nseg = Nsegments();
522
523 for(int iseg=0; iseg<nseg; iseg++) {
524 seg = (EdbSegment*)eSegments->At(iseg);
525
526 TPolyLine3D *line = new TPolyLine3D(2);
527 line->SetPoint(0, seg->GetX0(), seg->GetY0(), seg->GetZ0() );
528 line->SetPoint(1,
529 seg->GetX0() + seg->GetTx()*seg->GetDz(),
530 seg->GetY0() + seg->GetTy()*seg->GetDz(),
531 seg->GetZ0()+ seg->GetDz() );
532 line->Draw();
533 }
534 }
535
536}
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
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
TPolyMarker3D * DrawClustersSegments(Option_t *opt=0) const
Definition: EdbView.cxx:316
Int_t Ntracks() const
Definition: EdbView.h:217
TPolyMarker3D * DrawClustersFog(Option_t *opt=0) const
Definition: EdbView.cxx:300
Float_t GetZ1() const
Definition: EdbView.h:196
Float_t GetZ4() const
Definition: EdbView.h:199
void hist()
Definition: init.C:23

◆ DrawClustersFog()

TPolyMarker3D * EdbView::DrawClustersFog ( Option_t *  opt = 0) const
301{
302 EdbCluster *cluster=0;
303 int ncl = eClusters->GetLast()+1;
304
305 TPolyMarker3D *marks = new TPolyMarker3D(ncl,1);
306
307 for(int i=0; i<ncl; i++){
308 cluster = (EdbCluster*)(eClusters->At(i));
309 if( cluster->GetSegment() < 0 )
310 marks->SetNextPoint( cluster->GetX(), cluster->GetY(), cluster->GetZ() );
311 }
312 return marks;
313}
Float_t GetX() const
Definition: EdbCluster.h:51
Float_t GetY() const
Definition: EdbCluster.h:52
Float_t GetZ() const
Definition: EdbCluster.h:53

◆ DrawClustersSegments()

TPolyMarker3D * EdbView::DrawClustersSegments ( Option_t *  opt = 0) const
317{
318 EdbCluster *cluster=0;
319 int ncl = eClusters->GetLast()+1;
320
321 TPolyMarker3D *marks = new TPolyMarker3D(ncl,1);
322
323 for(int i=0; i<ncl; i++){
324 cluster = (EdbCluster*)(eClusters->At(i));
325 if( cluster->GetSegment() >= 0 )
326 marks->SetNextPoint( cluster->GetX(), cluster->GetY(), cluster->GetZ() );
327 }
328 return marks;
329}

◆ GenerateClustersFog()

void EdbView::GenerateClustersFog ( float  density)
617{
618 float x1=0, x2=512;
619 float y1=0, y2=512;
620
621 float x=0,y=0,z=0;
622 int area=0, volume=0, frame=0, side=0;
623
624 float z1 = GetZ1();
625 float z2 = GetZ2();
626 float z3 = GetZ3();
627 float z4 = GetZ4();
628
629 float activeVolume = TMath::Abs((x2-x1)*(y2-y1)*( z2-z1 + z4-z3 ));
630 int nClusters = (int)( density*activeVolume/(10*10*10) );
631
632 int ft = GetNframesTop();
633 int fb = GetNframesBot();
634 int nf = ft+fb;
635
636 for( int i=0; i<nClusters; i++ ){
637
638 frame = gRandom->Integer( nf );
639 if( frame < ft ) {
640 side = 0;
641 z = z1 + (z2-z1)/ft * frame;
642 }
643 else {
644 side = 1;
645 z = z3 + (z4-z3)/fb * (frame-ft);
646 }
647
648 x = x1 + gRandom->Rndm()*(x2-x1);
649 y = y1 + gRandom->Rndm()*(y2-y1);
650 area = (int)gRandom->Exp(5);
651 volume = gRandom->Poisson(50);
652
653 AddCluster(x,y,z, area,volume, frame, side);
654 }
655}
Int_t GetNframesTop() const
Definition: EdbView.h:207
Float_t GetZ3() const
Definition: EdbView.h:198
Int_t GetNframesBot() const
Definition: EdbView.h:208
EdbCluster * AddCluster(EdbCluster *c)
Definition: EdbView.h:225
Float_t GetZ2() const
Definition: EdbView.h:197

◆ GenerateClustersSegment()

void EdbView::GenerateClustersSegment ( EdbSegment segment,
int  n0 = 25,
float  sigma = .1 
)

n0 - mean number of grains/100 microns
sigma - XY dispersion in microns

568{
571
572 float xmin = 0;
573 float ymin = 0;
574 float xmax = 512;
575 float ymax = 512;
576
577 float length0 = 100./n0; // mean flight length [microns]
578
579 float x0 = segment->GetX0();
580 float y0 = segment->GetY0();
581 float z0 = segment->GetZ0();
582 float dz = segment->GetDz();
583 float x=0,y=0,z=0;
584 float path=0;
585 int area=0, volume=0, frame=-1;
586 int side=segment->GetSide();
587 int seg=segment->GetID();
588
589 float tx = segment->GetTx();
590 float ty = segment->GetTy();
591
592 length0 = length0/sqrt( 1 + tx*tx + ty*ty );
593
594 int ncl=0;
595
596 while( 1 ) {
597 path += gRandom->Exp(length0);
598 if( path>dz) break;
599 z = z0+path;
600 x = x0 + tx*path + gRandom->Gaus(0,sigma);
601 y = y0 + ty*path + gRandom->Gaus(0,sigma);
602
603 if( x < xmin ) continue;
604 if( x > xmax ) continue;
605 if( y < ymin ) continue;
606 if( y > ymax ) continue;
607
608 volume = gRandom->Poisson(50);
609 AddCluster(x,y,z, area, volume, frame, side, seg );
610 ncl++;
611 }
612 segment->SetPuls(ncl);
613}
brick z0
Definition: RecDispMC.C:106
void SetPuls(int puls)
Definition: EdbSegment.h:95
Int_t GetSide() const
Definition: EdbSegment.h:89

◆ GenerateClustersTrack()

void EdbView::GenerateClustersTrack ( EdbTrack track,
int  n0 = 25,
float  sigma = .1 
)
540{
541 EdbSegment segment;
542
543 int count = Nsegments();
544 float x0 = track->GetX0();
545 float y0 = track->GetY0();
546 float z1 = GetZ1();
547 float z2 = GetZ2();
548 float z3 = GetZ3();
549 float z4 = GetZ4();
550 float tx = track->GetTx();
551 float ty = track->GetTy();
552
553 segment.Set( x0,y0, z1, tx, ty, z2-z1, 0, -1, count );
554 AddSegment(&segment);
556
557 x0 = x0 + (z3-z1)*tx;
558 y0 = y0 + (z3-z1)*ty;
559 count = Nsegments();
560
561 segment.Set( x0,y0, z3, tx, ty, z4-z3, 1, -1, count );
562 AddSegment(&segment);
564}
void Set(float x, float y, float z, float tx, float ty, float dz=0, int side=0, int puls=0, int id=0)
Definition: EdbSegment.cxx:93
EdbSegment * GetSegment(int i) const
Definition: EdbView.h:219
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
void GenerateClustersSegment(EdbSegment *segment, int n0=25, float sigma=.1)
Definition: EdbView.cxx:567
Definition: bitview.h:14

◆ GenerateFrames()

void EdbView::GenerateFrames ( int  n = 32)
659{
660 if(eFrames) eFrames->Clear();
661
662 EdbFrame *frame=0;
663
664 for( int ifr=0; ifr<n; ifr++) {
665 frame = new EdbFrame( ifr,2,3,"123456" );
666 AddFrame(frame);
667 }
668}
void AddFrame(int id, float z, int ncl=0, int npix=0)
Definition: EdbView.cxx:268

◆ GetAreaID()

Int_t EdbView::GetAreaID ( ) const
inline
191{ if(eHeader) return GetHeader()->GetAreaID(); else return 0;}
Int_t GetAreaID() const
Definition: EdbView.h:91
EdbViewHeader * GetHeader() const
Definition: EdbView.h:163

◆ GetCluster()

EdbCluster * EdbView::GetCluster ( int  i) const
inline
218{ return (EdbCluster*)eClusters->At(i); }

◆ GetClusters()

TClonesArray * EdbView::GetClusters ( ) const
inline
164{ return eClusters; }

◆ GetClustersAddr()

void * EdbView::GetClustersAddr ( )
inline
170{ return &eClusters; }

◆ GetClustersFrame()

TList * EdbView::GetClustersFrame ( int  frame) const
286{
287 TList *list = new TList();
288
289 EdbCluster *cluster=0;
290 int ncl = eClusters->GetLast()+1;
291
292 for(int i=0; i<ncl; i++){
293 cluster = (EdbCluster*)(eClusters->UncheckedAt(i));
294 if( cluster->GetFrame() == frame ) list->Add(cluster);
295 }
296 return list;
297}
Int_t GetFrame() const
Definition: EdbCluster.h:56

◆ GetDZbot()

Float_t EdbView::GetDZbot ( ) const
inline
204{ if(eHeader) return GetHeader()->GetDZbot(); else return 0; }
Float_t GetDZbot() const
Definition: EdbView.h:118

◆ GetDZtop()

Float_t EdbView::GetDZtop ( ) const
inline
203{ if(eHeader) return GetHeader()->GetDZtop(); else return 0; }
Float_t GetDZtop() const
Definition: EdbView.h:117

◆ GetFrame()

EdbFrame * EdbView::GetFrame ( int  i) const
inline
221{ return (EdbFrame*)eFrames->At(i); }

◆ GetFrames()

TClonesArray * EdbView::GetFrames ( ) const
inline
167{ return eFrames; }

◆ GetFramesAddr()

void * EdbView::GetFramesAddr ( )
inline
173{ return &eFrames; }

◆ GetHeader()

EdbViewHeader * EdbView::GetHeader ( ) const
inline
163{ return eHeader; }

◆ GetHeaderAddr()

void * EdbView::GetHeaderAddr ( )
inline
169{ return &eHeader; }

◆ GetLastSystemTime()

Long_t EdbView::GetLastSystemTime ( ) const
inline
211{ return eLastSystemTime; }
Long_t eLastSystemTime
system time when view was saved
Definition: EdbView.h:145

◆ GetNframes()

Int_t EdbView::GetNframes ( ) const
inline
206{ return GetNframesTop()+GetNframesBot(); }

◆ GetNframesBot()

Int_t EdbView::GetNframesBot ( ) const
inline
208{ if(eHeader) return GetHeader()->GetNframesBot(); else return 0; }
Int_t GetNframesBot() const
Definition: EdbView.h:121

◆ GetNframesTop()

Int_t EdbView::GetNframesTop ( ) const
inline
207{ if(eHeader) return GetHeader()->GetNframesTop(); else return 0; }
Int_t GetNframesTop() const
Definition: EdbView.h:120

◆ GetSegment()

EdbSegment * EdbView::GetSegment ( int  i) const
inline
219{ return (EdbSegment*)eSegments->At(i); }

◆ GetSegments()

TClonesArray * EdbView::GetSegments ( ) const
inline
165{ return eSegments; }

◆ GetSegmentsAddr()

void * EdbView::GetSegmentsAddr ( )
inline
171{ return &eSegments; }

◆ GetSegmentsClusters()

TObjArray * EdbView::GetSegmentsClusters ( ) const
435{
436 int nseg = eSegments->GetLast()+1;
437 TObjArray *segments = new TObjArray(nseg);
438
439 EdbSegment *seg1=0, *seg2=0;
440
441 for(int iseg=0; iseg<nseg; iseg++) {
442 seg1 = (EdbSegment*)eSegments->At(iseg);
443 seg2 = new EdbSegment(*seg1);
444 segments->Add(seg2);
445 }
446
447 EdbCluster *cl1=0, *cl2=0;
448
449 int ncl = eClusters->GetLast()+1;
450 int segID=-1;
451
452 for(int i=0; i<ncl; i++) {
453 cl1 = (EdbCluster*)eClusters->At(i);
454 segID = cl1->GetSegment();
455 if( segID < 0) continue;
456 if( segID > nseg) continue;
457 cl2 = new EdbCluster(*cl1);
458 seg1 = (EdbSegment*)segments->At(segID);
459 seg1->AddElement(cl2);
460 }
461
462 return segments;
463}

◆ GetTime()

Int_t EdbView::GetTime ( ) const
inline
210{ if(eHeader) return GetHeader()->GetTime(); else return 0; }
Int_t GetTime()
Definition: EdbView.h:82

◆ GetTrack()

EdbTrack * EdbView::GetTrack ( int  i) const
inline
220{ return (EdbTrack*)eTracks->At(i); }

◆ GetTracks()

TClonesArray * EdbView::GetTracks ( ) const
inline
166{ return eTracks; }

◆ GetTracksAddr()

void * EdbView::GetTracksAddr ( )
inline
172{ return &eTracks; }

◆ GetViewID()

Int_t EdbView::GetViewID ( ) const
inline
190{ if(eHeader) return GetHeader()->GetViewID(); else return 0;}
Int_t GetViewID() const
Definition: EdbView.h:90

◆ GetXview()

Float_t EdbView::GetXview ( ) const
inline
193{ if(eHeader) return GetHeader()->GetXview(); else return 0; }
Float_t GetXview() const
Definition: EdbView.h:93

◆ GetYview()

Float_t EdbView::GetYview ( ) const
inline
194{ if(eHeader) return GetHeader()->GetYview(); else return 0; }
Float_t GetYview() const
Definition: EdbView.h:94

◆ GetZ0bot()

Float_t EdbView::GetZ0bot ( ) const
inline
202{ if(eHeader) return GetHeader()->GetZ0bot(); else return 0; }
Float_t GetZ0bot() const
Definition: EdbView.h:116

◆ GetZ0top()

Float_t EdbView::GetZ0top ( ) const
inline
201{ if(eHeader) return GetHeader()->GetZ0top(); else return 0; }
Float_t GetZ0top() const
Definition: EdbView.h:115

◆ GetZ1()

Float_t EdbView::GetZ1 ( void  ) const
inline
196{ if(eHeader) return GetHeader()->GetZ1(); else return 0; }
Float_t GetZ1() const
Definition: EdbView.h:110

◆ GetZ2()

Float_t EdbView::GetZ2 ( void  ) const
inline
197{ if(eHeader) return GetHeader()->GetZ2(); else return 0; }
Float_t GetZ2() const
Definition: EdbView.h:111

◆ GetZ3()

Float_t EdbView::GetZ3 ( void  ) const
inline
198{ if(eHeader) return GetHeader()->GetZ3(); else return 0; }
Float_t GetZ3() const
Definition: EdbView.h:112

◆ GetZ4()

Float_t EdbView::GetZ4 ( void  ) const
inline
199{ if(eHeader) return GetHeader()->GetZ4(); else return 0; }
Float_t GetZ4() const
Definition: EdbView.h:113

◆ Nclusters()

Int_t EdbView::Nclusters ( ) const
inline
215{ return eClusters->GetLast()+1; }

◆ Nsegments()

Int_t EdbView::Nsegments ( ) const
inline
216{ return eSegments->GetLast()+1; }

◆ Ntracks()

Int_t EdbView::Ntracks ( ) const
inline
217{ return eTracks->GetLast()+1; }

◆ Print()

void EdbView::Print ( Option_t *  opt = 0) const
186{
187 printf("\nEdbView: ");
188 if(eClusters) printf( "\t %d clusters\n", eClusters->GetLast()+1);
189 if(eSegments) printf("\t\t %d segments\n", eSegments->GetLast()+1);
190 if(eTracks) printf("\t\t %d tracks \n", eTracks->GetLast()+1);
191 if(eFrames) printf("\t\t %d frames \n", eFrames->GetLast()+1);
192 eHeader->Print();
193}
void Print() const
Definition: EdbView.cxx:749

◆ PrintClusters()

void EdbView::PrintClusters ( Option_t *  opt = 0) const
197{
198 EdbCluster *cluster =0;
199
200 int ncl = eClusters->GetLast()+1;
201 for( int icl=0; icl<ncl; icl++) {
202 cluster = (EdbCluster*)(eClusters->UncheckedAt(icl));
203 cluster->Print();
204 }
205}
void Print(Option_t *opt=0) const
Definition: EdbCluster.cxx:58

◆ ReadView()

int EdbView::ReadView ( char *  fname)
178{
179 FILE *fp = fopen(fname,"r");
180 fclose(fp);
181 return 0;
182}
const char * fname
Definition: mc2raw.cxx:41
fclose(pFile)

◆ Scale()

void EdbView::Scale ( float  zscale)

scale all z coordinates and angles to the given factor

114{
116 EdbCluster *c=0;
117 EdbSegment *s=0;
118 EdbTrack *t=0;
119 if(eClusters) {
120 int ncl=Nclusters();
121 for( int i=0; i<ncl; i++ ) {
122 c = GetCluster(i);
123 c->SetZ( c->GetZ()*zscale );
124 }
125 }
126 if(eSegments) {
127 int nseg=Nsegments();
128 for( int i=0; i<nseg; i++ ) {
129 s = GetSegment(i);
130 s->SetZ0( s->GetZ0()*zscale );
131 s->SetTx( s->GetTx()/zscale );
132 s->SetTy( s->GetTy()/zscale );
133 }
134 }
135 if(eTracks) {
136 int ntr=Ntracks();
137 for( int i=0; i<ntr; i++ ) {
138 t = GetTrack(i);
139 t->SetZ0( t->GetZ0()*zscale );
140 t->SetTx( t->GetTx()/zscale );
141 t->SetTy( t->GetTy()/zscale );
142 }
143 }
144}
virtual void SetZ(float z)
Definition: EdbCluster.h:83
EdbTrack * GetTrack(int i) const
Definition: EdbView.h:220
EdbCluster * GetCluster(int i) const
Definition: EdbView.h:218

◆ SetAreaID()

void EdbView::SetAreaID ( int  id)
inline
183{ GetHeader()->SetAreaID(id); }
void SetAreaID(int id)
Definition: EdbView.h:100

◆ SetCoordXY()

void EdbView::SetCoordXY ( float  x,
float  y 
)
inline
185{ GetHeader()->SetCoordXY(x,y); }
void SetCoordXY(float x, float y)
Definition: EdbView.h:101

◆ SetCoordZ()

void EdbView::SetCoordZ ( float  z1,
float  z2,
float  z3,
float  z4 
)
inline
187 { GetHeader()->SetCoordZ(z1,z2,z3,z4); }
void SetCoordZ(float z1, float z2, float z3, float z4)
Definition: EdbView.h:102

◆ SetLastSystemTime()

void EdbView::SetLastSystemTime ( Long_t  time)
inline
212{ eLastSystemTime=time; }

◆ SetNframes()

void EdbView::SetNframes ( int  top,
int  bot 
)
inline
184{ GetHeader()->SetNframes(top,bot); }
void SetNframes(int top, int bot)
Definition: EdbView.h:104

◆ Shift()

void EdbView::Shift ( float  zshift)

shift all z coordinates to the given offset

148{
150 EdbCluster *c=0;
151 EdbSegment *s=0;
152 EdbTrack *t=0;
153 if(eClusters) {
154 int ncl=Nclusters();
155 for( int i=0; i<ncl; i++ ) {
156 c = GetCluster(i);
157 c->SetZ( c->GetZ()+zshift );
158 }
159 }
160 if(eSegments) {
161 int nseg=Nsegments();
162 for( int i=0; i<nseg; i++ ) {
163 s = GetSegment(i);
164 s->SetZ0( s->GetZ0()+zshift );
165 }
166 }
167 if(eTracks) {
168 int ntr=Ntracks();
169 for( int i=0; i<ntr; i++ ) {
170 t = GetTrack(i);
171 t->SetZ0( t->GetZ0()+zshift );
172 }
173 }
174}

◆ Transform()

void EdbView::Transform ( const EdbAffine2D aff)
97{
98 if(eClusters) {
99 int ncl=Nclusters();
100 for( int i=0; i<ncl; i++ ) ((EdbPoint2D*)(GetCluster(i)))->Transform(aff);
101 }
102 if(eSegments) {
103 int nseg=Nsegments();
104 for( int i=0; i<nseg; i++ ) GetSegment(i)->Transform(aff);
105 }
106 if(eTracks) {
107 int ntr=Ntracks();
108 for( int i=0; i<ntr; i++ ) GetTrack(i)->Transform(aff);
109 }
110}
virtual 2D point
Definition: EdbVirtual.h:76
virtual void Transform(const EdbAffine2D *aff)
Definition: EdbSegment.cxx:59

◆ Xmin()

float EdbView::Xmin ( ) const
inline

View coordinates are in pixels.

◆ Ymin()

float EdbView::Ymin ( ) const
inline

starting from 0

◆ ZFrameMax()

float EdbView::ZFrameMax ( ) const
235{
236 float zmax=GetFrame(0)->GetZ();
237 for(int i=1; i<GetNframes(); i++)
238 {
239 float z = GetFrame(i)->GetZ();
240 if(z>zmax) zmax=z;
241 }
242 return zmax;
243}
float GetZ() const
Definition: EdbFrame.h:42
EdbFrame * GetFrame(int i) const
Definition: EdbView.h:221
Int_t GetNframes() const
Definition: EdbView.h:206

◆ ZFrameMin()

float EdbView::ZFrameMin ( ) const
223{
224 float zmin=GetFrame(0)->GetZ();
225 for(int i=1; i<GetNframes(); i++)
226 {
227 float z = GetFrame(i)->GetZ();
228 if(z<zmin) zmin=z;
229 }
230 return zmin;
231}

◆ Zmax()

float EdbView::Zmax ( ) const
inline
256{return TMath::Max( GetZ1(), GetZ4() );}

◆ Zmin()

float EdbView::Zmin ( ) const
inline
255{return TMath::Min( GetZ1(), GetZ4() );}

Member Data Documentation

◆ eClusters

TClonesArray* EdbView::eClusters
private

array of Clusters

◆ eFrames

TClonesArray* EdbView::eFrames
private

array of Frames (images)

◆ eHeader

EdbViewHeader* EdbView::eHeader
private

View header.

◆ eLastSystemTime

Long_t EdbView::eLastSystemTime
private

system time when view was saved

◆ eSegments

TClonesArray* EdbView::eSegments
private

array of Segments

◆ eTracks

TClonesArray* EdbView::eTracks
private

array of Tracks


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