FEDRA emulsion software from the OPERA Collaboration
EdbPlateAlignment Class Reference

plate-to-plate alignment More...

#include <EdbPlateAlignment.h>

Inheritance diagram for EdbPlateAlignment:
Collaboration diagram for EdbPlateAlignment:

Public Member Functions

void Align (EdbPattern &p1, EdbPattern &p2, float dz)
 
void CoarseAl (EdbPattern &p1, EdbPattern &p2)
 
void DoubletsFilterOut (EdbPattern &p1, EdbPattern &p2)
 
 EdbPlateAlignment ()
 
void FineAl (EdbPattern &p1, EdbPattern &p2)
 
void FineAlAff (EdbPattern &p1, EdbPattern &p2, EdbLayer &la1)
 
void ProduceReport ()
 
void RankCouples (TObjArray &arr1, TObjArray &arr2)
 
void SaveCouplesTree ()
 
void SetDoublets (float dx, float dy, float dtx, float dty)
 
void SetParCoarseAl (float zcorr, float dpos=300, float dang=0.015, float dz=122, float dphi=0.01)
 
void SetParFineAl ()
 
void SetParTestAl (float zcorr, float dz=500, float dphi=0.03)
 
void SetSigma (float spos, float sang)
 
void SlowAlignXY (EdbPattern &p1, EdbPattern &p2, EdbH2 &hxy, EdbH1 &hphi, const char *name="slowal")
 
void TestAl (EdbPattern &p1, EdbPattern &p2)
 
virtual ~EdbPlateAlignment ()
 
- Public Member Functions inherited from EdbAlignmentV
void AddSegCouple (EdbSegP *s1, EdbSegP *s2)
 
void ApplyLimitsOffset (float &xmin1, float &xmax1, float &xmin2, float &xmax2, float offsetMax)
 
Int_t CalcAffFull ()
 
Int_t CalcApplyFractMeanDiff ()
 
Int_t CalcApplyMeanDiff ()
 
Float_t CalcFractMeanDiff (int ivar, float fraction)
 
Float_t CalcMeanDiff (int ivar)
 
Float_t CalcMeanDZ (float tmin=0.1, float tmax=2.)
 
Float_t CalcMeanShr (float tmin=0.1, float tmax=2.)
 
Int_t CalculateAffTXTY (EdbAffine2D &aff)
 
Int_t CalculateAffTXTY (TObjArray &arr1, TObjArray &arr2, EdbAffine2D &aff)
 
Int_t CalculateAffTXTYTurn (TObjArray &arr1, TObjArray &arr2, EdbAffine2D &aff)
 
Int_t CalculateAffXY (EdbAffine2D &aff)
 
Int_t CalculateAffXY (TObjArray &arr1, TObjArray &arr2, EdbAffine2D &aff)
 
Int_t CalculateAffXYTurn (EdbAffine2D &aff)
 
Int_t CalculateAffXYTurn (TObjArray &arr1, TObjArray &arr2, EdbAffine2D &aff)
 
void CloseOutputFile ()
 
void Corr2Aff (EdbLayer &layer)
 
void Corr2Aff (EdbSegCorr &corr, EdbLayer &layer)
 
void CorrToCoG (int side, EdbPattern &p)
 functions alpplied to the individual patterns More...
 
void CorrToCoG (int side, TObjArray &p)
 
float CoupleQuality (EdbSegP &s1, EdbSegP &s2)
 
void DefineGuessCell (float xmin1, float xmax1, float ymin1, float ymax1, float xmin2, float xmax2, float ymin2, float ymax2, int np1, int np2, float binOK)
 Selector functions. More...
 
int DoubletsFilterOut (int checkview, TH2F *hxy=0, TH2F *htxty=0)
 
 EdbAlignmentV ()
 
void FillCell (int side, EdbPattern &pat)
 
void FillCell (int side, TObjArray &pat)
 
int FillCombinations ()
 
int FillCombinations (float dv[4], float dxMax, float dyMax, bool doFill)
 
void FillGuessCell (EdbPattern &p1, EdbPattern &p2, float binOK=1., float offsetMax=2000.)
 
void FillGuessCell (TObjArray &p1, TObjArray &p2, float binOK=1., float offsetMax=2000.)
 
void FillThetaHist (int side, EdbH2 &htxy)
 
Int_t FindCorrDiff (float dvsame[4], int side=0, int nlim=10)
 
Float_t FindDensityPeak (TArrayF &arr, float fraction)
 
Int_t FindDiff (TObjArray &arr1, TObjArray &arr2, float dvlim[4], float dvfound[4])
 
float FineCorrPhi (TObjArray &sel1, TObjArray &sel2)
 
float FineCorrZ ()
 
float FineCorrZ (TObjArray &sel1, TObjArray &sel2)
 
void HDistance (EdbPattern &p1, EdbPattern &p2, float dxMax, float dyMax)
 
void InitHphi (int n, float min, float max)
 
void InitHshr0 (int n, float min, float max)
 
void InitHshr1 (int n, float min, float max)
 
void InitHx (int n, float min, float max)
 
void InitHy (int n, float min, float max)
 
void InitHz (int n, float min, float max)
 
void InitOutputFile (const char *file="report_al.root", const char *option="RECREATE")
 IO, initialization and finalization functions. More...
 
void InitPatCellBin (int side, EdbPattern &pat, float binx, float biny)
 
void InitPatCellN (EdbCell2 &cell, EdbPattern &pat, int nx, int ny)
 
Bool_t IsInsideDVsame (EdbSegP &s1, EdbSegP &s2)
 
Int_t Ncoins (float dvlim[4], EdbH2 *hdxy=0, EdbH2 *hdtxty=0, TObjArray *sel1=0, TObjArray *sel2=0)
 
Int_t Ncp ()
 Functions applied to the selected parallel arrays. More...
 
int OptimiseVar1 (int side, int ivar, EdbH2 *hdxy=0, EdbH2 *hdtxy=0)
 
void OptimiseVar2 (int side1, int ivar1, int side2, int ivar2, EdbH2 &h12, EdbH2 *hdxy=0, EdbH2 *hdtxty=0)
 
int SelectBestCouple ()
 
int SelectIsolated ()
 
Bool_t SideOK (int side)
 
char * StrDVsame () const
 
float TX (int side, EdbSegP &s)
 
float TY (int side, EdbSegP &s)
 
Bool_t ValidCoinsidence (EdbSegP &s1, EdbSegP &s2, float dvlim[4], float dvfound[4])
 
float Var (int side, EdbSegP &s, int ivar)
 
float Var (int side, int iseg, int ivar)
 
float X (int side, EdbSegP &s)
 Correction parameters handling. More...
 
float Xmax (int side, EdbPattern &p)
 
float Xmax (int side, TObjArray &p)
 
float Xmin (int side, EdbPattern &p)
 
float Xmin (int side, TObjArray &p)
 
float Y (int side, EdbSegP &s)
 
float Ymax (int side, EdbPattern &p)
 
float Ymax (int side, TObjArray &p)
 
float Ymin (int side, EdbPattern &p)
 
float Ymin (int side, TObjArray &p)
 
virtual ~EdbAlignmentV ()
 

Public Attributes

Int_t eCoarseMin
 minimum coinsidences to accept alignment More...
 
Bool_t eCoarseOK
 
