FEDRA emulsion software from the OPERA Collaboration
EdbViewMap Class Reference

views map More...

#include <EdbViewMap.h>

Inheritance diagram for EdbViewMap:
Collaboration diagram for EdbViewMap:

Public Member Functions

void AddInverseAlign ()
 
void CheckView (AView *v)
 
void CheckView (int side, int area, int idview)
 
void CheckViewFrameAlignQuality (const char *file)
 
void ConvertRun (const char *fin, const char *fout)
 
void CorrectCols (TTree &tree)
 
void CorrectLines (TTree &tree)
 
void CorrectToStage ()
 
void DoCorrection ()
 
void DoCorrection1 (AArea &area)
 
void DoCorrectionAllNeib (AView &v)
 
void DoCorrectionBestNeib (AView &v)
 
 EdbViewMap ()
 
void FillAHT (AlignmentParView *apv, EdbViewHeader *vn, EdbViewHeader *vc, EdbAffine2D *aff, float w, TTree *&aht, const char *name)
 
void FillALcp ()
 
void FillAR ()
 
AAreaGetAArea (int side, int area)
 
AViewGetAView (EdbViewHeader &h)
 
AViewGetAView (int side, int area, int view)
 
EdbViewHeaderGetHeader (int i)
 
EdbViewHeaderGetPinHeader (int i)
 
EdbViewHeaderGetViewHeader (AlignmentParView &apv, int fs)
 
void InitAR ()
 
bool IsBug0 (AView &v1, AView &v2)
 
void MakeAHTnocorr ()
 
int PinPoint (int pinviewid)
 
int PinViewID (int pinpoint, int side)
 
float Quality (AApv &apv)
 
void ReadPinViewsHeaders (const char *file)
 
void ReadViewAlign (const char *file)
 
void ReadViewsHeaders (const char *file, TCut &cut)
 
void SaveColsCorrToRun (const char *fin)
 
void SaveCorrToRun (const char *fin)
 
void SaveLinesCorrToRun (const char *fin)
 
AViewSelectBestNotCorrected (AArea &area)
 
void t (const char *fin, const char *fout, int algorithm, const char *cut="1")
 
void Transform (AView &v, const AlignmentParView &apv)
 
float ViewCorrectability (AView &v)
 
float ViewCorrectability1 (AView &v)
 
float ViewCorrectabilityAll (AView &v)
 
int ViewSide (const EdbViewHeader &vh) const
 
virtual ~EdbViewMap ()
 

Static Public Member Functions

static bool IsApplied (const AlignmentParFrame &apf)
 
static bool IsApplied (const AlignmentParView &apv)
 
static bool IsFound (const AlignmentParFrame &apf)
 
static bool IsFound (const AlignmentParView &apv)
 
static bool IsPin (const AlignmentParView &apv)
 
static bool IsRecovered (const AlignmentParFrame &apf)
 

Public Attributes

TTree * eAHT
 aligned headers tree More...
 
AlignmentParView ** eALcp
 read tree with alignment parameters here More...
 
Int_t eAlgorithm
 algorithm used for correction: 0-bestNeigbour, 1-allNeigbours, 2,3 - hysteresis corr More...
 
ARun eAR
 structures with pointers to the above containers More...
 
TH2D * eHX0even
 even lines x-offset More...
 
TH2D * eHX0odd
 odd lines x-offset side=0 More...
 
TH2D * eHX1even
 even lines x-offset More...
 
TH2D * eHX1odd
 odd lines x-offset side=1 More...
 
TH2D * eHY0even
 even lines x-offset More...
 
TH2D * eHY0odd
 odd lines x-offset More...
 
TH2D * eHY1even
 even lines x-offset More...
 
TH2D * eHY1odd
 odd lines x-offset More...
 
Int_t eNcp
 number of objects in eALcp More...
 
Int_t eNpvh
 number of headers More...
 
Int_t eNvh
 number of headers More...
 
Int_t ePinPointID
 pin point ID to be used for this alignment pinID=(1000 + 2*pointID+side) More...
 
TObjArray ePinViewHeaders
 read headers for all pinned views here More...
 
TObjArray eViewHeaders
 

Detailed Description

views map

Constructor & Destructor Documentation

◆ EdbViewMap()

EdbViewMap::EdbViewMap ( )
22{
23 eNcp =0;
24 eALcp =0;
25 eNvh =0;
26 eAHT =0;
27 eAlgorithm=0;
28 eHX0odd=0;
29 eHX0even=0;
30 eHY0odd=0;
31 eHY0even=0;
32 eHX1odd=0;
33 eHX1even=0;
34 eHY1odd=0;
35 eHY1even=0;
36}
AlignmentParView ** eALcp
read tree with alignment parameters here
Definition: EdbViewMap.h:75
Int_t eAlgorithm
algorithm used for correction: 0-bestNeigbour, 1-allNeigbours, 2,3 - hysteresis corr
Definition: EdbViewMap.h:87
TH2D * eHX0odd
odd lines x-offset side=0
Definition: EdbViewMap.h:91
TTree * eAHT
aligned headers tree
Definition: EdbViewMap.h:85
TH2D * eHY1even
even lines x-offset
Definition: EdbViewMap.h:98
TH2D * eHY0odd
odd lines x-offset
Definition: EdbViewMap.h:93
Int_t eNcp
number of objects in eALcp
Definition: EdbViewMap.h:74
TH2D * eHY0even
even lines x-offset
Definition: EdbViewMap.h:94
TH2D * eHX1odd
odd lines x-offset side=1
Definition: EdbViewMap.h:95
TH2D * eHY1odd
odd lines x-offset
Definition: EdbViewMap.h:97
Int_t eNvh
number of headers
Definition: EdbViewMap.h:77
TH2D * eHX0even
even lines x-offset
Definition: EdbViewMap.h:92
TH2D * eHX1even
even lines x-offset
Definition: EdbViewMap.h:96

◆ ~EdbViewMap()

EdbViewMap::~EdbViewMap ( )
virtual
40{
41 if(eALcp) delete[] eALcp;
42}

Member Function Documentation

◆ AddInverseAlign()

void EdbViewMap::AddInverseAlign ( )
204{
205 for(int i=0; i<eNcp; i++)
206 {
208 AlignmentParView *vi = eALcp[i+eNcp] = new AlignmentParView(*v);
209 vi->view1 = v->view2;
210 vi->area1 = v->area2;
211 vi->side1 = v->side2;
212 vi->view2 = v->view1;
213 vi->area2 = v->area1;
214 vi->side2 = v->side1;
215 vi->dx = -v->dx;
216 vi->dy = -v->dy;
217 vi->dz = -v->dz;
218 }
219 eNcp *= 2;
220}
Definition: EdbRun.h:33
Float_t dy
found offset
Definition: EdbRun.h:41
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 view2
2-d view
Definition: EdbRun.h:35
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

