FEDRA emulsion software from the OPERA Collaboration
EdbRun Class Reference

#include <EdbRun.h>

Inheritance diagram for EdbRun:
Collaboration diagram for EdbRun:

Public Member Functions

int AddAsciiFile (const char *fname, const char *objname)
 
void AddFrameAlign (Int_t frame1, Int_t frame2, Int_t view, Int_t area, Int_t side, Float_t dxglobal, Float_t dyglobal, Float_t dx, Float_t dy, Float_t z1, Float_t z2, Int_t n1tot, Int_t n2tot, Int_t n1, Int_t n2, Int_t nsg, Int_t nbg, Int_t flag)
 
void AddPinViewHeader (EdbViewHeader &h)
 
void AddView ()
 
void AddView (EdbView *view)
 
void AddViewAlign (Int_t view1, Int_t view2, Int_t area1, Int_t area2, Int_t side1, Int_t side2, Float_t dx, Float_t dy, Float_t dz, Int_t n1tot, Int_t n2tot, Int_t n1, Int_t n2, Int_t nsg, Int_t nbg, Int_t flag, EdbViewHeader &vh1, EdbViewHeader &vh2)
 
void AddViewMerge (Int_t view1, Int_t view2, Int_t area1, Int_t area2, Int_t side1, Int_t side2, Float_t dx, Float_t dy, Float_t dz, Int_t n1tot, Int_t n2tot, Int_t n1, Int_t n2, Int_t nsg, Int_t nbg, Int_t flag)
 
void Close ()
 
void Create (const char *fname)
 
 EdbRun ()
 
 EdbRun (const char *fname, const char *status="READ")
 
 EdbRun (EdbRun &run, const char *fname)
 
 EdbRun (int id, const char *status="READ", const char *path=".")
 
int ExtractAsciiFile (const char *fname, const char *objname)
 
void GeneratePredictions (int n)
 
int GetEntries () const
 
EdbViewGetEntry (int entry, int ih=1, int icl=0, int iseg=1, int itr=0, int ifr=0)
 
TClonesArray * GetEntryClusters (int entry) const
 
TObjArray * GetEntryFrames (int entry) const
 
EdbViewHeaderGetEntryHeader (int entry) const
 
int GetEntryNumberWithIndex (int event, int view) const
 
TClonesArray * GetEntrySegments (int entry) const
 
TClonesArray * GetEntryTracks (int entry) const
 
TFile * GetFile ()
 
TDatime * GetFinishTime () const
 
EdbRunHeaderGetHeader () const
 
EdbMarksSetGetMarks () const
 
TTree * GetPinViews () const
 
EdbPredictionDCGetPrediction (int ip)
 
EdbPredictionsBoxGetPredictions () const
 
int GetRunID () const
 
TDatime * GetStartTime () const
 
TTree * GetTree () const
 
EdbViewGetView () const
 
void Init ()
 
int Npredictions () const
 
void Open (const char *fname)
 
void OpenUpdate (const char *fname)
 
void Print () const
 
void PrintBranchesStatus () const
 
void PrintLog (const char *fname) const
 
void Save ()
 
void SaveViews ()
 
void SelectOpenMode (const char *fname, const char *status="READ")
 
void SetComment (const char *com)
 
void SetMarks (EdbMarksSet *marks)
 
void SetRunID (int id)
 
void SetTitle (const char *tit)
 
void SetView ()
 
void SetView (EdbView *view)
 
void TransformDC ()
 
virtual ~EdbRun ()
 

Public Attributes

EdbMarksSeteMarks
 fiducial marks More...
 

Private Attributes

AlignmentParFrame eFA
 
TFile * eFile
 file associated with the Run More...
 
TTree * eFrameAlign
 frames alignment More...
 
EdbRunHeadereHeader
 run header More...
 
TString ePath
 runs directory path More...
 
TTree * ePinViews
 pinned views More...
 
EdbPredictionsBoxePredictions
 predictions to scan ($c) More...
 
EdbViewHeaderePVH
 to add pinned view header More...
 
TTree * eTree
 tree with View objects More...
 
AlignmentParView eVA
 
EdbViewHeadereVH1
 alignment headers More...
 
EdbViewHeadereVH2
 alignment headers More...
 
EdbVieweView
 view using for import and export data More...
 
TTree * eViewAlign
 view neighbours alignment More...
 
TTree * eViewMerge
 view merging alignment More...
 
AlignmentParView eVM
 

Constructor & Destructor Documentation

◆ EdbRun() [1/4]

EdbRun::EdbRun ( )

◆ EdbRun() [2/4]

