FEDRA emulsion software from the OPERA Collaboration
EdbPositionAlignment Class Reference

new alignment class developed mainly for compton search More...

#include <EdbPositionAlignment.h>

Inheritance diagram for EdbPositionAlignment:
Collaboration diagram for EdbPositionAlignment:

Public Member Functions

bool ActivatePosTree (const char *name="postree")
 
int Align ()
 
void AlignWithTracks (int nang=3, int npos=0)
 
int DoubletsFilterOut (bool checkview=true)
 
int DoubletsFilterOutSide (EdbPatCell2 &pc, bool checkview=true)
 
 EdbPositionAlignment ()
 
int FillArrays (EdbPattern &p1, EdbPattern &p2)
 
int FillArrays (TObjArray &arr1, TObjArray &arr2, float xymin[2], float xymax[2])
 
int FillCombinations ()
 
int FillCombinations (EdbPatCell2 &pc1, EdbPatCell2 &pc2, float dx, float dy, float dtx, float dty)
 
int FillPosTree (float dz1, float dz2, int flag)
 
void FindCorrections (EdbPattern &pat1, EdbPattern &pat2, float DZ, bool doShrink)
 
void FindDiff (TObjArray &arr1, TObjArray &arr2, EdbPatCell2 &pc1, EdbPatCell2 &pc2, float &dx, float &dy, float &dtx, float &dty)
 
void FindDiff10 (float &dx, float &dy, float &dtx, float &dty)
 
void FindDiff12 (float &dx, float &dy, float &dtx, float &dty)
 
void FindDiff20 (float &dx, float &dy, float &dtx, float &dty)
 
int Link ()
 
void PositionPlot (EdbH2 &hd)
 
void PositionPlot (float dz1, float dz2, EdbH2 &hd)
 
void PositionPlotA (EdbH2 &hd, float DR=0, TObjArray *arr1=0, TObjArray *arr2=0)
 
void PrintStat ()
 
void PrintSummary ()
 
void RankCouples ()
 
void RankCouples (TObjArray &arr1, TObjArray &arr2, TObjArray &arrtr)
 
void RankCouples0 ()
 
void RankCouplesFast (TObjArray &arrtr)
 
int ReadPosTree (TTree *tree, EdbPattern *p=0, EdbPattern *p1=0, EdbPattern *p2=0, TEventList *lst=0)
 
void ResetAngularCorr ()
 
void ResetPeak ()
 
void ResetPositionCorr ()
 
void SaveAsTree (EdbPattern &pat, const char *name)
 
int SelectNarrowPeakDXDY (float dr, EdbH2 &hxy)
 
int SelectPeak (EdbPattern &p1, EdbPattern &p2, float dx, float dy, float dz1, float dz2)
 
int SelectPeak (EdbPattern &p1, EdbPattern &p2, float dx=10, float dy=10)
 
int SelectZone (float min[2], float max[2], TObjArray &arr1, TObjArray &arr2, float maxDens)
 
int SelectZoneSide (EdbPatCell2 &pc, float min[2], float max[2], TObjArray &arr, int nmax=kMaxInt)
 
void ShrinkageSelection (float dzbase)
 
int ShrinkageSelectionSide (EdbPatCell2 &pc1, EdbPatCell2 &pc2, EdbH2 &hshr, int nstep, float deltaShr)
 
int SpotsFilterOut (int nmax)
 
int WideSearchXY (EdbPattern &pat1, EdbPattern &pat2, EdbH2 &hxy, float dz, float xmin, float xmax, float ymin, float ymax)
 
int WritePosTree ()
 
float Xcorr (EdbSegP &s, float dz)
 
float Xcorr (EdbSegP &s, float dz, float dx)
 
float Ycorr (EdbSegP &s, float dz)
 
float Ycorr (EdbSegP &s, float dz, float dy)
 
void Zselection (EdbH2 &hd)
 
virtual ~EdbPositionAlignment ()
 

Static Public Member Functions

static int SelectBestComptons (TObjArray &a, TObjArray &arr, int nmax)
 

Public Attributes

EdbAffine2DeAff
 the found affine transformation (when applied to pattern1 gives pattern2 ) More...
 
Float_t eBinX
 
Float_t eBinY
 bin size for the differential hist (for example 5 microns) More...
 
Float_t eChi2Max
 cut for the chi2 of the basetrack More...
 
TObjArray eComb1
 
TObjArray eComb2
 array with the selected combinations of segments More...
 
Bool_t eDoRot
 if true - perform the rotation selection More...
 
Float_t eDTXmax
 
Float_t eDTXmin
 
Float_t eDTYmax
 max angular acceptance (ex: 0.15) for the coinsidences More...
 
Float_t eDTYmin
 min angular difference for the dublets cutout More...
 
Float_t eDXmin
 
Float_t eDYmin
 min position difference for the dublets cutout More...
 
EdbH2 eHDZ
 histogram used for Z-selection More...
 
EdbH2 eHpeak
 histogram used for peak selection More...
 
EdbH2 eHShr1
 histogram used for the shrinkage-selection More...
 
EdbH2 eHShr2
 histogram used for the shrinkage-selection More...
 
Int_t eNpeak
 number of found combinations More...
 
Int_t eNShr1step
 
Int_t eNShr2step
 number of steps for the Shr-selection More...
 
Int_t eNZ1step
 
Int_t eNZ2step
 number of steps for the z-selection More...
 
EdbPatCell2 ePC
 cells with the pointers to tracks More...
 
EdbPatCell2 ePC1
 
EdbPatCell2 ePC2
 cells with the pointers to segments More...
 
TTree * ePosTree
 optional: tree with s1:s2 for the selected combinations More...
 
EdbH1 eRot
 definition of the rotation steps More...
 
Float_t eRpeak
 coordinate peak acceptance More...
 
Float_t eS0tx
 
Float_t eS0ty
 sigmas at 0 angle for the chi2 calculation More...
 
Float_t eS0x
 
Float_t eS0y
 
TObjArray eSegCouples
 segment couples objects to fill couples format tree More...
 
Float_t eShr1from
 
Float_t eShr1to
 limits in Shr for the peak search (ePat1) More...
 
Float_t eShr2from
 
Float_t eShr2to
 limits in Shr for the peak search (ePat2) More...
 
TString eSmoothKernel
 used to smooth histograms More...
 
TObjArray eTracks
 tracks created with segments of eComb1, eComb2 More...
 
Float_t eWbaseMin
 cut for the w of basetrack to accept it More...
 
float eX0
 
Float_t eXcell
 
Float_t eXpeak
 
float eY0
 coordinates of the center of the zone (for the ePeakNT only) More...
 
Float_t eYcell
 cell size (for example 50 microns) More...
 
Float_t eYpeak
 peak position in X,Y More...
 
Float_t eZ1from
 
Float_t eZ1peak
 
Float_t eZ1to
 limits in Z for the peak search (ePat1) More...
 
Float_t eZ2from
 
Float_t eZ2peak
 peak position in Z More...
 
Float_t eZ2to
 limits in Z for the peak search (ePat2) More...
 

Detailed Description

new alignment class developed mainly for compton search

Constructor & Destructor Documentation

◆ EdbPositionAlignment()

EdbPositionAlignment::EdbPositionAlignment ( )
32{
33 eXcell=eYcell=50; // bin size (for example 50 microns)
34 eDTXmax=eDTYmax=0.15; // max angular acceptance (ex: 0.15) for the coincidence
35 eDXmin=eDYmin=5; // min position difference for the doublets cutout
36 eDTXmin=eDTYmin=0.005; // min angular difference for the doublets cutout
37
38 eX0=eY0=0;
39
43 eNpeak = 0;
44 eAff=0;
45
46 eSmoothKernel = "0";
47 ePosTree=0;
48
49 // linking parameters
50 eS0x = eS0y = 2.;
51 eS0tx = eS0ty = 0.01; // sigmas at 0 angle for the chi2 calculation
52 eWbaseMin = 12; // cut for the w of basetrack to accept it
53 eChi2Max = 3.; // cut for the chi2 of the basetrack
54
55}
TTree * ePosTree
optional: tree with s1:s2 for the selected combinations
Definition: EdbPositionAlignment.h:49
Float_t eDTXmin
Definition: EdbPositionAlignment.h:31
Float_t eZ1to
limits in Z for the peak search (ePat1)
Definition: EdbPositionAlignment.h:35
Float_t eDTYmin
min angular difference for the dublets cutout
Definition: EdbPositionAlignment.h:31
Float_t eZ2to
limits in Z for the peak search (ePat2)
Definition: EdbPositionAlignment.h:36
Float_t eDTYmax
max angular acceptance (ex: 0.15) for the coinsidences
Definition: EdbPositionAlignment.h:27
Int_t eNpeak
number of found combinations
Definition: EdbPositionAlignment.h:42
Float_t eZ1from
Definition: EdbPositionAlignment.h:35
Float_t eWbaseMin
cut for the w of basetrack to accept it
Definition: EdbPositionAlignment.h:53
Float_t eS0tx
Definition: EdbPositionAlignment.h:52
Float_t eDXmin
Definition: EdbPositionAlignment.h:30
Float_t eZ2from
Definition: EdbPositionAlignment.h:36
float eY0
coordinates of the center of the zone (for the ePeakNT only)
Definition: EdbPositionAlignment.h:24
TString eSmoothKernel
used to smooth histograms
Definition: EdbPositionAlignment.h:47
Float_t eDTXmax
Definition: EdbPositionAlignment.h:27
Float_t eZ1peak
Definition: EdbPositionAlignment.h:40
float eX0
Definition: EdbPositionAlignment.h:24
Float_t eXcell
Definition: EdbPositionAlignment.h:26
Float_t eChi2Max
cut for the chi2 of the basetrack
Definition: EdbPositionAlignment.h:54
Float_t eS0x
Definition: EdbPositionAlignment.h:52
Float_t eYcell
cell size (for example 50 microns)
Definition: EdbPositionAlignment.h:26
Float_t eS0y
Definition: EdbPositionAlignment.h:52
Float_t eS0ty
sigmas at 0 angle for the chi2 calculation
Definition: EdbPositionAlignment.h:52
Int_t eNZ2step
number of steps for the z-selection
Definition: EdbPositionAlignment.h:37
EdbAffine2D * eAff
the found affine transformation (when applied to pattern1 gives pattern2 )
Definition: EdbPositionAlignment.h:43
Float_t eDYmin
min position difference for the dublets cutout
Definition: EdbPositionAlignment.h:30
Int_t eNZ1step
Definition: EdbPositionAlignment.h:37
Float_t eZ2peak
peak position in Z
Definition: EdbPositionAlignment.h:40