◆ CheckView() [1/2]

void EdbViewMap::CheckView ( AView v)
329{
330 v->header->Print();
331 float corr = ViewCorrectability(*v);
332 printf("%d corrected: %d nalv= %d correctability: %f \n", v->header->GetViewID(), v->isCorrected, v->nalv, corr);
333
334 for(int i=0; i<v->nalv; i++)
335 {
336 AlignmentParView *apv = v->alv[i]->apv;
337 if(apv->view1== v->header->GetViewID()|| apv->view2==v->header->GetViewID())
338 printf("%d %d %d\n",apv->view1,apv->view2,apv->nsg);
339 }
340}
void Print() const
Definition: EdbView.cxx:749
Int_t GetViewID() const
Definition: EdbView.h:90
float ViewCorrectability(AView &v)
Definition: EdbViewMap.cxx:419
AlignmentParView * apv
Definition: EdbViewMap.h:34
AApv * alv[16]
Definition: EdbViewMap.h:43
bool isCorrected
Definition: EdbViewMap.h:44
int nalv
Definition: EdbViewMap.h:42
EdbViewHeader * header
Definition: EdbViewMap.h:41
Int_t nsg
peak height
Definition: EdbRun.h:47

◆ CheckView() [2/2]

void EdbViewMap::CheckView ( int  side,
int  area,
int  idview 
)
323{
324 AView *v = eAR.side[side].area[area]->view[idview];
325 CheckView(v);
326}
void CheckView(int side, int area, int idview)
Definition: EdbViewMap.cxx:322
ARun eAR
structures with pointers to the above containers
Definition: EdbViewMap.h:83
AView ** view
Definition: EdbViewMap.h:54
ASide side[2]
Definition: EdbViewMap.h:67
AArea ** area
Definition: EdbViewMap.h:60
< view and neighbouring
Definition: EdbViewMap.h:40

◆ CheckViewFrameAlignQuality()

void EdbViewMap::CheckViewFrameAlignQuality ( const char *  file)

tag view as bad if there are notaligned frames inside

224{
226 TFile f(file);
227 TTree *t = (TTree*)(f.Get("FrameAlign"));
229 t->SetBranchAddress("frpar",&fr);
230 int n = t->GetEntries();
231 printf("%d entries in FrameAlign \n",n);
232 int view=-1, area=-1, side=-1;
233 int changeflag=0;
234 int flag=0;
235 for(int i=0; i<n; i++)
236 {
237 t->GetEntry(i);
238 if( fr.view != view || area != fr.area || side != fr.side) { // new view started here
239 if(changeflag>2) GetAView(side,area,view)->quality /= 100;
240 view = fr.view; area = fr.area; side = fr.side;
241 flag = 0;
242 changeflag=0;
243 }
244 if(fr.flag!=flag) { flag=fr.flag; changeflag++; }
245 }
246 if(changeflag>2) GetAView(side,area,view)->quality /= 100;
247
248 f.Close();
249 Log(2,"EdbViewMap::CheckViewFrameAlignQuality","Read %d frame alignments",n);
250}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
FILE * f
Definition: RecDispMC.C:150
AView * GetAView(EdbViewHeader &h)
Definition: EdbViewMap.cxx:253
void t(const char *fin, const char *fout, int algorithm, const char *cut="1")
Definition: EdbViewMap.cxx:862
TFile * file
Definition: write_pvr.C:3
float quality
Definition: EdbViewMap.h:46
Definition: EdbRun.h:53
Int_t view
view id
Definition: EdbRun.h:56
Int_t flag
in first 3 bits: {f_Applied = 0x01,f_Found = 0x02,f_Recovered = 0x04};
Definition: EdbRun.h:71
Int_t side
Definition: EdbRun.h:58
Int_t area
Definition: EdbRun.h:57

◆ ConvertRun()

void EdbViewMap::ConvertRun ( const char *  fin,
const char *  fout 
)
547{
548 EdbRun rin(fin);
549 EdbView *v = rin.GetView();
550 EdbRun rout(fout,"RECREATE");
551
552 int n=rin.GetEntries();
553 printf("entries=%d\n",n);
554
555 for(int i=0; i<n; i++)
556 {
557 rin.GetEntry(i,1,0,1,1,1);
558 //rin.GetEntry(i,1,0,0,0,0);
559
560 EdbViewHeader *h = v->GetHeader();
561 AView *av=GetAView(*h);
562 if(av) {
563 EdbAffine2D a( av->aff );
564 h->SetAffine( a.A11(),a.A12(),a.A21(),a.A22(),a.B1(),a.B2() );
565 }
566 rout.AddView(v);
567 }
568
569 if(eAHT) eAHT->Write();
570
571 //FillAHTapv();
572 //if(eAHTapv) eAHTapv->Write();
573
574 rout.Close();
575 //rin.Close();
576}
void a()
Definition: check_aligned.C:59
Definition: EdbAffine.h:17
Definition: EdbRun.h:75
view identification
Definition: EdbView.h:26
void SetAffine(float a11, float a12, float a21, float a22, float b1, float b2)
Definition: EdbView.h:67
Base scanning data object: entry into Run tree.
Definition: EdbView.h:134
EdbViewHeader * GetHeader() const
Definition: EdbView.h:163
EdbAffine2D aff
Definition: EdbViewMap.h:45

◆ CorrectCols()

void EdbViewMap::CorrectCols ( TTree &  tree)
718{
719 tree.SetAlias("xm","(vc.eXview+vn.eXview)/2.");
720 tree.SetAlias("ym","(vc.eYview+vn.eYview)/2.");
721
722 TCut cutx("cutx","flag!=0 && abs(vn.eXview-vc.eXview)<100");
723 TCut cuty("cuty","flag!=0 && abs(vn.eYview-vc.eYview)<100");
724
725 tree.Draw( "dx:ym:xm>>hx0even(30,0,120000,20,0,100000)", "!(vn.eCol%2)&&side1==0" && cutx ,"prof colz");
726 tree.Draw( "dx:ym:xm>>hx0odd(30,0,120000,20,0,100000)", "(vn.eCol%2)&&side1==0" && cutx ,"prof colz");
727 tree.Draw( "dy:ym:xm>>hy0even(30,0,120000,20,0,100000)", "!(vn.eCol%2)&&side1==0" && cuty ,"prof colz");
728 tree.Draw( "dy:ym:xm>>hy0odd(30,0,120000,20,0,100000)", "(vn.eCol%2)&&side1==0" && cuty ,"prof colz");
729 tree.Draw( "dx:ym:xm>>hx1even(30,0,120000,20,0,100000)", "!(vn.eCol%2)&&side1==1" && cutx ,"prof colz");
730 tree.Draw( "dx:ym:xm>>hx1odd(30,0,120000,20,0,100000)", "(vn.eCol%2)&&side1==1" && cutx ,"prof colz");
731 tree.Draw( "dy:ym:xm>>hy1even(30,0,120000,20,0,100000)", "!(vn.eCol%2)&&side1==1" && cuty ,"prof colz");
732 tree.Draw( "dy:ym:xm>>hy1odd(30,0,120000,20,0,100000)", "(vn.eCol%2)&&side1==1" && cuty ,"prof colz");
733 eHX0odd = (TH2D*)((TH2D*)gDirectory->Get("hx0odd"))->Clone("X0odd"); // positive
734 eHX0even = (TH2D*)((TH2D*)gDirectory->Get("hx0even"))->Clone("X0even"); // negative
735 eHY0odd = (TH2D*)((TH2D*)gDirectory->Get("hy0odd"))->Clone("Y0odd"); // positive
736 eHY0even = (TH2D*)((TH2D*)gDirectory->Get("hy0even"))->Clone("Y0even"); // negative
737 eHX1odd = (TH2D*)((TH2D*)gDirectory->Get("hx1odd"))->Clone("X1odd"); // positive
738 eHX1even = (TH2D*)((TH2D*)gDirectory->Get("hx1even"))->Clone("X1even"); // negative
739 eHY1odd = (TH2D*)((TH2D*)gDirectory->Get("hy1odd"))->Clone("Y1odd"); // positive
740 eHY1even = (TH2D*)((TH2D*)gDirectory->Get("hy1even"))->Clone("Y1even"); // negative
741}