EdbRun::EdbRun ( int  id,
const char *  status = "READ",
const char *  path = "." 
)
57{
58 Init();
59
60 eHeader = new EdbRunHeader(id);
61
62 ePath = path;
63
64 char fname[256]="";
65 sprintf(fname,"%s/run%4.4d.root",path,id);
66
67 SelectOpenMode( fname, status );
68}
Definition: EdbRunHeader.h:95
void SelectOpenMode(const char *fname, const char *status="READ")
Definition: EdbRun.cxx:164
TString ePath
runs directory path
Definition: EdbRun.h:85
void Init()
Definition: EdbRun.cxx:145
EdbRunHeader * eHeader
run header
Definition: EdbRun.h:79
const char * fname
Definition: mc2raw.cxx:41

◆ EdbRun() [3/4]

EdbRun::EdbRun ( const char *  fname,
const char *  status = "READ" 
)
72{
73 Init();
74
75 eHeader = new EdbRunHeader(0);
76
77 SelectOpenMode( fname, status );
78}

◆ EdbRun() [4/4]

EdbRun::EdbRun ( EdbRun run,
const char *  fname 
)
82{
83 Init();
84
85 if(run.GetHeader())
87 if(run.GetMarks())
88 eMarks = new EdbMarksSet(*(run.GetMarks()));
89 if(run.GetPredictions())
91 Log(3,"EdbRun::EdbRun","Run headers are copied");
92 SelectOpenMode( fname, "RECREATE" );
93}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
Definition: EdbFiducial.h:88
Definition: EdbPrediction.h:71
EdbMarksSet * eMarks
fiducial marks
Definition: EdbRun.h:101
EdbPredictionsBox * GetPredictions() const
Definition: EdbRun.h:119
EdbPredictionsBox * ePredictions
predictions to scan ($c)
Definition: EdbRun.h:87
EdbMarksSet * GetMarks() const
Definition: EdbRun.h:120
EdbRunHeader * GetHeader() const
Definition: EdbRun.h:138
EdbRun * run
Definition: check_raw.C:38

◆ ~EdbRun()

EdbRun::~EdbRun ( )
virtual
97{
98 //SafeDelete(eTree);
99
100 SafeDelete(eFile);
101 SafeDelete(eHeader);
102 SafeDelete(ePredictions);
103 SafeDelete(eMarks);
104
105 //SafeDelete(eView); // provocate crash in EdbRunAccess::CopyRawDataXY - to understand!
106/*
107 if(ePVH)
108 SafeDelete(ePVH);
109 if(eVH1)
110 SafeDelete(eVH1);
111 if(eVH2)
112 SafeDelete(eVH2);
113*/
114}
TFile * eFile
file associated with the Run
Definition: EdbRun.h:84

Member Function Documentation

◆ AddAsciiFile()

int EdbRun::AddAsciiFile ( const char *  fname,
const char *  objname 
)

can be used for reading relativly small ascii files from fname and save into the current run with name objname as TObjString

593{
596
597 FILE *fp=fopen(fname,"r");
598 if (fp==NULL) {
599 printf("ERROR open file: %s n", fname); return(-1);
600 }else
601 printf( "Read Ascii file: %s\n", fname );
602 TObjString str;
603 char buf[256];
604 while( fgets(buf,256,fp)!=NULL ) str.String().Append(buf);
605 fclose(fp);
606 return str.Write(objname);
607}
fclose(pFile)
#define NULL
Definition: nidaqmx.h:84

◆ AddFrameAlign()

void EdbRun::AddFrameAlign ( Int_t  frame1,
Int_t  frame2,
Int_t  view,
Int_t  area,
Int_t  side,
Float_t  dxglobal,
Float_t  dyglobal,
Float_t  dx,
Float_t  dy,
Float_t  z1,
Float_t  z2,
Int_t  n1tot,
Int_t  n2tot,
Int_t  n1,
Int_t  n2,
Int_t  nsg,
Int_t  nbg,
Int_t  flag 
)
388{
389 if(!eFrameAlign) {
390 eFrameAlign = new TTree("FrameAlign","Neighbour frames offsets");
391 eFrameAlign->Branch("frpar",&eFA,"frame1/I:frame2/I:view/I:area/I:side/I:dxglobal/F:dyglobal/F:dx/F:dy/F:z1/F:z2/F:n1tot/I:n2tot/I:n1/I:n2/I:nsg/I:nbg/I:flag/I");
392 }
393 if(eFrameAlign) {
394 eFA.frame1=frame1; eFA.frame2=frame2;
395 eFA.view=view;
396 eFA.area=area;
397 eFA.side=side;
398 eFA.dxglobal=dxglobal;
399 eFA.dyglobal=dyglobal;
400 eFA.dx=dx; eFA.dy=dy;
401 eFA.z1=z1; eFA.z2=z2;
402 eFA.n1tot=n1tot; eFA.n2tot=n2tot;
403 eFA.n1=n1; eFA.n2=n2;
404 eFA.nsg=nsg; eFA.nbg=nbg;
405 eFA.flag=flag;
406 eFrameAlign->Fill();
407 }
408}
AlignmentParFrame eFA
Definition: EdbRun.h:95
TTree * eFrameAlign
frames alignment
Definition: EdbRun.h:91
Int_t frame1
1-st frame id
Definition: EdbRun.h:54
Float_t dyglobal
global offset applied for all frames of this view
Definition: EdbRun.h:60
Int_t view
view id
Definition: EdbRun.h:56
Float_t dx
found offset
Definition: EdbRun.h:61
Int_t n2tot
Definition: EdbRun.h:66
Int_t nsg
peak height
Definition: EdbRun.h:69
Float_t z2
z of the second frame
Definition: EdbRun.h:64
Int_t n1
number of grains in intersection area
Definition: EdbRun.h:67
Float_t dy
found offset
Definition: EdbRun.h:62
Float_t dxglobal
global offset applied for all frames of this view (XfrStage=Xframecorrected-dxglobal-flag*dx)
Definition: EdbRun.h:59
Int_t flag
in first 3 bits: {f_Applied = 0x01,f_Found = 0x02,f_Recovered = 0x04};
Definition: EdbRun.h:71
Int_t nbg
bg level (random coinsidences)
Definition: EdbRun.h:70
Int_t side
Definition: EdbRun.h:58
Int_t frame2
2-d frame id
Definition: EdbRun.h:55
Float_t z1
z of the first frame
Definition: EdbRun.h:63
Int_t n2
Definition: EdbRun.h:68
Int_t n1tot
total number of clusters in the frame
Definition: EdbRun.h:65
Int_t area
Definition: EdbRun.h:57

