FEDRA emulsion software from the OPERA Collaboration
EdbRunAccess Class Reference

helper class for access to the run data More...

#include <EdbRunAccess.h>

Inheritance diagram for EdbRunAccess:
Collaboration diagram for EdbRunAccess:

Public Member Functions

bool AcceptRawSegment (EdbView *view, int ud, EdbSegP &segP, int side, int entry)
 
bool AddRWDToRun (char *rwdname, const char *options="")
 
void AddSegmentCut (int ud, int xi, float min[5], float max[5])
 
void AddSegmentCut (int ud, int xi, float var[10])
 
void AddSegmentCut (int xi, const char *cutline)
 
void ApplyCorrections (const EdbView &view, EdbSegment &s, const int rs)
 
float CalculateSegmentChi2 (EdbSegment *seg, float sx, float sy, float sz)
 
int Check0Views (EdbPattern &pat, int thres=1)
 
int CheckEmptyViews (EdbPattern &pat)
 
float CheckMeanSegsPerView (EdbPattern &pat)
 
void CheckRunLine ()
 
TH2F * CheckUpDownOffsets ()
 
void CheckViewSize ()
 
void CheckViewStep ()
 
void CheckViewStep (int ud)
 
void ClearCuts ()
 
bool CopyRawDataXY (float x0, float y0, float dR, const char *file)
 
 EdbRunAccess ()
 
 EdbRunAccess (const char *fname)
 
 EdbRunAccess (EdbRun *run)
 
void FillDZMaps ()
 
int FillVP ()
 
int FirstArea () const
 
EdbScanCondGetCond (int ud)
 
EdbSegmentCutGetCut (int ud, int i)
 
float GetCutGR () const
 
int GetEntryXY (int ud, float x, float y)
 
EdbLayerGetLayer (int id)
 
EdbScanCondGetMakeCond (int ud)
 
EdbLayerGetMakeLayer (int id)
 
Int_t GetNareas ()
 
Int_t GetNviewsPerArea ()
 
int GetPatternData (EdbPattern &pat, int side, int nviews, TArrayI &srt, int &nrej)
 
int GetPatternDataForPrediction (int id, int side, EdbPattern &pat)
 
int GetPatternXY (EdbSegP &s, int side, EdbPattern &pat, float rmin=200)
 
int GetPatternXYcut (EdbSegP &s, int side, EdbPattern &pat, float dr, float dt)
 
EdbSegmentGetRawSegment (EdbView &v, int sid, int rs=0)
 
EdbSegmentGetRawSegment (int vid, int sid, int rs=0)
 
EdbSegmentGetRawSegmentN (int vid, int sid, int rs=0)
 
float GetRawSegmentPix (EdbSegment *seg)
 
EdbRunGetRun () const
 
int GetViewsArea (int ud, TArrayI &entr, float xmin, float xmax, float ymin, float ymax)
 
int GetViewsArea (int ud, TArrayI &entr, int area, float &xmin, float &xmax, float &ymin, float &ymax)
 
int GetViewsAreaMarg (int ud, TArrayI &entr, int area, float xmarg, float ymarg)
 
int GetViewsXY (int ud, TArrayI &entr, float x, float y, float r=200.)
 
int GetVolumeArea (EdbPatternsVolume &vol, int area)
 
int GetVolumeData (EdbPatternsVolume &vol, int nviews, TArrayI &srt, int &nrej)
 
int GetVolumeXY (EdbSegP &s, EdbPatternsVolume &vol)
 
EdbPatternGetVP (int ud) const
 
void GuessNviewsPerArea ()
 
bool InitRun (const char *runfile=0, bool do_update=false)
 
bool InitRunFromRWC (char *rwcname, bool bAddRWD=true, const char *options="")
 
int LastArea () const
 
int NCuts (int ud)
 
float OverlapX (int ud)
 
float OverlapY (int ud)
 
bool PassCuts (int ud, EdbSegment &seg)
 
void Print ()
 
void PrintStat ()
 
void ReadVAfile ()
 
float SegmentWeight (const EdbSegment &s)
 
void Set0 ()
 
void SetCond (int ud, EdbScanCond &cond)
 
void SetCutBottom (int ud, float wmin)
 
void SetCutLeft (int ud, float wmin)
 
void SetCutRight (int ud, float wmin)
 
void SetCutTop (int ud, float wmin)
 
void SetPixelCorrection (const char *str)
 
bool SetSegmentAtExternalSurface (EdbSegment *seg, int side)
 
int ViewSide (const EdbView *view) const
 
float ViewX (const EdbView *view) const
 
float ViewY (const EdbView *view) const
 
virtual ~EdbRunAccess ()
 

Public Attributes

Int_t eAFID
 if =1 - use affine transformations of the fiducial marks More...
 
Int_t eCLUST
 1-use clusters, 0 - do not More...
 
Int_t eDoPixelCorr
 apply or not pix/mic correction when read data (default is 0) More...
 
Bool_t eDoViewAnalysis
 fill or not the histograms for optional view analysis More...
 
TGraph2D * eGraphDZ [3]
 keep the base/layer1/layer2 thickness calculated using eZ1/eZ2/eZ3/eZ4 for each view More...
 
TGraph2D * eGraphZ [4]
 keep z1/z2/z3/z4 surfaces using eZ1/eZ2/eZ3/eZ4 for each view More...
 
TCut eHeaderCut
 header cut to be applied in run initialization More...
 
EdbH2 eHViewXY [3]
 XY segments distribution in a view local coords. More...
 
Bool_t eInvertSides
 0 -do nothing, 1-invert sides More...
 
Float_t ePixelCorrX
 pixel/micron correction factor to be applied for data More...
 
Float_t ePixelCorrY
 
Int_t eTracking
 to test tracking alorithm: -1-ignored(def),0/1 - trackings to accept More...
 
Bool_t eUseDensityAsW
 in case of LASSO tracking possible to use eSigmaY as eW
More...
 
Bool_t eUseExternalSurface
 if true - set segment position corrisponding to the very external cluster More...
 
Int_t eWeightAlg
 0-puls, 1 - density(former eUseDensityAsW), 2-Likelyhood (eSigmaX) More...
 

Private Attributes

EdbAffine2DeAffStage2Abs
 affine transformation extracted from Marks (if AFID=11) More...
 
EdbScanCondeCond [3]
 
Float_t eCutGR
 grain cut (chi) More...
 
TObjArray * eCuts [3]
 arrays of cuts to be applied to segments More...
 
Int_t eFirstArea
 
Int_t eLastArea
 
EdbLayereLayers [3]
 base(0),up(1),down(2) layers More...
 
Int_t eNareas
 
Int_t eNviewsPerArea
 
EdbRuneRun
 pointer to the run to be accessed More...
 
TString eRunFileName
 
TObjArray * eViewCorr
 corrections obtained from the views alignment More...
 
Float_t eViewXmax [3]
 
Float_t eViewXmin [3]
 
Float_t eViewYmax [3]
 
Float_t eViewYmin [3]
 
EdbPatterneVP [3]
 
Float_t eXmax
 
Float_t eXmin
 
Float_t eXstep [3]
 
Float_t eYmax
 
Float_t eYmin
 
Float_t eYstep [3]
 

Detailed Description

helper class for access to the run data

//////////////////////////////////////////////////////////////////////// // EdbRunAccess // // OPERA data Run Access helper class // // ////////////////////////////////////////////////////////////////////////

Constructor & Destructor Documentation

◆ EdbRunAccess() [1/3]

EdbRunAccess::EdbRunAccess ( )

◆ EdbRunAccess() [2/3]

EdbRunAccess::EdbRunAccess ( EdbRun r)

56{
57 Set0();
58 eRun=r;
59}
EdbRun * eRun
pointer to the run to be accessed
Definition: EdbRunAccess.h:51
void Set0()
Definition: EdbRunAccess.cxx:62
void r(int rid=2)
Definition: test.C:201

◆ EdbRunAccess() [3/3]

EdbRunAccess::EdbRunAccess ( const char *  fname)

37{
38 Set0();
40}
TString eRunFileName
Definition: EdbRunAccess.h:50
const char * fname
Definition: mc2raw.cxx:41

◆ ~EdbRunAccess()

EdbRunAccess::~EdbRunAccess ( )
virtual

44{
45 if(eRun) { delete eRun; eRun=0; }
46 for(int i=0; i<3; i++) {
47 if(eVP[i]) { delete eVP[i]; eVP[i]=0; }
48 if(eLayers[i]) { delete eLayers[i]; eLayers[i]=0; }
49 if(eCuts[i]) { delete eCuts[i]; eCuts[i]=0; }
50 if(eCond[i]) { delete eCond[i]; eCond[i]=0; }
51 }
52}
EdbLayer * eLayers[3]
base(0),up(1),down(2) layers
Definition: EdbRunAccess.h:60
TObjArray * eCuts[3]
arrays of cuts to be applied to segments
Definition: EdbRunAccess.h:62
EdbPattern * eVP[3]
Definition: EdbRunAccess.h:63
EdbScanCond * eCond[3]
Definition: EdbRunAccess.h:61

Member Function Documentation

◆ AcceptRawSegment()

bool EdbRunAccess::AcceptRawSegment ( EdbView view,
int  id,
EdbSegP segP,
int  side,
int  entry 
)