◆ CorrectLines()

void EdbViewMap::CorrectLines ( TTree &  tree)
676{
677 tree.SetAlias("xm","(vc.eXview+vn.eXview)/2.");
678 tree.SetAlias("ym","(vc.eYview+vn.eYview)/2.");
679
680 TCut cutx("cutx","flag!=0 && abs(vn.eXview-vc.eXview)<100");
681 TCut cuty("cuty","flag!=0 && abs(vn.eYview-vc.eYview)<100");
682
683 tree.Draw( "dx:ym:xm>>hx0even(30,0,120000,20,0,100000)", "!(vn.eRow%2)&&side1==0" && cutx ,"prof colz");
684 tree.Draw( "dx:ym:xm>>hx0odd(30,0,120000,20,0,100000)", "(vn.eRow%2)&&side1==0" && cutx ,"prof colz");
685 tree.Draw( "dy:ym:xm>>hy0even(30,0,120000,20,0,100000)", "!(vn.eRow%2)&&side1==0" && cuty ,"prof colz");
686 tree.Draw( "dy:ym:xm>>hy0odd(30,0,120000,20,0,100000)", "(vn.eRow%2)&&side1==0" && cuty ,"prof colz");
687 tree.Draw( "dx:ym:xm>>hx1even(30,0,120000,20,0,100000)", "!(vn.eRow%2)&&side1==1" && cutx ,"prof colz");
688 tree.Draw( "dx:ym:xm>>hx1odd(30,0,120000,20,0,100000)", "(vn.eRow%2)&&side1==1" && cutx ,"prof colz");
689 tree.Draw( "dy:ym:xm>>hy1even(30,0,120000,20,0,100000)", "!(vn.eRow%2)&&side1==1" && cuty ,"prof colz");
690 tree.Draw( "dy:ym:xm>>hy1odd(30,0,120000,20,0,100000)", "(vn.eRow%2)&&side1==1" && cuty ,"prof colz");
691 eHX0odd = (TH2D*)((TH2D*)gDirectory->Get("hx0odd"))->Clone("X0odd"); // positive
692 eHX0even = (TH2D*)((TH2D*)gDirectory->Get("hx0even"))->Clone("X0even"); // negative
693 eHY0odd = (TH2D*)((TH2D*)gDirectory->Get("hy0odd"))->Clone("Y0odd"); // positive
694 eHY0even = (TH2D*)((TH2D*)gDirectory->Get("hy0even"))->Clone("Y0even"); // negative
695 eHX1odd = (TH2D*)((TH2D*)gDirectory->Get("hx1odd"))->Clone("X1odd"); // positive
696 eHX1even = (TH2D*)((TH2D*)gDirectory->Get("hx1even"))->Clone("X1even"); // negative
697 eHY1odd = (TH2D*)((TH2D*)gDirectory->Get("hy1odd"))->Clone("Y1odd"); // positive
698 eHY1even = (TH2D*)((TH2D*)gDirectory->Get("hy1even"))->Clone("Y1even"); // negative
699
700 /*
701 printf("eHXodd integral: %f eHXeven integral: %f \n",
702 eHXodd->Integral(),
703 eHXeven->Integral());
704 TH2D *hxdiff = (TH2D *)(eHXodd->Clone("hxdiff"));
705 hxdiff->Add(eHXeven, 1);
706 printf("hdiff integral: %f \n",hxdiff->Integral());
707 printf("eHYodd integral: %f eHYeven integral: %f \n",
708 eHYodd->Integral(),
709 eHYeven->Integral());
710 TH2D *hydiff = (TH2D *)(eHYodd->Clone("hydiff"));
711 hydiff->Add(eHYeven, 1);
712 printf("hydiff integral: %f \n",hydiff->Integral());
713 */
714}

◆ CorrectToStage()

void EdbViewMap::CorrectToStage ( )

for each area eliminate scaling&rotation in respect to the stage coord

511{
513 for( int side=0; side<2; side++)
514 for( int area=0; area<eAR.side[side].narea; area++)
515 {
516 AArea *a= eAR.side[side].area[area];
517 TArrayF xst(a->nview);
518 TArrayF yst(a->nview);
519 TArrayF x(a->nview);
520 TArrayF y(a->nview);
521 int cnt=0;
522 for(int i=0; i<a->nview; i++)
523 {
524 AView *v = GetAView(side,area,i);
525 if(!v) continue;
526 if(!v->header) continue;
527 xst[cnt]= v->header->GetXview();
528 yst[cnt]= v->header->GetYview();
529 v->aff.ShiftX( xst[cnt] );
530 v->aff.ShiftY( yst[cnt] );
531 x[cnt] = v->aff.B1();
532 y[cnt] = v->aff.B2();
533 cnt++;
534 }
535 EdbAffine2D aff;
536 aff.Calculate(cnt,x.GetArray(),y.GetArray(),xst.GetArray(),yst.GetArray());
537 aff.Print();
538 for(int i=0; i<a->nview; i++)
539 {
540 if(GetAView(side,area,i)) GetAView(side,area,i)->aff.Transform(aff);
541 }
542 }
543}
Float_t B2() const
Definition: EdbAffine.h:48
void Print(Option_t *opt="") const
Definition: EdbAffine.cxx:52
void ShiftX(float d)
Definition: EdbAffine.h:64
void ShiftY(float d)
Definition: EdbAffine.h:65
Int_t Calculate(EdbPointsBox2D *b1, EdbPointsBox2D *b2)
Definition: EdbAffine.cxx:231
Float_t B1() const
Definition: EdbAffine.h:47
void Transform(const EdbAffine2D *a)
Definition: EdbAffine.cxx:93
Float_t GetXview() const
Definition: EdbView.h:93
Float_t GetYview() const
Definition: EdbView.h:94
Definition: EdbViewMap.h:51
int narea
Definition: EdbViewMap.h:61