◆ AddPinViewHeader()

void EdbRun::AddPinViewHeader ( EdbViewHeader h)
412{
413 if(!ePinViews) {
414 ePinViews = new TTree("PinViews","Pinned Views");
415 ePinViews->Branch("headers","EdbViewHeader",&ePVH);
416 }
417 if(ePinViews) {
418 *ePVH=h;
419 ePinViews->Fill();
420 }
421}
EdbViewHeader * ePVH
to add pinned view header
Definition: EdbRun.h:96
TTree * ePinViews
pinned views
Definition: EdbRun.h:92

◆ AddView() [1/2]

void EdbRun::AddView ( )
306{
307 AddView(eView);
308}
EdbView * eView
view using for import and export data
Definition: EdbRun.h:81
void AddView()
Definition: EdbRun.cxx:305

◆ AddView() [2/2]

void EdbRun::AddView ( EdbView view)
312{
313
314 if( view != eView ) {
315 Log(3," EdbRun::AddView","WARNING: view!=eView - inefficient cycle!");
316 SetView(view);
317 }
318
319 EdbViewHeader *viewHeader = eView->GetHeader();
320
321 //Long_t oldtime = eView->GetLastSystemTime();
322 //Long_t newtime = gSystem->Now();
323 //viewHeader->SetTime( (Int_t)((ULong_t)newtime-(ULong_t)oldtime) );
324 //eView->SetLastSystemTime(newtime);
325
326 viewHeader->SetNclusters(eView->Nclusters());
327 viewHeader->SetNsegments(eView->Nsegments());
328 viewHeader->SetEvent( GetHeader()->GetFlag(8) );
329 viewHeader->SetTrack( GetHeader()->GetFlag(9) );
330
331 eTree->Fill();
332 eView->Clear();
333}
void SetView()
Definition: EdbRun.cxx:285
TTree * eTree
tree with View objects
Definition: EdbRun.h:83
view identification
Definition: EdbView.h:26
void SetTrack(int track)
Definition: EdbView.h:107
void SetNsegments(int nseg)
Definition: EdbView.h:86
void SetEvent(int event)
Definition: EdbView.h:108
void SetNclusters(int nclu)
Definition: EdbView.h:85
EdbViewHeader * GetHeader() const
Definition: EdbView.h:163
Int_t Nsegments() const
Definition: EdbView.h:216
Int_t Nclusters() const
Definition: EdbView.h:215
void Clear()
Definition: EdbView.cxx:79

◆ AddViewAlign()