◆ ~EdbPositionAlignment()

EdbPositionAlignment::~EdbPositionAlignment ( )
virtual
59{
60 eComb1.Clear();
61 eComb2.Clear();
62 SafeDelete(eAff);
63 eTracks.Delete();
64}
TObjArray eComb1
Definition: EdbPositionAlignment.h:18
TObjArray eComb2
array with the selected combinations of segments
Definition: EdbPositionAlignment.h:18
TObjArray eTracks
tracks created with segments of eComb1, eComb2
Definition: EdbPositionAlignment.h:19

Member Function Documentation

◆ ActivatePosTree()

bool EdbPositionAlignment::ActivatePosTree ( const char *  name = "postree")
106{
107 if(ePosTree) SafeDelete(ePosTree);
108 ePosTree = new TTree(name,"postree (s1:s2)");
109 EdbSegP *s1 =0, *s2 =0, *s=0;
110 EdbSegCouple *cp=0;
111 Float_t dz1,dz2;
112 Int_t flag;
113 Float_t chi2;
114 int pid1=0,pid2=0;
115 float xv=0,yv=0;
116
117 ePosTree->Branch("pid1", &pid1,"pid1/I");
118 ePosTree->Branch("pid2", &pid2,"pid2/I");
119 ePosTree->Branch("xv", &xv,"xv/F");
120 ePosTree->Branch("yv", &yv,"yv/F");
121 ePosTree->Branch("cp", "EdbSegCouple",&cp,32000,99);
122 ePosTree->Branch("s1.", "EdbSegP",&s1,32000,99);
123 ePosTree->Branch("s2.", "EdbSegP",&s2,32000,99);
124 ePosTree->Branch("s." , "EdbSegP",&s ,32000,99);
125
126 ePosTree->Branch("dz1", &dz1,"dz1/F");
127 ePosTree->Branch("dz2", &dz2,"dz2/F");
128 ePosTree->Branch("flag", &flag,"flag/I");
129 ePosTree->Branch("chi2", &chi2,"chi2/F");
130 return true;
131}
Definition: EdbSegCouple.h:17
Definition: EdbSegP.h:21
s
Definition: check_shower.C:55
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30
EdbSegCouple * cp
Definition: tlg2couples.C:28
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24
Float_t chi2
Definition: testBGReduction_By_ANN.C:14

◆ Align()

int EdbPositionAlignment::Align ( )

find the best alignment of the given zone

696{
698
699 int n[2] = { (int)(2*(eXcell/eBinX)+1), (int)(2*(eYcell/eBinY)+1) };
700 float min[2] = {-eXcell,-eYcell};
701 float max[2] = { eXcell, eYcell};
702 EdbH2 hd; // define the differential hist for the peaks search
703 hd.InitH2(n, min, max);
704 Zselection(hd); // roughly define eZ1peak, eZ2peak, eXpeak,eYpeak
705
706 if(ePosTree) FillPosTree( eZ1peak, eZ2peak, 0); // flag=0
707
708 EdbPattern p1, p2;
709 int npeak = 0;
710 npeak = SelectPeak(p1,p2, eRpeak, eRpeak );
711 p1.Clear(); p2.Clear();
712 npeak = SelectPeak(p1,p2, eRpeak, eRpeak );
713 //if(ePosTree) FillPosTree( p1, p2, 0, 0, 1); // flag=1
714
715 return npeak;
716}
float min(TClonesArray *t)
Definition: bitview.cxx:275
fast 2-dim histogram class (used as a basis for EdbCell2)
Definition: EdbCell2.h:19
int InitH2(const EdbH2 &h)
Definition: EdbCell2.cpp:78
Definition: EdbPattern.h:273
Float_t eRpeak
coordinate peak acceptance
Definition: EdbPositionAlignment.h:32
int SelectPeak(EdbPattern &p1, EdbPattern &p2, float dx, float dy, float dz1, float dz2)
Definition: EdbPositionAlignment.cxx:719
Float_t eBinY
bin size for the differential hist (for example 5 microns)
Definition: EdbPositionAlignment.h:28
Float_t eBinX
Definition: EdbPositionAlignment.h:28
int FillPosTree(float dz1, float dz2, int flag)
Definition: EdbPositionAlignment.cxx:189
void Zselection(EdbH2 &hd)
Definition: EdbPositionAlignment.cxx:624
int max
Definition: check_shower.C:41

◆ AlignWithTracks()

void EdbPositionAlignment::AlignWithTracks ( int  nang = 3,
int  npos = 0 
)

assuming that the combinations are found and tracks are reconstructed
correct DTX,DTY of the sides in respect of the tracks

920{
923
924 for(int i=0; i<nang; i++) {
925 float dx, dy, dtx, dty;
926 FindDiff10(dx, dy, dtx, dty);
927 Log(3,"AlignWithTracks", "%d angle side 1 diff: %f %f %f %f with %d ", i,dx, dy, dtx, dty, eComb1.GetEntries());
928 ePC1.eDTX += dtx;
929 ePC1.eDTY += dty;
930 FindDiff20(dx, dy, dtx, dty);
931 Log(3,"AlignWithTracks", "%d angle side 2 diff: %f %f %f %f with %d ", i,dx, dy, dtx, dty, eComb2.GetEntries());
932 ePC2.eDTX += dtx;
933 ePC2.eDTY += dty;
934 }
935
936 for(int i=0; i<npos; i++) {
937 float dx, dy, dtx, dty;
938 FindDiff12(dx, dy, dtx, dty);
939 Log(3,"AlignWithTracks", "%d pos side 1 diff: %f %f %f %f with %d ", i, dx, dy, dtx, dty, eComb1.GetEntries());
940 ePC1.eDX += dx/2.;
941 ePC1.eDY += dy/2.;
942 ePC2.eDX -= dx/2.;
943 ePC2.eDY -= dy/2.;
944 }
945
946}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
float eDX
Definition: EdbPatCell2.h:18
float eDTX
Definition: EdbPatCell2.h:20
float eDY
Definition: EdbPatCell2.h:18
float eDTY
corrections to be applied if eApplyCorr==true
Definition: EdbPatCell2.h:20
void FindDiff10(float &dx, float &dy, float &dtx, float &dty)
Definition: EdbPositionAlignment.cxx:955
void FindDiff12(float &dx, float &dy, float &dtx, float &dty)
Definition: EdbPositionAlignment.cxx:949
void FindDiff20(float &dx, float &dy, float &dtx, float &dty)
Definition: EdbPositionAlignment.cxx:961
EdbPatCell2 ePC2
cells with the pointers to segments
Definition: EdbPositionAlignment.h:17
EdbPatCell2 ePC1
Definition: EdbPositionAlignment.h:17

◆ DoubletsFilterOut()

int EdbPositionAlignment::DoubletsFilterOut ( bool  checkview = true)
766{
767 return DoubletsFilterOutSide(ePC1,checkview) + DoubletsFilterOutSide(ePC2,checkview);
768}
int DoubletsFilterOutSide(EdbPatCell2 &pc, bool checkview=true)
Definition: EdbPositionAlignment.cxx:771

◆ DoubletsFilterOutSide()

int EdbPositionAlignment::DoubletsFilterOutSide ( EdbPatCell2 pc,
bool  checkview = true 
)

assign flag = -10 for the duplicated segments
if(checkview) - do not remove close segments if them are in the same view

772{
775
776 int nout =0;
777 int n = FillCombinations(pc, pc, eDXmin, eDYmin, eDTXmin, eDTYmin);
778 EdbSegP *s1,*s2;
779 for(int i=0; i<n; i++) {
780 s1 = (EdbSegP*)eComb1.UncheckedAt(i);
781 s2 = (EdbSegP*)eComb2.UncheckedAt(i);
782 if(checkview) if( s2->Aid(0)==s1->Aid(0) && s2->Aid(1)==s1->Aid(1) && s2->Side()==s1->Side() ) continue;
783 if( s2->W()>s1->W() ) s1->SetFlag(-10);
784 else s2->SetFlag(-10);
785 nout++;
786 }
787
788 if(ePosTree) FillPosTree( 0., 0., -10);
789
790 Log(2,"DubletsFilterOut","%d segments discarded with DX,DY,DTX,DTY: (%5.1f %5.1f %5.3f %5.3f) checkview =%d",
791 nout,eDXmin, eDYmin, eDTXmin, eDTYmin, checkview );
792 return nout;
793}
int FillCombinations()
Definition: EdbPositionAlignment.cxx:534
Int_t Side() const
mandatory virtual functions:
Definition: EdbSegP.h:170
Float_t W() const
Definition: EdbSegP.h:151
Int_t Aid(int i) const
Definition: EdbSegP.h:169
void SetFlag(int flag)
Definition: EdbSegP.h:130