Bool_t eDoCoarse
 
Bool_t eDoCorrectBeforeSaving
 apply corrections before saving the couples tree in al.root More...
 
Bool_t eDoFine
 
Bool_t eDoTestAl
 
Float_t eDoublets [4]
 dx,dy,dtx,dty for the dublets cutout More...
 
Float_t eDPHI
 the range +- dphi will be scanned by coarce align More...
 
Float_t eDZ
 the range +- dz will be scanned by coarce align More...
 
Int_t eFineMin
 minimum coinsidences to accept alignment More...
 
Bool_t eFineOK
 
EdbPeak2 eH_xy_coarse
 
EdbPeak2 eH_xy_final
 the final alignment peak More...
 
EdbPeak2 eH_zphi_coarse
 the results of the coarse alignment More...
 
Int_t eNcoins
 final number of coinsidence used for affine calc More...
 
Bool_t eNoScale
 if 1 - disable scaling for calculation of affine transformations More...
 
Float_t eOffsetMax
 the maximal offset to be looked for More...
 
Bool_t eRankCouples
 rank couples More...
 
Bool_t eSaveCouples
 save couples tree with the report file More...
 
TObjArray eSegCouples
 segment couples objects to fill couples format tree More...
 
Float_t eSigma [2]
 sigma of the bt useful for the fine alignment ie:(10,0.01) More...
 
Bool_t eStatus
 overall alignment result (true - OK) More...
 
Bool_t eTestAlOK
 
- Public Attributes inherited from EdbAlignmentV
EdbSegCorr eCorr [2]
 corrections for side 1 and 2 (v[7]) - the result of the alignment More...
 
EdbLayer eCorrL [2]
 corrections in form of affine transformations - the final output More...
 
TH1I * eDoubletsRate
 can be filled in FillCombinations() More...
 
Float_t eDVsame [4]
 (dx,dy,dtx,dty) condition for the coinsidence More...
 
EdbH1 eH [2][7]
 
EdbH2 eHxy
 position 2d histo to be used in OptimiseVar2 More...
 
TFile * eOutputFile
 
EdbCell2 ePC [2]
 2-d position cells with patterns segments More...
 
TObjArray eS [2]
 "parallel" arrays with the selected combinations of segments More...
 
Bool_t eUseAffCorr
 if "true" - use eCorrL for corrections More...
 
Float_t eXmarg
 
Float_t eYmarg
 margins for the cell definition More...
 

Additional Inherited Members

- Static Public Member Functions inherited from EdbAlignmentV
static Int_t CheckEqualArr (TObjArray &arr1, TObjArray &arr2)
 

Detailed Description

plate-to-plate alignment

Constructor & Destructor Documentation

◆ EdbPlateAlignment()

EdbPlateAlignment::EdbPlateAlignment ( )
32{
33 SetSigma( 10, 0.005 );
34 eOffsetMax = 500.;
35 eDZ = 120;
36 eDPHI = 0.009;
37 eDoublets[0]=eDoublets[1]=5;
38 eDoublets[2]=eDoublets[3]=0.003;
39
40 eStatus =false;
41 eDoTestAl =false; eTestAlOK =true;
42 eDoCoarse =true; eCoarseOK =false;
43 eDoFine =true; eFineOK =false;
44 eSaveCouples=false;
45 eRankCouples=false;
47 eNoScale=0;
48
49 eNcoins = 0;
50 eCoarseMin = 5;
51 eFineMin = 5;
52}
Float_t eOffsetMax
the maximal offset to be looked for
Definition: EdbPlateAlignment.h:12
Bool_t eCoarseOK
Definition: EdbPlateAlignment.h:20
Bool_t eSaveCouples
save couples tree with the report file
Definition: EdbPlateAlignment.h:23
Bool_t eFineOK
Definition: EdbPlateAlignment.h:21
void SetSigma(float spos, float sang)
Definition: EdbPlateAlignment.h:53
Bool_t eDoFine
Definition: EdbPlateAlignment.h:21
Int_t eFineMin
minimum coinsidences to accept alignment
Definition: EdbPlateAlignment.h:27
Float_t eDZ
the range +- dz will be scanned by coarce align
Definition: EdbPlateAlignment.h:14
Int_t eCoarseMin
minimum coinsidences to accept alignment
Definition: EdbPlateAlignment.h:28
Bool_t eStatus
overall alignment result (true - OK)
Definition: EdbPlateAlignment.h:24
Float_t eDoublets[4]
dx,dy,dtx,dty for the dublets cutout
Definition: EdbPlateAlignment.h:16
Float_t eDPHI
the range +- dphi will be scanned by coarce align
Definition: EdbPlateAlignment.h:15
Int_t eNcoins
final number of coinsidence used for affine calc
Definition: EdbPlateAlignment.h:25
Bool_t eDoTestAl
Definition: EdbPlateAlignment.h:19
Bool_t eTestAlOK
Definition: EdbPlateAlignment.h:19
Bool_t eDoCoarse
Definition: EdbPlateAlignment.h:20
Bool_t eRankCouples
rank couples
Definition: EdbPlateAlignment.h:22
Bool_t eDoCorrectBeforeSaving
apply corrections before saving the couples tree in al.root
Definition: EdbPlateAlignment.h:35
Bool_t eNoScale
if 1 - disable scaling for calculation of affine transformations
Definition: EdbPlateAlignment.h:18

◆ ~EdbPlateAlignment()

EdbPlateAlignment::~EdbPlateAlignment ( )
virtual
56{
58}
void CloseOutputFile()
Definition: EdbAlignmentV.cxx:48

Member Function Documentation

◆ Align()