void EdbRun::AddViewAlign ( Int_t  view1,
Int_t  view2,
Int_t  area1,
Int_t  area2,
Int_t  side1,
Int_t  side2,
Float_t  dx,
Float_t  dy,
Float_t  dz,
Int_t  n1tot,
Int_t  n2tot,
Int_t  n1,
Int_t  n2,
Int_t  nsg,
Int_t  nbg,
Int_t  flag,
EdbViewHeader vh1,
EdbViewHeader vh2 
)
362{
363 if(!eViewAlign) {
364 eViewAlign = new TTree("ViewAlign","Alignment offsets");
365 eViewAlign->Branch("alpar",&eVA,"view1/I:view2/I:area1/I:area2/I:side1/I:side2/I:dx/F:dy/F:dz/F:n1tot/I:n2tot/I:n1/I:n2/I:nsg/I:nbg/I:flag/I");
366 eViewAlign->Branch("vh1","EdbViewHeader",&eVH1);
367 eViewAlign->Branch("vh2","EdbViewHeader",&eVH2);
368 }
369 if(eViewAlign) {
370 eVA.view1=view1; eVA.view2=view2;
371 eVA.area1=area1; eVA.area2=area2;
372 eVA.side1=side1; eVA.side2=side2;
373 eVA.dx=dx; eVA.dy=dy; eVA.dz=dz;
374 eVA.n1tot=n1tot; eVA.n2tot=n2tot;
375 eVA.n1=n1; eVA.n2=n2;
376 eVA.nsg=nsg; eVA.nbg=nbg;
377 eVA.flag=flag;
378 *eVH1=vh1;
379 *eVH2=vh2;
380 eViewAlign->Fill();
381 }
382}
brick dz
Definition: RecDispMC.C:107
EdbViewHeader * eVH1
alignment headers
Definition: EdbRun.h:97
AlignmentParView eVA
Definition: EdbRun.h:94
EdbViewHeader * eVH2
alignment headers
Definition: EdbRun.h:98
TTree * eViewAlign
view neighbours alignment
Definition: EdbRun.h:90
Float_t dy
found offset
Definition: EdbRun.h:41
Int_t n1
number of grains in intersection area
Definition: EdbRun.h:45
Int_t n2tot
Definition: EdbRun.h:44
Int_t side1
Definition: EdbRun.h:38
Float_t dz
found offset
Definition: EdbRun.h:42
Float_t dx
found offset
Definition: EdbRun.h:40
Int_t n1tot
total number of grains in the view
Definition: EdbRun.h:43
Int_t nsg
peak height
Definition: EdbRun.h:47
Int_t view2
2-d view
Definition: EdbRun.h:35
Int_t n2
Definition: EdbRun.h:46
Int_t flag
enum Flags {f_Applied = 0x01,f_Pin = 0x02,f_Found = 0x04}; f_Pin: view2 is a pinned view
Definition: EdbRun.h:49
Int_t area2
Definition: EdbRun.h:37
Int_t side2
Definition: EdbRun.h:39
Int_t view1
1-st view
Definition: EdbRun.h:34
Int_t area1
Definition: EdbRun.h:36
Int_t nbg
bg level (random coinsidences)
Definition: EdbRun.h:48

◆ AddViewMerge()

void EdbRun::AddViewMerge ( Int_t  view1,
Int_t  view2,
Int_t  area1,
Int_t  area2,
Int_t  side1,
Int_t  side2,
Float_t  dx,
Float_t  dy,
Float_t  dz,
Int_t  n1tot,
Int_t  n2tot,
Int_t  n1,
Int_t  n2,
Int_t  nsg,
Int_t  nbg,
Int_t  flag 
)
339{
340 if(!eViewMerge) {
341 eViewMerge = new TTree("ViewMerge","Merging offsets");
342 eViewMerge->Branch("alpar",&eVM,"view1/I:view2/I:area1/I:area2/I:side1/I:side2/I:dx/F:dy/F:dz/F:n1tot/I:n2tot/I:n1/I:n2/I:nsg/I:nbg/I:flag/I");
343 }
344 if(eViewMerge) {
345 eVM.view1=view1; eVM.view2=view2;
346 eVM.area1=area1; eVM.area2=area2;
347 eVM.side1=side1; eVM.side2=side2;
348 eVM.dx=dx; eVM.dy=dy; eVM.dz=dz;
349 eVM.n1tot=n1tot; eVM.n2tot=n2tot;
350 eVM.n1=n1; eVM.n2=n2;
351 eVM.nsg=nsg; eVM.nbg=nbg;
352 eVM.flag=flag;
353 eViewMerge->Fill();
354 }
355}
TTree * eViewMerge
view merging alignment
Definition: EdbRun.h:89
AlignmentParView eVM
Definition: EdbRun.h:93

◆ Close()

void EdbRun::Close ( )
440{
441 eFile = eTree->GetCurrentFile();
442 const char *status = eFile->GetOption();
443
444 GetFile()->cd();
445 eHeader->GetFinishTime()->Set();
446
447 if( strcmp(status,"READ") ) { // file is in "write" mode
448 if( !strcmp(status,"CREATE") ) { // save header in CREATE mode only
449 Save();
450 }
451 if(eTree) eTree->Write();
452 if(eViewMerge) eViewMerge->Write();
453 if(eViewAlign) eViewAlign->Write();
454 if(eFrameAlign) eFrameAlign->Write();
455 if(ePinViews) ePinViews->Write();
456 }
457 eFile->Purge();
458 eFile->Close();
459}
TDatime * GetFinishTime()
Definition: EdbRunHeader.h:145
void Save()
Definition: EdbRun.cxx:424
TFile * GetFile()
Definition: EdbRun.h:149