◆ FillArrays() [1/2]

int EdbPositionAlignment::FillArrays ( EdbPattern p1,
EdbPattern p2 
)
417{
418 ePC1.Reset();
419 ePC2.Reset();
420
421 float marg = eXcell/2.;
422 float min1[2] = { Min(p1.Xmin(),p2.Xmin())-marg, Min(p1.Ymin(),p2.Ymin())-marg };
423 float max1[2] = { Max(p1.Xmax(),p2.Xmax())+marg, Max(p1.Ymax(),p2.Ymax())+marg };
424 int n1[2] = { (int)((max1[0]-min1[0])/eXcell+1), (int)((max1[1]-min1[1])/eYcell+1) };
425 //eXcell = Max((max1[0]-min1[0])/n1[0],mincell);
426 //eYcell = Max((max1[1]-min1[1])/n1[1],mincell);
427 float min2[2] = { min1[0]-eXcell, min1[1]-eYcell };
428 float max2[2] = { max1[0]+eXcell, max1[1]+eYcell };
429 int n2[2] = { n1[0]+2 , n1[1]+2 };
430
431 if( n1[0]<1 || n1[1]<1 ) return 0;
432 int meancell1 = (int)(p1.N()/n1[0]/n1[1]);
433 int meancell2 = (int)(p2.N()/n2[0]/n2[1]);
434 int maxcell1 = (int)(10 + meancell1 + 5*Sqrt(meancell1));
435 int maxcell2 = (int)(10 + meancell2 + 5*Sqrt(meancell2));
436
437 ePC1.InitCell(maxcell1, n1, min1, max1);
438 ePC2.InitCell(maxcell2, n2, min2, max2);
439
440 int nc1 = ePC1.FillCell(p1);
441 int nc2 = ePC2.FillCell(p2);
442
443 Log(3,"FillArrays","filled objects/cells: %d / %d and %d / %d in cells of (%6.1f x %6.1f)",
444 p1.N(), nc1, p2.N(), nc2, eXcell, eYcell);
445 return nc1+nc2;
446}
void Reset()
Definition: EdbCell2.h:171
int InitCell(EdbCell2 &c)
Definition: EdbCell2.h:167
int FillCell(EdbPattern &pat)
limits should be already defined
Definition: EdbPatCell2.cxx:96
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

◆ FillArrays() [2/2]

int EdbPositionAlignment::FillArrays ( TObjArray &  arr1,
TObjArray &  arr2,
float  xymin[2],
float  xymax[2] 
)
450{
451 ePC1.Reset();
452 ePC2.Reset();
453
454 float marg = 0.0001;
455 float min1[2] = { xymin[0]-marg, xymin[1]-marg };
456 float max1[2] = { xymax[0]+marg, xymax[1]+marg };
457 int n1[2] = { (int)((max1[0]-min1[0])/eXcell+1), (int)((max1[1]-min1[1])/eYcell+1) };
458 eXcell = (max1[0]-min1[0])/n1[0];
459 eYcell = (max1[1]-min1[1])/n1[1];
460 float min2[2] = { min1[0]-eXcell, min1[1]-eYcell };
461 float max2[2] = { max1[0]+eXcell, max1[1]+eYcell };
462 int n2[2] = { n1[0]+2 , n1[1]+2 };
463
464 if( n1[0]<1 || n1[1]<1 ) return 0;
465 int entr1 = arr1.GetEntriesFast();
466 int entr2 = arr2.GetEntriesFast();
467 int meancell1 = (int)(entr1/n1[0]/n1[1]);
468 int meancell2 = (int)(entr2/n2[0]/n2[1]);
469 int maxcell = (int)Max(meancell1,meancell2)+10;
470 maxcell += (int)(5*Sqrt(maxcell)); // define safe value for max/cell limit
471
472 ePC1.InitCell(maxcell, n1, min1, max1);
473 ePC2.InitCell(maxcell, n2, min2, max2);
474
475 int nc1 = ePC1.FillCell(arr1);
476 int nc2 = ePC2.FillCell(arr2);
477
478 Log(3,"FillArrays","%d and %d cells filled", nc1,nc2);
479 return nc1;
480}

◆ FillCombinations() [1/2]

int EdbPositionAlignment::FillCombinations ( )
535{
537}

◆ FillCombinations() [2/2]

int EdbPositionAlignment::FillCombinations ( EdbPatCell2 pc1,
EdbPatCell2 pc2,
float  dx,
float  dy,
float  dtx,
float  dty 
)

the cells must be already filled

541{
543 pc1.eDXlim = dx; pc1.eDYlim = dy;
544 pc1.eDTXlim = dtx; pc1.eDTYlim = dty;
545 eComb1.Clear();
546 eComb2.Clear();
547 int ir[2] = { int((dx-0.0001)/eXcell)+1, int((dy-0.0001)/eYcell)+1 };
548 int ncomb = pc1.FillCombinations(pc2, ir, eComb1, eComb2);
549 Log(3,"FillCombinations","%d selected with the acceptance: %7.2f %7.2f ( %d x %4.1f, %d x %4.1f ) and %7.4f %7.4f ",
550 ncomb, dx,dy, ir[0], eXcell, ir[1], eYcell, dtx, dty );
551 return ncomb;
552}
int FillCombinations(EdbPatCell2 &cell, int ir2[2], TObjArray &arrC1, TObjArray &arrC2, bool doFill=1)
Definition: EdbPatCell2.cxx:122
float eDYlim
acceptance limits for the combinations selection
Definition: EdbPatCell2.h:12
float eDXlim
Definition: EdbPatCell2.h:12
float eDTYlim
Definition: EdbPatCell2.h:13
float eDTXlim
Definition: EdbPatCell2.h:13

◆ FillPosTree()

int EdbPositionAlignment::FillPosTree ( float  dz1,
float  dz2,
int  flag 
)
190{
191 if(!ePosTree) return 0;
192 EdbSegP *sb1 = new EdbSegP();
193 EdbSegP *sb2 = new EdbSegP();
194 EdbSegP *sb = new EdbSegP();
196
197 ePosTree->SetBranchAddress("cp" , &cp );
198 ePosTree->SetBranchAddress("s1." , &sb1);
199 ePosTree->SetBranchAddress("s2." , &sb2);
200 ePosTree->SetBranchAddress("s." , &sb);
201 ePosTree->SetBranchAddress("dz1" , &dz1);
202 ePosTree->SetBranchAddress("dz2" , &dz2);
203 ePosTree->SetBranchAddress("flag", &flag);
204 Float_t chi2;
205 ePosTree->SetBranchAddress("chi2", &chi2);
206
207 EdbSegCouple *cpt=0;
208 EdbSegP *s1, *s2;
209 EdbTrackP *t;
210 int nsaved = 0;
211
212 int n = eComb1.GetEntries();
213 int nt = eTracks.GetEntries();
214 int ncp = eSegCouples.GetEntries();
215
216 for(int i=0; i<n; i++) {
217 s1 = (EdbSegP*)eComb1.UncheckedAt(i);
218 s2 = (EdbSegP*)eComb2.UncheckedAt(i);
219 if(nt==n) {
220 t = (EdbTrackP*)eTracks.UncheckedAt(i);
221 if(t->Chi2()>eChi2Max) continue;
222 sb->Copy( *((EdbSegP*)t) );
223 }
224 if(ncp==n) {
225 cpt = (EdbSegCouple*)eSegCouples.UncheckedAt(i);
226 *cp = *cpt;
227 }
228 sb1->Copy(*s1);
229 sb2->Copy(*s2);
230 ePosTree->Fill();
231 nsaved++;
232 }
233 Log(2,"FillPosTree","%s %d couples saved with chi2 <= %f", ePosTree->GetName(), nsaved, eChi2Max);
234 return n;
235}
TObjArray eSegCouples
segment couples objects to fill couples format tree
Definition: EdbPositionAlignment.h:20
void Copy(const EdbSegP &s)
Definition: EdbSegP.cxx:105
Definition: EdbPattern.h:113
TTree * t
Definition: check_shower.C:4

◆ FindCorrections()

void EdbPositionAlignment::FindCorrections ( EdbPattern pat1,
EdbPattern pat2,
float  DZ,
bool  doShrink 
)

find corrections assuming the patterns already closely aligned (linking case)
DZ is the distance between patterns(fixed)