◆ DoCorrection()

void EdbViewMap::DoCorrection ( )
503{
504 for( int iside=0; iside<2; iside++)
505 for( int iarea=0; iarea<eAR.side[iside].narea; iarea++)
506 DoCorrection1( *(eAR.side[iside].area[iarea]) );
507}
void DoCorrection1(AArea &area)
Definition: EdbViewMap.cxx:486

◆ DoCorrection1()

void EdbViewMap::DoCorrection1 ( AArea area)
487{
488 Log(3,"DoCorrection","area with %d views, alg=%d",area.nview,eAlgorithm);
489 for(int i=0; i<area.nview; i++)
490 {
491 AView *v = SelectBestNotCorrected(area);
492 if(v)
493 {
494 if ( eAlgorithm==0 ) DoCorrectionBestNeib(*v);
495 else if( eAlgorithm==1 ) DoCorrectionAllNeib(*v);
496 }
497 }
498 Log(3,"DoCorrection","ok");
499}
void DoCorrectionAllNeib(AView &v)
Definition: EdbViewMap.cxx:376
AView * SelectBestNotCorrected(AArea &area)
Definition: EdbViewMap.cxx:465
void DoCorrectionBestNeib(AView &v)
Definition: EdbViewMap.cxx:343
int nview
Definition: EdbViewMap.h:53

◆ DoCorrectionAllNeib()

void EdbViewMap::DoCorrectionAllNeib ( AView v)

allNeigbours algorithm
convention for apv: view1+offset => view2

TODO: insert sigma approach

377{
380
382
383 int ncorr=0;
384 EdbAffine2D aff[10]; // corrections obtained from different neighbours
385 float weight[10];
386 for( int inb=0; inb < v.nalv; inb++ )
387 {
388 AApv *apv = v.alv[inb];
389 AView *vc=GetAView( apv->apv->side2, apv->apv->area2, apv->apv->view2);
390 if(!vc) continue;
391 if(vc->isCorrected)
392 {
393 aff[ncorr].Set( vc->aff );
394 aff[ncorr].ShiftX(apv->apv->dx);
395 aff[ncorr].ShiftY(apv->apv->dy);
396 weight[ncorr] = Quality(*apv);
397 //printf("%5d %5d %9.1f ",apv->apv->view1, apv->apv->view2, weight[ncorr]); aff[ncorr].Print();
398 FillAHT( apv->apv, v.header, vc->header, &(aff[ncorr]), weight[ncorr], eAHT, "AHT" );
399 ncorr++;
400 }
401 }
402
403 float xmean=0, ymean=0;
404 float wmean=0;
405 for( int i=0; i < ncorr; i++ )
406 {
407 xmean += aff[i].B1()*weight[i];
408 ymean += aff[i].B2()*weight[i];
409 wmean += weight[i];
410 }
411 xmean /= wmean;
412 ymean /= wmean;
413 v.aff.ShiftX(xmean);
414 v.aff.ShiftY(ymean);
415 v.isCorrected=true;
416}
void Set(EdbAffine2D &a)
Definition: EdbAffine.h:36
float Quality(AApv &apv)
Definition: EdbViewMap.cxx:301
void FillAHT(AlignmentParView *apv, EdbViewHeader *vn, EdbViewHeader *vc, EdbAffine2D *aff, float w, TTree *&aht, const char *name)
Definition: EdbViewMap.cxx:604
< side-by-side alignment
Definition: EdbViewMap.h:33

◆ DoCorrectionBestNeib()

void EdbViewMap::DoCorrectionBestNeib ( AView v)

bestNeigbours algorithm
convention for apv: view1+offset => view2

344{
347 Log(3,"DoCorrectionBestNeib","%d alv",v.nalv);
348 AView *vcbest=0;
349 AApv *apv0=0;
350 float wmax=0;
351 for( int inb=0; inb < v.nalv; inb++ )
352 {
353 AApv *apv = v.alv[inb];
354 AView *vc=GetAView( apv->apv->side2, apv->apv->area2, apv->apv->view2);
355 //printf("i=%d\n",inb);
356 if(!vc) continue;
357 if(vc->isCorrected)
358 {
359 float w=Quality(*apv);
360 //printf("w=%f\n",w);
361 if( w > wmax) { wmax=w; vcbest=vc; apv0=apv;}
362 }
363 }
364 if(vcbest&&apv0) {
365 v.aff.Transform( vcbest->aff );
366
367 Transform( v, *(apv0->apv) );
368 //printf( "%d %d %d ", apv0->apv->side1, apv0->apv->area1, v.header->GetViewID()); v.aff.Print();
369 float w=Quality(*apv0);
370 v.isCorrected=true;
371 FillAHT( apv0->apv, v.header, vcbest->header, &(v.aff), w, eAHT, "AHT" );
372 }
373}
void Transform(AView &v, const AlignmentParView &apv)
Definition: EdbViewMap.cxx:315
void w(int rid=2, int nviews=2)
Definition: test.C:27

◆ FillAHT()

void EdbViewMap::FillAHT ( AlignmentParView apv,
EdbViewHeader vn,
EdbViewHeader vc,
EdbAffine2D aff,
float  w,
TTree *&  aht,
const char *  name 
)
605{
606 if(!aht) {
607 aht=new TTree( name , "Aligned Headers");
608 aht->Branch("alpar",apv,"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");
609 aht->Branch("vn",vn);
610 aht->Branch("vc",vc);
611 aht->Branch("aff",aff);
612 aht->Branch("w",&w,"w/F");
613 } else {
614 aht->SetBranchAddress("alpar",&(apv->view1) );
615 aht->SetBranchAddress("vn",&vn);
616 aht->SetBranchAddress("vc",&vc);
617 aht->SetBranchAddress("aff",&aff);
618 aht->SetBranchAddress("w",&w);
619 }
620 aht->Fill();
621}
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24

◆ FillALcp()