◆ Create()

void EdbRun::Create ( const char *  fname)
239{
241 if(!eMarks) eMarks = new EdbMarksSet();
242
243 if( !gSystem->AccessPathName(fname, kFileExists) )
244 Log(2,"EdbRun::Create","WARNING : file %s already exists!", fname);
245
246 Log(2,"EdbRun::Create","Open new file %s", fname);
247 eFile = new TFile(fname,"RECREATE");
248
249 eFile->SetCompressionLevel(2);
250
251 eTree = new TTree("Views","Scanning Viewes data");
252 eTree->SetAutoSave(100000000); // autosave each 100Mb
253 eTree->SetMaxTreeSize(32000000000LL); //set 32 Gb file size limit
254 Log(3,"EdbRun::Create","set maxtreesize as: %lld",eTree->GetMaxTreeSize());
255
256 TClonesArray *eClusters = eView->GetClusters();
257 TClonesArray *eSegments = eView->GetSegments();
258 TClonesArray *eTracks = eView->GetTracks();
259 TClonesArray *eFrames = eView->GetFrames();
260 int savelevel = gErrorIgnoreLevel;
261 gErrorIgnoreLevel = kError;
262 eTree->Branch( "clusters", &eClusters);
263 eTree->Branch( "segments", &eSegments);
264 eTree->Branch( "tracks", &eTracks);
265 eTree->Branch( "frames", &eFrames);
266 eTree->Branch("headers", "EdbViewHeader", eView->GetHeaderAddr());
267 gErrorIgnoreLevel = savelevel;
268 SetView();
269}
TClonesArray * GetSegments() const
Definition: EdbView.h:165
TClonesArray * GetClusters() const
Definition: EdbView.h:164
TClonesArray * GetFrames() const
Definition: EdbView.h:167
TClonesArray * GetTracks() const
Definition: EdbView.h:166
void * GetHeaderAddr()
Definition: EdbView.h:169

◆ ExtractAsciiFile()

int EdbRun::ExtractAsciiFile ( const char *  fname,
const char *  objname 
)
611{
612 FILE *fp=fopen(fname,"w");
613 if (fp==NULL) {
614 printf("ERROR open file: %s n", fname); return(-1);
615 }else
616 printf( "Write to Ascii file: %s\n", fname );
617 TObjString *str=(TObjString *)(eFile->Get(objname));
618 fprintf(fp,"%s",str->String().Data());
619 fclose(fp);
620 return str->String().Sizeof();
621}

◆ GeneratePredictions()

void EdbRun::GeneratePredictions ( int  n)
640{
642
644}
void Generate(int n)
Definition: EdbPrediction.cxx:193

◆ GetEntries()

int EdbRun::GetEntries ( ) const
inline
136{ return (Int_t)eTree->GetEntries(); }

◆ GetEntry()

EdbView * EdbRun::GetEntry ( int  entry,
int  ih = 1,
int  icl = 0,
int  iseg = 1,
int  itr = 0,
int  ifr = 0 
)
463{
464 if(eView) eView->Clear();
465 else SetView();
466
467 if( !eView->GetClusters() ) SetView();
468
469 if(ih) GetEntryHeader( entry );
470 if(icl) GetEntryClusters( entry );
471 if(iseg) GetEntrySegments( entry );
472 if(itr) GetEntryTracks( entry );
473 if(ifr) GetEntryFrames( entry );
474
475 return eView;
476}
TLegendEntry * entry
Definition: Canv_SYSTEMATICS_ALLCOMBINED__RMSEnergy__vs__Energy__ELECTRON.C:130
TClonesArray * GetEntryTracks(int entry) const
Definition: EdbRun.cxx:479
TClonesArray * GetEntryClusters(int entry) const
Definition: EdbRun.cxx:507
TObjArray * GetEntryFrames(int entry) const
Definition: EdbRun.cxx:534
EdbViewHeader * GetEntryHeader(int entry) const
Definition: EdbRun.cxx:521
TClonesArray * GetEntrySegments(int entry) const
Definition: EdbRun.cxx:493

◆ GetEntryClusters()

TClonesArray * EdbRun::GetEntryClusters ( int  entry) const
508{
509 TClonesArray *clusters = eView->GetClusters();
510 TBranchClones *clustersBranch = (TBranchClones*)(eTree->GetBranch("clusters"));
511 if(!clustersBranch) return 0;
512 //clustersBranch->SetAutoDelete( kTRUE );
513 clustersBranch->ResetReadEntry();
514 clustersBranch->SetAddress(&clusters);
515 clustersBranch->GetEntry(entry);
516
517 return clusters;
518}

◆ GetEntryFrames()