void EdbPlateAlignment::Align ( EdbPattern p1,
EdbPattern p2,
float  dz 
)
62{
63 eStatus = false;
64 eHxy.Delete();
68
69 Log(2,"EdbPlateAlignment::Align","pattern 1 with %d in limits x:(%12.2f %12.2f) y:(%12.2f %12.2f)",
70 p1.N(), p1.Xmin(),p1.Xmax(),p1.Ymin(),p1.Ymax() );
71 Log(2,"EdbPlateAlignment::Align","pattern 2 with %d in limits x:(%12.2f %12.2f) y:(%12.2f %12.2f)",
72 p2.N(), p2.Xmin(),p2.Xmax(),p2.Ymin(),p2.Ymax() );
73
74 if(p1.N()<1||p2.N()<1) return;
75 DoubletsFilterOut(p1,p2);
76
78 if(eDoTestAl) {
79 TestAl(p1,p2);
80 if(!eTestAlOK) goto END;
81 }
82
83 CorrToCoG(1,p2);
84
85 if(eDoCoarse) {
86 SetParCoarseAl( eCorr[0].V(2) , eOffsetMax, 3*eSigma[1], eDZ, eDPHI );
87 //SetParCoarseAl( eCorr[0].V(2) , 3*eSigma[0], 3*eSigma[1], eDZ, eDPHI );
88 CoarseAl(p1,p2);
89 if(!eCoarseOK) goto END;
90 }
91
92 if(eDoFine) {
94 for(int i=0; i<5; i++) FineAl(p1,p2);
95
96 eUseAffCorr=true;
97 Corr2Aff( eCorrL[0] );
98
99 for(int i=0; i<5; i++) FineAlAff(p1,p2,eCorrL[0] );
101 if(!eFineOK) goto END;
102 }
103
104 END:
105
106 if(eRankCouples) {
107// eCorr[0].SetV(2,0); eCorr[1].SetV(2, 0); ??
108 RankCouples( eS[0], eS[1] );
109 }
112}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
brick dz
Definition: RecDispMC.C:107
EdbH2 eHxy
position 2d histo to be used in OptimiseVar2
Definition: EdbAlignmentV.h:31
void CorrToCoG(int side, EdbPattern &p)
functions alpplied to the individual patterns
Definition: EdbAlignmentV.cxx:430
void Corr2Aff(EdbLayer &layer)
Definition: EdbAlignmentV.cxx:1027
Bool_t eUseAffCorr
if "true" - use eCorrL for corrections
Definition: EdbAlignmentV.h:23
EdbLayer eCorrL[2]
corrections in form of affine transformations - the final output
Definition: EdbAlignmentV.h:25
EdbSegCorr eCorr[2]
corrections for side 1 and 2 (v[7]) - the result of the alignment
Definition: EdbAlignmentV.h:24
TObjArray eS[2]
"parallel" arrays with the selected combinations of segments
Definition: EdbAlignmentV.h:21
void Delete()
Definition: EdbCell2.cpp:68
void Init(const EdbH2 &h, int npeaks=10)
Definition: EdbCell2.cpp:250
void Delete()
Definition: EdbCell2.cpp:243
EdbPeak2 eH_zphi_coarse
the results of the coarse alignment
Definition: EdbPlateAlignment.h:30
void FineAlAff(EdbPattern &p1, EdbPattern &p2, EdbLayer &la1)
Definition: EdbPlateAlignment.cxx:343
Float_t eSigma[2]
sigma of the bt useful for the fine alignment ie:(10,0.01)
Definition: EdbPlateAlignment.h:11
void TestAl(EdbPattern &p1, EdbPattern &p2)
Definition: EdbPlateAlignment.cxx:479
void SetParFineAl()
Definition: EdbPlateAlignment.cxx:532
void FineAl(EdbPattern &p1, EdbPattern &p2)
Definition: EdbPlateAlignment.cxx:382
void SetParCoarseAl(float zcorr, float dpos=300, float dang=0.015, float dz=122, float dphi=0.01)
Definition: EdbPlateAlignment.cxx:541
EdbPeak2 eH_xy_coarse
Definition: EdbPlateAlignment.h:31
void ProduceReport()
Definition: EdbPlateAlignment.cxx:133
void SetParTestAl(float zcorr, float dz=500, float dphi=0.03)
Definition: EdbPlateAlignment.cxx:567
void DoubletsFilterOut(EdbPattern &p1, EdbPattern &p2)
Definition: EdbPlateAlignment.cxx:328
void CoarseAl(EdbPattern &p1, EdbPattern &p2)
Definition: EdbPlateAlignment.cxx:429
void SaveCouplesTree()
Definition: EdbPlateAlignment.cxx:115
EdbPeak2 eH_xy_final
the final alignment peak
Definition: EdbPlateAlignment.h:32
void RankCouples(TObjArray &arr1, TObjArray &arr2)
Definition: EdbPlateAlignment.cxx:634
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
Int_t N() const
Definition: EdbPattern.h:86

◆ CoarseAl()

void EdbPlateAlignment::CoarseAl ( EdbPattern p1,
EdbPattern p2 
)

all pars should be already defined

430{
432
433 if( Log(3, "CoarseAl","with patterns %d %d", p1.N(), p2.N() ) ) {
434 eCorr[0].Print(); eCorr[1].Print();
435 }
436 FillGuessCell(p1,p2,1.,eOffsetMax);
438 int n[2] = { eH[0][0].N() , eH[0][1].N() };
439 float min[2] = { eH[0][0].Xmin(), eH[0][1].Xmin() };
440 float max[2] = { eH[0][0].Xmax(), eH[0][1].Xmax() };
441 eHxy.InitH2(n, min, max);
442
443 //OptimiseVar1(0, 6, &eHxy);
444 //OptimiseVar1(0, 2, &eHxy);
445
446 EdbH2 h_z_phi;
447 OptimiseVar2( 0, 6, 0, 2, h_z_phi, &eHxy ); // variate rotation and z of the p1
448
449 eH_zphi_coarse.Init(h_z_phi);
450
451 EdbPeak2 pk2z(h_z_phi);
452 float dphi, dz;
453 pk2z.FindPeak( dphi, dz );
454 pk2z.ProbPeak();
455 eCorr[0].SetV(2, dz );
456 eCorr[0].SetV(6, dphi);
457
459 Log(2,"CoarseAl","DVsame: %s",StrDVsame());
461
462 EdbPeak2 pk2(eHxy);
463 float dx, dy;
464 int npk = (int)pk2.FindPeak9( dx, dy );
465 eCorr[0].AddV(0,dx); eCorr[0].AddV(1,dy);
466 pk2.ProbPeak();
467 //pk2.Print();
468
469 Log(2,"CoarseAl","%d coins with peak of %d/%6.3f with offsets (dx:dy:dz:dphi): %7.2f %7.2f %8.2f %7.3f with patterns of: %d %d",
470 npk, (int)pk2.ePeak[0], pk2.eMean[0],
471 eCorr[0].V(0)-eCorr[1].V(0), eCorr[0].V(1)-eCorr[1].V(1), eCorr[0].V(2), eCorr[0].V(6),
472 p1.N(),p2.N() );
473
474 if(npk>3) eCoarseOK = true;
475 else eCoarseOK = false;
476}
float min(TClonesArray *t)
Definition: bitview.cxx:275
Int_t Ncoins(float dvlim[4], EdbH2 *hdxy=0, EdbH2 *hdtxty=0, TObjArray *sel1=0, TObjArray *sel2=0)
Definition: EdbAlignmentV.cxx:722
Float_t eDVsame[4]
(dx,dy,dtx,dty) condition for the coinsidence
Definition: EdbAlignmentV.h:16
EdbH1 eH[2][7]
Definition: EdbAlignmentV.h:27
void OptimiseVar2(int side1, int ivar1, int side2, int ivar2, EdbH2 &h12, EdbH2 *hdxy=0, EdbH2 *hdtxty=0)
Definition: EdbAlignmentV.cxx:323
char * StrDVsame() const
Definition: EdbAlignmentV.h:50
void FillGuessCell(EdbPattern &p1, EdbPattern &p2, float binOK=1., float offsetMax=2000.)
Definition: EdbAlignmentV.cxx:880
int FillCombinations()
Definition: EdbAlignmentV.cxx:208
int N() const
Definition: EdbCell1.h:48
float Xmax() const
Definition: EdbCell1.h:55
float Xmin() const
Definition: EdbCell1.h:54
fast 2-dim histogram class (used as a basis for EdbCell2)
Definition: EdbCell2.h:19
int InitH2(const EdbH2 &h)
Definition: EdbCell2.cpp:78
peak analyser for EdbH2
Definition: EdbCell2.h:105
void Print()
float V(int i)
Definition: EdbSegCorr.h:23
void SetV(int i, float x)
Definition: EdbSegCorr.h:21
void AddV(int i, float x)
Definition: EdbSegCorr.h:22
int max
Definition: check_shower.C:41