1009{
1012
1013 FillArrays(pat1, pat2); // fill at nominal quote - important for doublets and for shrinkage check
1014 //DoubletsFilterOut(true);
1015 if(doShrink) ShrinkageSelection(DZ);
1016
1017 ePC1.eDZ = DZ/2;
1018 ePC2.eDZ = -DZ/2;
1019 ePC1.eApplyCorr = ePC2.eApplyCorr = true;
1020 FillArrays(pat1, pat2); // fill at the same quote - to increase the coinsidences
1022 //Link();
1023
1024 int nxy[2] = {51,51}; // TODO
1025 float min[2] = {-2*eXcell, -2*eXcell };
1026 float max[2] = { 2*eXcell, 2*eXcell };
1027 EdbH2 hxy;
1028 hxy.InitH2(nxy, min, max);
1029
1030 SelectNarrowPeakDXDY(10., hxy); // only peak couples remaining in eComb1/2
1031 RankCouples0(); // construct tracks (without correction)
1032
1033 if(ActivatePosTree("postreeDXDY")) { FillPosTree(0, 0, 9); WritePosTree(); }
1034
1037
1038
1039 AlignWithTracks(3,0); // only angular correction of segments to keep original basetrack angle
1040
1041 if(ActivatePosTree("postreeDXDYalign")) { FillPosTree(0, 0, 8); WritePosTree(); }
1042
1043 Log(3,"FindCorrections","side1 (shr,dx,dy,dtx,dty): %f %f %f %f %f", ePC1.eShr, ePC1.eDX, ePC1.eDY, ePC1.eDTX, ePC1.eDTY);
1044 Log(3,"FindCorrections","side2 (shr,dx,dy,dtx,dty): %f %f %f %f %f", ePC2.eShr, ePC2.eDX, ePC2.eDY, ePC2.eDTX, ePC2.eDTY);
1045}
float eShr
corrections to be applied if eApplyCorr==true
Definition: EdbPatCell2.h:19
float eDZ
corrections to be applied if eApplyCorr==true
Definition: EdbPatCell2.h:18
bool eApplyCorr
Definition: EdbPatCell2.h:17
void ResetPositionCorr()
Definition: EdbPositionAlignment.h:130
void AlignWithTracks(int nang=3, int npos=0)
Definition: EdbPositionAlignment.cxx:919
void ShrinkageSelection(float dzbase)
Definition: EdbPositionAlignment.cxx:806
bool ActivatePosTree(const char *name="postree")
Definition: EdbPositionAlignment.cxx:105
int FillArrays(EdbPattern &p1, EdbPattern &p2)
Definition: EdbPositionAlignment.cxx:416
void ResetAngularCorr()
Definition: EdbPositionAlignment.h:131
void RankCouples0()
Definition: EdbPositionAlignment.cxx:294
int WritePosTree()
Definition: EdbPositionAlignment.cxx:134
int SelectNarrowPeakDXDY(float dr, EdbH2 &hxy)
Definition: EdbPositionAlignment.cxx:878
float DZ
Definition: hwinit.C:66

◆ FindDiff()

void EdbPositionAlignment::FindDiff ( TObjArray &  arr1,
TObjArray &  arr2,
EdbPatCell2 pc1,
EdbPatCell2 pc2,
float &  dx,
float &  dy,
float &  dtx,
float &  dty 
)
969{
970 Log(3,"FindDiff","pc1: %f %f %f %f %f %f %f",pc1.eDX, pc1.eDY, pc1.eDZ, pc1.eDTX, pc1.eDTY, pc1.ePhi, pc1.eShr);
971 Log(3,"FindDiff","pc2: %f %f %f %f %f %f %f",pc2.eDX, pc2.eDY, pc2.eDZ, pc2.eDTX, pc2.eDTY, pc2.ePhi, pc2.eShr);
972
973 int n = arr1.GetEntries();
974 EdbSegP *s1,*s2;
975 Double_t sdx=0, sdy=0, sdtx=0, sdty=0;
976 for(int i=0; i<n; i++) {
977 s1 = (EdbSegP*)arr1.UncheckedAt(i);
978 s2 = (EdbSegP*)arr2.UncheckedAt(i);
979 sdx += pc2.Xs(*s2)-pc1.Xs(*s1);
980 sdy += pc2.Ys(*s2)-pc1.Ys(*s1);
981 sdtx += pc2.TXs(*s2)-pc1.TXs(*s1);
982 sdty += pc2.TYs(*s2)-pc1.TYs(*s1);
983 }
984 if(Abs(sdtx)>10) Log(1,"FindDiff","%f %f %f %f",sdx,sdy,sdtx,sdty);
985 dx = sdx/n;
986 dy = sdy/n;
987 dtx = sdtx/n;
988 dty = sdty/n;
989}
float TYs(EdbSegP &s)
Definition: EdbPatCell2.h:50
float Xs(EdbSegP &s)
Definition: EdbPatCell2.h:40
float TXs(EdbSegP &s)
Definition: EdbPatCell2.h:49
float ePhi
rotation angle (using dx,dy one can set-up the center of rotation)
Definition: EdbPatCell2.h:21
float Ys(EdbSegP &s)
Definition: EdbPatCell2.h:41

◆ FindDiff10()

void EdbPositionAlignment::FindDiff10 ( float &  dx,
float &  dy,
float &  dtx,
float &  dty 
)
956{
957 FindDiff(eComb1, eTracks, ePC1, ePC, dx,dy,dtx,dty);
958}
EdbPatCell2 ePC
cells with the pointers to tracks
Definition: EdbPositionAlignment.h:21
void FindDiff(TObjArray &arr1, TObjArray &arr2, EdbPatCell2 &pc1, EdbPatCell2 &pc2, float &dx, float &dy, float &dtx, float &dty)
Definition: EdbPositionAlignment.cxx:967

◆ FindDiff12()

void EdbPositionAlignment::FindDiff12 ( float &  dx,
float &  dy,
float &  dtx,
float &  dty 
)
950{
951 FindDiff(eComb1, eComb2, ePC1, ePC2, dx,dy,dtx,dty);
952}

◆ FindDiff20()

void EdbPositionAlignment::FindDiff20 ( float &  dx,
float &  dy,
float &  dtx,
float &  dty 
)
962{
963 FindDiff(eComb2, eTracks, ePC2, ePC, dx,dy,dtx,dty);
964}

◆ Link()

int EdbPositionAlignment::Link ( )
993{
994 //FillArrays(p1, p2);
995 eComb1.Clear();
996 eComb2.Clear();
998 ePC1.eDTXlim = 0.12; ePC1.eDTYlim = 0.12;
999 int ir[2] = {2,2};
1000 int ncomb = ePC1.FillCombinations(ePC2, ir, eComb1, eComb2);
1001 //RankCouples();
1002
1003 Log(1,"WARNING: deprecated function: EdbPositionAlignment::Link"," %8d combinations found", ncomb );
1004 return ncomb;
1005}

◆ PositionPlot() [1/2]

void EdbPositionAlignment::PositionPlot ( EdbH2 hd)

output: hd - differential hist, should be already defined
use correction defined in ePC1, ePC2 if requested

600{
603
604 int n1 = eComb1.GetEntries();
605 EdbSegP *s1,*s2;
606 for(int i=0; i<n1; i++) {
607 s1 = (EdbSegP*)eComb1.UncheckedAt(i);
608 s2 = (EdbSegP*)eComb2.UncheckedAt(i);
609 hd.Fill( ePC2.Xs(*s2) - ePC1.Xs(*s1),
610 ePC2.Ys(*s2) - ePC1.Ys(*s1) );
611 }
612}
int Fill(float x, float y)
Definition: EdbCell2.h:88

◆ PositionPlot() [2/2]

void EdbPositionAlignment::PositionPlot ( float  dz1,
float  dz2,
EdbH2 hd 
)

dz1,dz2 offsets in Z,
output: hd - differential hist, should be already defined

585{
588
589 int n1 = eComb1.GetEntries();
590 EdbSegP *s1,*s2;
591 for(int i=0; i<n1; i++) {
592 s1 = (EdbSegP*)eComb1.UncheckedAt(i);
593 s2 = (EdbSegP*)eComb2.UncheckedAt(i);
594 hd.Fill( Xcorr(*s2,dz2)-Xcorr(*s1,dz1), Ycorr(*s2,dz2)-Ycorr(*s1,dz1) );
595 }
596}
float Ycorr(EdbSegP &s, float dz)
Definition: EdbPositionAlignment.h:115
float Xcorr(EdbSegP &s, float dz)
Definition: EdbPositionAlignment.h:114

◆ PositionPlotA()

void EdbPositionAlignment::PositionPlotA ( EdbH2 hd,
float  DR = 0,
TObjArray *  arr1 = 0,
TObjArray *  arr2 = 0 
)

dz1,dz2 offsets in Z should be defined in ePC1, ePC2
output: hd - differential hist, should be already defined
fill arr1 and arr2 if requested, using dr

556{
561
562 hd.CleanCells();
563 int n1 = eComb1.GetEntries();
564 EdbSegP *s1,*s2;
565 bool fillarr=false;
566 if(arr1&&arr2&&DR>0) fillarr=true;
567 float dx, dy, dr;
568 for(int i=0; i<n1; i++) {
569 s1 = (EdbSegP*)eComb1.UncheckedAt(i);
570 s2 = (EdbSegP*)eComb2.UncheckedAt(i);
571 dx = ePC2.Xs(*s2)-ePC1.Xs(*s1);
572 dy = ePC2.Ys(*s2)-ePC1.Ys(*s1);
573 hd.Fill( dx, dy );
574 if(fillarr) {
575 dr = Sqrt(dx*dx+dy*dy);
576 if(dr>DR) continue;
577 arr1->Add(s1);
578 arr2->Add(s2);
579 }
580 }
581}
void CleanCells()
Definition: EdbCell2.cpp:114

◆ PrintStat()