void EdbViewMap::FillALcp ( )
267{
268 for(int i=0; i<eNcp; i++) {
269 AlignmentParView *apv = eALcp[i];
270 if(!apv) continue;
271 AView *v1 = GetAView( apv->side1, apv->area1, apv->view1 );
272 if(!v1) continue;
273 AView *v2 = GetAView( apv->side2, apv->area2, apv->view2 );
274 if(!v2) continue;
275
276 float quality=1;
277 if( IsBug0(*v1,*v2) ) continue;
278 //printf("%d %d %d - %d %d %d\n",apv->side1, apv->area1, apv->view1, apv->side2, apv->area2, apv->view2);
279
280 if(v1) {
281 if(v1->nalv>15) { Log(1,"EdbViewMap::FillALcp","ERROR: nalv out of limits: %d",v1->nalv); continue;}
282 v1->alv[v1->nalv] = new AApv;
283 v1->alv[v1->nalv]->apv = apv;
284 v1->alv[v1->nalv]->quality=quality;
285 v1->nalv++;
286 }
287 /*
288 if(v2) {
289 if(v2->nalv>15) { Log(1,"EdbViewMap::FillALcp","ERROR: nalv out of limits: %d",v2->nalv); continue;}
290 v2->alv[v2->nalv] = new AApv;
291 v2->alv[v2->nalv]->apv = apv;
292 v2->alv[v2->nalv]->quality=quality;
293 v2->nalv++;
294 }
295 */
296 }
297 Log(3,"FillALcp","filled with %d cp",eNcp);
298}
bool IsBug0(AView &v1, AView &v2)
Definition: EdbViewMap.cxx:174
float quality
Definition: EdbViewMap.h:35

◆ FillAR()

void EdbViewMap::FillAR ( )
93{
94 for(int i=0; i<eNvh; i++)
95 {
96 EdbViewHeader *vh = GetHeader(i);
97 int side=ViewSide(*vh);
98 int area=vh->GetAreaID();
99 int view = vh->GetViewID();
100 //printf("%d %d (%d) %d (%d)\n",side,area, eAR.side[side].narea ,view, eAR.side[side].area[area]->nview );
101 eAR.side[side].area[area]->view[view]->header=vh;
102 if( view > eAR.side[side].area[area]->ilast) eAR.side[side].area[area]->ilast=view;
103 }
104 Log(2,"EdbViewMap::FillARun","with %d headers",eNvh );
105
106 for(int i=0; i<eNpvh; i++) //pinned
107 {
109 int side=ViewSide(*vh);
110 int area=vh->GetAreaID();
111 int oldid = vh->GetViewID();
112 int newid = (eAR.side[side].area[area]->ilast++);
113 printf("oldid %d newid %d\n",oldid,newid);
114 eAR.side[side].area[area]->view[newid]->header=vh;
115 eAR.side[side].area[area]->view[ newid ]->isCorrected=true;
116 eAR.side[side].area[area]->view[ newid ]->quality=1;
117 vh->SetViewID(newid); // update header with new id
118 for(int ia=0; ia<eNcp; ia++)
119 {
120 AlignmentParView *apv = eALcp[ia];
121 if( apv->side1==side && apv->area1==area )
122 {
123 if( apv->view1==oldid) apv->view1=newid; // update APV vith newid
124 if( apv->view2==oldid) apv->view2=newid;
125 }
126 }
127 }
128 Log(2,"EdbViewMap::FillARun","with %d pinned headers",eNpvh );
129}
void SetViewID(int id)
Definition: EdbView.h:99
Int_t GetAreaID() const
Definition: EdbView.h:91
int ViewSide(const EdbViewHeader &vh) const
Definition: EdbViewMap.h:115
EdbViewHeader * GetPinHeader(int i)
Definition: EdbViewMap.h:114
EdbViewHeader * GetHeader(int i)
Definition: EdbViewMap.h:113
Int_t eNpvh
number of headers
Definition: EdbViewMap.h:80
int ilast
Definition: EdbViewMap.h:52

◆ GetAArea()

AArea * EdbViewMap::GetAArea ( int  side,
int  area 
)
inline
121{ return eAR.side[side].area[area]; }

◆ GetAView() [1/2]

AView * EdbViewMap::GetAView ( EdbViewHeader h)
254{
255 return GetAView( ViewSide(h), h.GetAreaID(), h.GetViewID() );
256}

◆ GetAView() [2/2]

AView * EdbViewMap::GetAView ( int  side,
int  area,
int  view 
)
258{
259 if( side<0||side>1) return 0;
260 if( area<0||area>=eAR.side[side].narea) return 0;
261 if( view<0||view>=eAR.side[side].area[area]->nview) return 0;
262 return eAR.side[side].area[area]->view[view];
263}

◆ GetHeader()

EdbViewHeader * EdbViewMap::GetHeader ( int  i)
inline
113{ return (EdbViewHeader*)(eViewHeaders.At(i)); }
TObjArray eViewHeaders
Definition: EdbViewMap.h:78

◆ GetPinHeader()

EdbViewHeader * EdbViewMap::GetPinHeader ( int  i)
inline
114{ return (EdbViewHeader*)(ePinViewHeaders.At(i)); }
TObjArray ePinViewHeaders
read headers for all pinned views here
Definition: EdbViewMap.h:81

◆ GetViewHeader()

EdbViewHeader * EdbViewMap::GetViewHeader ( AlignmentParView apv,
int  fs 
)
625{
626 EdbViewHeader *h=0;
627 int side=-1;
628 int vid =-1;
629 int area = -1;
630 if(fs==1) {
631 side = apv.side1;
632 area = apv.area1;
633 vid = apv.view1;
634 } else if(fs==2) {
635 side = apv.side2;
636 area = apv.area2;
637 vid = apv.view2;
638 }
639 else return h;
640
641 for(int i=0; i<eNvh; i++)
642 {
643 h = GetHeader(i);
644 if(ViewSide(*h) == side)
645 if( h->GetViewID()==vid )
646 if(h->GetAreaID()==area)
647 return h;
648 }
649 return 0;
650}

◆ InitAR()

void EdbViewMap::InitAR ( )
46{
47 int narea=0;
48 for(int i=0; i<eNvh; i++)
49 {
50 EdbViewHeader *vh = GetHeader(i);
51 if( vh) if(vh->GetAreaID()>narea) narea=vh->GetAreaID();
52 }
53 narea++;
54 for(int iside=0; iside<2; iside++) {
55 eAR.side[iside].narea = narea;
56 eAR.side[iside].area = new AArea*[narea];
57 for(int iarea=0; iarea<narea; iarea++) {
58 eAR.side[iside].area[iarea] = new AArea;
59 eAR.side[iside].area[iarea]->nview=0;
60 //eAR.side[iside].area[iarea]->nalv=0;
61 }
62 }
63 printf("%d areas\n", narea);
64
65 for(int i=0; i<eNvh; i++)
66 {
67 EdbViewHeader *vh = GetHeader(i);
68 //eAR.side[ViewSide(*vh)].area[vh->GetAreaID()]->nview++;
69 if( eAR.side[ViewSide(*vh)].area[vh->GetAreaID()]->nview < vh->GetViewID())
70 eAR.side[ViewSide(*vh)].area[vh->GetAreaID()]->nview = vh->GetViewID();
71 }
72
73 for(int iside=0; iside<2; iside++)
74 for(int iarea=0; iarea<narea; iarea++)
75 {
76 eAR.side[iside].area[iarea]->nview += (eNpvh +1); // keep space for pinned views at the end of the buffer
77 eAR.side[iside].area[iarea]->ilast = 0;
78 printf( "%d %d %d\n", iside, iarea, eAR.side[iside].area[iarea]->nview );
79 eAR.side[iside].area[iarea]->view = new AView*[eAR.side[iside].area[iarea]->nview];
80 for( int iv=0; iv<eAR.side[iside].area[iarea]->nview; iv++ )
81 {
82 AView *v = eAR.side[iside].area[iarea]->view[iv] = new AView;
83 v->header=0;
84 v->quality=1;
85 v->nalv=0;
86 v->isCorrected=0;
87 }
88 }
89}