TObjArray * EdbRun::GetEntryFrames ( int  entry) const
535{
536 TObjArray *frames = eView->GetFrames();
537 TBranch *framesBranch = eTree->GetBranch("frames");
538 if(!framesBranch) return 0;
539 framesBranch->SetAddress(&frames);
540 framesBranch->GetEntry(entry);
541
542 return frames;
543}

◆ GetEntryHeader()

EdbViewHeader * EdbRun::GetEntryHeader ( int  entry) const
522{
523 EdbViewHeader *header = eView->GetHeader();
524 TBranch *headerBranch = eTree->GetBranch("headers");
525 if(!headerBranch) return 0;
526 headerBranch->SetAddress(&header);
527 headerBranch->GetEntry(entry);
528 // eTree->GetEntry(entry);
529
530 return header;
531}

◆ GetEntryNumberWithIndex()

int EdbRun::GetEntryNumberWithIndex ( int  event,
int  view 
) const
inline
152 { return eTree->GetEntryNumberWithIndex( event,view ); }
int event
Definition: shower_tr.C:25

◆ GetEntrySegments()

TClonesArray * EdbRun::GetEntrySegments ( int  entry) const
494{
495 TClonesArray *segments = eView->GetSegments();
496 TBranchClones *segmentsBranch = (TBranchClones*)(eTree->GetBranch("segments"));
497 if(!segmentsBranch) return 0;
498
499 segmentsBranch->ResetReadEntry();
500 segmentsBranch->SetAddress(&segments);
501 segmentsBranch->GetEntry(entry);
502
503 return segments;
504}

◆ GetEntryTracks()

TClonesArray * EdbRun::GetEntryTracks ( int  entry) const
480{
481 TClonesArray *tracks=eView->GetTracks();
482 TBranchClones *tracksBranch = (TBranchClones*)(eTree->GetBranch("tracks"));
483 if(!tracksBranch) return 0;
484
485 tracksBranch->ResetReadEntry();
486 tracksBranch->SetAddress(&tracks);
487 tracksBranch->GetEntry(entry);
488
489 return tracks;
490}
TTree * tracks
Definition: check_tr.C:19

◆ GetFile()

TFile * EdbRun::GetFile ( )
inline
149{ return eFile; }

◆ GetFinishTime()

TDatime * EdbRun::GetFinishTime ( ) const
inline
129{ return eHeader->GetFinishTime(); }

◆ GetHeader()

EdbRunHeader * EdbRun::GetHeader ( ) const
inline
138{ return eHeader; }

◆ GetMarks()

EdbMarksSet * EdbRun::GetMarks ( ) const
inline
120{ return eMarks; }

◆ GetPinViews()

TTree * EdbRun::GetPinViews ( ) const
inline
114{return ePinViews;}

◆ GetPrediction()

EdbPredictionDC * EdbRun::GetPrediction ( int  ip)
inline
116{ return ePredictions->GetPrediction(ip); }
EdbPredictionDC * GetPrediction(int i) const
Definition: EdbPrediction.cxx:281

◆ GetPredictions()

EdbPredictionsBox * EdbRun::GetPredictions ( ) const
inline
119{ return ePredictions; }

◆ GetRunID()

int EdbRun::GetRunID ( ) const
inline
127{ return eHeader->GetRunID(); }
Int_t GetRunID() const
Definition: EdbRunHeader.h:142

◆ GetStartTime()

TDatime * EdbRun::GetStartTime ( ) const
inline
128{ return eHeader->GetStartTime(); }
TDatime * GetStartTime()
Definition: EdbRunHeader.h:144

◆ GetTree()

TTree * EdbRun::GetTree ( ) const
inline
113{return eTree; }

◆ GetView()

EdbView * EdbRun::GetView ( ) const
inline
110{ return eView; }

◆ Init()

void EdbRun::Init ( void  )
146{
147 eView = new EdbView(); //each run has his own view
148 eHeader = 0;
149 eTree = 0;
150 eFile = 0;
151 ePath = "";
152 ePredictions = 0;
153 eMarks = 0;
154 eViewMerge = 0;
155 eViewAlign = 0;
156 eFrameAlign = 0;
157 ePinViews = 0;
158 ePVH=new EdbViewHeader();
159 eVH1=new EdbViewHeader();
160 eVH2=new EdbViewHeader();
161}
Base scanning data object: entry into Run tree.
Definition: EdbView.h:134

◆ Npredictions()

int EdbRun::Npredictions ( ) const
inline
117{ return ePredictions->N(); }
Int_t N() const
mandatory virtual functions:
Definition: EdbPrediction.h:92

◆ Open()