944{
945 EdbSegment *seg = view->GetSegment(id);
946
947 if( !PassCuts(side,*seg) ) return false;
948
949 if(eTracking>=0)
950 if(((Int_t)(seg->GetUniqueID()>>16))!=eTracking) return false;
951
952 float pix, chi2;
953 if(eCLUST) {
956 } else {
957 pix = GetRawSegmentPix(seg);
958 segP.SetVolume( pix );
960 GetCond(1)->SigmaXgr(), //TODO: side logic
961 GetCond(1)->SigmaYgr(),
962 GetCond(1)->SigmaZgr());
963 if(chi2>GetCutGR()) return false;
964 segP.SetChi2( chi2 );
965 }
966 }
967
968 if(eDoViewAnalysis) {
969 eHViewXY[side].Fill( seg->GetX0(), seg->GetY0());
970 }
971
972 EdbLayer *layer=GetLayer(side);
973 const EdbAffine2D *affview = 0;
974 if (eAFID==1) affview = view->GetHeader()->GetAffine();
975 else if (eAFID==11) affview = eAffStage2Abs;
976 else if(eAFID==2) {
977 affview = (EdbAffine2D*)eViewCorr->At(entry);
978 if(id==0) {
979 view->GetHeader()->GetAffine()->Print();
980 affview->Print();
981 }
982 }
983
984 float x,y,z,tx,ty,puls;
985 tx = seg->GetTx()/layer->Shr();
986 ty = seg->GetTy()/layer->Shr();
987 x = seg->GetX0();
988 y = seg->GetY0();
989
990 if(eDoPixelCorr) {
991 x *= ePixelCorrX;
992 y *= ePixelCorrY;
993 //printf("%f %f \n",ePixelCorrX,ePixelCorrY);
994 }
995
996 x += layer->Zcorr()*tx;
997 y += layer->Zcorr()*ty;
998 z = layer->Z() + layer->Zcorr();
999
1000 float xx, yy, txr, tyr;
1001 if( eAFID==0 ) {
1002 xx = x+view->GetXview();
1003 yy = y+view->GetYview();
1004 txr = tx;
1005 tyr = ty;
1006 }
1007 else if(eAFID==11) {
1008 x += view->GetXview();
1009 y += view->GetYview();
1010 xx = affview->A11()*x+affview->A12()*y+affview->B1();
1011 yy = affview->A21()*x+affview->A22()*y+affview->B2();
1012 txr = affview->A11()*tx+affview->A12()*ty; // rotate also angles
1013 tyr = affview->A21()*tx+affview->A22()*ty;
1014 }
1015 else {
1016 //seg->Transform( affview );
1017 xx = affview->A11()*x+affview->A12()*y+affview->B1();
1018 yy = affview->A21()*x+affview->A22()*y+affview->B2();
1019 txr = affview->A11()*tx+affview->A12()*ty; // rotate also angles
1020 tyr = affview->A21()*tx+affview->A22()*ty;
1021 }
1022
1023 if(eWeightAlg==2) puls = seg->GetPuls(); // use SigmaX for cuts only
1024 else puls = SegmentWeight(*seg);
1025
1026 EdbAffine2D *aff = layer->GetAffineTXTY();
1027 float txx = aff->A11()*txr+aff->A12()*tyr+aff->B1();
1028 float tyy = aff->A21()*txr+aff->A22()*tyr+aff->B2();
1029
1030 segP.Set( seg->GetID(),xx,yy,txx,tyy,puls,0);
1031 segP.SetZ( z );
1032 segP.SetDZ( seg->GetDz()*layer->Shr() );
1033 segP.SetW( puls );
1034 segP.SetVolume( seg->GetVolume() );
1035 segP.SetChi2( seg->GetSigmaX() );
1036 segP.SetProb( seg->GetSigmaY() ); //test: grains density after LASSO
1037 segP.SetVid(entry,id);
1038 segP.SetAid(view->GetAreaID(),view->GetViewID(), side);
1039
1040 return true;
1041}
TLegendEntry * entry
Definition: Canv_SYSTEMATICS_ALLCOMBINED__RMSEnergy__vs__Energy__ELECTRON.C:130
Definition: EdbAffine.h:17
Float_t B2() const
Definition: EdbAffine.h:48
Float_t A22() const
Definition: EdbAffine.h:46
void Print(Option_t *opt="") const
Definition: EdbAffine.cxx:52
Float_t A21() const
Definition: EdbAffine.h:45
Float_t A12() const
Definition: EdbAffine.h:44
Float_t B1() const
Definition: EdbAffine.h:47
Float_t A11() const
Definition: EdbAffine.h:43
int Fill(float x, float y)
Definition: EdbCell2.h:88
Definition: EdbLayer.h:39
float Zcorr() const
Definition: EdbLayer.h:90
float Shr() const
Definition: EdbLayer.h:89
EdbAffine2D * GetAffineTXTY()
Definition: EdbLayer.h:120
float Z() const
Definition: EdbLayer.h:77
bool SetSegmentAtExternalSurface(EdbSegment *seg, int side)
Definition: EdbRunAccess.cxx:1167
EdbH2 eHViewXY[3]
XY segments distribution in a view local coords.
Definition: EdbRunAccess.h:31
float CalculateSegmentChi2(EdbSegment *seg, float sx, float sy, float sz)
Definition: EdbRunAccess.cxx:1061
EdbLayer * GetLayer(int id)
Definition: EdbRunAccess.h:113
TObjArray * eViewCorr
corrections obtained from the views alignment
Definition: EdbRunAccess.h:75
Float_t ePixelCorrY
Definition: EdbRunAccess.h:41
bool PassCuts(int ud, EdbSegment &seg)
Definition: EdbRunAccess.cxx:1044
Int_t eCLUST
1-use clusters, 0 - do not
Definition: EdbRunAccess.h:27
float GetCutGR() const
Definition: EdbRunAccess.h:165
Int_t eAFID
if =1 - use affine transformations of the fiducial marks
Definition: EdbRunAccess.h:26
Int_t eTracking
to test tracking alorithm: -1-ignored(def),0/1 - trackings to accept
Definition: EdbRunAccess.h:44
Bool_t eDoViewAnalysis
fill or not the histograms for optional view analysis
Definition: EdbRunAccess.h:30
Bool_t eUseExternalSurface
if true - set segment position corrisponding to the very external cluster
Definition: EdbRunAccess.h:28
float SegmentWeight(const EdbSegment &s)
Definition: EdbRunAccess.cxx:935
EdbScanCond * GetCond(int ud)
Definition: EdbRunAccess.h:175
Int_t eDoPixelCorr
apply or not pix/mic correction when read data (default is 0)
Definition: EdbRunAccess.h:39
Int_t eWeightAlg
0-puls, 1 - density(former eUseDensityAsW), 2-Likelyhood (eSigmaX)
Definition: EdbRunAccess.h:33
float GetRawSegmentPix(EdbSegment *seg)
Definition: EdbRunAccess.cxx:1151
EdbAffine2D * eAffStage2Abs
affine transformation extracted from Marks (if AFID=11)
Definition: EdbRunAccess.h:53
Float_t ePixelCorrX
pixel/micron correction factor to be applied for data
Definition: EdbRunAccess.h:40
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 GetY0() const
Definition: EdbSegment.h:39
virtual Float_t GetTy() const
Definition: EdbSegment.h:42
void SetVolume(float w)
Definition: EdbSegP.h:136
void SetProb(float prob)
Definition: EdbSegP.h:134
void SetZ(float z)
Definition: EdbSegP.h:125
void SetW(float w)
Definition: EdbSegP.h:132
void SetChi2(float chi2)
Definition: EdbSegP.h:135
void SetDZ(float dz)
Definition: EdbSegP.h:126
void SetVid(int vid, int sid)
Definition: EdbSegP.h:137
void SetAid(int a, int v, int side=0)
Definition: EdbSegP.h:138
void Set(int id, float x, float y, float tx, float ty, float w, int flag)
Definition: EdbSegP.h:87
segment of the track
Definition: EdbSegment.h:63
Int_t GetVolume() const
Definition: EdbSegment.h:91
Int_t GetPuls() const
Definition: EdbSegment.h:90
float GetSigmaX() const
Definition: EdbSegment.h:86
Int_t GetID() const
Definition: EdbSegment.h:92
float GetSigmaY() const
Definition: EdbSegment.h:87
EdbAffine2D const * GetAffine() const
Definition: EdbView.h:72
Float_t GetXview() const
Definition: EdbView.h:193
EdbViewHeader * GetHeader() const
Definition: EdbView.h:163
EdbSegment * GetSegment(int i) const
Definition: EdbView.h:219
Int_t GetViewID() const
Definition: EdbView.h:190
Int_t GetAreaID() const
Definition: EdbView.h:191
Float_t GetYview() const
Definition: EdbView.h:194
Float_t chi2
Definition: testBGReduction_By_ANN.C:14

◆ AddRWDToRun()

bool EdbRunAccess::AddRWDToRun ( char *  rwdname,
const char *  options = "" 
)

258{
259 Log(2,"EdbRunAccess::AddRWDToRun"," %s\n",rwdname);
260 if(!eRun) return false;
261 if( gSystem->AccessPathName(rwdname, kFileExists) ) {
262 Log(1,"EdbRunAccess::AddRWDToRun","ERROR open file: %s\n",rwdname);
263 return false;
264 }
265
266 int f;//fragment index
268 f++;
269 if( !AddRWD( eRun, rwdname, f, options) ) return false;
271 return true;
272}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
FILE * f
Definition: RecDispMC.C:150
void SetNareas(int n)
Definition: EdbRunHeader.h:134
Int_t GetNareas() const
Definition: EdbRunHeader.h:148
EdbRunHeader * GetHeader() const
Definition: EdbRun.h:138
int AddRWD(EdbRun *run, char *rwdname, int fragID, const char *options)
Definition: libDataConversion.cpp:176

◆ AddSegmentCut() [1/3]

void EdbRunAccess::AddSegmentCut ( int  layer,
int  xi,
float  min[5],
float  max[5] 
)

1207{
1208 float var[10]={ min[0],max[0], min[1],max[1], min[2],max[2], min[3],max[3], min[4],max[4] };
1209 AddSegmentCut(layer,xi,var);
1210}
float min(TClonesArray *t)
Definition: bitview.cxx:275
void AddSegmentCut(int xi, const char *cutline)
Definition: EdbRunAccess.cxx:1187
int max
Definition: check_shower.C:41

◆ AddSegmentCut() [2/3]

void EdbRunAccess::AddSegmentCut ( int  layer,
int  xi,
float  var[10] 
)

1200{
1201 if(!eCuts[layer]) eCuts[layer] = new TObjArray();
1202 eCuts[layer]->Add( new EdbSegmentCut(xi,var) );
1203}
Definition: EdbSegmentCut.h:6

◆ AddSegmentCut() [3/3]

void EdbRunAccess::AddSegmentCut ( int  xi,
const char *  cutline 
)

1188{
1189 float var[10];
1190 int onoff=-1;
1191 if(sscanf(cutline,"%d %f %f %f %f %f %f %f %f %f %f",&onoff,
1192 &var[0],&var[1], &var[2],&var[3], &var[4],&var[5], &var[6],&var[7], &var[8], &var[9] )==11)
1193 if(onoff>-1) {
1194 AddSegmentCut(1,xi,var);
1195 AddSegmentCut(2,xi,var);
1196 }
1197}

◆ ApplyCorrections()

void EdbRunAccess::ApplyCorrections ( const EdbView view,
EdbSegment s,
const int  rs 
)