void EdbPositionAlignment::PrintStat ( )
68{
71 printf("Z variate as: pat1: %d (%6.1f %6.1f) pat2: %d (%6.1f %6.1f) \n",
73}
void PrintStat()
Definition: EdbCell2.cpp:674

◆ PrintSummary()

void EdbPositionAlignment::PrintSummary ( )
77{
78 printf("Alignment summary: ncomb = %7d\n",eComb1.GetEntries());
79 printf("%5d coinsidence in the peak with dXY: (%6.2f %6.2f) dz12: (%6.1f %6.1f) at XY: (%10.1f %10.1f)\n",
81}
Float_t eXpeak
Definition: EdbPositionAlignment.h:41
Float_t eYpeak
peak position in X,Y
Definition: EdbPositionAlignment.h:41

◆ RankCouples() [1/2]

void EdbPositionAlignment::RankCouples ( )
288{
289 eTracks.Clear();
291}
void RankCouples()
Definition: EdbPositionAlignment.cxx:287

◆ RankCouples() [2/2]

void EdbPositionAlignment::RankCouples ( TObjArray &  arr1,
TObjArray &  arr2,
TObjArray &  arrtr 
)
324{
325 int n = arr1.GetEntries();
326 Log(3,"RankCouples","%d" ,n);
327
328 EdbSegP *s1, *s2;
329 EdbTrackP *t;
331 EdbScanCond cond1,cond2;
332 cond1.SetSigma0(1.,1., 0.010,0.010); cond1.SetDegrad(4); cond1.SetPulsRamp0(6,9); cond1.SetPulsRamp04(6,9);
333 cond2.SetSigma0(1.,1., 0.010,0.010); cond2.SetDegrad(4); cond2.SetPulsRamp0(6,9); cond2.SetPulsRamp04(6,9);
334 EdbSegP seg;
335 for(int i=0; i<n; i++) {
336 s1 = (EdbSegP*)arr1.At(i);
337 s2 = (EdbSegP*)arr2.At(i);
338 seg.Copy(*s1); // to set correctly vid, aid, etc
339 //tf.Chi2SegM( *s1, *s2, seg, cond1, cond2);
340 tf.Chi2ASeg( *s1, *s2, seg, cond1, cond2);
341 //seg.SetChi2( tf.Chi2ACP( *s1, *s2, cond1) ); //TODO test!!
342
343 t = new EdbTrackP(seg);
344 t->AddSegment(s1);
345 t->AddSegment(s2);
346 // s2->SetMC(s.MCEvt(),s.MCTrack());
347 arrtr.Add( t );
348 }
349
351 Log(2,"RankCouples","%d tracks ok", arrtr.GetEntriesFast() );
352}
TObjArray * arrtr
Definition: RecDispMC.C:129
void RankCouplesFast(TObjArray &arrtr)
Definition: EdbPositionAlignment.cxx:355
Definition: EdbScanCond.h:10
void SetPulsRamp0(float p1, float p2)
Definition: EdbScanCond.h:74
void SetDegrad(float d)
Definition: EdbScanCond.h:71
void SetSigma0(float x, float y, float tx, float ty)
Definition: EdbScanCond.h:62
void SetPulsRamp04(float p1, float p2)
Definition: EdbScanCond.h:75
Definition: EdbTrackFitter.h:17
static float Chi2ASeg(EdbSegP &s1, EdbSegP &s2, EdbSegP &s, EdbScanCond &cond1, EdbScanCond &cond2)
Definition: EdbTrackFitter.cxx:217

◆ RankCouples0()

void EdbPositionAlignment::RankCouples0 ( )

construct basetracks only using the position information

295{
297
298 eTracks.Clear();
299 int n = eComb1.GetEntries();
300 Log(3,"RankCouples0","%d" ,n);
301 EdbSegP *s1, *s2, *s;
302 float dz;
303 for(int i=0; i<n; i++) {
304 s1 = (EdbSegP*)eComb1.At(i);
305 s2 = (EdbSegP*)eComb2.At(i);
306 dz = s2->Z()-s1->Z();
307 s = new EdbSegP( i,
308 0.5*(s1->X() + s2->X()),
309 0.5*(s1->Y() + s2->Y()),
310 (s2->X() - s1->X())/dz,
311 (s2->Y() - s1->Y())/dz,
312 s1->W() + s2->W()
313 );
314 s->SetZ(0.5*(s1->Z() + s2->Z()));
315 eTracks.Add( s );
316 }
317
318 Log(2,"RankCouples0","%d tracks ok", eTracks.GetEntriesFast() );
319
320}
brick dz
Definition: RecDispMC.C:107
Float_t X() const
Definition: EdbSegP.h:173
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174

◆ RankCouplesFast()

void EdbPositionAlignment::RankCouplesFast ( TObjArray &  arrtr)
356{
357 int ntr = arrtr.GetEntries();
358 Log(1,"RankCouplesFast","%d tracks", ntr );
359 TArrayF chi2arr(ntr);
360 TArrayI ind(ntr);
361 eSegCouples.Delete();
362 eSegCouples.Expand(ntr);
363 EdbTrackP *t=0;
364 for(int i=0; i<ntr; i++) {
365 t = ((EdbTrackP*)arrtr.At(i));
366 chi2arr[i]= t->Chi2();
367 t->GetSegment(0)->SetFlag(0); //use flag as a counter
368 t->GetSegment(1)->SetFlag(0);
369 }
370 Sort(ntr,chi2arr.GetArray(),ind.GetArray(),0);
371
372 EdbSegP *s0=0, *s1=0;
373 for(int i=0; i<ntr; i++) {
374 t = ((EdbTrackP*)arrtr.At(ind[i]));
375 s0 = t->GetSegment(0);
376 s1 = t->GetSegment(1);
377 s0->SetFlag( s0->Flag()+1 );
378 s1->SetFlag( s1->Flag()+1 );
379 EdbSegCouple *sc = new EdbSegCouple();
380 sc->SetN1(s0->Flag());
381 sc->SetN2(s1->Flag());
382 eSegCouples.Add(sc);
383 }
384
385 for(int i=0; i<ntr; i++) {
386 t = ((EdbTrackP*)arrtr.At(ind[i]));
387 s0 = t->GetSegment(0);
388 s1 = t->GetSegment(1);
389 EdbSegCouple *sc = (EdbSegCouple *)eSegCouples.UncheckedAt(i);
390 sc->SetCHI2P( t->Chi2() );
391 sc->SetN1tot(s0->Flag());
392 sc->SetN2tot(s1->Flag());
393 s0->SetFlag(0);
394 s1->SetFlag(0);
395 }
396
397}
void SetN1(int n1)
Definition: EdbSegCouple.h:43
void SetN2(int n2)
Definition: EdbSegCouple.h:44
void SetN1tot(int n)
Definition: EdbSegCouple.h:45
void SetN2tot(int n)
Definition: EdbSegCouple.h:46
void SetCHI2P(float chi2)
Definition: EdbSegCouple.h:48
Int_t Flag() const
Definition: EdbSegP.h:149

◆ ReadPosTree()

int EdbPositionAlignment::ReadPosTree ( TTree *  tree,
EdbPattern p = 0,
EdbPattern p1 = 0,
EdbPattern p2 = 0,
TEventList *  lst = 0 
)
239{
240 if(!tree) return 0;
241 EdbSegP *sb1 = 0;//new EdbSegP();
242 EdbSegP *sb2 = 0;//new EdbSegP();
243 EdbSegP *sb = 0;//new EdbSegP();
244 Float_t dz1,dz2;
245 Int_t flag;
246 tree->SetBranchAddress("s1." , &sb1);
247 tree->SetBranchAddress("s2." , &sb2);
248 tree->SetBranchAddress("s." , &sb);
249 tree->SetBranchAddress("dz1" , &dz1);
250 tree->SetBranchAddress("dz2" , &dz2);
251 tree->SetBranchAddress("flag", &flag);
252 Float_t chi2;
253 tree->SetBranchAddress("chi2", &chi2);
254
255 int n = tree->GetEntries();
256 if(lst) n = lst->GetN();
257
258 for(int i=0; i<n; i++) {
259 if(lst) tree->GetEntry(lst->GetEntry(i));
260 else tree->GetEntry(i);
261 if(p) p->AddSegment( *sb );
262 if(p1) p1->AddSegment( *sb1 );
263 if(p2) p2->AddSegment( *sb2 );
264 }
265 Log(2,"ReadPosTree","%d basetracks read from %s", n,tree->GetName());
266 return n;
267}
EdbSegP * AddSegment(int i, EdbSegP &s)
Definition: EdbPattern.cxx:72
p
Definition: testBGReduction_AllMethods.C:8

◆ ResetAngularCorr()

void EdbPositionAlignment::ResetAngularCorr ( )
inline
131{ ePC1.eDTX = ePC1.eDTY = ePC2.eDTX = ePC2.eDTY = 0;}

◆ ResetPeak()

void EdbPositionAlignment::ResetPeak ( )
inline

◆ ResetPositionCorr()

void EdbPositionAlignment::ResetPositionCorr ( )
inline
130{ ePC1.eDX = ePC1.eDY = ePC2.eDX = ePC2.eDY = 0;}

◆ SaveAsTree()

void EdbPositionAlignment::SaveAsTree ( EdbPattern pat,
const char *  name 
)
271{
272 TTree *tree = new TTree(name,"pattern tree");
273 EdbSegP *s=new EdbSegP();
274 tree->Branch("s" , "EdbSegP",&s ,32000,99);
275 int n = pat.N(), nsaved=0;
276 for(int i=0; i<n; i++) {
277 s->Copy( *(pat.GetSegment(i)) );
278 tree->Fill();
279 nsaved++;
280 }
281 tree->Write();
282 SafeDelete(tree);
283 Log(2,"SaveAsTree","%d segments saved", nsaved);
284}
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66

◆ SelectBestComptons()

int EdbPositionAlignment::SelectBestComptons ( TObjArray &  a,
TObjArray &  arr,
int  nmax 
)
static
506{
507 int ic=0;
508 int n = a.GetEntriesFast();
509 EdbSegP *s;
510 if(n>nmax) {
511 TArrayF comptonity(n); // should be defined as: the bigger - the better
512 TArrayI ind(n);
513 for(int i=0; i<n; i++) {
514 s = (EdbSegP*)a.At(i);
515 if(s->Flag()<0) comptonity[i] =0;
516 else comptonity[i] = s->W()*100 - s->Theta(); // proportional to ngrains and low theta
517 }
518 Sort(n,comptonity.GetArray(),ind.GetArray(),1); // sort in decreasing order
519 for(int i=0; i<nmax; i++) {
520 s = (EdbSegP*)a.At(ind[i]);
521 arr.Add(s); ic++;
522 }
523 }
524 else {
525 for(int i=0; i<n; i++) {
526 s = (EdbSegP*)a.UncheckedAt(i);
527 if(s->Flag()>-1) { arr.Add(s); ic++; }
528 }
529 }
530 return ic;
531}
void a()
Definition: check_aligned.C:59

◆ SelectNarrowPeakDXDY()

int EdbPositionAlignment::SelectNarrowPeakDXDY ( float  dr,
EdbH2 hxy 
)

assuming that the combinations are found
make the position plot and find the maximum,
select only the combinations belongs to the maximum cutting out the tails

879{
883
884 float vpeak[2]={0,0};
885 PositionPlotA(hxy);
886 EdbPeak2 p2d(hxy);
887 p2d.FindPeak(vpeak);
888 ePC1.eDX += vpeak[0];
889 ePC1.eDY += vpeak[1];
890
891 for(int i=0; i<6; i++) {
892 TObjArray comb1,comb2;
893 PositionPlotA(hxy, dr, &comb1, &comb2);
894 float dx, dy, dtx, dty;
895 FindDiff(comb1,comb2,ePC1,ePC2, dx, dy, dtx, dty);
896 Log(3,"SelectNarrowPeakDXDY", "diff: %f %f %f %f with %d eDX,eDY = %f %f", dx, dy, dtx, dty, comb1.GetEntries(), ePC1.eDX,ePC1.eDY);
897 ePC1.eDX += dx;
898 ePC1.eDY += dy;
899 }
900
901 TObjArray comb1,comb2;
902 PositionPlotA(hxy, dr, &comb1, &comb2);
903 hxy.DrawH2("dxdy")->Write("dxdy");
904
905 if(ActivatePosTree("postreeDXDYwide")) { FillPosTree(0, 0, 9); WritePosTree(); }
906
907 eComb1.Clear(); eComb2.Clear();
908 int nsel = comb1.GetEntries();
909 for(int i=0; i<nsel; i++) {
910 eComb1.Add(comb1.UncheckedAt(i));
911 eComb2.Add(comb2.UncheckedAt(i));
912 }
913
914 Log(2,"SelectNarrowPeakDXDY", "peak of %d found at %f %f with dr = %f", eComb1.GetEntries(), ePC1.eDX,ePC1.eDY, dr);
915 return comb1.GetEntries();
916}
TH2F * DrawH2(const char *name="plot2d", const char *title="EdbH2plot2D")
Definition: EdbCell2.cpp:187
peak analyser for EdbH2
Definition: EdbCell2.h:105
void PositionPlotA(EdbH2 &hd, float DR=0, TObjArray *arr1=0, TObjArray *arr2=0)
Definition: EdbPositionAlignment.cxx:555

◆ SelectPeak() [1/2]

int EdbPositionAlignment::SelectPeak ( EdbPattern p1,
EdbPattern p2,
float  dx,
float  dy,
float  dz1,
float  dz2 
)

Select the couples formed the peak
assumed that the peak eXpeak, eYpeak is approximately found
Input: dz1,dz2 - for projection
dxMax, dYmax - for selection
Output p1,p2 with found segments

721{
727
728 int n1 = eComb1.GetEntries();
729
730 TH2F *hpeak = new TH2F("hpeak","Selected peak",100,-dxMax+eXpeak,dxMax+eXpeak,100,-dyMax+eYpeak,dyMax+eYpeak);
731 int ic=0;
732 EdbSegP *s1,*s2;
733 EdbSegP s;
734 float dx,dy;
735 for(int i=0; i<n1; i++) {
736 s1 = (EdbSegP*)eComb1.At(i);
737 s2 = (EdbSegP*)eComb2.At(i);
738 dx = Xcorr(*s2,dz2) - Xcorr(*s1,dz1);
739 dy = Ycorr(*s2,dz2) - Ycorr(*s1,dz1);
740 if(Abs(dx-eXpeak)>dxMax) continue;
741 if(Abs(dy-eYpeak)>dyMax) continue;
742 s.Copy(*s1); s.SetX(Xcorr(*s1,dz1)); s.SetY(Ycorr(*s1,dz1)); s.SetZ(0); p1.AddSegment(s);
743 s.Copy(*s2); s.SetX(Xcorr(*s2,dz2)); s.SetY(Ycorr(*s2,dz2)); s.SetZ(0); p2.AddSegment(s);
744 hpeak->Fill(dx,dy);
745 ic++;
746 }
747 eNpeak = ic;
748 eXpeak = hpeak->GetMean(1);
749 eYpeak = hpeak->GetMean(2);
750 Log(2,"SelectPeak","%4d coincidences selected at XY:(%8.2f %8.2f) +-%6.2f with RMS:(%5.2f %5.2f) Uniformity:(%5.2f %5.2f)",
751 ic, hpeak->GetMean(1), hpeak->GetMean(2), dxMax, hpeak->GetRMS(1), hpeak->GetRMS(2),
752 hpeak->GetRMS(1)/(2*dxMax/Sqrt(12)), hpeak->GetRMS(2)/(2*dyMax/Sqrt(12)) );
753
754 SafeDelete(eAff);
755 eAff = new EdbAffine2D();
756 eAff->ShiftX(hpeak->GetMean(1));
757 eAff->ShiftY(hpeak->GetMean(2));
758
759 SafeDelete(hpeak);
760
761 return ic;
762}
Definition: EdbAffine.h:17
void ShiftX(float d)
Definition: EdbAffine.h:64
void ShiftY(float d)
Definition: EdbAffine.h:65

◆ SelectPeak() [2/2]

int EdbPositionAlignment::SelectPeak ( EdbPattern p1,
EdbPattern p2,
float  dx = 10,
float  dy = 10 
)
inline
121 { return SelectPeak(p1,p2,dx,dy,eZ1peak,eZ2peak); }

◆ SelectZone()

int EdbPositionAlignment::SelectZone ( float  min[2],
float  max[2],
TObjArray &  arr1,
TObjArray &  arr2,
float  maxDens 
)

density unit is N/100/100 microns

484{
486 int ic=0;
487 int nmax = (int)(maxDens*(max[0]-min[0])*(max[1]-min[1]) /100/100);
488 ic += SelectZoneSide(ePC1, min, max, arr1, nmax);
489 float min2[2] = {min[0]-eXcell, min[1]-eYcell };
490 float max2[2] = {max[0]+eXcell, max[1]+eYcell };
491 nmax = (int)(maxDens*(max2[0]-min2[0])*(max2[1]-min2[1]) /100/100);
492 ic += SelectZoneSide(ePC2, min2, max2, arr2, nmax);
493 return ic;
494}
int SelectZoneSide(EdbPatCell2 &pc, float min[2], float max[2], TObjArray &arr, int nmax=kMaxInt)
Definition: EdbPositionAlignment.cxx:497

◆ SelectZoneSide()

int EdbPositionAlignment::SelectZoneSide ( EdbPatCell2 pc,
float  min[2],
float  max[2],
TObjArray &  arr,
int  nmax = kMaxInt 
)
498{
499 TObjArray a(arr.GetSize());
501 return SelectBestComptons(a,arr,nmax);
502}
int SelectObjects(TObjArray &arr)
Definition: EdbCell2.cpp:738
static int SelectBestComptons(TObjArray &a, TObjArray &arr, int nmax)
Definition: EdbPositionAlignment.cxx:505

◆ ShrinkageSelection()

void EdbPositionAlignment::ShrinkageSelection ( float  dzbase)

select the shrinkage

807{
809 ePC1.eShr=ePC2.eShr=1;
810
811 ePC1.eDZ = dzbase; ePC2.eDZ = 0;
812 ShrinkageSelectionSide( *(&ePC1), *(&ePC2), eHShr1, 10, 0.5);
813 ePC2.eDZ = -dzbase; ePC1.eDZ = 0;
814 ShrinkageSelectionSide( *(&ePC2), *(&ePC1), eHShr2, 10, 0.5);
815
816 ePC1.eDZ = dzbase; ePC2.eDZ = 0;
817 ShrinkageSelectionSide( *(&ePC1), *(&ePC2), eHShr1, 50, 0.5);
818 ePC2.eDZ = -dzbase; ePC1.eDZ = 0;
819 ShrinkageSelectionSide( *(&ePC2), *(&ePC1), eHShr2, 50, 0.5);
820
821 Log(2,"ShrinkageSelection","found %7.3f %7.3f", ePC1.eShr, ePC2.eShr );
822}
EdbH2 eHShr2
histogram used for the shrinkage-selection
Definition: EdbPositionAlignment.h:62
EdbH2 eHShr1
histogram used for the shrinkage-selection
Definition: EdbPositionAlignment.h:61
int ShrinkageSelectionSide(EdbPatCell2 &pc1, EdbPatCell2 &pc2, EdbH2 &hshr, int nstep, float deltaShr)
Definition: EdbPositionAlignment.cxx:825

◆ ShrinkageSelectionSide()

int EdbPositionAlignment::ShrinkageSelectionSide ( EdbPatCell2 pc1,
EdbPatCell2 pc2,
EdbH2 hshr,
int  nstep,
float  deltaShr 
)

select the shrinkage for one side

826{
828
829 int npeak=0;
830
831 int n[2] = { nstep, 1 };
832 float min[2] = {pc1.eShr-deltaShr, pc2.eShr-0.001};
833 float max[2] = {pc1.eShr+deltaShr, pc2.eShr+0.001};
834 EdbH2 h12;
835
836 h12.InitH2(n, min, max);
837
838 //pc1.eDXlim = eXcell; pc1.eDYlim = eYcell;
839 //pc1.eDTXlim = 0.1; pc1.eDTYlim = 0.1;
840 //int ir[2] = {1,1};
841
842 pc1.eDXlim = pc1.eDYlim = 15;
843 pc1.eDTXlim = 0.1; pc1.eDTYlim = 0.1;
844 int ir[2] = {0,0};
845
846 pc1.eApplyCorr = pc2.eApplyCorr = true;
847
848 eComb1.Clear();
849 eComb2.Clear();
850
851 float x1,x2;
852 for(int i1=0; i1<n[0]; i1++)
853 {
854 for(int i2=0; i2<n[1]; i2++)
855 {
856 x1 = pc1.eShr = h12.X(i1);
857 x2 = pc2.eShr = h12.Y(i2);
858 int ncomb = pc1.FillCombinations(pc2, ir, eComb1, eComb2, 0);
859 h12.Fill(x1,x2,ncomb);
860 Log(3,"ShrinkageSelection:","ncomb = %6d at (%7.2f %7.2f)", ncomb, x1, x2);
861 }
862 }
863
864 if(gEDBDEBUGLEVEL>2) h12.PrintStat();
865
866 hshr.Copy(h12);
867
868 float vpeak[2];
869 EdbPeak2 p2d(h12);
870 npeak = p2d.FindPeak(vpeak);
871 pc1.eShr = vpeak[0];
872
873 Log(3,"ShrSelection","maxcomb: %6d (%7.2f %7.2f)", npeak,vpeak[0],vpeak[1]);
874 return npeak;
875}
float X(int ix) const
Definition: EdbCell2.h:60
float Y(int iy) const
Definition: EdbCell2.h:61
void Copy(const EdbH2 &h)
Definition: EdbCell2.cpp:45
void PrintStat()
Definition: EdbCell2.cpp:199
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ SpotsFilterOut()

int EdbPositionAlignment::SpotsFilterOut ( int  nmax)

discard cells with nenries > nmax

797{
799 int nout1 = ePC1.DiscardHighCells(nmax);
800 int nout2 = ePC2.DiscardHighCells(nmax);
801 Log(3,"SpotsFilterOut","%d and %d cells discarded", nout1, nout2);
802 return nout1+nout2;
803}
int DiscardHighCells(int nmax)
Definition: EdbCell2.cpp:132

◆ WideSearchXY()

int EdbPositionAlignment::WideSearchXY ( EdbPattern pat1,
EdbPattern pat2,
EdbH2 hxy,
float  dz,
float  xmin,
float  xmax,
float  ymin,
float  ymax 
)

wide prealignment (in case of the basetracks)

1049{
1051
1052 int npeak=0;
1053
1054 ePC1.eDZ = dz/2;
1055 ePC2.eDZ = -dz/2;
1056 ePC1.eApplyCorr = ePC.eApplyCorr = true;
1057 FillArrays(pat1, pat2); // fill at central quote
1058 DoubletsFilterOut(false);
1059
1060 int n[2] = { (int)((xmax-xmin)/eXcell) +1 , (int)((ymax-ymin)/eYcell) +1 };
1061 float min[2] = { -((int)(n[0]/2)+0.5)*eXcell , -((int)(n[1]/2)+0.5)*eYcell };
1062 float max[2] = { ((int)(n[0]/2)+0.5)*eXcell , ((int)(n[1]/2)+0.5)*eYcell };
1063 EdbH2 h12;
1064 h12.InitH2(n, min, max);
1065
1067 ePC1.eDTXlim = 0.01; ePC1.eDTYlim = 0.01;
1068
1069 int ir[2] = {1,1};
1070 float x1,x2;
1071 for(int i1=0; i1<n[0]; i1++)
1072 {
1073 for(int i2=0; i2<n[1]; i2++)
1074 {
1075 x1 = ePC1.eDX = h12.X(i1);
1076 x2 = ePC1.eDY = h12.Y(i2);
1077 eComb1.Clear();
1078 eComb2.Clear();
1079 int ncomb = ePC1.FillCombinations(ePC2, ir, eComb1, eComb2);
1080 h12.Fill(x1,x2,ncomb);
1081 Log(3,"WideSearchXY:","ncomb = %6d at (%7.2f %7.2f)", ncomb, x1, x2);
1082 }
1083 }
1084
1085 h12.PrintStat();
1086 hxy.Copy(h12);
1087
1088 float vpeak[2];
1089 EdbPeak2 p2d(h12);
1090 npeak = p2d.FindPeak(vpeak);
1091 ePC1.eDX = vpeak[0];
1092 ePC1.eDY = vpeak[1];
1093 eComb1.Clear();
1094 eComb2.Clear();
1095 int ncomb = ePC1.FillCombinations(ePC2, ir, eComb1, eComb2);
1096 Log(3,"WideSearchXY","ncomb = %d",ncomb);
1097
1098 int nxy[2] = {21,21}; // TODO
1099 float minxy[2] = {-1000, -1000 };
1100 float maxxy[2] = { 1000, 1000 };
1101 EdbH2 hxyn;
1102 hxyn.InitH2(nxy, minxy, maxxy);
1103 SelectNarrowPeakDXDY(1000, hxyn);
1104
1105 Log(2,"WideSearchXY","PC1: %f %f %f %f PC2:%f %f %f %f ", ePC1.eDX, ePC1.eDY, ePC1.eDZ,ePC1.eShr, ePC2.eDX, ePC2.eDY, ePC2.eDZ, ePC2.eShr);
1106
1107 hxyn.DrawH2("hxyn")->Write("hxyn");
1108
1109 ActivatePosTree("align");
1110 FillPosTree(dz/2, -dz/2, 0);
1111 WritePosTree();
1112
1113 Log(2,"WideSearchXY","maxcomb: %6d (%7.2f %7.2f)", npeak,vpeak[0],vpeak[1]);
1114 return npeak;
1115}
int DoubletsFilterOut(bool checkview=true)
Definition: EdbPositionAlignment.cxx:765

◆ WritePosTree()

Int_t EdbPositionAlignment::WritePosTree ( )
135{
136 if(!ePosTree) return 0;
137 ePosTree->SetAlias("x1","s1.eX+dz1*s1.eTX");
138 ePosTree->SetAlias("x2","s2.eX+dz2*s2.eTX");
139 ePosTree->SetAlias("y1","s1.eY+dz1*s1.eTY");
140 ePosTree->SetAlias("y2","s2.eY+dz2*s2.eTY");
141
142 ePosTree->SetAlias( "tx1c" , Form("s1.eTX/(%f)+(%f)",ePC1.eShr,ePC1.eDTX) );
143 ePosTree->SetAlias( "ty1c" , Form("s1.eTY/(%f)+(%f)",ePC1.eShr,ePC1.eDTY) );
144 ePosTree->SetAlias( "tx2c" , Form("s2.eTX/(%f)+(%f)",ePC2.eShr,ePC2.eDTX) );
145 ePosTree->SetAlias( "ty2c" , Form("s2.eTY/(%f)+(%f)",ePC2.eShr,ePC2.eDTY) );
146
147 ePosTree->SetAlias("x1c" , Form("(%f) + s1.eX + tx1c*(%f)",ePC1.eDX,ePC1.eDZ) );
148 ePosTree->SetAlias("y1c" , Form("(%f) + s1.eY + ty1c*(%f)",ePC1.eDY,ePC1.eDZ) );
149 ePosTree->SetAlias("x2c" , Form("(%f) + s2.eX + tx2c*(%f)",ePC2.eDX,ePC2.eDZ) );
150 ePosTree->SetAlias("y2c" , Form("(%f) + s2.eY + ty2c*(%f)",ePC2.eDY,ePC2.eDZ) );
151
152 ePosTree->SetAlias("dx" , "x2c-x1c" );
153 ePosTree->SetAlias("dy" , "y2c-y1c" );
154 ePosTree->SetAlias("dtx" , "tx2c-tx1c" );
155 ePosTree->SetAlias("dty" , "ty2c-ty1c" );
156
157 ePosTree->SetAlias("dr" , "sqrt(dx*dx+dy*dy)" );
158 ePosTree->SetAlias("dt" , "sqrt(dtx*dtx+dty*dty)" );
159
160 ePosTree->Write();
161 SafeDelete(ePosTree);
162 return 1;
163}

◆ Xcorr() [1/2]

float EdbPositionAlignment::Xcorr ( EdbSegP s,
float  dz 
)
inline
114{return s.X() + dz*s.TX();}

◆ Xcorr() [2/2]

float EdbPositionAlignment::Xcorr ( EdbSegP s,
float  dz,
float  dx 
)
inline
116{return dx + s.X() + dz*s.TX();}

◆ Ycorr() [1/2]

float EdbPositionAlignment::Ycorr ( EdbSegP s,
float  dz 
)
inline
115{return s.Y() + dz*s.TY();}

◆ Ycorr() [2/2]

float EdbPositionAlignment::Ycorr ( EdbSegP s,
float  dz,
float  dy 
)
inline
117{return dy + s.Y() + dz*s.TY();}

◆ Zselection()

void EdbPositionAlignment::Zselection ( EdbH2 hd)
625{
626 if(eNZ1step<1 || eNZ2step<1) {
627 Log(1,"Zselection","ERROR: eNZ1step = %d, eNZ2step = %d",eNZ1step,eNZ2step);
628 return;
629 }
630
631 float maxMax=0, meanMax=0, meanMean=0, meanbin=0;
632 int maxbin=0,ic=0;
633 //float maxPeak = 0;
634
635 int nz[2] = {eNZ1step,eNZ2step};
636 float minz[2] = {Min(eZ1from,eZ1to), Min(eZ2from,eZ2to)};
637 float maxz[2] = {Max(eZ1from,eZ1to), Max(eZ2from,eZ2to)};
638 EdbH2 hz12;
639 hz12.InitH2(nz, minz, maxz);
640 //hz12.PrintStat();
641
642 float dz1,dz2;
643 for(int i1=0; i1<nz[0]; i1++)
644 {
645 for(int i2=0; i2<nz[1]; i2++)
646 {
647 dz1 = hz12.X(i1);
648 dz2 = hz12.Y(i2);
649 hd.CleanCells();
650 PositionPlot(dz1,dz2,hd);
651
652 EdbPeak2 p2di(hd);
653 p2di.Smooth(eSmoothKernel);
654 //maxbin = hd.MaxBin();
655 maxbin = p2di.MaxBin();
656 if(maxbin>maxMax) maxMax=maxbin;
657 hz12.Fill(dz1,dz2,maxbin);
658 //meanbin = hd.Mean();
659 meanbin = p2di.Mean();
660 meanMax += maxbin;
661 meanMean += meanbin;
662 ic++;
663 Log(3,"Zselection","max: %5d mean: %8.4f at dZ(%7.2f %7.2f)", maxbin,meanbin, dz1, dz2);
664 }
665 }
666 meanMax /= ic;
667 meanMean /= ic;
668
669 EdbPeak2 p2z(hz12);
670 p2z.Smooth("k5a");
671 eHDZ.Copy(*((EdbH2*)(&p2z)));
672 p2z.FindGlobalPeak(eZ1peak, eZ2peak, 0.05);
673
674 hd.CleanCells();
676
677 EdbPeak2 p2d(hd);
678 p2d.Smooth(eSmoothKernel);
679 eHpeak.Copy(*((EdbH2*)(&p2d)));
680 p2d.FindPeak(eXpeak, eYpeak);
681
682 // eHpeak.Copy(hd);
683 //float vpeak[2];
684 //EdbPeak2 p2d(eHpeak);
685 //p2d.FindPeak(vpeak);
686 //eXpeak=vpeak[0];
687 //eYpeak=vpeak[1];
688
689 Log(2,"Zselection","maxMax: %8.2f meanMax: %8.4f meanMean: %8.4f at dZ(%7.2f %7.2f) XY(%7.2f %7.2f)",
690 maxMax, meanMax, meanMean, eZ1peak, eZ2peak, eXpeak,eYpeak);
691}
EdbH2 eHDZ
histogram used for Z-selection
Definition: EdbPositionAlignment.h:45
void PositionPlot(float dz1, float dz2, EdbH2 &hd)
Definition: EdbPositionAlignment.cxx:584
EdbH2 eHpeak
histogram used for peak selection
Definition: EdbPositionAlignment.h:44

Member Data Documentation

◆ eAff

EdbAffine2D* EdbPositionAlignment::eAff

the found affine transformation (when applied to pattern1 gives pattern2 )

◆ eBinX

Float_t EdbPositionAlignment::eBinX

◆ eBinY

Float_t EdbPositionAlignment::eBinY

bin size for the differential hist (for example 5 microns)

◆ eChi2Max

Float_t EdbPositionAlignment::eChi2Max

cut for the chi2 of the basetrack

◆ eComb1

TObjArray EdbPositionAlignment::eComb1

◆ eComb2

TObjArray EdbPositionAlignment::eComb2

array with the selected combinations of segments

◆ eDoRot

Bool_t EdbPositionAlignment::eDoRot

if true - perform the rotation selection

◆ eDTXmax

Float_t EdbPositionAlignment::eDTXmax

◆ eDTXmin

Float_t EdbPositionAlignment::eDTXmin

◆ eDTYmax

Float_t EdbPositionAlignment::eDTYmax

max angular acceptance (ex: 0.15) for the coinsidences

◆ eDTYmin

Float_t EdbPositionAlignment::eDTYmin

min angular difference for the dublets cutout

◆ eDXmin

Float_t EdbPositionAlignment::eDXmin

◆ eDYmin

Float_t EdbPositionAlignment::eDYmin

min position difference for the dublets cutout

◆ eHDZ

EdbH2 EdbPositionAlignment::eHDZ

histogram used for Z-selection

◆ eHpeak

EdbH2 EdbPositionAlignment::eHpeak

histogram used for peak selection

◆ eHShr1

EdbH2 EdbPositionAlignment::eHShr1

histogram used for the shrinkage-selection

◆ eHShr2

EdbH2 EdbPositionAlignment::eHShr2

histogram used for the shrinkage-selection

◆ eNpeak

Int_t EdbPositionAlignment::eNpeak

number of found combinations

◆ eNShr1step

Int_t EdbPositionAlignment::eNShr1step

◆ eNShr2step

Int_t EdbPositionAlignment::eNShr2step

number of steps for the Shr-selection

◆ eNZ1step

Int_t EdbPositionAlignment::eNZ1step

◆ eNZ2step

Int_t EdbPositionAlignment::eNZ2step

number of steps for the z-selection

◆ ePC

EdbPatCell2 EdbPositionAlignment::ePC

cells with the pointers to tracks

◆ ePC1

EdbPatCell2 EdbPositionAlignment::ePC1

◆ ePC2

EdbPatCell2 EdbPositionAlignment::ePC2

cells with the pointers to segments

◆ ePosTree

TTree* EdbPositionAlignment::ePosTree

optional: tree with s1:s2 for the selected combinations

◆ eRot

EdbH1 EdbPositionAlignment::eRot

definition of the rotation steps

◆ eRpeak

Float_t EdbPositionAlignment::eRpeak

coordinate peak acceptance

◆ eS0tx

Float_t EdbPositionAlignment::eS0tx

◆ eS0ty

Float_t EdbPositionAlignment::eS0ty

sigmas at 0 angle for the chi2 calculation

◆ eS0x

Float_t EdbPositionAlignment::eS0x

◆ eS0y

Float_t EdbPositionAlignment::eS0y

◆ eSegCouples

TObjArray EdbPositionAlignment::eSegCouples

segment couples objects to fill couples format tree

◆ eShr1from

Float_t EdbPositionAlignment::eShr1from

◆ eShr1to

Float_t EdbPositionAlignment::eShr1to

limits in Shr for the peak search (ePat1)

◆ eShr2from

Float_t EdbPositionAlignment::eShr2from

◆ eShr2to

Float_t EdbPositionAlignment::eShr2to

limits in Shr for the peak search (ePat2)

◆ eSmoothKernel

TString EdbPositionAlignment::eSmoothKernel

used to smooth histograms

◆ eTracks

TObjArray EdbPositionAlignment::eTracks

tracks created with segments of eComb1, eComb2

◆ eWbaseMin

Float_t EdbPositionAlignment::eWbaseMin

cut for the w of basetrack to accept it

◆ eX0

float EdbPositionAlignment::eX0

◆ eXcell

Float_t EdbPositionAlignment::eXcell

◆ eXpeak

Float_t EdbPositionAlignment::eXpeak

◆ eY0

float EdbPositionAlignment::eY0

coordinates of the center of the zone (for the ePeakNT only)

◆ eYcell

Float_t EdbPositionAlignment::eYcell

cell size (for example 50 microns)

◆ eYpeak

Float_t EdbPositionAlignment::eYpeak

peak position in X,Y

◆ eZ1from

Float_t EdbPositionAlignment::eZ1from

◆ eZ1peak

Float_t EdbPositionAlignment::eZ1peak

◆ eZ1to

Float_t EdbPositionAlignment::eZ1to

limits in Z for the peak search (ePat1)

◆ eZ2from

Float_t EdbPositionAlignment::eZ2from

◆ eZ2peak

Float_t EdbPositionAlignment::eZ2peak

peak position in Z

◆ eZ2to

Float_t EdbPositionAlignment::eZ2to

limits in Z for the peak search (ePat2)


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