◆ DoubletsFilterOut()

void EdbPlateAlignment::DoubletsFilterOut ( EdbPattern p1,
EdbPattern p2 
)
329{
330 EdbAlignmentV adup;
331 for(int i=0; i<4; i++) adup.eDVsame[i]=eDoublets[i];
332
333 adup.FillGuessCell(p1,p1,1.,eOffsetMax);
334 adup.FillCombinations();
335 adup.DoubletsFilterOut(0); // assign flag -10 to the duplicated segments
336
337 adup.FillGuessCell(p2,p2,1.,eOffsetMax);
338 adup.FillCombinations();
339 adup.DoubletsFilterOut(0); // assign flag -10 to the duplicated segments
340}
universal basic alignment class
Definition: EdbAlignmentV.h:13
int DoubletsFilterOut(int checkview, TH2F *hxy=0, TH2F *htxty=0)
Definition: EdbAlignmentV.cxx:67

◆ FineAl()

void EdbPlateAlignment::FineAl ( EdbPattern p1,
EdbPattern p2 
)
383{
384 Log(3, "FineAl","with patterns %d %d", p1.N(), p2.N() );
385 FillGuessCell(p1,p2,1.,eOffsetMax);
386 float dxlim = eH[0][0].Xmax()-eH[0][0].Xmin();
387 float dylim = eH[0][1].Xmax()-eH[0][1].Xmin();
388
389 eDVsame[0] = dxlim/2.; eDVsame[1] = dylim/2.;
390 //eDVsame[0] = dxlim; eDVsame[1] = dylim;
391
392 //eDVsame[2]=eDVsame[3]=0.01;
393 eDVsame[2]=eDVsame[3]=3.*eSigma[1];
394 FillCombinations( eDVsame, dxlim, dylim, 1);
395
396 int n[2] = { 9,9 };
397 float min[2] = { eH[0][0].Xmin(), eH[0][1].Xmin() };
398 float max[2] = { eH[0][0].Xmax(), eH[0][1].Xmax() };
399 eHxy.InitH2(n, min, max);
400
401 for(int i=0; i<3; i++) FindCorrDiff(eDVsame);
402
403 TObjArray sel1,sel2;
404 Ncoins(eDVsame, &eHxy, 0, &sel1, &sel2);
405
406 for(int i=0; i<3; i++) {
407 float dz = FineCorrZ(sel1,sel2);
408 eCorr[0].AddV(2,dz);
409 }
410
411// for(int i=0; i<3; i++) {
412// float dphi = FineCorrPhi(sel1,sel2);
413// eCorr[0].AddV(6,dphi);
414// }
415
416 int npk=Ncoins(eDVsame, &eHxy);
417
418 Log(2,"FineAl","peak of %d with offsets (dx:dy:dz:dphi): %7.2f %7.2f %8.2f %7.3f with patterns of: %d %d",
419 npk,
420 eCorr[0].V(0)-eCorr[1].V(0), eCorr[0].V(1)-eCorr[1].V(1), eCorr[0].V(2), eCorr[0].V(6),
421 p1.N(),p2.N() );
422
423
424 if(npk>3) eFineOK = true;
425 else eFineOK = false;
426}
float FineCorrZ()
Definition: EdbAlignmentV.h:89
Int_t FindCorrDiff(float dvsame[4], int side=0, int nlim=10)
Definition: EdbAlignmentV.cxx:665

◆ FineAlAff()

void EdbPlateAlignment::FineAlAff ( EdbPattern p1,
EdbPattern p2,
EdbLayer la1 
)

assuming that the patterns are already aligned with the accuracy of a few
microns or the corrections are setted.
find affine transfornations

344{
348 Log(3, "FineAlAff","with patterns %d %d", p1.N(), p2.N() );
349 FillGuessCell(p1,p2,1.,eOffsetMax);
350 float dxlim = eH[0][0].Xmax()-eH[0][0].Xmin();
351 float dylim = eH[0][1].Xmax()-eH[0][1].Xmin();
352
353 eDVsame[0] = dxlim/2.; eDVsame[1] = dylim/2.;
354 //eDVsame[0] = dxlim; eDVsame[1] = dylim;
355 eDVsame[2]=eDVsame[3]=3.*eSigma[1];
356 FillCombinations( eDVsame, dxlim, dylim, 1);
357
358 int n[2] = { 9,9 };
359 float min[2] = { eH[0][0].Xmin(), eH[0][1].Xmin() };
360 float max[2] = { eH[0][0].Xmax(), eH[0][1].Xmax() };
361 eHxy.InitH2(n, min, max);
362
363 TObjArray sel1, sel2;
364 int npk= Ncoins(eDVsame, &eHxy, 0, &sel1, &sel2);
365 EdbAffine2D aff;
366 if(eNoScale) CalculateAffXYTurn(sel1,sel2,aff);
367 else CalculateAffXY(sel1,sel2,aff);
368 la1.GetAffineXY()->Transform(&aff);
369
370 EdbAffine2D afftxty;
371 CalculateAffTXTYTurn(sel1,sel2,afftxty);
372 la1.GetAffineTXTY()->Transform(&afftxty);
373
374 Log(2,"FineAlAff","peak of %d with patterns of: %d %d", npk, p1.N(),p2.N() );
375 if(gEDBDEBUGLEVEL>2) la1.GetAffineXY()->Print();
376
377 if(npk>3) eFineOK = true;
378 else eFineOK = false;
379}
Definition: EdbAffine.h:17
void Print(Option_t *opt="") const
Definition: EdbAffine.cxx:52
void Transform(const EdbAffine2D *a)
Definition: EdbAffine.cxx:93
Int_t CalculateAffXY(EdbAffine2D &aff)
Definition: EdbAlignmentV.h:95
Int_t CalculateAffTXTYTurn(TObjArray &arr1, TObjArray &arr2, EdbAffine2D &aff)
Definition: EdbAlignmentV.cxx:998
Int_t CalculateAffXYTurn(EdbAffine2D &aff)
Definition: EdbAlignmentV.h:94
EdbAffine2D * GetAffineTXTY()
Definition: EdbLayer.h:120
EdbAffine2D * GetAffineXY()
Definition: EdbLayer.h:119
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ ProduceReport()