◆ IsApplied() [1/2]

static bool EdbViewMap::IsApplied ( const AlignmentParFrame apf)
inlinestatic
148{return test_ith_bit(0,apf.flag);}
bool test_ith_bit(int i_, int flags_)
Definition: EdbViewMap.h:26

◆ IsApplied() [2/2]

static bool EdbViewMap::IsApplied ( const AlignmentParView apv)
inlinestatic
152{return test_ith_bit(0,apv.flag);}
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

◆ IsBug0()

bool EdbViewMap::IsBug0 ( AView v1,
AView v2 
)

protection for bug in LASSO 11/06/2015

175{
177 if( Abs( v1.header->GetXview() - v2.header->GetXview()) > 1000. ||
178 Abs( v1.header->GetYview() - v2.header->GetYview()) > 1000. ) return 1;
179 return 0;
180}

◆ IsFound() [1/2]

static bool EdbViewMap::IsFound ( const AlignmentParFrame apf)
inlinestatic
149{return test_ith_bit(1,apf.flag);}

◆ IsFound() [2/2]

static bool EdbViewMap::IsFound ( const AlignmentParView apv)
inlinestatic
154{return test_ith_bit(2,apv.flag);}

◆ IsPin()

static bool EdbViewMap::IsPin ( const AlignmentParView apv)
inlinestatic
153{return test_ith_bit(1,apv.flag);}

◆ IsRecovered()

static bool EdbViewMap::IsRecovered ( const AlignmentParFrame apf)
inlinestatic
150{return test_ith_bit(2,apf.flag);}

◆ MakeAHTnocorr()

void EdbViewMap::MakeAHTnocorr ( )
654{
655 TTree *tree=0;
656 EdbAffine2D *aff = new EdbAffine2D();
657 for(int i=0; i<eNcp; i++)
658 {
659 AlignmentParView *apv=eALcp[i];
660 EdbViewHeader *v1 = GetViewHeader(*apv,1);
661 EdbViewHeader *v2 = GetViewHeader(*apv,2);
662 if(v1&&v2) FillAHT(apv, v1, v2, aff, 0, tree, "AHTnocorr");
663 }
664
665 if (eAlgorithm==2) CorrectLines( *tree );
666 else if(eAlgorithm==3) CorrectCols( *tree );
667
668 TFile rout( "aht.root","RECREATE" );
669 tree->Write();
670 rout.Close();
671
672}
EdbViewHeader * GetViewHeader(AlignmentParView &apv, int fs)
Definition: EdbViewMap.cxx:624
void CorrectLines(TTree &tree)
Definition: EdbViewMap.cxx:675
void CorrectCols(TTree &tree)
Definition: EdbViewMap.cxx:717

◆ PinPoint()

int EdbViewMap::PinPoint ( int  pinviewid)
inline
128{return (-1000-pinviewid)/2; }

◆ PinViewID()

int EdbViewMap::PinViewID ( int  pinpoint,
int  side 
)
inline
127{return -(1000+2*pinpoint+side); }

◆ Quality()

float EdbViewMap::Quality ( AApv apv)

the bigger the better

302{
303
304 float q = 0;
305 if(apv.apv->view1>=0&&apv.apv->view2>=0)
306 {
307 AView *v1 = GetAView( apv.apv->side1, apv.apv->area1, apv.apv->view1 );
308 AView *v2 = GetAView( apv.apv->side2, apv.apv->area2, apv.apv->view2 );
309 if(v1&&v2) q = apv.apv->nsg * apv.quality * v1->quality * v2->quality;
310 }
311 return q;
312}
q
Definition: testBGReduction_AllMethods.C:55

◆ ReadPinViewsHeaders()

void EdbViewMap::ReadPinViewsHeaders ( const char *  file)
155{
156 EdbRun r(file);
157 TTree *t = r.GetPinViews();
158 if(t) {
160 t->SetBranchAddress("headers",&h);
161 int n = t->GetEntries();
162 for(int i=0; i<n; i++)
163 {
164 t->GetEntry(i);
165 if(h) ePinViewHeaders.Add( new EdbViewHeader( *h ) );
166 }
167 eNpvh = n;
168 }
169 r.Close();
170 Log(2,"EdbViewMap::ReadPinViewsHeaders","Read %d headers",eNpvh);
171}
void r(int rid=2)
Definition: test.C:201

◆ ReadViewAlign()

void EdbViewMap::ReadViewAlign ( const char *  file)
184{
185 eNcp=0;
186 TFile f(file);
187 TTree *t = (TTree*)(f.Get("ViewAlign"));
189 t->SetBranchAddress("alpar",&v);
190 int n = t->GetEntries();
191 eALcp = new AlignmentParView*[2*n]; // doubled space for inverse align
192 for(int i=0; i<n; i++)
193 {
194 t->GetEntry(i);
195 eALcp[i] = new AlignmentParView(v);
196 //eNcp++;
197 }
198 eNcp=n;
199 f.Close();
200 Log(2,"EdbViewMap::ReadViewAlign","Read %d aligned views couples",eNcp);
201}

◆ ReadViewsHeaders()

void EdbViewMap::ReadViewsHeaders ( const char *  file,
TCut &  cut 
)

read headers from runfile, fill eViewHeaders

133{
135 EdbRun r(file);
136 EdbView *view = r.GetView();
137 r.GetTree()->Draw(">>lst", cut );
138 TEventList *lst = (TEventList*)(gDirectory->GetList()->FindObject("lst"));
139 if(!lst) {Log(1,"EdbViewMap::ReadViewsHeaders","ERROR!: events list (lst) did not found!"); return;}
140 int n = lst->GetN();
141 for(int i=0; i<n; i++)
142 {
143 view = r.GetEntry(lst->GetEntry(i),1,0,0,0);
144 if(view)
145 if(view->GetHeader())
146 eViewHeaders.Add( new EdbViewHeader( *(view->GetHeader()) ) );
147 }
148 eNvh = n;
149 r.Close();
150 Log(2,"EdbViewMap::ReadViewsHeaders","Read %d headers",eNvh);
151}
TCut cut
Definition: check_shower.C:6