void EdbRun::Open ( const char *  fname)
179{
180 Log(3,"EdbRun::Open","Open an existing file %s", fname);
181 eFile = new TFile(fname);
182 if(!eFile) return;
183 eTree = (TTree*)eFile->Get("Views");
184 eViewMerge = (TTree*)eFile->Get("ViewMerge");
185 eViewAlign = (TTree*)eFile->Get("ViewAlign");
186 eFrameAlign = (TTree*)eFile->Get("FrameAlign");
187 ePinViews = (TTree*)eFile->Get("PinViews");
188 if(!eTree) {
189 Log(1,"EdbRun::Open","ERROR: %s has no Views tree",fname);
190 return;
191 }
192
193 SetView();
194
195 SafeDelete(eHeader);
196 eHeader = (EdbRunHeader*)eFile->Get("RunHeader");
197
198 SafeDelete(ePredictions);
199 ePredictions = (EdbPredictionsBox*)eFile->Get("Predictions");
200
201 SafeDelete(eMarks);
202 eMarks = (EdbMarksSet*)eFile->Get("Marks");
203
204 if(gEDBDEBUGLEVEL>2) Print();
205}
void Print() const
Definition: EdbRun.cxx:580
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ OpenUpdate()

void EdbRun::OpenUpdate ( const char *  fname)
209{
210 Log(3,"EdbRun::OpenUpdate","Open an existing file for update %s", fname);
211 eFile = new TFile(fname,"UPDATE");
212
213 eTree = (TTree*)eFile->Get("Views");
214 SetView();
215
216 eViewMerge = (TTree*)eFile->Get("ViewMerge");
217 eViewAlign = (TTree*)eFile->Get("ViewAlign");
218 eFrameAlign = (TTree*)eFile->Get("FrameAlign");
219 ePinViews = (TTree*)eFile->Get("PinViews");
220 eViewAlign->SetBranchAddress("alpar",&eVA);
221 eViewMerge->SetBranchAddress("alpar",&eVM);
222 eFrameAlign->SetBranchAddress("alpar",&eFA);
223 ePinViews->SetBranchAddress("pin",&ePVH);
224
225 SafeDelete(eHeader);
226 eHeader = (EdbRunHeader*)eFile->Get("RunHeader");
227
228 SafeDelete(ePredictions);
229 ePredictions = (EdbPredictionsBox*)eFile->Get("Predictions");
230
231 SafeDelete(eMarks);
232 eMarks = (EdbMarksSet*)eFile->Get("Marks");
233
234 if(gEDBDEBUGLEVEL>2) Print();
235}

◆ Print()

void EdbRun::Print ( ) const
581{
582 if(eHeader) eHeader->Print();
583
584 int nentr=0;
585 if(eTree) nentr = (Int_t)eTree->GetEntries();
586 printf("Number of entries: \t %d \n", nentr);
587 printf("--------------------------------------------------------------------\n\n");
588 //PrintBranchesStatus();
589}
void Print()
other routines
Definition: EdbRunHeader.cxx:73

◆ PrintBranchesStatus()

void EdbRun::PrintBranchesStatus ( ) const

not supported...

547{
549
550 EdbViewHeader *header = GetEntryHeader( 0 );
551 if(header) printf("header branch exists\n");
552 TClonesArray *clusters = GetEntryClusters( 0 );
553 if(clusters) printf("clusters branch exists\n");
554 TObjArray *frames = GetEntryFrames( 0 );
555 if(frames) printf("frames branch exists\n");
556}

◆ PrintLog()

void EdbRun::PrintLog ( const char *  fname) const
560{
561 char buf[256];
562
563 sprintf( buf,"run%4.4d %s \t%s \t%s",
564 GetRunID(),
565 eHeader->GetName(),
566 eHeader->GetPlate()->GetName(),
567 eHeader->GetTitle() );
568
569 printf("**************************************************************************\n");
570 printf("\n%s\n", buf);
571 printf("**************************************************************************\n");
572
573 FILE *fp = fopen( fname,"a" );
574 fprintf(fp,"%s", buf);
575 fclose(fp);
576
577}
EdbPlate * GetPlate() const
Definition: EdbRunHeader.h:141
int GetRunID() const
Definition: EdbRun.h:127

◆ Save()

void EdbRun::Save ( )
425{
426 if(eHeader) eHeader->Write("RunHeader");
427 if(ePredictions) ePredictions->Write("Predictions");
428 if(eMarks) eMarks->Write("Marks");
429
430 //if(eViewMerge) eViewMerge->Write();
431 //if(eViewAlign) eViewAlign->Write();
432 //if(eFrameAlign) eViewMerge->Write();
433
434 //SaveViews();
435 //eFile->Purge();
436}

◆ SaveViews()

void EdbRun::SaveViews ( )
inline
168{ eTree->AutoSave(); }

◆ SelectOpenMode()