void EdbPlateAlignment::ProduceReport ( )
134{
135 TH2F *zphiCoarse=0, *xyCoarse=0;
136 TH2F *xyFine=0;
137 gStyle->SetPalette(1);
138 if( eDoTestAl&&eTestAlOK ) {
139 Log(1,"Alignment Report","TestAl OK");
140 }
141 if( eDoCoarse&&eCoarseOK ) {
142 zphiCoarse = eH_zphi_coarse.DrawH2("zphi_coarse", "Z vs Phi after coarse align");
143 xyCoarse = eH_xy_coarse.DrawH2("xy_coarse", "dY vs dX after coarse align");
144 float dx,dy;
145 int npk = (int)eH_xy_coarse.FindPeak9( dx, dy );
147
148 if (npk<eCoarseMin) eStatus = false; // too few coinsidences
149 else if(npk < 3*Sqrt(eH_xy_coarse.eMean[0]) ) eStatus = false; // too small signal/noise
150 else eStatus = true; // good alignment
151
152 if(eStatus)
153 Log(2,"Alignment Report","Coarse %s: %d coins with peak of %d/%6.3f","OK",
154 npk, (int)eH_xy_coarse.ePeak[0], eH_xy_coarse.eMean[0]);
155 else {
156 Log(1,"Alignment Report","Coarse %s: %d coins with peak of %d/%6.3f","FAILED!",
157 npk, (int)eH_xy_coarse.ePeak[0], eH_xy_coarse.eMean[0]);
158 return;
159 }
160 }
161
162
163 if( eDoFine&&eFineOK ) {
164 xyFine = eH_xy_final.DrawH2("xy_final","dY vs dX after final alignment");
166 if( eNcoins < eFineMin ) eStatus = false; // too few coinsidences
167 else eStatus = true; // good alignment
168
169 if(eStatus)
170 Log(2,"Alignment Report","Fine %s: %d coins inside (%6.1f %6.1f) with angular acc: (%6.3f %6.3f)", "OK",
171 (int)(eH_xy_final.Integral()),
173 eDVsame[2], eDVsame[3] );
174 else {
175 Log(2,"Alignment Report","Fine %s: %d coins inside (%6.1f %6.1f) with angular acc: (%6.3f %6.3f)", "FAILED!",
176 (int)(eH_xy_final.Integral()),
178 eDVsame[2], eDVsame[3] );
179 return;
180 }
181 }
182
183
184 if(eOutputFile) {
185 Log(2,"Alignment Report","Save to file %s", eOutputFile->GetName());
186 gStyle->SetPalette(1);
187 bool batch = gROOT->IsBatch();
188 gROOT->SetBatch();
189
190 TCanvas *cc = new TCanvas("cc","Alignment report",800,800);
191
192 EdbAffine2D *aXY = eCorrL[0].GetAffineXY();
193 float xcenter1 = (ePC[0].Xmin()+ePC[0].Xmax())/2.;
194 float ycenter1 = (ePC[0].Ymin()+ePC[0].Ymax())/2.;
195 float xoffset= aXY->A11()*xcenter1 + aXY->A12()*ycenter1 + aXY->B1() - xcenter1;
196 float yoffset= aXY->A21()*xcenter1 + aXY->A22()*ycenter1 + aXY->B2() - ycenter1;
197
198 eCorrL[0].SetXY( xcenter1, ycenter1 );
199
200 const char *str = Form( "Nfinal= %5d Peak: %5d/%6.3f dx,dy,dz = %7.3f %7.3f %7.3f",
202 xoffset, yoffset, eCorrL[0].Zcorr() );
203 Log(1,"Alignment Report","%s", str);
204
205 TPaveText *ctit = new TPaveText(0.01,0.943,0.99,0.998);
206 ctit->AddText( Form("Alignment of %s",eOutputFile->GetName()) );
207 ctit->AddText( str );
208 ctit->Draw();
209
210 TPad *c = new TPad("c","plots",0.01,0.05,0.99,0.94);
211 c->Divide(3,3);
212 c->Draw();
213 TH2F *hp1 = ePC[0].DrawH2("hp1","Pattern1"); c->cd(1); hp1->SetStats(0); hp1->Draw("colz");
214 TH2F *hp2 = ePC[1].DrawH2("hp2","Pattern2"); c->cd(2); hp2->SetStats(0); hp2->Draw("colz");
215 if(zphiCoarse) { c->cd(3); zphiCoarse->SetStats(0); zphiCoarse->Draw("colz"); }
216 if(xyCoarse) { c->cd(4); xyCoarse->SetStats(0); xyCoarse->Draw("colz"); }
217 if(xyFine) { c->cd(5); xyFine->SetStats(0); xyFine->Draw("colz"); }
218
219
220 float xmin = Min( ePC[0].Xmin(), ePC[1].Xmin() );
221 float xmax = Max( ePC[0].Xmax(), ePC[1].Xmax() );
222 float ymin = Min( ePC[0].Ymin(), ePC[1].Ymin() );
223 float ymax = Max( ePC[0].Ymax(), ePC[1].Ymax() );
224
225 Log(2,"Alignment Report","xmax-xmin: %f (%f-%f) ymax-ymin: %f (%f-%f)",
226 (xmax-xmin),xmax,xmin, (ymax-ymin),ymax,ymin);
227
228
229 int n = CheckEqualArr(eS[0],eS[1]);
230
231 c->cd(6);
232
233 TH2I *hdtxy = new TH2I("dTXdTY","dTY vs dTX after final alignment",25,-3*eSigma[1], 3*eSigma[1], 25, -3*eSigma[1], 3*eSigma[1]);
234 for(int i=0; i<n; i++) {
235 EdbSegP *s1 = (EdbSegP*)eS[0].UncheckedAt(i);
236 EdbSegP *s2 = (EdbSegP*)eS[1].UncheckedAt(i);
237 float dtx = TX(1,*s2) - TX(0,*s1);
238 float dty = TY(1,*s2) - TY(0,*s1);
239 hdtxy->Fill(dtx,dty);
240 }
241 hdtxy->SetStats(0); hdtxy->Draw("colz");
242
243 c->cd(8);
244
245 TH2I *h2 = new TH2I("plateXY","plateXY",100,xmin-10, xmax+10, 100, ymin-10,ymax+10);
246 h2->SetStats(0);
247 h2->Draw();
248
249 if(n<1000) {
250 for(int i=0; i<n; i++) {
251 EdbSegP *s1 = (EdbSegP*)eS[0].UncheckedAt(i);
252 EdbSegP *s2 = (EdbSegP*)eS[1].UncheckedAt(i);
253 float x1 = X(0,*s1), y1 = Y(0,*s1);
254 float x2 = X(1,*s2), y2 = Y(1,*s2);
255 TMarker *m1 = new TMarker(x1,y1,22);
256 m1->SetMarkerColor(kRed);
257 m1->SetMarkerSize(0.6);
258 m1->Draw();
259 TMarker *m2 = new TMarker(x2,y2,23);
260 m2->SetMarkerColor(kBlue);
261 m2->SetMarkerSize(0.6);
262 m2->Draw();
263 }
264 }
265
266 c->cd(9);
267 TH2I *h2t = new TH2I("plateTXTY","plateTXTY",100,-1, 1, 100, -1,1 );
268 h2t->SetStats(0);
269 h2t->Draw();
270 if(n<1000) {
271 for(int i=0; i<n; i++) {
272 EdbSegP *s1 = (EdbSegP*)eS[0].UncheckedAt(i);
273 EdbSegP *s2 = (EdbSegP*)eS[1].UncheckedAt(i);
274 float tx1 = TX(0, *s1), ty1 = TY(0, *s1);
275 float tx2 = TX(1, *s2), ty2 = TY(1, *s2);
276 TMarker *m1 = new TMarker(tx1,ty1,22);
277 m1->SetMarkerColor(kRed);
278 m1->SetMarkerSize(0.6);
279 m1->Draw();
280 TMarker *m2 = new TMarker(tx2,ty2,23);
281 m2->SetMarkerColor(kBlue);
282 m2->SetMarkerSize(0.6);
283 m2->Draw();
284 }
285 }
286
287
288 c->cd(7);
289 int nbin=25;
290 TArrayF xbins(nbin+1);
291 xbins[0]= 0.15;
292 float s0 = xbins[0]*xbins[0];
293 for(int i=1; i<nbin+1; i++) xbins[i] = Sqrt( xbins[i-1]*xbins[i-1] + s0 );
294
295 // TH1F *thetadens = new TH1F("thetadens","theta density",21,0, 0.7 );
296 TH1F *thetadens = new TH1F("thetadens","theta density",nbin,xbins.GetArray() );
297 for(int i=0; i<n; i++) {
298 //EdbSegP *s1 = (EdbSegP*)eS[0].UncheckedAt(i);
299 EdbSegP *s2 = (EdbSegP*)eS[1].UncheckedAt(i);
300 float tx = TX(1,*s2), ty = TY(1,*s2);
301 float theta = Sqrt(tx*tx+ty*ty);
302 thetadens->Fill(theta);
303 }
304 thetadens->Draw();
305
306 //--
307
308 eH_xy_coarse.Write("peak2c");
309 if(zphiCoarse) zphiCoarse->Write();
310 if(xyCoarse) xyCoarse->Write();
311 if(xyFine) xyFine->Write();
312 cc->Write("report_al");
313 eCorrL[0].Write("corr_layer1");
314 eCorrL[1].Write("corr_layer2");
315 SafeDelete(c);
316
317 SafeDelete(cc);
318
319 gROOT->SetBatch(batch);
320 }
321
322 SafeDelete(zphiCoarse);
323 SafeDelete(xyCoarse);
324 SafeDelete(xyFine);
325}
cout<< tr-> GetEntries()<< endl
Float_t B2() const
Definition: EdbAffine.h:48
Float_t A22() const
Definition: EdbAffine.h:46
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
float Y(int side, EdbSegP &s)
Definition: EdbAlignmentV.h:119
static Int_t CheckEqualArr(TObjArray &arr1, TObjArray &arr2)
Definition: EdbAlignmentV.cxx:389
float TY(int side, EdbSegP &s)
Definition: EdbAlignmentV.h:121
float Ymin(int side, EdbPattern &p)
Definition: EdbAlignmentV.cxx:1065
TFile * eOutputFile
Definition: EdbAlignmentV.h:35
EdbCell2 ePC[2]
2-d position cells with patterns segments
Definition: EdbAlignmentV.h:18
float Xmax(int side, EdbPattern &p)
Definition: EdbAlignmentV.cxx:1054
float X(int side, EdbSegP &s)
Correction parameters handling.
Definition: EdbAlignmentV.h:118
float Ymax(int side, EdbPattern &p)
Definition: EdbAlignmentV.cxx:1076
float TX(int side, EdbSegP &s)
Definition: EdbAlignmentV.h:120
float Xmin(int side, EdbPattern &p)
Definition: EdbAlignmentV.cxx:1043
float Xmax() const
Definition: EdbCell2.h:65
float Ymin() const
Definition: EdbCell2.h:66
Long_t Integral()
Definition: EdbCell2.cpp:216
TH2F * DrawH2(const char *name="plot2d", const char *title="EdbH2plot2D")
Definition: EdbCell2.cpp:187
float Xmin() const
Definition: EdbCell2.h:64
float Ymax() const
Definition: EdbCell2.h:67
float Zcorr() const
Definition: EdbLayer.h:90
void SetXY(float x, float y)
Definition: EdbLayer.h:96
TArrayF eMean
Definition: EdbCell2.h:111
TArrayF ePeak
Definition: EdbCell2.h:109
float ProbPeak()
Definition: EdbCell2.cpp:287
float FindPeak9(float &x, float &y)
Definition: EdbCell2.cpp:345
Definition: EdbSegP.h:21
TH1F * h2
Definition: energy.C:19
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30
new TCanvas()