◆ SaveColsCorrToRun()

void EdbViewMap::SaveColsCorrToRun ( const char *  fin)
804{
805 EdbRun rin(fin);
806 EdbView *v = rin.GetView();
807 int n=rin.GetEntries();
808 printf("entries=%d\n",n);
809 TObjArray viewcorr(n);
810
811 for(int i=0; i<n; i++)
812 {
813 rin.GetEntry(i,1,0,0,0,0);
814 EdbViewHeader *h = v->GetHeader();
815 float x = h->GetXview();
816 float y = h->GetYview();
817 EdbAffine2D *aff = new EdbAffine2D(*h->GetAffine());
818
819 if(h->GetCol()%2) //odd
820 {
821 if(ViewSide(*h)) {
822 aff->ShiftX( 0.5*eHX1odd->Interpolate(x,y) );
823 aff->ShiftY( -0.5*eHY1odd->Interpolate(x,y) );
824 }
825 else {
826 aff->ShiftX( 0.5*eHX0odd->Interpolate(x,y) );
827 aff->ShiftY( -0.5*eHY0odd->Interpolate(x,y) );
828 }
829 }
830 else if(!(h->GetCol()%2)) //even
831 {
832 if(ViewSide(*h)) {
833 aff->ShiftX( 0.5*eHX1even->Interpolate(x,y) );
834 aff->ShiftY( -0.5*eHY1even->Interpolate(x,y) );
835 }
836 else {
837 aff->ShiftX( 0.5*eHX0even->Interpolate(x,y) );
838 aff->ShiftY( -0.5*eHY0even->Interpolate(x,y) );
839 }
840 }
841 viewcorr.Add( new EdbAffine2D(*aff) );
842 }
843 int k = strlen(fin);
844 TString s(fin,k-4); s+="va.root";
845 TFile rout( s.Data(),"RECREATE" );
846
847 viewcorr.Write("viewcorr",1);
848 if(eAHT) eAHT->Write();
849 eHX0odd->Write();
850 eHY0odd->Write();
851 eHX0even->Write();
852 eHY0even->Write();
853 eHX1odd->Write();
854 eHY1odd->Write();
855 eHX1even->Write();
856 eHY1even->Write();
857 rout.Close();
858 //rin.Close();
859}
Int_t GetCol() const
Definition: EdbView.h:123
EdbAffine2D const * GetAffine() const
Definition: EdbView.h:72
s
Definition: check_shower.C:55

◆ SaveCorrToRun()

void EdbViewMap::SaveCorrToRun ( const char *  fin)
580{
581 EdbRun rin(fin);
582 EdbView *v = rin.GetView();
583 int n=rin.GetEntries();
584 printf("entries=%d\n",n);
585 TObjArray viewcorr(n);
586
587 for(int i=0; i<n; i++)
588 {
589 rin.GetEntry(i,1,0,0,0,0);
590 EdbViewHeader *h = v->GetHeader();
591 viewcorr.Add( new EdbAffine2D(GetAView(*h)->aff) );
592 }
593 int k = strlen(fin);
594 TString s(fin,k-4); s+="va.root";
595 TFile rout( s.Data(),"RECREATE" );
596
597 viewcorr.Write("viewcorr",1);
598 if(eAHT) eAHT->Write();
599 rout.Close();
600 //rin.Close();
601}

◆ SaveLinesCorrToRun()

void EdbViewMap::SaveLinesCorrToRun ( const char *  fin)
745{
746 EdbRun rin(fin);
747 EdbView *v = rin.GetView();
748 int n=rin.GetEntries();
749 printf("entries=%d\n",n);
750 TObjArray viewcorr(n);
751
752 for(int i=0; i<n; i++)
753 {
754 rin.GetEntry(i,1,0,0,0,0);
755 EdbViewHeader *h = v->GetHeader();
756 float x = h->GetXview();
757 float y = h->GetYview();
758 EdbAffine2D *aff = new EdbAffine2D(*h->GetAffine());
759
760 if(h->GetRow()%2) //odd
761 {
762 if(ViewSide(*h)) {
763 aff->ShiftX( 0.5*eHX1odd->Interpolate(x,y) );
764 aff->ShiftY( -0.5*eHY1odd->Interpolate(x,y) );
765 }
766 else {
767 aff->ShiftX( 0.5*eHX0odd->Interpolate(x,y) );
768 aff->ShiftY( -0.5*eHY0odd->Interpolate(x,y) );
769 }
770 }
771 else if(!(h->GetRow()%2)) //even
772 {
773 if(ViewSide(*h)) {
774 aff->ShiftX( 0.5*eHX1even->Interpolate(x,y) );
775 aff->ShiftY( -0.5*eHY1even->Interpolate(x,y) );
776 }
777 else {
778 aff->ShiftX( 0.5*eHX0even->Interpolate(x,y) );
779 aff->ShiftY( -0.5*eHY0even->Interpolate(x,y) );
780 }
781 }
782 viewcorr.Add( new EdbAffine2D(*aff) );
783 }
784 int k = strlen(fin);
785 TString s(fin,k-4); s+="va.root";
786 TFile rout( s.Data(),"RECREATE" );
787
788 viewcorr.Write("viewcorr",1);
789 if(eAHT) eAHT->Write();
790 eHX0odd->Write();
791 eHY0odd->Write();
792 eHX0even->Write();
793 eHY0even->Write();
794 eHX1odd->Write();
795 eHY1odd->Write();
796 eHX1even->Write();
797 eHY1even->Write();
798 rout.Close();
799 //rin.Close();
800}
Int_t GetRow() const
Definition: EdbView.h:124

◆ SelectBestNotCorrected()

AView * EdbViewMap::SelectBestNotCorrected ( AArea area)
466{
467 AView *vbest=0;
468 float maxcorr=0;
469 for(int i=0; i<area.nview; i++)
470 {
471 AView *v=area.view[i];
472 if(!v) Log(1,"EdbViewMap::SelectBestNotCorrected","ERROR; zero pointer");
473 float corr = ViewCorrectability(*v);
474 if( corr > maxcorr)
475 {
476 //printf("%d nalv=%d correatability=%f\n",i, v->nalv, corr);
477 vbest=v;
478 maxcorr = corr;
479 }
480
481 }
482 return vbest;
483}

◆ t()