void EdbRun::SelectOpenMode ( const char *  fname,
const char *  status = "READ" 
)
165{
166 if ( !strcmp(status,"RECREATE") ) Create(fname);
167 else if( !strcmp(status,"READ") ) Open(fname);
168 else if( !strcmp(status,"UPDATE") ) OpenUpdate(fname);
169 else{
170 if( !gSystem->AccessPathName(fname, kFileExists) ) Open(fname);
171 else Create(fname);
172 }
173
174 Log(3,"EdbRun::SelectOpenMode","root file version: %d",eFile->GetVersion());
175}
void Open(const char *fname)
Definition: EdbRun.cxx:178
void OpenUpdate(const char *fname)
Definition: EdbRun.cxx:208
void Create(const char *fname)
Definition: EdbRun.cxx:238

◆ SetComment()

void EdbRun::SetComment ( const char *  com)
inline
132{ eHeader->SetComment(com); }
void SetComment(const char *com)
Definition: EdbRunHeader.h:135

◆ SetMarks()

void EdbRun::SetMarks ( EdbMarksSet marks)
inline
121{ eMarks = marks; }

◆ SetRunID()

void EdbRun::SetRunID ( int  id)
inline
131{ eHeader->SetRunID(id); }
void SetRunID(int id)
Definition: EdbRunHeader.h:130

◆ SetTitle()

void EdbRun::SetTitle ( const char *  tit)
inline
133{ eHeader->SetTitle(tit); }

◆ SetView() [1/2]

void EdbRun::SetView ( )
286{
287 if(!eView) {
288 Log(2,"EdbRun::SetView","WARNING: *eView=0");
289 eView = new EdbView();
290 }
291 else if( !eView->GetClusters() ) {
292 Log(1,"EdbRun::SetView","ERROR: eView was deleted!");
293 eView = new EdbView();
294 }
295
296 Log(3,"EdbRun::SetView","Note: EdbRun::SetView");
297 if(eView->GetHeader()) eTree->SetBranchAddress("headers", eView->GetHeaderAddr() );
298 if(eView->GetClusters()) eTree->SetBranchAddress("clusters", eView->GetClustersAddr() );
299 if(eView->GetSegments()) eTree->SetBranchAddress("segments", eView->GetSegmentsAddr() );
300 if(eView->GetTracks()) eTree->SetBranchAddress("tracks", eView->GetTracksAddr() );
301 if(eView->GetFrames()) eTree->SetBranchAddress("frames", eView->GetFramesAddr() );
302}
void * GetSegmentsAddr()
Definition: EdbView.h:171
void * GetClustersAddr()
Definition: EdbView.h:170
void * GetTracksAddr()
Definition: EdbView.h:172
void * GetFramesAddr()
Definition: EdbView.h:173

◆ SetView() [2/2]

void EdbRun::SetView ( EdbView view)
273{
274 if(!view) { Log(1,"EdbRun::SetView","ERROR: EdbRun::SetView: *view=0"); return; }
275
276 if(eView != view) {
277 if(eView) { eView->Clear(); SafeDelete(eView); }
278 eView = view;
279 }
280
281 SetView();
282}

◆ TransformDC()

void EdbRun::TransformDC ( )
626{
627 if(ePredictions){
628 if(eMarks){
629 EdbAffine2D *aff = eMarks->Abs2Stage();
631 printf("transfrom predictions according to fiducial marks:\n");
632 aff->Print();
633 SafeDelete(aff);
634 }
635 }
636}
Definition: EdbAffine.h:17
void Print(Option_t *opt="") const
Definition: EdbAffine.cxx:52
EdbAffine2D * Abs2Stage() const
Definition: EdbFiducial.cxx:271
virtual void Transform(const EdbAffine2D *a)
Definition: EdbVirtual.cxx:155

Member Data Documentation

◆ eFA

AlignmentParFrame EdbRun::eFA
private

◆ eFile

TFile* EdbRun::eFile
private

file associated with the Run

◆ eFrameAlign

TTree* EdbRun::eFrameAlign
private

frames alignment

◆ eHeader

EdbRunHeader* EdbRun::eHeader
private

run header

◆ eMarks

EdbMarksSet* EdbRun::eMarks

fiducial marks

◆ ePath

TString EdbRun::ePath
private

runs directory path

◆ ePinViews

TTree* EdbRun::ePinViews
private

pinned views

◆ ePredictions

EdbPredictionsBox* EdbRun::ePredictions
private

predictions to scan ($c)

◆ ePVH

EdbViewHeader* EdbRun::ePVH
private

to add pinned view header

◆ eTree

TTree* EdbRun::eTree
private

tree with View objects

◆ eVA

AlignmentParView EdbRun::eVA
private

◆ eVH1

EdbViewHeader* EdbRun::eVH1
private

alignment headers

◆ eVH2

EdbViewHeader* EdbRun::eVH2
private

alignment headers

◆ eView

EdbView* EdbRun::eView
private

view using for import and export data

◆ eViewAlign

TTree* EdbRun::eViewAlign
private

view neighbours alignment

◆ eViewMerge

TTree* EdbRun::eViewMerge
private

view merging alignment

◆ eVM

AlignmentParView EdbRun::eVM
private

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