◆ RankCouples()

void EdbPlateAlignment::RankCouples ( TObjArray &  arr1,
TObjArray &  arr2 
)
635{
636 int n = arr1.GetEntries();
637 Log(3,"RankCouples","%d" ,n);
638
640 EdbScanCond cond;
641 cond.SetSigma0(eSigma[0],eSigma[0],eSigma[1],eSigma[1]);
642 eSegCouples.Delete();
643 EdbSegP seg, seg1, seg2;
644 for(int i=0; i<n; i++) {
645 EdbSegP *s1 = ((EdbSegP*)arr1.UncheckedAt(i));
646 EdbSegP *s2 = ((EdbSegP*)arr2.UncheckedAt(i));
647
648 seg.Copy(*s1); // to set correctly vid, aid, etc
649 seg1.Copy(*s1);
650 seg2.Copy(*s2);
651
652 eCorr[0].ApplyCorrections(seg1);
653 eCorr[1].ApplyCorrections(seg2);
654
655 tf.Chi2PSeg( seg1, seg2, seg, cond, cond);
656
657 s1->SetFlag(0);
658 s2->SetFlag(0);
659 seg.SetFlag(0);
660 seg.SetSide(0);
661 seg.SetVolume(seg1.Volume()+seg2.Volume());
662 seg.SetDZ( Abs(seg1.eZ - seg2.eZ) );
663
664 seg.SetDZem( seg1.Chi2() + seg2.Chi2() ); // HACK: use eDZem variable to keep microtracking Likelihood
665
666 EdbSegCouple *sc=new EdbSegCouple();
667 sc->eS1=s1;
668 sc->eS2=s2;
669 sc->eS = new EdbSegP(seg);
670 sc->SetCHI2P( seg.Chi2() );
671 eSegCouples.Add(sc);
672 }
673
674 EdbSegCouple::SetSortFlag(0); // sort by CHI2P
675 eSegCouples.UnSort();
676 eSegCouples.Sort();
677
678 int ncp = eSegCouples.GetEntries();
679
680 for(int i=0; i<ncp; i++) {
681 EdbSegCouple *sc = (EdbSegCouple*)(eSegCouples.UncheckedAt(i));
682 sc->eS1->SetFlag( sc->eS1->Flag()+1 );
683 sc->eS2->SetFlag( sc->eS2->Flag()+1 );
684 sc->SetN1(sc->eS1->Flag());
685 sc->SetN2(sc->eS2->Flag());
686 }
687
688 for(int i=0; i<ncp; i++) {
689 EdbSegCouple *sc = (EdbSegCouple*)(eSegCouples.UncheckedAt(i));
690 sc->SetN1tot(sc->eS1->Flag());
691 sc->SetN2tot(sc->eS2->Flag());
692 }
693
694 Log(2,"RankCouples","%d couples ok", ncp );
695}
TObjArray eSegCouples
segment couples objects to fill couples format tree
Definition: EdbPlateAlignment.h:34
Definition: EdbScanCond.h:10
void SetSigma0(float x, float y, float tx, float ty)
Definition: EdbScanCond.h:62
void ApplyCorrections(EdbSegP &s)
Definition: EdbSegCorr.cxx:13
Definition: EdbSegCouple.h:17
void SetN1(int n1)
Definition: EdbSegCouple.h:43
void SetN2(int n2)
Definition: EdbSegCouple.h:44
EdbSegP * eS
the result of the fit
Definition: EdbSegCouple.h:27
void SetN1tot(int n)
Definition: EdbSegCouple.h:45
EdbSegP * eS2
Definition: EdbSegCouple.h:29
void SetN2tot(int n)
Definition: EdbSegCouple.h:46
void SetCHI2P(float chi2)
Definition: EdbSegCouple.h:48
static void SetSortFlag(int s=0)
Definition: EdbSegCouple.cxx:62
EdbSegP * eS1
pointers - useful when all segments are in memory
Definition: EdbSegCouple.h:28
void SetVolume(float w)
Definition: EdbSegP.h:136
void SetSide(int side=0)
Definition: EdbSegP.h:139
void SetDZem(float dz)
Definition: EdbSegP.h:127
Float_t eZ
coordinates
Definition: EdbSegP.h:31
Float_t Volume() const
Definition: EdbSegP.h:158
Float_t Chi2() const
Definition: EdbSegP.h:157
void SetDZ(float dz)
Definition: EdbSegP.h:126
void Copy(const EdbSegP &s)
Definition: EdbSegP.cxx:105
void SetFlag(int flag)
Definition: EdbSegP.h:130
Int_t Flag() const
Definition: EdbSegP.h:149
Definition: EdbTrackFitter.h:17
static float Chi2PSeg(EdbSegP &s1, EdbSegP &s2, EdbSegP &seg, EdbScanCond &cond1, EdbScanCond &cond2)
Definition: EdbTrackFitter.cxx:263