rs=0 local view RS
rs==1 plate RS (use layer 1&2 information Note that the last emlink
should be done with all corrections off to get the consistent data here
rs==1 brick RS - use layer 0 info

1106{
1111
1112 if(rs==0) return; // do nothing - remain in local view RS
1113 if(rs>0) // use layers information
1114 {
1115 if(eAFID) s.Transform( view.GetHeader()->GetAffine() );
1116 else {
1117 s.SetX0( s.GetX0()+view.GetXview());
1118 s.SetY0( s.GetY0()+view.GetYview());
1119 }
1120 EdbSegP sp(0,s.GetX0(),s.GetY0(),s.GetTx(),s.GetTy());
1121 sp.SetZ( GetLayer( ViewSide(&view) )->Zxy() ); //forget original Z here...
1122 GetLayer( ViewSide(&view) )->CorrectSeg(sp); //get to plate RS
1123 if(rs==2) {
1124 GetLayer(0)->CorrectSeg(sp); //get to brick RS
1125 sp.SetZ( sp.Z() + GetLayer(0)->Zxy() );
1126 }
1127 s.SetTx( sp.TX() );
1128 s.SetTy( sp.TY() );
1129 s.SetX0( sp.X() );
1130 s.SetY0( sp.Y() );
1131 s.SetZ0( sp.Z() );
1132 }
1133}
void CorrectSeg(EdbSegP &s)
Definition: EdbLayer.cxx:69
int ViewSide(const EdbView *view) const
Definition: EdbRunAccess.cxx:443
Definition: EdbSegP.h:21
s
Definition: check_shower.C:55

◆ CalculateSegmentChi2()

float EdbRunAccess::CalculateSegmentChi2 ( EdbSegment seg,
float  sx,
float  sy,
float  sz 
)

assumed that clusters are attached to segments

1062{
1064 double chi2=0;
1065 EdbCluster *cl=0;
1066 TObjArray *clusters = seg->GetElements();
1067 if(!clusters) return 0;
1068 int ncl = clusters->GetLast()+1;
1069 if(ncl<=0) return 0;
1070
1071 float xyz1[3], xyz2[3]; // segment line parametrized as 2 points
1072 float xyz[3];
1073 bool inside=true;
1074
1075 xyz1[0] = seg->GetX0() /sx;
1076 xyz1[1] = seg->GetY0() /sy;
1077 xyz1[2] = seg->GetZ0() /sz;
1078 xyz2[0] = (seg->GetX0() + seg->GetDz()*seg->GetTx()) /sx;
1079 xyz2[1] = (seg->GetY0() + seg->GetDz()*seg->GetTy()) /sy;
1080 xyz2[2] = (seg->GetZ0() + seg->GetDz()) /sz;
1081
1082 double d;
1083 for(int i=0; i<ncl; i++ ) {
1084 cl = (EdbCluster*)clusters->At(i);
1085 xyz[0] = cl->GetX()/sx;
1086 xyz[1] = cl->GetY()/sy;
1087 xyz[2] = cl->GetZ()/sz;
1088 d = EdbMath::DistancePointLine3(xyz,xyz1,xyz2, inside);
1089 chi2 += d*d;
1090 }
1091
1092 return TMath::Sqrt(chi2/ncl);
1093}
void d()
Definition: RecDispEX.C:381
Definition: EdbCluster.h:19
Float_t GetX() const
Definition: EdbCluster.h:51
Float_t GetY() const
Definition: EdbCluster.h:52
Float_t GetZ() const
Definition: EdbCluster.h:53
static double DistancePointLine3(float Point[3], float LineStart[3], float LineEnd[3], bool inside)
Definition: EdbMath.cxx:61
virtual Float_t GetZ0() const
Definition: EdbSegment.h:40
TObjArray * GetElements() const
Definition: EdbSegment.h:117

◆ Check0Views()

int EdbRunAccess::Check0Views ( EdbPattern pat,
int  thres = 1 
)

774{
775 Long_t n0=0, nview=pat.N();
776 if(!nview) return 0;
777 for(int i=0; i<nview; i++) if( pat.GetSegment(i)->Flag() < thres) n0++;
778 return n0;
779}
Int_t Flag() const
Definition: EdbSegP.h:149
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66

◆ CheckEmptyViews()

int EdbRunAccess::CheckEmptyViews ( EdbPattern pat)

794{
795 int nempty=0;
796 EdbSegP *s=0;
797 for(int i=0; i<pat.N(); i++) {
798 s = pat.GetSegment(i);
799 if( s->Vid(0)>0 ) continue;
800 if( s->Vid(1)>0 ) continue;
801 nempty++;
802 }
803 return nempty;
804}

◆ CheckMeanSegsPerView()

float EdbRunAccess::CheckMeanSegsPerView ( EdbPattern pat)

766{
767 Long_t nseg=0, nview=pat.N();
768 if(!nview) return 0;
769 for(int i=0; i<nview; i++) nseg += pat.GetSegment(i)->Flag();
770 return 1.*nseg/nview;
771}

◆ CheckRunLine()

void EdbRunAccess::CheckRunLine ( )

710{
711 if(!eRun) return;
712
713 cerr << endl << eRunFileName.Data();
714 int nviews=eRun->GetEntries();
715 cerr<<"\tv: "<< nviews;
716 eNareas=0;
717 if(eRun->GetHeader()) {
719 } else
720 cerr <<" RUN HEADER is missing!";
721
722 cerr<<" a: "<<eNareas;
723 int var=0;
724 if(eNareas>0) {
725 var = nviews/eNareas/2;
726 cerr<<" v/a/s: "<<var;
727 }
728 cerr<<endl;
729
730 if(eNareas>0)
731 if(var*2*eNareas != nviews) {
732 cerr<<"\tWARNING: areas and views numbers are mismatching!\n";
733 for(int ud=0; ud<3; ud++) {
734 EdbPattern *pat=GetVP(ud);
735 if(pat) {
736 int nv = pat->N();
737 cerr<<"\t\t views("<<ud<<"): "<<nv;
738 if(nv>0) {
739 cerr<<" Areas in range: ";
740 cerr<<pat->GetSegment(0)->Aid(0)<<" : "<<pat->GetSegment(nv-1)->Aid(0);
741 cerr<<" Views in range: ";
742 cerr<<pat->GetSegment(0)->Aid(1)<<" : "<<pat->GetSegment(nv-1)->Aid(1);
743 }
744 cerr<<endl;
745 }
746 }
747 }
748
749 for(int ud=0; ud<3; ud++) {
750 EdbPattern *pat=GetVP(ud);
751 int nempty=0;
752 if(pat) {
753 nempty = CheckEmptyViews(*pat);
754 if(nempty>0) cerr<<"\tWARNING: "<<nempty<<" empty views in layer "<<ud<<endl;
755 float meanseg = CheckMeanSegsPerView(*pat);
756 printf("Mean Segments/view = %10.1f\n",meanseg);
757 int thres=1;
758 int n0v = Check0Views(*pat, thres);
759 if(n0v) printf("%d views with <%d segments!\n",n0v,thres);
760 }
761 }
762}
Definition: EdbPattern.h:273
int CheckEmptyViews(EdbPattern &pat)
Definition: EdbRunAccess.cxx:793
EdbPattern * GetVP(int ud) const
Definition: EdbRunAccess.h:117
int Check0Views(EdbPattern &pat, int thres=1)
Definition: EdbRunAccess.cxx:773
float CheckMeanSegsPerView(EdbPattern &pat)
Definition: EdbRunAccess.cxx:765
Int_t eNareas
Definition: EdbRunAccess.h:57
int GetEntries() const
Definition: EdbRun.h:136
Int_t Aid(int i) const
Definition: EdbSegP.h:169

◆ CheckUpDownOffsets()

TH2F * EdbRunAccess::CheckUpDownOffsets ( )

1247{
1248 TH2F *h = new TH2F("UpDownOffsets","dXdY offsets between Up and correspondent Down views",100,-10,10,100,-10,10 );
1249 EdbPattern *p1 = GetVP(1); if(!p1) { Log(1,"EdbRunAccess::CheckUpDownOffsets","Error: GetVP(1) empty"); return h;}
1250 EdbPattern *p2 = GetVP(2); if(!p1) { Log(1,"EdbRunAccess::CheckUpDownOffsets","Error: GetVP(2) empty"); return h;}
1251 int N1 = p1->N(), N2 = p2->N();
1252 if(N1 != N2) Log(1,"EdbRunAccess::CheckUpDownOffsets","Warning: different number of views: %d -up, %d - down", p1->N(), p2->N() );
1253 int n = TMath::Min( N1, N2 );
1254
1255//segP.SetAid( view->GetAreaID() , view->GetViewID() );
1256
1257 TArrayL viewid1(N1), viewid2(N2); //view/area id to sincronize the views
1258 TArrayI ind1(N1), ind2(N2);
1259 for(int i=0; i<N1; i++) viewid1[i] = p1->GetSegment(i)->Aid(0) * 1000000 + p1->GetSegment(i)->Aid(1);
1260 for(int i=0; i<N2; i++) viewid2[i] = p2->GetSegment(i)->Aid(0) * 1000000 + p2->GetSegment(i)->Aid(1);
1261 TMath::Sort(N1,viewid1.GetArray(),ind1.GetArray(),0);
1262 TMath::Sort(N2,viewid2.GetArray(),ind2.GetArray(),0);
1263
1264 for(int i=0; i<n; i++) {
1265 EdbSegP *s1=p1->GetSegment(ind1[i]);
1266 EdbSegP *s2=p2->GetSegment(ind2[i]);
1267 h->Fill( s2->X()-s1->X(), s2->Y()-s1->Y() );
1268 }
1269 float sx = h->GetRMS(1);
1270 float sy = h->GetRMS(2);
1271 if(sx>0.5 || sy>0.5)
1272 Log(1,"EdbRunAccess::CheckUpDownOffsets","WARNING!!! XY stage possible malfunction (rms>0.5 micron) rmsX = %7.4f rmsY = %7.4f (by %d views)",
1273 h->GetRMS(1), h->GetRMS(2), n );
1274 else Log(2,"EdbRunAccess::CheckUpDownOffsets","rmsX = %7.4f rmsY = %7.4f (by %d views)", sx, sy, n );
1275 return h;
1276}
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ CheckViewSize()

void EdbRunAccess::CheckViewSize ( )

1234{
1235 for(int ud=1; ud<3; ud++ ) {
1236 eViewXmin[ud] = eHViewXY[ud].XminA(0.1);
1237 eViewXmax[ud] = eHViewXY[ud].XmaxA(0.1);
1238 eViewYmin[ud] = eHViewXY[ud].YminA(0.1);
1239 eViewYmax[ud] = eHViewXY[ud].YmaxA(0.1);
1240 Log(2,"EdbRunAccess::CheckViewSize"," side%d X:[%f %f] Y:[%f %f]", ud,
1241 eViewXmin[ud],eViewXmax[ud],eViewYmin[ud],eViewYmax[ud] );
1242 }
1243}
float XmaxA(float level=0)
Definition: EdbCell2.cpp:162
float XminA(float level=0)
Definition: EdbCell2.cpp:157
float YmaxA(float level=0)
Definition: EdbCell2.cpp:172
float YminA(float level=0)
Definition: EdbCell2.cpp:167
Float_t eViewXmin[3]
Definition: EdbRunAccess.h:72
Float_t eViewYmin[3]
Definition: EdbRunAccess.h:73
Float_t eViewXmax[3]
Definition: EdbRunAccess.h:72
Float_t eViewYmax[3]
Definition: EdbRunAccess.h:73

◆ CheckViewStep() [1/2]

void EdbRunAccess::CheckViewStep ( )
808{
809 CheckViewStep(1);
810 CheckViewStep(2);
811}
void CheckViewStep()
Definition: EdbRunAccess.cxx:807

◆ CheckViewStep() [2/2]

void EdbRunAccess::CheckViewStep ( int  ud)

assumed that in the tree (and in eVP) views are consecutive in the order of acquisition

815{
817 if(ud<0 || ud>2) return;
818 EdbPattern *pat=GetVP(ud); if(!pat) return;
819
820 //TFile f("test.root","RECREATE");
821 TH1F hx("viewXstep","viewXstep",3000,200,500);
822 TH1F hy("viewYstep","viewYstep",3000,200,500);
823 int n=pat->N(); if(n<2) return;
824 for( int i=0; i<n-1; i++ ) {
825 EdbSegP *s1 = pat->GetSegment(i);
826 EdbSegP *s2 = pat->GetSegment(i+1);
827 if( s1->Aid(0) != s1->Aid(0) ) continue; // skip different areas
828 if( s1->Aid(1) == s2->Aid(1) ) continue; // skip same view
829 float dx = Abs(s2->X()-s1->X());
830 float dy = Abs(s2->Y()-s1->Y());
831 if(dx>dy)
832 if(dx>50&&dx<500&&dy<50) //check that steps are in reasonable limits
833 hx.Fill(dx);
834 if(dy>dx)
835 if(dy>50&&dy<500&&dx<50) //check that steps are in reasonable limits
836 hy.Fill(dy);
837 }
838 //hx.Write(); hy.Write(); f.Close();
839 eXstep[ud] = hx.GetMean();
840 eYstep[ud] = hy.GetMean();
841 Log(2,"EdbRunAccess::CheckViewStep","side %d: stepX(%d) = %f +- %f stepY(%d) = %f +- %f",
842 ud, (int)(hx.GetEntries()), hx.GetMean(),hx.GetRMS(), (int)(hy.GetEntries()), hy.GetMean(),hy.GetRMS() );
843}
Float_t eXstep[3]
Definition: EdbRunAccess.h:70
Float_t eYstep[3]
Definition: EdbRunAccess.h:71

◆ ClearCuts()

void EdbRunAccess::ClearCuts ( )

96{
97 for(int i=0; i<3; i++) {
98 if(eCuts[i]) eCuts[i]->Clear();
99 }
100}

◆ CopyRawDataXY()

bool EdbRunAccess::CopyRawDataXY ( float  x0,
float  y0,
float  dR,
const char *  file 
)

1214{
1215 Log(2,"CopyRawDataXY", "dR = %f around (%f %f) into file %s",dR, x0,y0, file);
1216 if(!eRun) return false;
1217 EdbRun r(*eRun, file);
1218 EdbView *view = eRun->GetView();
1219 for(int side=1; side<=2; side++) {
1220 TArrayI entr(100000);
1221 int nv = GetViewsXY(side,entr,x0,y0, dR);
1222 for(int i=0; i<nv; i++) {
1223 view = eRun->GetEntry(entr[i]);
1224 r.AddView(view);
1225 }
1226 Log(2,"CopyRawDataXY", "%d views copied for side %d",nv, side);
1227 }
1228 r.Close();
1229 return true;
1230}
int GetViewsXY(int ud, TArrayI &entr, float x, float y, float r=200.)
Definition: EdbRunAccess.cxx:283
Definition: EdbRun.h:75
EdbView * GetView() const
Definition: EdbRun.h:110
EdbView * GetEntry(int entry, int ih=1, int icl=0, int iseg=1, int itr=0, int ifr=0)
Definition: EdbRun.cxx:462
Base scanning data object: entry into Run tree.
Definition: EdbView.h:134
TFile * file
Definition: write_pvr.C:3

◆ FillDZMaps()

void EdbRunAccess::FillDZMaps ( )

fill DZ maps

1281{
1283
1284 int eNsegMin=100;
1285
1286 Log(2,"EdbRunAccess::FillDZmaps","Get data from %s with AFID=%d", eRunFileName.Data(),eAFID);
1287 if(!eRun) { Log(1,"EdbRunAccess::FillDZMaps","ERROR: run is not opened\n"); return; }
1288 EdbView *view = eRun->GetView();
1289 EdbViewHeader *head = view->GetHeader();
1290 int nentries = eRun->GetEntries();
1291
1292 Double_t *x1 = new Double_t[nentries];
1293 Double_t *y1 = new Double_t[nentries];
1294 Double_t *x2 = new Double_t[nentries];
1295 Double_t *y2 = new Double_t[nentries];
1296 Double_t *z1 = new Double_t[nentries];
1297 Double_t *z2 = new Double_t[nentries];
1298 Double_t *z3 = new Double_t[nentries];
1299 Double_t *z4 = new Double_t[nentries];
1300 Int_t cnt_good=0;
1301
1302 int cnt1=0, cnt2=0;
1303 for(int ie=0; ie<nentries; ie++ ) {
1304 view = eRun->GetEntry(ie,1,0,0,0,0);
1305 if( head->GetNsegments() < eNsegMin ) continue;
1306 int side = ViewSide(view);
1307 if( side == 1 ) {
1308 x1[cnt1] = ViewX(view);
1309 y1[cnt1] = ViewY(view);
1310 z1[cnt1] = view->GetZ1();
1311 z2[cnt1] = view->GetZ2();
1312 cnt1++;
1313 }
1314 else if( side == 2 ) {
1315 x2[cnt2] = ViewX(view);
1316 y2[cnt2] = ViewY(view);
1317 z3[cnt2] = view->GetZ3();
1318 z4[cnt2] = view->GetZ4();
1319 cnt2++;
1320 }
1321 }
1322
1323 for(int i=0; i<4; i++) if( eGraphZ[i] ) delete eGraphZ[i];
1324 for(int i=0; i<3; i++) if( eGraphDZ[i] ) delete eGraphDZ[i];
1325
1326 eGraphZ[0] = new TGraph2D("graphZ1","graphZ1",cnt1,x1,y1,z1);
1327 eGraphZ[1] = new TGraph2D("graphZ2","graphZ2",cnt1,x1,y1,z2);
1328 eGraphZ[2] = new TGraph2D("graphZ3","graphZ3",cnt2,x2,y2,z3);
1329 eGraphZ[3] = new TGraph2D("graphZ4","graphZ4",cnt2,x2,y2,z4);
1330
1331 TFile f("dz.root","RECREATE");
1332 for(int i=0; i<4; i++) if( eGraphZ[i] ) eGraphZ[i]->Write();
1333 f.Close();
1334}
float ViewX(const EdbView *view) const
Definition: EdbRunAccess.h:145
TGraph2D * eGraphDZ[3]
keep the base/layer1/layer2 thickness calculated using eZ1/eZ2/eZ3/eZ4 for each view
Definition: EdbRunAccess.h:47
TGraph2D * eGraphZ[4]
keep z1/z2/z3/z4 surfaces using eZ1/eZ2/eZ3/eZ4 for each view
Definition: EdbRunAccess.h:46
float ViewY(const EdbView *view) const
Definition: EdbRunAccess.h:146
view identification
Definition: EdbView.h:26
Int_t GetNsegments() const
Definition: EdbView.h:88
Float_t GetZ3() const
Definition: EdbView.h:198
Float_t GetZ1() const
Definition: EdbView.h:196
Float_t GetZ4() const
Definition: EdbView.h:199
Float_t GetZ2() const
Definition: EdbView.h:197
int nentries
Definition: check_shower.C:40

◆ FillVP()

int EdbRunAccess::FillVP ( )

fill patterns up/down with dummy segments with center of view coordinates
to profit of the pattern tools for the fast views selection
the convention is:
segp->ID() - tree entry number
segp->X() - coordinates in the external (brick) RS (after transformation)
segp->Y() - /
segp->Aid(areaID,viewID);

847{
855
856
857 Log(2,"EdbRunAccess::FillVP","Get data from %s with AFID=%d",eRunFileName.Data(),eAFID);
858
859 if(!eRun) { Log(1,"EdbRunAccess::FillVP","ERROR: run is not opened\n"); return 0; }
860 EdbView *view = eRun->GetView();
861 EdbViewHeader *head = view->GetHeader();
862
863 int nentr = eRun->GetEntries();
864
865 eRun->GetTree()->Draw(">>lst", eHeaderCut );
866 TEventList *lst = (TEventList*)(gDirectory->GetList()->FindObject("lst"));
867 if(!lst) {Log(1,"EdbRunAccess::FillVP","ERROR!: events list (lst) did not found!"); return 0;}
868 nentr =lst->GetN();
869
870
871 if( eVP[1] ) delete eVP[1];
872 if( eVP[2] ) delete eVP[2];
873 eVP[1] = new EdbPattern( 0.,0., 0., nentr);
874 eVP[2] = new EdbPattern( 0.,0., 200., nentr); //TODO
875
876 EdbSegP segP;
877
878 eFirstArea = 9999999;
879 eLastArea = 0;
880
881 int side;
882 int nseg,ncl;
883 for(int ie=0; ie<nentr; ie++ ) {
884 int iv = lst->GetEntry(ie);
885 view = eRun->GetEntry(iv,1,0,0,0,0);
886
887 nseg = head->GetNsegments();
888 ncl = head->GetNclusters();
889 if(eAFID==0){
890 segP.Set( iv,view->GetXview(),view->GetYview(), 0,0,ncl,nseg);
891 }
892 else if(eAFID==1 || eAFID==2) {
893 segP.Set( iv,0,0, 0,0,ncl,nseg);
894 segP.Transform( head->GetAffine() );
895 }
896 else if(eAFID==11) {
897 segP.Set( iv,0,0, 0,0,ncl,nseg);
898 if(!eAffStage2Abs) Log(1,"","eAffStage2Abs =0");
899 segP.Transform( eAffStage2Abs );
900 }
901 else if(eAFID==20) { //to remove
902 segP.Set( iv,0,0, 0,0,ncl,nseg);
903 EdbAffine2D *aff = (EdbAffine2D*)eViewCorr->At(iv);
904 if(aff) {
905 segP.Transform( aff );
906 head->GetAffine()->Print();
907 aff->Print();
908 } else Log(1,"EdbRunAccess::FillVP", "WARNING! aff is missing for entry %d",iv);
909 }
910 if( view->GetAreaID() < eFirstArea ) eFirstArea= view->GetAreaID();
911 if( view->GetAreaID() > eLastArea ) eLastArea = view->GetAreaID();
912
913 segP.SetAid( view->GetAreaID() , view->GetViewID() );
914 segP.SetVid( head->GetNframesTop(), head->GetNframesBot() );
915
916 side = ViewSide(view);
917 if(side<1||side>2) continue;
918 eVP[side]->AddSegment(segP);
919 }
920 Log(2,"EdbRunAccess::FillVP"," %d entries in limits X=(%f : %f) Y=(%f %f)\n", nentr,
921 eVP[1]->Xmin(), eVP[1]->Xmax(),
922 eVP[1]->Ymin(), eVP[1]->Ymax() );
923 return nentr;
924}
Int_t eFirstArea
Definition: EdbRunAccess.h:55
Int_t eLastArea
Definition: EdbRunAccess.h:56
TCut eHeaderCut
header cut to be applied in run initialization
Definition: EdbRunAccess.h:43
TTree * GetTree() const
Definition: EdbRun.h:113
EdbSegP * AddSegment(int i, EdbSegP &s)
Definition: EdbPattern.cxx:72
virtual void Transform(const EdbAffine2D *a)
Int_t GetNclusters() const
Definition: EdbView.h:87
Int_t GetNframesTop() const
Definition: EdbView.h:120
Int_t GetNframesBot() const
Definition: EdbView.h:121

◆ FirstArea()

int EdbRunAccess::FirstArea ( ) const
inline
119{ return eFirstArea;}

◆ GetCond()

EdbScanCond * EdbRunAccess::GetCond ( int  ud)
inline
176 { if(eCond[ud]) return (EdbScanCond *)eCond[ud]; else return 0; }
Definition: EdbScanCond.h:10

◆ GetCut()

EdbSegmentCut * EdbRunAccess::GetCut ( int  ud,
int  i 
)
inline
164 { return (EdbSegmentCut *)(eCuts[ud]->UncheckedAt(i)); }

◆ GetCutGR()

float EdbRunAccess::GetCutGR ( ) const
inline
165{return eCutGR;}
Float_t eCutGR
grain cut (chi)
Definition: EdbRunAccess.h:66

◆ GetEntryXY()

int EdbRunAccess::GetEntryXY ( int  ud,
float  x,
float  y 
)

305{
306 // find view with the center closest to x,y and return it's entry
307
308 int entry = -1;
309 EdbPattern *pat=GetVP(ud);
310 if(!pat) return -1;
311 EdbSegP *seg=0;
312 float r,rmin=9999999.;
313 int imin=-1;
314 for( int i=0; i<pat->N(); i++ ) {
315 seg = pat->GetSegment(i);
316 r = Sqrt( (seg->X()-x)*(seg->X()-x)+(seg->Y()-y)*(seg->Y()-y));
317 if( r>rmin ) continue;
318 rmin = r;
319 entry = seg->ID();
320 imin = i;
321 }
322 //if(imin>0) pat->GetSegment(imin)->Print();
323 return entry;
324}
Int_t ID() const
Definition: EdbSegP.h:147

◆ GetLayer()

EdbLayer * EdbRunAccess::GetLayer ( int  id)
inline
114 { if(eLayers[id]) return (EdbLayer *)eLayers[id]; else return 0; }
UInt_t id
Definition: tlg2couples.C:117

◆ GetMakeCond()

EdbScanCond * EdbRunAccess::GetMakeCond ( int  ud)
inline
174{return 0;} //TODO?

◆ GetMakeLayer()

EdbLayer * EdbRunAccess::GetMakeLayer ( int  id)

276{
277 if(id<0) return 0; if(id>2) return 0;
278 if(!GetLayer(id)) eLayers[id] = new EdbLayer();
279 return GetLayer(id);
280}

◆ GetNareas()

Int_t EdbRunAccess::GetNareas ( )
inline
99{ return eNareas; }

◆ GetNviewsPerArea()

Int_t EdbRunAccess::GetNviewsPerArea ( )
inline
100{ return eNviewsPerArea; }
Int_t eNviewsPerArea
Definition: EdbRunAccess.h:58

◆ GetPatternData()

int EdbRunAccess::GetPatternData ( EdbPattern pat,
int  side,
int  nviews,
TArrayI &  srt,
int &  nrej 
)

506{
507 // get raw segments as a pattern for the given array of entries and given side
508
509 EdbSegP segP;
510 EdbView *view = eRun->GetView();
511 int nseg=0;
512 nrej=0;
513 int entry;
514 int nsegV;
515
516 for(int iu=0; iu<nviews; iu++) {
517 entry = srt[iu];
518 if(eCLUST) {
519 view = eRun->GetEntry(entry,1,1,1);
521 }
522 else view = eRun->GetEntry(entry);
523
524 nsegV = view->Nsegments();
525 if( ViewSide(view) != side ) continue;
526
527 for(int j=0;j<nsegV;j++) {
528 if(!AcceptRawSegment(view,j,segP,side,entry)) {
529 nrej++;
530 continue;
531 }
532 nseg++;
533 pat.AddSegment( segP);
534 }
535 }
536 return nseg;
537}
bool AcceptRawSegment(EdbView *view, int ud, EdbSegP &segP, int side, int entry)
Definition: EdbRunAccess.cxx:943
int AttachClustersToSegments()
Definition: EdbView.cxx:360
Int_t Nsegments() const
Definition: EdbView.h:216

◆ GetPatternDataForPrediction()

int EdbRunAccess::GetPatternDataForPrediction ( int  id,
int  side,
EdbPattern pat 
)

get raw segments belonging to views having a given id (view->GetHeader()->GetTrack()) and a given side
the result is stored into an EdbPattern object
this routine is called by EdbRunTracking::FindCandidates and to select microtracks from a raw.root file
acquired for a specific prediction
Special case: if id=-1 - accept all views for the given side

541{
547
548 EdbSegP segP;
549 EdbView *view = eRun->GetView();
550 //int nviews = eRun->GetEntries();
551 int nviews = eVP[side]->N();
552 Log(2,"EdbRunAccess::GetPatternDataForPrediction","%d views are selected for side %d",nviews,side);
553 int nseg=0;
554 int nsegV;
555
556 for(int ie=0; ie<nviews; ie++) {
557 int entry=eVP[side]->GetSegment(ie)->ID();
558 if(eCLUST) {
559 view = eRun->GetEntry(entry,1,1,1);
561 }
562 else view = eRun->GetEntry(entry);
563
564 if( id>=0 && view->GetHeader()->GetTrack()!=id) continue;
565
566 nsegV = view->Nsegments();
567 if( ViewSide(view) != side ) continue;
568 //Log(2,"EdbRunAccess::GetPatternDataForPrediction","ie = %d nsegV = %d",ie,nsegV);
569
570 for(int j=0;j<nsegV;j++) {
571 if(!AcceptRawSegment(view,j,segP,side,entry)) continue;
572 nseg++;
573 segP.SetScanID(pat.ScanID());
574 segP.SetSide(pat.Side());
575 pat.AddSegment( segP);
576 }
577 }
578 Log(2,"GetPatternDataForPrediction","%d segments for pred %d, side %d from %s",nseg,id,side, eRunFileName.Data());
579 return nseg;
580}
EdbID ScanID() const
Definition: EdbPattern.h:326
Int_t Side() const
Definition: EdbPattern.h:328
void SetSide(int side=0)
Definition: EdbSegP.h:139
void SetScanID(EdbID id)
Definition: EdbSegP.h:143
Int_t GetTrack() const
Definition: EdbView.h:126

◆ GetPatternXY()

int EdbRunAccess::GetPatternXY ( EdbSegP s,
int  side,
EdbPattern pat,
float  rmin = 200 
)

378{
379 // assuming s as the prediction in plate RS (also Z), fill vol with the segments
380 // of closest up and down views
381 // todo: Z of layers must be setted!
382
383 if(side<1||side>2) return 0;
384 TArrayI entr(100000);
385 float x = s.X()+ (eLayers[side]->Z()-s.Z())*s.TX();
386 float y = s.Y()+ (eLayers[side]->Z()-s.Z())*s.TY();
387 int nv = GetViewsXY(side,entr,x,y, Max(rmin,(float)300.) );
388 pat.SetZ(eLayers[side]->Z());
389 int nrej=0;
390 return GetPatternData(pat,side,nv,entr,nrej);
391}
int GetPatternData(EdbPattern &pat, int side, int nviews, TArrayI &srt, int &nrej)
Definition: EdbRunAccess.cxx:502
void SetZ(float z)
Definition: EdbPattern.h:41
Double_t Z
Definition: tlg2couples.C:104

◆ GetPatternXYcut()

int EdbRunAccess::GetPatternXYcut ( EdbSegP s0,
int  side,
EdbPattern pat,
float  dr,
float  dt 
)

352{
353 Log(2,"EdbRunAccess::GetPatternXYcut","side %d [%f %f %f %f] +-%5.1f +-%5.3f:", side, s0.X(), s0.Y(), s0.TX(),s0.TY(),dr,dt );
354
355 EdbPattern p0;
356 int nseg = GetPatternXY(s0, side, p0, dr);
357 Log(2,"EdbRunAccess::GetPatternXYcut","nseg = %d",nseg);
358 Log(2,"EdbRunAccess::GetPatternXYcut","side = %d zlayer=%f zseg= %f",side,eLayers[side]->Z(),s0.Z() );
359 float x = s0.X()+ (eLayers[side]->Z()-s0.Z())*s0.TX();
360 float y = s0.Y()+ (eLayers[side]->Z()-s0.Z())*s0.TY();
361 int ic=0;
362 for(int i=0; i<nseg; i++)
363 {
364 EdbSegP *s = p0.GetSegment(i);
365 if( Abs(s->X()-x) > dr) continue;
366 if( Abs(s->Y()-y) > dr) continue;
367 if( Abs(s->TX()-s0.TX()) > dt) continue;
368 if( Abs(s->TY()-s0.TY()) > dt) continue;
369 pat.AddSegment(*s);
370 ic++;
371 }
372 Log(2,"EdbRunAccess::GetPatternXYcut","selected = %d",ic);
373 return ic;
374}
int GetPatternXY(EdbSegP &s, int side, EdbPattern &pat, float rmin=200)
Definition: EdbRunAccess.cxx:377
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t Z() const
Definition: EdbSegP.h:153
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176

◆ GetRawSegment() [1/2]

EdbSegment * EdbRunAccess::GetRawSegment ( EdbView view,
int  sid,
int  rs = 0 
)

1137{
1138 EdbSegment *s = view.GetSegment(sid);
1139 ApplyCorrections( view, *s, rs );
1140 return s;
1141}
void ApplyCorrections(const EdbView &view, EdbSegment &s, const int rs)
Definition: EdbRunAccess.cxx:1105

◆ GetRawSegment() [2/2]

EdbSegment * EdbRunAccess::GetRawSegment ( int  vid,
int  sid,
int  rs = 0 
)

1145{
1146 EdbView *view = eRun->GetEntry(vid);
1147 return GetRawSegment( *view, sid, rs );
1148}
EdbSegment * GetRawSegment(int vid, int sid, int rs=0)
Definition: EdbRunAccess.cxx:1144

◆ GetRawSegmentN()

EdbSegment * EdbRunAccess::GetRawSegmentN ( int  vid,
int  sid,
int  rs = 0 
)

1097{
1098 EdbView *view = eRun->GetEntry(vid);
1099 EdbSegment *s = view->GetSegment(sid);
1100// printf("(%d %d) %d %f %f\n",vid,sid, s->GetPuls(), view->GetXview(),view->GetYview());
1101 return s;
1102}

◆ GetRawSegmentPix()

float EdbRunAccess::GetRawSegmentPix ( EdbSegment seg)

assumed that clusters are attached to segments

1152{
1154 float pix=0;
1155 EdbCluster *cl=0;
1156 TObjArray *clusters = seg->GetElements();
1157 if(!clusters) return 0;
1158 int ncl = clusters->GetLast()+1;
1159 for(int i=0; i<ncl; i++ ) {
1160 cl = (EdbCluster*)clusters->At(i);
1161 pix += cl->GetArea();
1162 }
1163 return pix;
1164}
Float_t GetArea() const
Definition: EdbCluster.h:54

◆ GetRun()

EdbRun * EdbRunAccess::GetRun ( ) const
inline
111{return eRun;}

◆ GetViewsArea() [1/2]

int EdbRunAccess::GetViewsArea ( int  ud,
TArrayI &  entr,
float  xmin,
float  xmax,
float  ymin,
float  ymax 
)

get all views in a given limits

614{
616
617 int nv=0;
618 if(ud<0) return nv; if(ud>2) return nv;
619
620 EdbPattern *pat=GetVP(ud);
621 if(!pat) return nv;
622 EdbSegP *seg=0;
623 for( int i=0; i<pat->N(); i++ ) {
624 seg = pat->GetSegment(i);
625 if(seg->X()<xmin) continue;
626 if(seg->X()>xmax) continue;
627 if(seg->Y()<ymin) continue;
628 if(seg->Y()>ymax) continue;
629 entr[nv++]=seg->ID();
630 }
631 return nv;
632}

◆ GetViewsArea() [2/2]

int EdbRunAccess::GetViewsArea ( int  ud,
TArrayI &  entr,
int  area,
float &  xmin,
float &  xmax,
float &  ymin,
float &  ymax 
)

get all views of "area" for the side "ud"; fill entr with the entries, calculate area limits

587{
590
591 int nv=0;
592 if(ud<0) return nv; if(ud>2) return nv;
593
594 EdbPattern *pat=GetVP(ud);
595 if(!pat) return nv;
596 xmin=eXmax; xmax=eXmin;
597 ymin=eYmax; ymax=eYmin;
598 EdbSegP *seg=0;
599 for( int i=0; i<pat->N(); i++ ) {
600 seg = pat->GetSegment(i);
601 if( seg->Aid(0) != area ) continue;
602 if(xmin>seg->X()) xmin=seg->X();
603 if(xmax<seg->X()) xmax=seg->X();
604 if(ymin>seg->Y()) ymin=seg->Y();
605 if(ymax<seg->Y()) ymax=seg->Y();
606 entr[nv++]=seg->ID();
607 }
608 return nv;
609}
Float_t eYmin
Definition: EdbRunAccess.h:68
Float_t eXmin
Definition: EdbRunAccess.h:68
Float_t eYmax
Definition: EdbRunAccess.h:68
Float_t eXmax
Definition: EdbRunAccess.h:68
Double_t X
Definition: tlg2couples.C:76
Double_t Y
Definition: tlg2couples.C:76

◆ GetViewsAreaMarg()

int EdbRunAccess::GetViewsAreaMarg ( int  ud,
TArrayI &  entr,
int  area,
float  xmarg,
float  ymarg 
)

return area views with margins
ud: 0-base, 1-down, 2-up //TODO - to check this convention!

637{
640
641 int nv=0;
642 if(ud<0 || ud>2) return nv;
643 EdbPattern *pat=GetVP(ud);
644 if(!pat) return nv;
645
646 EdbPattern ptmp;
647 EdbSegP *seg=0;
648 int N=pat->N();
649 for( int i=0; i<N; i++ ) {
650 seg = pat->GetSegment(i);
651 if( seg->Aid(0) == area ) ptmp.AddSegment( *seg );
652 }
653
654 float xmin = ptmp.Xmin();
655 float xmax = ptmp.Xmax();
656 float ymin = ptmp.Ymin();
657 float ymax = ptmp.Ymax();
658
659 for( int i=0; i<N; i++ ) {
660 seg = pat->GetSegment(i);
661 if( seg->X() < xmin-xmarg ) continue;
662 if( seg->X() > xmax+xmarg ) continue;
663 if( seg->Y() < ymin-ymarg ) continue;
664 if( seg->Y() > ymax+ymarg ) continue;
665 entr[nv++] = seg->ID();
666 }
667
668 return nv;
669}
virtual Float_t Xmax() const
Definition: EdbVirtual.cxx:196
virtual Float_t Ymin() const
Definition: EdbVirtual.cxx:206
virtual Float_t Xmin() const
Definition: EdbVirtual.cxx:186
virtual Float_t Ymax() const
Definition: EdbVirtual.cxx:216

◆ GetViewsXY()

int EdbRunAccess::GetViewsXY ( int  ud,
TArrayI &  entr,
float  x,
float  y,
float  rv = 200. 
)

285{
286 // return all views in rv-surrounds of the x-y point
287 int nv=0;
288 if(ud<0||ud>2) return nv;
289 EdbPattern *pat=GetVP(ud);
290 if(!pat) return nv;
291 EdbSegP *seg=0;
292 float r;
293 for( int i=0; i<pat->N(); i++ ) {
294 seg = pat->GetSegment(i);
295 r = Sqrt( (seg->X()-x)*(seg->X()-x)+(seg->Y()-y)*(seg->Y()-y));
296 if( r>rv ) continue;
297 entr[nv++] = seg->ID();
298 }
299 return nv;
300}

◆ GetVolumeArea()

int EdbRunAccess::GetVolumeArea ( EdbPatternsVolume vol,
int  area 
)

396{
397 // get raw segments as a volume for the given "area"
398
399 int maxA=10000; //TODO
400 TArrayI entr1(maxA);
401 float xmin,xmax,ymin,ymax;
402 int n1 = GetViewsArea( 1, entr1, area, xmin,xmax,ymin,ymax);
403 TArrayI entr2(maxA);
404 xmin -= eXstep[1]+eXstep[1]/10.;
405 xmax += eXstep[1]+eXstep[1]/10.;
406 ymin -= eYstep[1]+eYstep[1]/10.;
407 ymax += eYstep[1]+eYstep[1]/10.;
408 int n2 = GetViewsArea( 2, entr2, xmin,xmax,ymin,ymax );
409 // int n2 = GetViewsAreaMarg( 2, entr2, area, eXstep+eXstep/10.,eYstep+eYstep/10. );
410
411 printf("n1=%d n2=%d\n",n1,n2);
412
413 TArrayI eall(n1+n2);
414 int ic=0;
415 for(int i=0; i<n1; i++) eall[ic++]=entr1[i];
416 for(int i=0; i<n2; i++) eall[ic++]=entr2[i];
417
418 TArrayI ind(ic);
419 TArrayI srt(ic);
420 TMath::Sort(ic,eall.GetArray(),ind.GetArray(),0);
421 for (int i=0; i<ic; i++) srt[i] = eall[ind[i]];
422
423 int nrej=-1;
424 int nseg = GetVolumeData(vol, ic, srt, nrej);
425
426#ifdef _WINDOWS
427 Log(2,"GetVolumeArea","Area:%3d (%3.0f %%) views:%4d/%4d %6d +%6d segments are accepted; %6d - rejected"
428#else
429 Log(2,"GetVolumeArea","Area:%3d (%3.0f \%%) views:%4d/%4d %6d +%6d segments are accepted; %6d - rejected"
430#endif
431 ,area, 100.*area/eNareas
432 ,n1,n2
433 ,vol.GetPattern(0)->N()
434 ,vol.GetPattern(1)->N()
435 ,nrej );
436
437 if(nseg != vol.GetPattern(0)->N()+vol.GetPattern(1)->N())
438 Log(2,"GetVolumeArea","WARNING!!! nseg = %d != %d\n", nseg,vol.GetPattern(0)->N()+vol.GetPattern(1)->N());
439 return ic;
440}
EdbPatternsVolume * vol
Definition: RecDispNU.C:116
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
int GetViewsArea(int ud, TArrayI &entr, int area, float &xmin, float &xmax, float &ymin, float &ymax)
Definition: EdbRunAccess.cxx:585
int GetVolumeData(EdbPatternsVolume &vol, int nviews, TArrayI &srt, int &nrej)
Definition: EdbRunAccess.cxx:454

◆ GetVolumeData()

int EdbRunAccess::GetVolumeData ( EdbPatternsVolume vol,
int  nviews,
TArrayI &  srt,
int &  nrej 
)

458{
459 // get raw segments as a volume for the given array of entries
460
461 EdbSegP segP;
462 EdbView *view = eRun->GetView();
463 int side;
464 int nseg=0;
465 nrej=0;
466 int entry;
467 int nsegV;
468
469 for(int iu=0; iu<nviews; iu++) {
470 entry = srt[iu];
471 if(eCLUST) {
472 view = eRun->GetEntry(entry,1,1,1);
474 }
475 else view = eRun->GetEntry(entry);
476
477 nsegV = view->Nsegments();
478 side = ViewSide(view);
479
480 //printf(" EdbRunAccess::GetVolumeData: %d nsegV = %d side=%d\n",iu, nsegV,side);
481 if(side<1||side>2) continue;
482
483 //vol.Print();
484 for(int j=0;j<nsegV;j++) {
485 if(!AcceptRawSegment(view,j,segP,side,entry)) {
486 nrej++;
487 continue;
488 }
489 nseg++;
490 vol.GetPattern(side-1)->AddSegment( segP);
491 }
492 }
493
494 //printf("nseg: %d \t rejected: %d\n", nseg, nrej );
495
498 return nseg;
499}
void SetSegmentsZ()
Definition: EdbPattern.cxx:250

◆ GetVolumeXY()

int EdbRunAccess::GetVolumeXY ( EdbSegP s,
EdbPatternsVolume vol 
)

328{
329 // assuming s as the prediction in plate RS (also Z), fill vol with the segments
330 // of closest up and down views
331 // todo: Z of layers must be setted!
332
333 TArrayI entr(2);
334 float x,y;
335 x = s.X()+ (eLayers[1]->Z()-s.Z())*s.TX();
336 y = s.Y()+ (eLayers[1]->Z()-s.Z())*s.TY();
337 entr[0] = GetEntryXY(1,x,y);
338 //printf(":GetVolumeXY: 1 %f %f \n",x,y);
339
340 x = s.X()+ (eLayers[2]->Z()-s.Z())*s.TX();
341 y = s.Y()+ (eLayers[2]->Z()-s.Z())*s.TY();
342 entr[1] = GetEntryXY(2,x,y);
343 //printf(":GetVolumeXY: 2 %f %f \n",x,y);
344
345 //printf("EdbRunAccess::GetVolumeXY: entr = %d %d\n",entr[0],entr[1]);
346 int nrej=0;
347 return GetVolumeData(vol,2,entr,nrej);
348}
int GetEntryXY(int ud, float x, float y)
Definition: EdbRunAccess.cxx:304

◆ GetVP()

EdbPattern * EdbRunAccess::GetVP ( int  ud) const
inline
117{ if(ud>-1&&ud<3) return eVP[ud]; else return 0; }

◆ GuessNviewsPerArea()

void EdbRunAccess::GuessNviewsPerArea ( )

166{
167 eNviewsPerArea = 0;
168 int nviews=eRun->GetEntries();
169 if(!nviews) return;
170
171 int nareas=0;
172 int idmin = kMaxInt;
173 int idmax = kMinInt;
174 for(int iview=0; iview<nviews; iview++) {
175 int aid = eRun->GetEntry(iview,1,0,0,0,0)->GetAreaID();
176 if(aid>idmax)
177 {
178 idmax=aid;
179 nareas++;
180 }
181 if(aid<idmin) idmin=aid;
182 }
183 eNareas=nareas;
184 eFirstArea=idmin;
185 eLastArea=idmax;
186 eNviewsPerArea=nviews/eNareas/2;
187 if( eNviewsPerArea*eNareas*2 != nviews )
188 { // try if value in run header is correct
189 int nva=0;
190 if( eRun->GetHeader() && eRun->GetHeader()->GetArea() )
191 {
192 nva = eRun->GetHeader()->GetArea()->GetN();
193 if( nva>0 && nva*eNareas*2==nviews ) eNviewsPerArea=nva;
194 }
195 }
196 if(eNviewsPerArea*eNareas*2 != nviews)
197 Log(1,"EdbRunAccess::GuessNviewsPerArea","WARNING: eNviewsPerArea*eNareas*2 != nviews: %d*%d*2 != %d",
198 eNviewsPerArea, eNareas, nviews );
199}
Int_t GetN() const
Definition: EdbFiducial.h:160
EdbArea * GetArea() const
Definition: EdbRunHeader.h:139

◆ InitRun()

bool EdbRunAccess::InitRun ( const char *  runfile = 0,
bool  do_update = false 
)

113{
114 if(eRun) SafeDelete(eRun);
115 if(runfile) eRunFileName=runfile;
116 if(eAFID==2) ReadVAfile();
117 if( gSystem->AccessPathName(eRunFileName.Data(), kFileExists) ) {
118 Log(1,"EdbRunAccess::InitRun","ERROR! open file: %s\n",eRunFileName.Data());
119 return false;
120 }
121 if(do_update) eRun=new EdbRun(eRunFileName.Data(),"UPDATE");
122 else eRun=new EdbRun(eRunFileName.Data());
123 if(!eRun) {Log(1,"EdbRunAccess::InitRun","ERROR! can't open file %s\n",eRunFileName.Data()); return false; }
124 if(!(eRun->GetTree())) {Log(1,"EdbRunAccess::InitRun","ERROR! Vews tree is missed %s\n",eRunFileName.Data());return false;}
125 for(int i=0; i<3; i++) if(eLayers[i]==0) GetMakeLayer(i);
126 if(eAFID==11) {
127 eRun->GetMarks()->Print();
129 //if(eRun->GetMarks()) eAffStage2Abs = eRun->GetMarks()->Abs2Stage();
131 else Log(1,"EdbRunAccess::InitRun","ERROR! Stage2Abs do not available (eAFID=11)");
132 }
133 if(do_update) {
134 if(FillVP()<1) return false;
135 eVP[1]->Write("views_s1");
136 eVP[2]->Write("views_s2");
137 eRun->GetTree()->GetCurrentFile()->Purge();
138 } else {
139 if( eVP[1] ) delete eVP[1];
140 if( eVP[2] ) delete eVP[2];
141 eVP[1] = dynamic_cast<EdbPattern*>(gDirectory->Get("views_s1"));
142 eVP[2] = dynamic_cast<EdbPattern*>(gDirectory->Get("views_s2"));
143 }
144
145 int nn=0;
146 if( eVP[1] && eVP[2] ) nn = eVP[1]->N()+eVP[2]->N();
147 if(nn != eRun->GetEntries()) if(FillVP()<1) return false;
148
149 eVP[1]->Transform(eLayers[1]->GetAffineXY());
150 eVP[2]->Transform(eLayers[2]->GetAffineXY());
151 eXmin = eVP[1]->Xmin();
152 eXmax = eVP[1]->Xmax();
153 eYmin = eVP[1]->Ymin();
154 eYmax = eVP[1]->Ymax();
155
158 for(int i=0; i<3; i++) { eXstep[i]=400; eYstep[i]=400; }
159
160 Log(2,"EdbRunAccess::InitRun","%s AFID=%d",runfile, eAFID);
161 return true;
162}
EdbAffine2D * Stage2Abs() const
Definition: EdbFiducial.cxx:279
void Print(Option_t *opt="") const
Definition: EdbFiducial.cxx:593
virtual void Transform(const EdbAffine2D *a)
Definition: EdbVirtual.cxx:155
void ReadVAfile()
Definition: EdbRunAccess.cxx:202
int FillVP()
Definition: EdbRunAccess.cxx:846
EdbLayer * GetMakeLayer(int id)
Definition: EdbRunAccess.cxx:275
void GuessNviewsPerArea()
Definition: EdbRunAccess.cxx:165
void PrintStat()
Definition: EdbRunAccess.cxx:672
EdbMarksSet * GetMarks() const
Definition: EdbRun.h:120
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ InitRunFromRWC()

bool EdbRunAccess::InitRunFromRWC ( char *  rwcname,
bool  bAddRWD = true,
const char *  options = "" 
)

243{
244 if(!eRun)
245 if( gSystem->AccessPathName(rwcname, kFileExists) ) {
246 Log(1,"EdbRunAccess::InitRunFromRWC","ERROR open file: %s\n",rwcname);
247 return false;
248 }
249 eRun=new EdbRun(eRunFileName.Data(),"RECREATE");
250 if(!eRun) return false;
251
252 AddRWC(eRun, rwcname, bAddRWD, options);
253
254 return true;
255}
int AddRWC(EdbRun *run, char *rwcname, int bAddRWD, const char *options)
Definition: libDataConversion.cpp:64

◆ LastArea()

int EdbRunAccess::LastArea ( ) const
inline
120{ return eLastArea;}

◆ NCuts()

int EdbRunAccess::NCuts ( int  ud)
inline
151 {
152 if(!eCuts[ud]) return 0;
153 return eCuts[ud]->GetEntriesFast();
154 }

◆ OverlapX()

float EdbRunAccess::OverlapX ( int  ud)
inline
178{return (ud>0&&ud<3)? (eViewXmax[ud]-eViewXmin[ud] - eXstep[ud]): 0.; }

◆ OverlapY()

float EdbRunAccess::OverlapY ( int  ud)
inline
179{return (ud>0&&ud<3)? (eViewYmax[ud]-eViewYmin[ud] - eYstep[ud]): 0.; }

◆ PassCuts()

bool EdbRunAccess::PassCuts ( int  id,
EdbSegment seg 
)

1045{
1046 float var[5];
1047 var[0] = seg.GetX0();
1048 var[1] = seg.GetY0();
1049 var[2] = seg.GetTx();
1050 var[3] = seg.GetTy();
1051 var[4] = SegmentWeight(seg);
1052
1053 int nc = NCuts(id);
1054 for(int i=0; i<nc; i++) {
1055 if( !(GetCut(id,i)->PassCut(var)) ) return false;
1056 }
1057 return true;
1058}
int NCuts(int ud)
Definition: EdbRunAccess.h:151
EdbSegmentCut * GetCut(int ud, int i)
Definition: EdbRunAccess.h:163

◆ Print()

void EdbRunAccess::Print ( )

104{
105 printf(" EdbRunAccess: eAFID=%d eCLUST=%d\n",eAFID, eCLUST );
106 printf("file: %s\n", eRunFileName.Data());
107 for(int i=0; i<3; i++)
108 if(eLayers[i]) eLayers[i]->Print();
109 }
void Print()
Definition: EdbLayer.cxx:149

◆ PrintStat()

void EdbRunAccess::PrintStat ( )

673{
674 if(!eRun) return;
675
676 printf("\nSome run statistics:\n");
677 printf( "--------------------\n");
678 int nviews=eRun->GetEntries();
680 printf("entries : %d\n", nviews );
681 printf("areas : %d in the range: ( %d : %d )\n",
683 int var=0;
684 if(eNareas>0) {
685 var = nviews/eNareas;
686 printf("views/area : %d (%d/side)\n", var, var/2 );
687 }
688 if(var*eNareas != nviews) printf("\t WARNING: areas are not equal\n");
689
690 for(int ud=0; ud<3; ud++) {
691 EdbPattern *pat=GetVP(ud);
692 int nempty=0;
693 if(pat) {
694 printf( "side:%d xmin,xmax = %f %f ymin,ymax = %f %f\n",
695 ud,
696 pat->Xmin(),
697 pat->Xmax(),
698 pat->Ymin(),
699 pat->Ymax()
700 );
701 nempty = CheckEmptyViews(*pat);
702 if(nempty>0) printf("\t WARNING: %d empty views in layer %d\n",nempty,ud);
703 }
704 }
705 printf("\n");
706}

◆ ReadVAfile()

void EdbRunAccess::ReadVAfile ( )

203{
204 int k = eRunFileName.Length();
205 TString s(eRunFileName.Data(),k-4); s+="va.root";
206 Log( 2,"EdbRunAccess::ReadVAfile","Get view corrections from %s",s.Data() );
207 TFile f( s.Data() );
208 eViewCorr = (TObjArray*)f.Get("viewcorr");
209 f.Close();
210}

◆ SegmentWeight()

float EdbRunAccess::SegmentWeight ( const EdbSegment s)

936{
937 if (eWeightAlg==1) return Sqrt(s.GetPuls()*s.GetSigmaY()/2.);
938 else if(eWeightAlg==2) return s.GetSigmaX();
939 return s.GetPuls();
940}

◆ Set0()

void EdbRunAccess::Set0 ( )

63{
64 eRun=0;
65 for(int i=0; i<3; i++) {
66 eVP[i]=0;
67 eLayers[i]=0;
68 eCuts[i]=0;
69 eCond[i]=0;
70 }
71 eAFID=0;
72 eCLUST=0;
73 eFirstArea = 9999999;
74 eLastArea = 0;
75 eNareas=0;
79 eDoViewAnalysis=true;
80 eHViewXY[1].InitH2(500,-250,250, 500, -250, 250);
81 eHViewXY[2].InitH2(500,-250,250, 500, -250, 250);
84 eHeaderCut="1";
85 eTracking=-1;
86 eWeightAlg=0;
87 eViewCorr = 0;
90 for(int i=0; i<4; i++) eGraphZ[i] = 0;
91 for(int i=0; i<3; i++) eGraphDZ[i] = 0;
92}
int InitH2(const EdbH2 &h)
Definition: EdbCell2.cpp:78
Bool_t eInvertSides
0 -do nothing, 1-invert sides
Definition: EdbRunAccess.h:32

◆ SetCond()

void EdbRunAccess::SetCond ( int  ud,
EdbScanCond cond 
)
inline
167 {
168 if(ud<0) return;
169 if(ud>2) return;
170 if(eCond[ud]) delete (eCond[ud]);
171 (eCond[ud]) = new EdbScanCond(cond);
172 }

◆ SetCutBottom()

void EdbRunAccess::SetCutBottom ( int  side,
float  wmin 
)

235{
236 float xmin[5]={ eViewXmin[side] + OverlapX(side), eViewYmin[side] , -1, -1, wmin};
237 float xmax[5]={ eViewXmax[side] - OverlapX(side), eViewYmin[side]+ OverlapY(side), 1, 1, 50};
238 ClearCuts(); AddSegmentCut(side,1,xmin,xmax);
239}
void ClearCuts()
Definition: EdbRunAccess.cxx:95
float OverlapX(int ud)
Definition: EdbRunAccess.h:178
float OverlapY(int ud)
Definition: EdbRunAccess.h:179

◆ SetCutLeft()

void EdbRunAccess::SetCutLeft ( int  side,
float  wmin 
)

214{
215 float xmin[5]={ eViewXmin[side] , eViewYmin[side], -1, -1, wmin};
216 float xmax[5]={ eViewXmin[side] + OverlapX(side), eViewYmax[side], 1, 1, 50};
217 ClearCuts(); AddSegmentCut(side,1,xmin,xmax);
218}

◆ SetCutRight()

void EdbRunAccess::SetCutRight ( int  side,
float  wmin 
)

221{
222 float xmin[5]={ eViewXmax[side] - OverlapX(side), eViewYmin[side], -1, -1, wmin};
223 float xmax[5]={ eViewXmax[side] , eViewYmax[side], 1, 1, 50};
224 ClearCuts(); AddSegmentCut(side,1,xmin,xmax);
225}

◆ SetCutTop()

void EdbRunAccess::SetCutTop ( int  side,
float  wmin 
)

228{
229 float xmin[5]={ eViewXmin[side] + OverlapX(side), eViewYmax[side] - OverlapY(side), -1, -1, wmin};
230 float xmax[5]={ eViewXmax[side] - OverlapX(side), eViewYmax[side] , 1, 1, 50};
231 ClearCuts(); AddSegmentCut(side,1,xmin,xmax);
232}

◆ SetPixelCorrection()

void EdbRunAccess::SetPixelCorrection ( const char *  str)

928{
929 if( 3 != sscanf(str,"%d %f %f", &eDoPixelCorr, &ePixelCorrX, &ePixelCorrY) )
931 if(eDoPixelCorr) Log(2,"EdbRunAccess::SetPixelCorrection","%s",str);
932}

◆ SetSegmentAtExternalSurface()

bool EdbRunAccess::SetSegmentAtExternalSurface ( EdbSegment seg,
int  side 
)

set the segment coordinates correspondent to the cluster coord with maximal Z
assumed that clusters are attached to segments

1168{
1171 TObjArray *clusters = seg->GetElements();
1172 if(!clusters) return false;
1173 int ncl = clusters->GetLast()+1;
1174 EdbCluster *cl=0, *clZ=(EdbCluster*)clusters->At(0);
1175 for(int i=0; i<ncl; i++ ) {
1176 cl = (EdbCluster*)clusters->At(i);
1177 if (side==1) { if( clZ->Z() < cl->Z() ) clZ=cl; } //side=1: top use point of max Z
1178 else if(side==2) { if( clZ->Z() > cl->Z() ) clZ=cl; } //side=2: bottom use point of min Z
1179 }
1180 seg->SetX0(clZ->X());
1181 seg->SetY0(clZ->Y());
1182 seg->SetZ0(clZ->Z());
1183 return true;
1184}
virtual Float_t Z() const
Definition: EdbCluster.h:80
virtual void SetX0(float x)
Definition: EdbSegment.h:45
virtual void SetZ0(float z)
Definition: EdbSegment.h:47
virtual void SetY0(float y)
Definition: EdbSegment.h:46

◆ ViewSide()

int EdbRunAccess::ViewSide ( const EdbView view) const

444{
445 int side=0;
446 if( view->GetNframesTop()==0) side=2; // 2- bottom
447 else if(view->GetNframesBot()==0) side=1; // 1- top
448 else return 0;
449 if(eInvertSides) side = 3-side;
450 return side;
451}
Int_t GetNframesTop() const
Definition: EdbView.h:207
Int_t GetNframesBot() const
Definition: EdbView.h:208

◆ ViewX()

float EdbRunAccess::ViewX ( const EdbView view) const
inline
145{return view->GetXview();} // todo

◆ ViewY()

float EdbRunAccess::ViewY ( const EdbView view) const
inline
146{return view->GetYview();} // todo

Member Data Documentation

◆ eAffStage2Abs

EdbAffine2D* EdbRunAccess::eAffStage2Abs
private

affine transformation extracted from Marks (if AFID=11)

◆ eAFID

Int_t EdbRunAccess::eAFID

if =1 - use affine transformations of the fiducial marks

◆ eCLUST

Int_t EdbRunAccess::eCLUST

1-use clusters, 0 - do not

◆ eCond

EdbScanCond* EdbRunAccess::eCond[3]
private

◆ eCutGR

Float_t EdbRunAccess::eCutGR
private

grain cut (chi)

◆ eCuts

TObjArray* EdbRunAccess::eCuts[3]
private

arrays of cuts to be applied to segments

◆ eDoPixelCorr

Int_t EdbRunAccess::eDoPixelCorr

apply or not pix/mic correction when read data (default is 0)

◆ eDoViewAnalysis

Bool_t EdbRunAccess::eDoViewAnalysis

fill or not the histograms for optional view analysis

◆ eFirstArea

Int_t EdbRunAccess::eFirstArea
private

◆ eGraphDZ

TGraph2D* EdbRunAccess::eGraphDZ[3]

keep the base/layer1/layer2 thickness calculated using eZ1/eZ2/eZ3/eZ4 for each view

◆ eGraphZ

TGraph2D* EdbRunAccess::eGraphZ[4]

keep z1/z2/z3/z4 surfaces using eZ1/eZ2/eZ3/eZ4 for each view

◆ eHeaderCut

TCut EdbRunAccess::eHeaderCut

header cut to be applied in run initialization

◆ eHViewXY

EdbH2 EdbRunAccess::eHViewXY[3]

XY segments distribution in a view local coords.

◆ eInvertSides

Bool_t EdbRunAccess::eInvertSides

0 -do nothing, 1-invert sides

◆ eLastArea

Int_t EdbRunAccess::eLastArea
private

◆ eLayers

EdbLayer* EdbRunAccess::eLayers[3]
private

base(0),up(1),down(2) layers

◆ eNareas

Int_t EdbRunAccess::eNareas
private

◆ eNviewsPerArea

Int_t EdbRunAccess::eNviewsPerArea
private

◆ ePixelCorrX

Float_t EdbRunAccess::ePixelCorrX

pixel/micron correction factor to be applied for data

◆ ePixelCorrY

Float_t EdbRunAccess::ePixelCorrY

◆ eRun

EdbRun* EdbRunAccess::eRun
private

pointer to the run to be accessed

◆ eRunFileName

TString EdbRunAccess::eRunFileName
private

◆ eTracking

Int_t EdbRunAccess::eTracking

to test tracking alorithm: -1-ignored(def),0/1 - trackings to accept

◆ eUseDensityAsW

Bool_t EdbRunAccess::eUseDensityAsW

in case of LASSO tracking possible to use eSigmaY as eW

◆ eUseExternalSurface

Bool_t EdbRunAccess::eUseExternalSurface

if true - set segment position corrisponding to the very external cluster

◆ eViewCorr

TObjArray* EdbRunAccess::eViewCorr
private

corrections obtained from the views alignment

◆ eViewXmax

Float_t EdbRunAccess::eViewXmax[3]
private

◆ eViewXmin

Float_t EdbRunAccess::eViewXmin[3]
private

◆ eViewYmax

Float_t EdbRunAccess::eViewYmax[3]
private

◆ eViewYmin

Float_t EdbRunAccess::eViewYmin[3]
private

◆ eVP

EdbPattern* EdbRunAccess::eVP[3]
private

base/up/down side patterns (0,1,2) with the views coordinates inside

◆ eWeightAlg

Int_t EdbRunAccess::eWeightAlg

0-puls, 1 - density(former eUseDensityAsW), 2-Likelyhood (eSigmaX)

◆ eXmax

Float_t EdbRunAccess::eXmax
private

◆ eXmin

Float_t EdbRunAccess::eXmin
private

◆ eXstep

Float_t EdbRunAccess::eXstep[3]
private

◆ eYmax

Float_t EdbRunAccess::eYmax
private

◆ eYmin

Float_t EdbRunAccess::eYmin
private

◆ eYstep

Float_t EdbRunAccess::eYstep[3]
private

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