void EdbViewMap::t ( const char *  fin,
const char *  fout,
int  algorithm,
const char *  cut = "1" 
)
863{
864 TCut cut("cut",c);
865 EdbViewMap vm;
866 vm.ReadViewsHeaders(fin, cut); // read headers from runfile, fill eViewHeaders
867 vm.ReadPinViewsHeaders(fin); // read headers from runfile, fill ePinViewHeaders
868 vm.ReadViewAlign(fin); // read ViewAlign from runfile, fill eALcp
869 vm.eAlgorithm=algorithm;
870
871 if(algorithm==2) {
872 vm.MakeAHTnocorr(); // create tree with alignment parameters and views headers
873 vm.SaveLinesCorrToRun(fin); // correct hysteresis only
874 }
875 else if(algorithm==3) {
876 vm.MakeAHTnocorr(); // create tree with alignment parameters and views headers
877 vm.SaveColsCorrToRun(fin); // correct hysteresis only (colomns)
878 }
879 else {
880 vm.AddInverseAlign();
881
882 vm.InitAR(); // init eAR structure
883 vm.FillAR(); // fill eAR structure (assumed that view entry==h.GetViewID())
885
886 vm.FillALcp(); //fill neighbouring
887 //vm.CheckView(0,0,225);
888 //vm.CheckView(0,0,112);
889
890 vm.DoCorrection();
891 vm.CorrectToStage();
892 //vm.ConvertRun(fin,fout);
893 vm.SaveCorrToRun(fin);
894 }
895}
views map
Definition: EdbViewMap.h:71
void MakeAHTnocorr()
Definition: EdbViewMap.cxx:653
void DoCorrection()
Definition: EdbViewMap.cxx:502
void ReadViewsHeaders(const char *file, TCut &cut)
Definition: EdbViewMap.cxx:132
void FillALcp()
Definition: EdbViewMap.cxx:266
void CorrectToStage()
Definition: EdbViewMap.cxx:510
void SaveCorrToRun(const char *fin)
Definition: EdbViewMap.cxx:579
void CheckViewFrameAlignQuality(const char *file)
Definition: EdbViewMap.cxx:223
void ReadPinViewsHeaders(const char *file)
Definition: EdbViewMap.cxx:154
void InitAR()
Definition: EdbViewMap.cxx:45
void SaveLinesCorrToRun(const char *fin)
Definition: EdbViewMap.cxx:744
void AddInverseAlign()
Definition: EdbViewMap.cxx:203
void SaveColsCorrToRun(const char *fin)
Definition: EdbViewMap.cxx:803
void FillAR()
Definition: EdbViewMap.cxx:92
void ReadViewAlign(const char *file)
Definition: EdbViewMap.cxx:183

◆ Transform()

void EdbViewMap::Transform ( AView v,
const AlignmentParView apv 
)
316{
317 v.aff.ShiftX(apv.dx);
318 v.aff.ShiftY(apv.dy);
319}

◆ ViewCorrectability()

float EdbViewMap::ViewCorrectability ( AView v)
420{
421 if (eAlgorithm==0) return ViewCorrectability1(v);
422 else if(eAlgorithm==1) return ViewCorrectabilityAll(v);
423 return 0;
424}
float ViewCorrectability1(AView &v)
Definition: EdbViewMap.cxx:427
float ViewCorrectabilityAll(AView &v)
Definition: EdbViewMap.cxx:446

◆ ViewCorrectability1()

float EdbViewMap::ViewCorrectability1 ( AView v)
428{
429 float correctability=0;
430 if(!v.isCorrected)
431 {
432 for(int ia=0; ia<v.nalv; ia++ )
433 {
434 AApv *apv=v.alv[ia];
435 if(!apv) Log(1,"EdbViewMap::ViewCorrectability1","ERROR! apv = zero pointer");
436 if(apv->apv->view2<0) continue;
437 AView *v = GetAView( apv->apv->side2, apv->apv->area2, apv->apv->view2);
438 if(!v) { Log(1,"EdbViewMap::ViewCorrectability","ERROR! v(%d,%d,%d) = zero pointer",apv->apv->side2, apv->apv->area2, apv->apv->view2); continue;}
439 if( v->isCorrected ) if(Quality(*apv)>correctability) correctability=Quality(*apv);
440 }
441 }
442 return correctability;
443}

◆ ViewCorrectabilityAll()

float EdbViewMap::ViewCorrectabilityAll ( AView v)
447{
448 float correctability=0;
449 if(!v.isCorrected)
450 {
451 for(int ia=0; ia<v.nalv; ia++ )
452 {
453 AApv *apv=v.alv[ia];
454 if(!apv) Log(1,"EdbViewMap::ViewCorrectabilityAll","ERROR! apv = zero pointer");
455 if(apv->apv->view2<0) continue;
456 AView *v = GetAView( apv->apv->side2, apv->apv->area2, apv->apv->view2);
457 if(!v) { Log(1,"EdbViewMap::ViewCorrectability","ERROR! v(%d,%d,%d) = zero pointer",apv->apv->side2, apv->apv->area2, apv->apv->view2); continue;}
458 if( v->isCorrected ) correctability+=Quality(*apv);
459 }
460 }
461 return correctability;
462}

◆ ViewSide()

int EdbViewMap::ViewSide ( const EdbViewHeader vh) const
inline
115{ return vh.GetNframesTop()==0?0:1; }
Int_t GetNframesTop() const
Definition: EdbView.h:120

Member Data Documentation

◆ eAHT

TTree* EdbViewMap::eAHT

aligned headers tree

◆ eALcp

AlignmentParView** EdbViewMap::eALcp

read tree with alignment parameters here

◆ eAlgorithm

Int_t EdbViewMap::eAlgorithm

algorithm used for correction: 0-bestNeigbour, 1-allNeigbours, 2,3 - hysteresis corr

◆ eAR

ARun EdbViewMap::eAR

structures with pointers to the above containers

◆ eHX0even

TH2D* EdbViewMap::eHX0even

even lines x-offset

◆ eHX0odd

TH2D* EdbViewMap::eHX0odd

odd lines x-offset side=0

◆ eHX1even

TH2D* EdbViewMap::eHX1even

even lines x-offset

◆ eHX1odd

TH2D* EdbViewMap::eHX1odd

odd lines x-offset side=1

◆ eHY0even

TH2D* EdbViewMap::eHY0even

even lines x-offset

◆ eHY0odd

TH2D* EdbViewMap::eHY0odd

odd lines x-offset

◆ eHY1even

TH2D* EdbViewMap::eHY1even

even lines x-offset

◆ eHY1odd

TH2D* EdbViewMap::eHY1odd

odd lines x-offset

◆ eNcp

Int_t EdbViewMap::eNcp

number of objects in eALcp

◆ eNpvh

Int_t EdbViewMap::eNpvh

number of headers

◆ eNvh

Int_t EdbViewMap::eNvh

number of headers

◆ ePinPointID

Int_t EdbViewMap::ePinPointID

pin point ID to be used for this alignment pinID=(1000 + 2*pointID+side)

◆ ePinViewHeaders

TObjArray EdbViewMap::ePinViewHeaders

read headers for all pinned views here

◆ eViewHeaders

TObjArray EdbViewMap::eViewHeaders

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