◆ SaveCouplesTree()

void EdbPlateAlignment::SaveCouplesTree ( )
116{
117 EdbCouplesTree ect;
118 ect.InitCouplesTree("couples",0,"NEW");
119 int nseg = CheckEqualArr(eS[0],eS[1]);
120 for(int i=0; i<nseg; i++)
121 if(eDoCorrectBeforeSaving) eCorrL[0].CorrectSeg(*(EdbSegP*)eS[0].UncheckedAt(i));
122 for(int i=0; i<nseg; i++) {
123 if( eRankCouples ) {
125 ect.Fill( sc->eS1, sc->eS2, sc->eS, sc );
126 }
127 else ect.Fill( (EdbSegP*)eS[0].UncheckedAt(i), (EdbSegP*)eS[1].UncheckedAt(i) );
128 }
129 ect.eTree->AutoSave();
130}
Definition: EdbCouplesTree.h:17
bool InitCouplesTree(const char *name="couples", const char *fname=0, Option_t *mode="READ")
Definition: EdbCouplesTree.cxx:79
TTree * eTree
couples tree
Definition: EdbCouplesTree.h:25
Int_t Fill(EdbSegP *s1, EdbSegP *s2, EdbSegP *s=0, EdbSegCouple *cp=0, float xv=0, float yv=0, int pid1=0, int pid2=0)
Definition: EdbCouplesTree.cxx:172
void CorrectSeg(EdbSegP &s)
Definition: EdbLayer.cxx:69

◆ SetDoublets()

void EdbPlateAlignment::SetDoublets ( float  dx,
float  dy,
float  dtx,
float  dty 
)
inline
55 { eDoublets[0]=dx; eDoublets[1]=dy; eDoublets[2]=dtx; eDoublets[3]=dty; }

◆ SetParCoarseAl()

void EdbPlateAlignment::SetParCoarseAl ( float  zcorr,
float  dpos = 300,
float  dang = 0.015,
float  dz = 122,
float  dphi = 0.01 
)

medium-wide parameters setting for the coarse alignment

542{
544 Log(3,"EdbPlateAlignment::SetParCoarseAl","zcorr = %f dpos = %f dang=%f dz=%f dphi=%f",zcorr, dpos, dang, dz, dphi);
545 eDVsame[0] = eDVsame[1] = dpos;
546 //eDVsame[0] = eDVsame[1] = 3*eSigma[0];
547 eDVsame[2] = eDVsame[3] = dang;
548 int n = (int)( 2*dpos/ (4*eSigma[0]) );
549 n = n%2 ? n : n+1; // make n - odd
550 Log(3, "SetParCoarseAl"," n = %d dpos = %f ", n, dpos);
551 InitHx( n, -dpos, dpos );
552 InitHy( n, -dpos, dpos );
553
554 float zstepmin=0.5; // 0.5 microns is the minimum step size
555 int nz = Min( (int)( 2*dz/ zstepmin ), 51);
556 nz = nz%2 ? nz : nz+1; // make n - odd
557 InitHz( nz, zcorr-dz, zcorr+dz );
558 eCorr[0].SetV(2,zcorr);
559
560 float phimin=0.00001; // minimum step size
561 int nphi = Min( (int)( 2*dphi/ phimin ), 51);
562 nphi = nphi%2 ? nphi : nphi+1; // make n - odd
563 InitHphi( nphi, eCorr[0].V(6)-dphi, eCorr[0].V(6)+dphi );
564}
void InitHz(int n, float min, float max)
Definition: EdbAlignmentV.h:46
void InitHphi(int n, float min, float max)
Definition: EdbAlignmentV.h:47
void InitHx(int n, float min, float max)
Definition: EdbAlignmentV.h:44
void InitHy(int n, float min, float max)
Definition: EdbAlignmentV.h:45

◆ SetParFineAl()

void EdbPlateAlignment::SetParFineAl ( )

narrow parameters setting for the fine alignment

533{
535 InitHx( 1, -4.*eSigma[0], 4.*eSigma[0] );
536 InitHy( 1, -4.*eSigma[0], 4.*eSigma[0] );
537}

◆ SetParTestAl()

void EdbPlateAlignment::SetParTestAl ( float  zcorr,
float  dz = 500,
float  dphi = 0.03 
)

wide parameters setting for the test pre-alignment

568{
570 InitHx( 81, -2000, 2000 );
571 InitHy( 81, -2000, 2000 );
572 int nbinz = Min(101, Max(21, (int)(dz/50.) ));
573 InitHz( nbinz, zcorr-dz, zcorr+dz );
574 int nbinphi = Min(101, Max(21, (int)(dphi/0.0001) ));
575 InitHphi( nbinphi, -dphi, dphi );
576 eCorr[0].SetV(2,zcorr);
577}

◆ SetSigma()

void EdbPlateAlignment::SetSigma ( float  spos,
float  sang 
)
inline
53{ eSigma[0]=spos; eSigma[1]=sang; }

◆ SlowAlignXY()

void EdbPlateAlignment::SlowAlignXY ( EdbPattern p1,
EdbPattern p2,
EdbH2 hxy,
EdbH1 hphi,
const char *  name = "slowal" 
)
581{
582 int n1= pat1.N();
583 int n2= pat2.N();
584 float x0 = (Min(pat1.Xmin(),pat2.Xmin())+Max(pat1.Xmax(),pat2.Xmax()))/2;
585 float y0 = (Min(pat1.Ymin(),pat2.Ymin())+Max(pat1.Ymax(),pat2.Ymax()))/2;
586 struct point{ int i; float x; float y; };
587 point **p1 = new point*[n1];
588 point **p2 = new point*[n2];
589
590
591 TFile f( Form("%s.root",name),"RECREATE");
592 TNtuple *ntpeak = new TNtuple("ntpeak","","phi:dx:dy:peak");
593
594 for(int i=0; i<n1; i++) p1[i] = new point;
595 for(int i=0; i<n2; i++) p2[i] = new point;
596
597 for(int iphi=0; iphi<hphi.N(); iphi++) {
598 float phi=hphi.X(iphi);
599 for(int i=0; i<n1; i++) {
600 EdbSegP *s= pat1.GetSegment(i);
601 p1[i]->x = (s->X()-x0)*Cos(phi) - (s->Y()-y0)*Sin(phi);
602 p1[i]->y = (s->X()-x0)*Sin(phi) + (s->Y()-y0)*Cos(phi);
603 p1[i]->i = i;
604 }
605 for(int i=0; i<n2; i++) {
606 EdbSegP *s= pat2.GetSegment(i);
607 p2[i]->x = (s->X()-x0);
608 p2[i]->y = (s->Y()-y0);
609 p2[i]->i = i;
610 }
611
612 for(int i1=0; i1<n1; i1++) {
613 for(int i2=0; i2<n2; i2++) {
614 float dx= p2[i2]->x-p1[i1]->x;
615 float dy= p2[i2]->y-p1[i1]->y;
616 hxy.Fill(dx,dy);
617 }
618 }
619 EdbPeak2 pk2(hxy);
620 float dx, dy;
621 int npk = pk2.FindPeak( dx, dy );
622 printf("phi = %f dx = %f dy = %f npk = %d\n",phi,dx,dy,npk);
623 hphi.Fill(phi,npk);
624 ntpeak->Fill(phi,dx,dy,npk);
625 hxy.Write(Form("hxy%d",iphi));
626 hxy.CleanCells();
627 }
628 hphi.Write("hphi");
629 ntpeak->Write();
630 f.Close();
631}
FILE * f
Definition: RecDispMC.C:150
float X(int i) const
Definition: EdbCell1.h:53
int Fill(float x)
Definition: EdbCell1.h:63
void CleanCells()
Definition: EdbCell2.cpp:114
int Fill(float x, float y)
Definition: EdbCell2.h:88
s
Definition: check_shower.C:55
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24

◆ TestAl()

void EdbPlateAlignment::TestAl ( EdbPattern p1,
EdbPattern p2 
)

alignment of 2 plates with testal-like preselection.
Application: basetracks patterns with: well defined angles, big possible offsets,

big possible rotations; z-variation
Result: find best dz and rotation and position offset and set them as:

eCorr[0].SetV(2, dz), eCorr[0].SetV(6, dphi); eCorr[0].AddV(0, dx); eCorr[0].AddV(1, dy);

480{
486
487 float dxlim = eH[0][0].Xmax()-eH[0][0].Xmin();
488 float dylim = eH[0][1].Xmax()-eH[0][1].Xmin();
489 HDistance(p1,p2, dxlim, dylim );
490
491 float x0 = p2.Xmean(), y0 = p2.Ymean(); // centralize patterns
492 eCorr[0].SetV(0,-x0); eCorr[0].SetV(1,-y0);
493 eCorr[1].SetV(0,-x0); eCorr[1].SetV(1,-y0);
494
495 int n[2] = { eH[0][0].N() , eH[0][1].N() };
496 float min[2] = { eH[0][0].Xmin(), eH[0][1].Xmin() };
497 float max[2] = { eH[0][0].Xmax(), eH[0][1].Xmax() };
498 eHxy.InitH2(n, min, max);
499
500 eDVsame[0] = dxlim/2.;
501 eDVsame[1] = dylim/2.;
502 eDVsame[2]=eDVsame[3]=eSigma[1];
503
504 EdbH2 h_z_phi;
505 OptimiseVar2( 0, 6, 0, 2, h_z_phi, &eHxy ); // variate rotation and z of the p1
506
507 h_z_phi.DrawH2()->Draw("colz");
508
509 EdbPeak2 pk2z(h_z_phi);
510 float dphi, dz;
511 pk2z.FindPeak( dphi, dz );
512 pk2z.ProbPeak();
513
514 eCorr[0].SetV(2, dz );
515 eCorr[0].SetV(6, dphi);
516
518 EdbPeak2 pk2(eHxy);
519 float dx, dy;
520 pk2.FindPeak( dx, dy );
521 pk2.ProbPeak();
522 eCorr[0].AddV(0,dx); // set x-y corrections
523 eCorr[0].AddV(1,dy);
524
525 Log(2,"TestAl","peak of %d/%6.3f with offsets (dx:dy:dz:dphi): %7.2f %7.2f %8.2f %7.3f with patterns of: %d %d",
526 (int)pk2.ePeak[0], pk2.eMean[0],
527 eCorr[0].V(0)-eCorr[1].V(0), eCorr[0].V(1)-eCorr[1].V(1), eCorr[0].V(2), eCorr[0].V(6),
528 p1.N(),p2.N() );
529}
void HDistance(EdbPattern &p1, EdbPattern &p2, float dxMax, float dyMax)
Definition: EdbAlignmentV.cxx:752
float Xmean()
Definition: EdbPattern.cxx:1517
float Ymean()
Definition: EdbPattern.cxx:1527

Member Data Documentation

◆ eCoarseMin

Int_t EdbPlateAlignment::eCoarseMin

minimum coinsidences to accept alignment

◆ eCoarseOK

Bool_t EdbPlateAlignment::eCoarseOK

◆ eDoCoarse

Bool_t EdbPlateAlignment::eDoCoarse

◆ eDoCorrectBeforeSaving

Bool_t EdbPlateAlignment::eDoCorrectBeforeSaving

apply corrections before saving the couples tree in al.root

◆ eDoFine

Bool_t EdbPlateAlignment::eDoFine

◆ eDoTestAl

Bool_t EdbPlateAlignment::eDoTestAl

◆ eDoublets

Float_t EdbPlateAlignment::eDoublets[4]

dx,dy,dtx,dty for the dublets cutout

◆ eDPHI

Float_t EdbPlateAlignment::eDPHI

the range +- dphi will be scanned by coarce align

◆ eDZ

Float_t EdbPlateAlignment::eDZ

the range +- dz will be scanned by coarce align

◆ eFineMin

Int_t EdbPlateAlignment::eFineMin

minimum coinsidences to accept alignment

◆ eFineOK

Bool_t EdbPlateAlignment::eFineOK

◆ eH_xy_coarse

EdbPeak2 EdbPlateAlignment::eH_xy_coarse

◆ eH_xy_final

EdbPeak2 EdbPlateAlignment::eH_xy_final

the final alignment peak

◆ eH_zphi_coarse

EdbPeak2 EdbPlateAlignment::eH_zphi_coarse

the results of the coarse alignment

◆ eNcoins

Int_t EdbPlateAlignment::eNcoins

final number of coinsidence used for affine calc

◆ eNoScale

Bool_t EdbPlateAlignment::eNoScale

if 1 - disable scaling for calculation of affine transformations

◆ eOffsetMax

Float_t EdbPlateAlignment::eOffsetMax

the maximal offset to be looked for

◆ eRankCouples

Bool_t EdbPlateAlignment::eRankCouples

rank couples

◆ eSaveCouples

Bool_t EdbPlateAlignment::eSaveCouples

save couples tree with the report file

◆ eSegCouples

TObjArray EdbPlateAlignment::eSegCouples

segment couples objects to fill couples format tree

◆ eSigma

Float_t EdbPlateAlignment::eSigma[2]

sigma of the bt useful for the fine alignment ie:(10,0.01)

◆ eStatus

Bool_t EdbPlateAlignment::eStatus

overall alignment result (true - OK)

◆ eTestAlOK

Bool_t EdbPlateAlignment::eTestAlOK

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