FEDRA emulsion software from the OPERA Collaboration
EdbAlignmentMap Class Reference

2-d alignment map finder More...

#include <EdbAlignmentMap.h>

Inheritance diagram for EdbAlignmentMap:
Collaboration diagram for EdbAlignmentMap:

Public Member Functions

void AlignMap ()
 
int ApplyMap (EdbPattern &pat)
 
int CheckDZbase (TEnv &cenv, EdbPattern &p1, EdbPattern &p2)
 
void CheckPattern (EdbPattern &p, const char *suffix)
 
int CheckXY (TEnv &cenv, EdbPattern &p1all, EdbPattern &p2all)
 
 EdbAlignmentMap (const char *file=0, const char *mode=0)
 
void ExtractMapFromTree ()
 
int FillMapTree (EdbPositionAlignment &pol, int izone)
 
void get_run_patterns (const char *runfile, TEnv *cenv, EdbPattern &p1, EdbPattern &p2)
 
int InitFile (const char *file=0, const char *mode=0)
 
int Link ()
 
void Link (const char *file, EdbPlateP &plate)
 
int Link (TEnv &cenv, EdbPattern &p1, EdbPattern &p2, EdbPlateP &plate)
 
void SaveAll ()
 
void SaveMap (const char *file)
 
int SelectSampleForShrinkageCorr (EdbPattern &p, EdbPattern &psel, float wmin, float wmax, float tmin, float tmax)
 
virtual ~EdbAlignmentMap ()
 

Static Public Member Functions

static void GetPostreeAsPat (EdbPattern &pat, const char *filename)
 

Public Attributes

float eDensityMax
 the max segments density for patterns selection (in N/100/100 microns) More...
 
TEnv * eEnv
 environment used to pass the parameters More...
 
EdbPositionAlignment eGlobal
 service object for the zones selection More...
 
TGraph2D * eGraphDX
 
TGraph2D * eGraphDY
 
TGraph2D * eGraphDZ1
 
TGraph2D * eGraphDZ2
 
TVector3 eGV1
 
TVector3 eGV2
 global offsets found by the test alignment (dx:dy:dz) More...
 
EdbCell2 eMap
 in each cell TArratF with (x,y, dx,dy,dz1,dz2, n) More...
 
TTree * eMapTree
 
int eNx
 
int eNy
 number of divisions calculated using eXcell, eYcell More...
 
TFile * eOutputFile
 
EdbPatternePat1
 big patterns to be splitted and aligned piece by piece More...
 
EdbPatternePat2
 
EdbPlatePePlate
 
float eXcell
 
float eYcell
 approximate zones size (for example 10000 microns) More...
 

Detailed Description

2-d alignment map finder

Constructor & Destructor Documentation

◆ EdbAlignmentMap()

EdbAlignmentMap::EdbAlignmentMap ( const char *  file = 0,
const char *  mode = 0 
)
26{
27 eEnv = 0;
28 ePat1 = ePat2 = 0;
29 eOutputFile = 0;
30 eMapTree = 0;
31 eNx=eNy=0;
33 InitFile(file, mode);
34}
EdbPattern * ePat2
Definition: EdbAlignmentMap.h:20
TEnv * eEnv
environment used to pass the parameters
Definition: EdbAlignmentMap.h:18
int InitFile(const char *file=0, const char *mode=0)
Definition: EdbAlignmentMap.cxx:47
int eNy
number of divisions calculated using eXcell, eYcell
Definition: EdbAlignmentMap.h:30
EdbPattern * ePat1
big patterns to be splitted and aligned piece by piece
Definition: EdbAlignmentMap.h:19
int eNx
Definition: EdbAlignmentMap.h:30
TFile * eOutputFile
Definition: EdbAlignmentMap.h:39
TGraph2D * eGraphDY
Definition: EdbAlignmentMap.h:35
TGraph2D * eGraphDZ1
Definition: EdbAlignmentMap.h:36
TGraph2D * eGraphDX
Definition: EdbAlignmentMap.h:34
TGraph2D * eGraphDZ2
Definition: EdbAlignmentMap.h:37
TTree * eMapTree
Definition: EdbAlignmentMap.h:40
TFile * file
Definition: write_pvr.C:3

◆ ~EdbAlignmentMap()

EdbAlignmentMap::~EdbAlignmentMap ( )
virtual
38{
39 SaveAll();
40 if(eOutputFile) {
41 eOutputFile->Close();
42 SafeDelete(eOutputFile);
43 }
44}
void SaveAll()
Definition: EdbAlignmentMap.cxx:86

Member Function Documentation

◆ AlignMap()

void EdbAlignmentMap::AlignMap ( )
101{
102 eDensityMax = eEnv->GetValue("comptonmap.DensityMax" , 1.);
103
104 for(int i=0; i<ePat1->N(); i++) ePat1->GetSegment(i)->SetFlag(0);
105 for(int i=0; i<ePat2->N(); i++) ePat2->GetSegment(i)->SetFlag(0);
106 eYcell = eXcell = eEnv->GetValue("comptonmap.zonesize" , 9000);
107
108 eGlobal.eXcell = 200; // set approximate bin size in microns
109 eGlobal.eYcell = 200;
110 eGlobal.eDXmin = eGlobal.eDYmin = eEnv->GetValue("comptonmap.DoubletsDX" , 5.); // acceptance to remove doublets
111 eGlobal.eDTXmin = eGlobal.eDTYmin = eEnv->GetValue("comptonmap.DoubletsDTX" , 0.005);
113
114 TH2F *hpat1 = eGlobal.ePC1.DrawH2("hpat1");
115 TH2F *hpat2 = eGlobal.ePC2.DrawH2("hpat2");
116 hpat1->Write("hpat1");
117 hpat2->Write("hpat2");
118
119 TH1I *spectr1 = eGlobal.ePC1.DrawSpectrum("spectr1");
120 spectr1->Write("spectrum1");
121 TH1I *spectr2 = eGlobal.ePC2.DrawSpectrum("spectr2");
122 spectr2->Write("spectrum2");
123
124 //eGlobal.DoubletsFilterOut(true);
125
126 eGlobal.ActivatePosTree("doublets");
129
130 eGlobal.SpotsFilterOut(30); // remove overoccupated cells
131
132 if(gEDBDEBUGLEVEL>1) {
135 EdbH2 *pat1sel = eGlobal.ePC1.FillSelectedH2();
136 EdbH2 *pat2sel = eGlobal.ePC2.FillSelectedH2();
137 hpat1 = pat1sel->DrawH2("hpat1sel"); hpat1->Write("hpat1sel");
138 hpat2 = pat2sel->DrawH2("hpat2sel"); hpat2->Write("hpat2sel");
139 }
140
142 local.eYcell = local.eXcell = eEnv->GetValue("comptonmap.Xcell" , 100); // set approximate bin size in microns
143 local.eDTXmax = local.eDTYmax = eEnv->GetValue("comptonmap.DTmax" , 0.15); // acceptance for coincidence selections
144 local.eNZ1step = eEnv->GetValue("comptonmap.NZ1step" , 100);
145 local.eZ1from = eEnv->GetValue("comptonmap.Z1from" , -150);
146 local.eZ1to = eEnv->GetValue("comptonmap.Z1to" , 50);
147 local.eNZ2step = eEnv->GetValue("comptonmap.NZ2step" , 100);
148 local.eZ2from = eEnv->GetValue("comptonmap.Z2from" , -50);
149 local.eZ2to = eEnv->GetValue("comptonmap.Z2to" , 150);
150 //local.eNZ1step = local.eNZ2step = 15; // set big step in z for the global rough alignment
151 local.eBinY = local.eBinX = eEnv->GetValue("comptonmap.BinXY" , 10); // bin size for the diff hist
152 local.eRpeak = eEnv->GetValue("comptonmap.Rpeak" , 4.); // position peak acceptance
153 float testzone = eEnv->GetValue("comptonmap.TestZone" , 30000);
154 //eEnv->Print();
155
156 float center[2] = { (eGlobal.ePC1.Xmax()+eGlobal.ePC1.Xmin())/2, (eGlobal.ePC1.Ymax()+eGlobal.ePC1.Ymin())/2};
157 float min[2] = { Max( eGlobal.ePC1.Xmin(), center[0]-testzone/2),
158 Max( eGlobal.ePC1.Ymin(), center[1]-testzone/2) };
159 float max[2] = { Min( eGlobal.ePC1.Xmax(), center[0]+testzone/2),
160 Min( eGlobal.ePC1.Ymax(), center[1]+testzone/2) };
161 Log(1,"Global alignment","select area in limits: X(%f %f), Y(%f %f) with max dens: %f ",min[0],max[0],min[1],max[1],eDensityMax);
162
163 TObjArray arr1(50000), arr2(50000);
165 float areaGlobal = (max[0]-min[0])*(max[1]-min[1]);
166
167 printf("******* Global alignment: n1 = %d n2 = %d *******\n",arr1.GetEntriesFast(), arr2.GetEntriesFast() );
168 local.FillArrays(arr1, arr2, min, max);
169 local.FillCombinations();
170 local.eSmoothKernel = "k5a";
171 float normKernel = 25;
172 local.Align(); // global peak selection
173 //local.eAff->Print();
174 TH2F *hdxy = local.eHpeak.DrawH2("global_hdxy");
175 TH2F *hdz12 = local.eHDZ.DrawH2("global_hdz12");
176 hdxy->Write("global_hdxy");
177 hdz12->Write("global_hdz12");
178
179 local.PrintSummary();
180 eGV1.SetXYZ(local.eXpeak,local.eYpeak,local.eZ1peak);
181 eGV2.SetXYZ( 0, 0,local.eZ2peak);
182 EdbPeak2 p2d(local.eHpeak);
183 p2d.ProbPeaks(2); // prob for 2 most significant peaks
184 float peakVolume = p2d.EstimatePeakVolumeSafe(0) / normKernel;
185 printf("******* End of Global alignment *******\n\n");
186
187 // do global correction:
188 ePat1->ProjectTo(local.eZ1peak);
189 ePat2->ProjectTo(local.eZ2peak);
190 ePat1->Transform(local.eAff);
193
194 // define the best zone size using the results of the global alignment
195 float requestPeakVolume = eEnv->GetValue("comptonmap.PeakVolume" , 150.);
196 float zonesize = Sqrt( areaGlobal * requestPeakVolume/peakVolume );
197 Log(1,"EdbAlignmentMap::AlignMap",
198 "zonesize selected as %10.2f considering global peakVolume = %10.2f and requested zonePeakVolume as %10.2f",
199 zonesize, peakVolume, requestPeakVolume);
200
201 // fill initial map
202 min[0] = eGlobal.ePC1.Xmin(); min[1] = eGlobal.ePC1.Ymin();
203 max[0] = eGlobal.ePC1.Xmax(); max[1] = eGlobal.ePC1.Ymax();
204 //int n[2] = { (int)((max[0]-min[0])/eXcell+1), (int)((max[1]-min[1])/eYcell+1) };
205 int n[2] = { (int)((max[0]-min[0])/zonesize+1), (int)((max[1]-min[1])/zonesize+1) };
206 eNx = n[0];
207 eNy = n[1];
208 eEnv->SetValue("comptonmap.eNx",eNx);
209 eEnv->SetValue("comptonmap.eNy",eNy);
210 eMap.InitCell(1,n,min,max);
211 local.ResetPeak();
212
213 local.eNZ1step = 25;
214 local.eZ1from = -50;
215 local.eZ1to = 50;
216 local.eNZ2step = 25;
217 local.eZ2from = -50;
218 local.eZ2to = 50;
219 //local.eXcell = local.eYcell = 100;
220
221 eMap.PrintStat();
222
223 int izone = 0;
224 for(int iy=0; iy<eMap.NY(); iy++)
225 for(int ix=0; ix<eMap.NX(); ix++) {
226 min[0] = eMap.X(ix) - eMap.Xbin()/2;
227 min[1] = eMap.Y(iy) - eMap.Ybin()/2;
228 max[0] = eMap.X(ix) + eMap.Xbin()/2;
229 max[1] = eMap.Y(iy) + eMap.Ybin()/2;
230 local.eX0 = eMap.X(ix);
231 local.eY0 = eMap.Y(iy);
232 printf("\n(%d:%d)\n", ix,iy);
233 arr1.Clear();
234 arr2.Clear();
236 local.eXcell = local.eYcell = eEnv->GetValue("comptonmap.Xcell" , 100);
237 local.FillArrays(arr1, arr2, min, max);
238 printf("align with n1 = %d n2 = %d\n",arr1.GetEntriesFast(), arr2.GetEntriesFast());
239 local.FillCombinations();
240 local.Align(); // local peak selection
241
242 FillMapTree( local, izone );
243 local.PrintSummary();
244 local.ResetPeak();
245 izone++;
246 }
247}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
float min(TClonesArray *t)
Definition: bitview.cxx:275
float eDensityMax
the max segments density for patterns selection (in N/100/100 microns)
Definition: EdbAlignmentMap.h:21
TVector3 eGV1
Definition: EdbAlignmentMap.h:27
TVector3 eGV2
global offsets found by the test alignment (dx:dy:dz)
Definition: EdbAlignmentMap.h:27
int FillMapTree(EdbPositionAlignment &pol, int izone)
Definition: EdbAlignmentMap.cxx:250
float eYcell
approximate zones size (for example 10000 microns)
Definition: EdbAlignmentMap.h:29
float eXcell
Definition: EdbAlignmentMap.h:29
EdbPositionAlignment eGlobal
service object for the zones selection
Definition: EdbAlignmentMap.h:25
EdbCell2 eMap
in each cell TArratF with (x,y, dx,dy,dz1,dz2, n)
Definition: EdbAlignmentMap.h:32
int InitCell(EdbCell2 &c)
Definition: EdbCell2.h:167
void PrintStat()
Definition: EdbCell2.cpp:674
fast 2-dim histogram class (used as a basis for EdbCell2)
Definition: EdbCell2.h:19
float X(int ix) const
Definition: EdbCell2.h:60
float Xmax() const
Definition: EdbCell2.h:65
int NX() const
Definition: EdbCell2.h:50
float Ymin() const
Definition: EdbCell2.h:66
float Ybin() const
Definition: EdbCell2.h:78
float Y(int iy) const
Definition: EdbCell2.h:61
float Xbin() const
Definition: EdbCell2.h:77
TH1I * DrawSpectrum(const char *name="plot1d", const char *title="EdbH2 DrawSpectrun")
Definition: EdbCell2.cpp:178
int NY() const
Definition: EdbCell2.h:51
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
EdbH2 * FillSelectedH2()
Definition: EdbPatCell2.cxx:171
peak analyser for EdbH2
Definition: EdbCell2.h:105
virtual void Transform(const EdbAffine2D *a)
Definition: EdbVirtual.cxx:155
new alignment class developed mainly for compton search
Definition: EdbPositionAlignment.h:14
Float_t eRpeak
coordinate peak acceptance
Definition: EdbPositionAlignment.h:32
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
int SelectZone(float min[2], float max[2], TObjArray &arr1, TObjArray &arr2, float maxDens)
Definition: EdbPositionAlignment.cxx:483
Float_t eDTYmax
max angular acceptance (ex: 0.15) for the coinsidences
Definition: EdbPositionAlignment.h:27
void PrintSummary()
Definition: EdbPositionAlignment.cxx:76
Float_t eZ1from
Definition: EdbPositionAlignment.h:35
Float_t eXpeak
Definition: EdbPositionAlignment.h:41
bool ActivatePosTree(const char *name="postree")
Definition: EdbPositionAlignment.cxx:105
int DoubletsFilterOut(bool checkview=true)
Definition: EdbPositionAlignment.cxx:765
Float_t eYpeak
peak position in X,Y
Definition: EdbPositionAlignment.h:41
EdbH2 eHDZ
histogram used for Z-selection
Definition: EdbPositionAlignment.h:45
int Align()
Definition: EdbPositionAlignment.cxx:695
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
int FillArrays(EdbPattern &p1, EdbPattern &p2)
Definition: EdbPositionAlignment.cxx:416
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 eBinY
bin size for the differential hist (for example 5 microns)
Definition: EdbPositionAlignment.h:28
Float_t eBinX
Definition: EdbPositionAlignment.h:28
void ResetPeak()
Definition: EdbPositionAlignment.h:123
Float_t eYcell
cell size (for example 50 microns)
Definition: EdbPositionAlignment.h:26
int WritePosTree()
Definition: EdbPositionAlignment.cxx:134
int SpotsFilterOut(int nmax)
Definition: EdbPositionAlignment.cxx:796
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
EdbPatCell2 ePC2
cells with the pointers to segments
Definition: EdbPositionAlignment.h:17
Float_t eDYmin
min position difference for the dublets cutout
Definition: EdbPositionAlignment.h:30
EdbH2 eHpeak
histogram used for peak selection
Definition: EdbPositionAlignment.h:44
Int_t eNZ1step
Definition: EdbPositionAlignment.h:37
Float_t eZ2peak
peak position in Z
Definition: EdbPositionAlignment.h:40
EdbPatCell2 ePC1
Definition: EdbPositionAlignment.h:17
int FillCombinations()
Definition: EdbPositionAlignment.cxx:534
void SetFlag(int flag)
Definition: EdbSegP.h:130
void SetZ(float z)
Definition: EdbPattern.h:41
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
void ProjectTo(const float dz)
Definition: EdbPattern.cxx:268
void SetSegmentsZ()
Definition: EdbPattern.cxx:250
int max
Definition: check_shower.C:41
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ ApplyMap()

int EdbAlignmentMap::ApplyMap ( EdbPattern pat)

applied to "pat1" gives "pat2"

319{
321
322 Log(1,"ApplyMap","to pattern with %d segments",pat.N());
323 if(!eGraphDX) return 0;
324 if(!eGraphDY) return 0;
325 EdbSegP *s;
326 int n = pat.N();
327 for(int i=0; i<n; i++) {
328 s = pat.GetSegment(i);
329 float dx = eGraphDX->Interpolate(s->X(),s->Y());
330 float dy = eGraphDY->Interpolate(s->X(),s->Y());
331 //float dz = eGraphDZ1->Interpolate(s->X(),s->Y());
332 float dz = 0;
333 s->SetZ( dz + s->Z());
334 s->SetX( dx + (s->X() + s->TX()*dz));
335 s->SetY( dy + (s->Y() + s->TY()*dz));
336 }
337
338 /*
339 n = pat2.N();
340 for(int i=0; i<n; i++) {
341 s = pat2.GetSegment(i);
342 float dz = eGraphDZ2->Interpolate(s->X(),s->Y());
343 float dx=0, dy=0;
344 s->SetZ( dz + s->Z());
345 s->SetX( dx + (s->X() + s->TX()*dz));
346 s->SetY( dy + (s->Y() + s->TY()*dz));
347 }
348 */
349 Log(1,"ApplyMap","%d segments corrected",n);
350
351 return n;
352}
brick dz
Definition: RecDispMC.C:107
Definition: EdbSegP.h:21
s
Definition: check_shower.C:55

◆ CheckDZbase()

int EdbAlignmentMap::CheckDZbase ( TEnv &  cenv,
EdbPattern p1,
EdbPattern p2 
)

test "alignment" to check dz

572{
574
575 float wmin = cenv.GetValue("fedra.link.shr.Wmin", 9 );
576 float wmax = cenv.GetValue("fedra.link.shr.Wmax", 100 );
577 float tmin = cenv.GetValue("fedra.link.shr.Tmin", 0.01 );
578 float tmax = cenv.GetValue("fedra.link.shr.Tmax", 0.5 );
579
580 EdbPattern p1sel,p2sel;
581 SelectSampleForShrinkageCorr(p1all,p1sel, wmin, wmax, tmin, tmax);
582 SelectSampleForShrinkageCorr(p2all,p2sel, wmin, wmax, tmin, tmax);
583
585 //atest.eXcell = atest.eYcell = cenv.GetValue("fedra.link.CellXY", 10 );
586 atest.eXcell = atest.eYcell = cenv.GetValue("fedra.link.CellXY", 100 ); // enlarge
587 atest.eBinX = atest.eBinY = cenv.GetValue("fedra.link.BinXY", 4 );
588 atest.eDTXmax = atest.eDTYmax = cenv.GetValue("fedra.link.DTMax", 0.02 );
589
590 atest.eNZ1step = 35;
591 atest.eZ1from = -350;
592 atest.eZ1to = 350;
593 atest.eNZ2step = 35;
594 atest.eZ2from = -350;
595 atest.eZ2to = 350;
596
597
598 float DZbase = cenv.GetValue("fedra.link.DZbase", 214.);
599 //p1sel.SetZ(0); p1sel.SetSegmentsZ();
600 //p2sel.SetZ(DZbase); p2sel.SetSegmentsZ();
601 p1sel.ProjectTo(DZbase/2.);
602 p2sel.ProjectTo(-DZbase/2.);
603
604 atest.SaveAsTree(p1sel,"p1sel");
605 atest.SaveAsTree(p2sel,"p2sel");
606
607 atest.FillArrays(p1sel,p2sel);
608 printf( "align with n1 = %d n2 = %d\n",p1sel.N(), p2sel.N() );
609 atest.FillCombinations();
610 atest.Align(); // local peak selection
611
612 FillMapTree( atest, 0 );
613
614 return 0;
615}
int SelectSampleForShrinkageCorr(EdbPattern &p, EdbPattern &psel, float wmin, float wmax, float tmin, float tmax)
Definition: EdbAlignmentMap.cxx:618
Definition: EdbPattern.h:273
void SaveAsTree(EdbPattern &pat, const char *name)
Definition: EdbPositionAlignment.cxx:270
TEnv cenv("emrec")

◆ CheckPattern()

void EdbAlignmentMap::CheckPattern ( EdbPattern p,
const char *  suffix 
)

produce control plots for this pattern

641{
643
644 TH1F *hw = new TH1F( Form("w_%s",suffix), "W", 100, 0, 100);
645 TH2F *htxty = new TH2F( Form("txty_%s",suffix), "W",
646 100, -1.5, 1.5, 100, -1.5,1.5 );
647
648 int n = p.N();
649 EdbSegP *s;
650 for(int i=0; i<n; i++) {
651 s = p.GetSegment(i);
652 hw->Fill( s->W() );
653 htxty->Fill( s->TX(), s->TY());
654 }
655
656 hw->Write();
657 htxty->Write();
658
659 Log(2,"CheckPattern","with %d segments",n);
660}
p
Definition: testBGReduction_AllMethods.C:8

◆ CheckXY()

int EdbAlignmentMap::CheckXY ( TEnv &  cenv,
EdbPattern p1all,
EdbPattern p2all 
)

test "alignment" to check DX and DY

743{
745
746 float wmin = cenv.GetValue("fedra.align.Wmin", 16 );
747 float wmax = cenv.GetValue("fedra.align.Wmax", 100 );
748 float tmin = cenv.GetValue("fedra.align.Tmin", 0.03 );
749 float tmax = cenv.GetValue("fedra.align.Tmax", 0.5 );
750
751 EdbPattern p1sel,p2sel;
752 SelectSampleForShrinkageCorr(p1all,p1sel, wmin, wmax, tmin, tmax);
753 SelectSampleForShrinkageCorr(p2all,p2sel, wmin, wmax, tmin, tmax);
754
756 atest.eXcell = atest.eYcell = cenv.GetValue("fedra.align.CellXY", 100 ); // enlarge
757 atest.eBinX = atest.eBinY = cenv.GetValue("fedra.align.BinXY", 4 );
758 atest.eDTXmax = atest.eDTYmax = cenv.GetValue("fedra.align.DTMax", 0.01 );
759
760 float DZ = cenv.GetValue("fedra.align.DZ", 1300);
761
762 EdbH2 hxy;
763 float xmin = -2000, ymin = -2000, xmax = 2000, ymax = 2000;
764
765 //for(int i=0; i<=11; i++) {
766 //float dz = DZ+(i-5)*500;
767 atest.WideSearchXY(p1sel, p2sel, hxy, DZ, xmin, xmax, ymin, ymax);
768 hxy.DrawH2("hxy")->Write(Form("hxy"));
769 //}
770
771 return 0;
772}
int WideSearchXY(EdbPattern &pat1, EdbPattern &pat2, EdbH2 &hxy, float dz, float xmin, float xmax, float ymin, float ymax)
Definition: EdbPositionAlignment.cxx:1048
float DZ
Definition: hwinit.C:66

◆ ExtractMapFromTree()

void EdbAlignmentMap::ExtractMapFromTree ( )

assume that in eMap is defined the map structure

356{
358 int n[2] = {eMap.NX()+2, eMap.NY()+2}; // define map with margins of one bin
359 float min[2] = {eMap.Xmin()-eMap.Xbin(), eMap.Ymin()-eMap.Ybin()};
360 float max[2] = {eMap.Xmax()+eMap.Xbin(), eMap.Ymax()+eMap.Ybin()};
361
362 EdbCell2 map;
363 map.InitCell(1,n,min,max); // in each cell: TarrayF of (x,y, dx,dy,dz1,dz2, n)
364
365 int N,IZONE;
366 float X, Y, DX, DY, DZ1, DZ2;
367 eMapTree->ResetBranchAddresses();
368 eMapTree->SetBranchAddress("izone",&IZONE);
369 eMapTree->SetBranchAddress("peak" ,&N);
370 eMapTree->SetBranchAddress("x0" ,&X);
371 eMapTree->SetBranchAddress("y0" ,&Y);
372 eMapTree->SetBranchAddress("dx" ,&DX);
373 eMapTree->SetBranchAddress("dy" ,&DY);
374 eMapTree->SetBranchAddress("dz1" ,&DZ1);
375 eMapTree->SetBranchAddress("dz2" ,&DZ2);
376 int nentr = eMapTree->GetEntries();
377 for(int i=0; i<nentr; i++) {
378 eMapTree->GetEntry(i);
379 TArrayF *arr = new TArrayF(7);
380 arr->AddAt( X , 0);
381 arr->AddAt( Y , 1);
382 arr->AddAt( DX +eGV1.X(), 2);
383 arr->AddAt( DY +eGV1.Y(), 3);
384 arr->AddAt( DZ1+eGV1.Z(), 4);
385 arr->AddAt( DZ2 , 5);
386 arr->AddAt( N , 6);
387 map.AddObject(X,Y,(TObject*)arr);
388 }
389
390 int ix0,iy0;
391 for(int ix=0; ix<map.NX(); ix++)
392 for(int iy=0; iy<map.NY(); iy++) {
393 ix0 = ix; iy0 = iy;
394 if(ix==0) ix0=1;
395 if(ix==map.NX()-1 ) ix0=map.NX()-2;
396 if(iy==0) iy0=1;
397 if(iy==map.NY()-1 ) iy0=map.NY()-2;
398 if( ix0!=ix || iy0!=iy ) {
399 TArrayF *arr0 = (TArrayF*)map.GetObject(ix0,iy0,0);
400 TArrayF *arr = new TArrayF(*arr0);
401 arr->AddAt( map.X(ix), 0 );
402 arr->AddAt( map.Y(iy), 1 );
403 map.AddObject(ix,iy,(TObject*)arr);
404 }
405 }
406
407 Log(2, "EdbAlignmentMap::ExtractMapFromTree","with %d modes", map.Ncell() );
408 if(gEDBDEBUGLEVEL>=2) map.PrintStat();
409
410 int nentries = map.Ncell();
411 Double_t *vx = new Double_t[nentries];
412 Double_t *vy = new Double_t[nentries];
413 Double_t *vdx = new Double_t[nentries];
414 Double_t *vdy = new Double_t[nentries];
415 Double_t *vdz1 = new Double_t[nentries];
416 Double_t *vdz2 = new Double_t[nentries];
417
418 for(int i=0; i<nentries; i++) {
419 TArrayF *arr = (TArrayF*)map.GetObject(i,0);
420 if(!arr) break;
421 vx[i] = arr->At(0);
422 vy[i] = arr->At(1);
423 vdx[i] = arr->At(2);
424 vdy[i] = arr->At(3);
425 vdz1[i] = arr->At(4);
426 vdz2[i] = arr->At(5);
427 }
428
429 eGraphDX = new TGraph2D("graphDX","graphDX",nentries,vx,vy,vdx);
430 eGraphDY = new TGraph2D("graphDY","graphDY",nentries,vx,vy,vdy);
431 eGraphDZ1 = new TGraph2D("graphDZ1","graphDZ1",nentries,vx,vy,vdz1);
432 eGraphDZ2 = new TGraph2D("graphDZ2","graphDZ2",nentries,vx,vy,vdz2);
433 //eGraphDZ = new TGraph2D("graphDZ","graphDZ",nentries,vx,vy,vdz);
434
435 if(eOutputFile) {
436 if(eOutputFile->IsWritable()) {
437 eGraphDX->Write( "graphDX");
438 eGraphDY->Write( "graphDY");
439 eGraphDZ1->Write("graphDZ1");
440 eGraphDZ2->Write("graphDZ2");
441 eOutputFile->Purge();
442 }
443 }
444}
class to group 2-dim objects
Definition: EdbCell2.h:148
TObject * GetObject(float x, float y, int ientr) const
Definition: EdbCell2.h:195
bool AddObject(float v[2], TObject *obj)
Definition: EdbCell2.h:173
int Ncell() const
Definition: EdbCell2.h:49
int nentries
Definition: check_shower.C:40
Double_t X
Definition: tlg2couples.C:76
Double_t Y
Definition: tlg2couples.C:76

◆ FillMapTree()

int EdbAlignmentMap::FillMapTree ( EdbPositionAlignment pol,
int  izone 
)
251{
252 if(!eMapTree) return 0;
253
254 float peakbin=0, mean3=0, mean=0;
255 if(pol.eHpeak.Integral()>0) {
256 EdbPeak2 p2d(pol.eHpeak);
257 p2d.ProbPeaks(2); // prob for 2 most significant peaks
258 peakbin = p2d.ePeak[0];
259 mean3 = p2d.eMean3[0];
260 mean = p2d.eMean[0];
261 }
262
263 TH2F *hdxy = pol.eHpeak.DrawH2("hdxy");
264 TH2F *hdz12 = pol.eHDZ.DrawH2("hdz12");
265
266 eMapTree->SetBranchAddress("izone", &izone);
267 eMapTree->SetBranchAddress("x0", &pol.eX0);
268 eMapTree->SetBranchAddress("y0", &pol.eY0);
269 eMapTree->SetBranchAddress("peak", &pol.eNpeak);
270 eMapTree->SetBranchAddress("peakbin", &peakbin);
271 eMapTree->SetBranchAddress("mean3", &mean3);
272 eMapTree->SetBranchAddress("mean", &mean);
273 eMapTree->SetBranchAddress("dx", &pol.eXpeak);
274 eMapTree->SetBranchAddress("dy", &pol.eYpeak);
275 eMapTree->SetBranchAddress("dz1", &pol.eZ1peak);
276 eMapTree->SetBranchAddress("dz2", &pol.eZ2peak);
277 eMapTree->SetBranchAddress("Hdxy", &hdxy);
278 eMapTree->SetBranchAddress("Hdz12", &hdz12);
279 eMapTree->Fill();
280
281 SafeDelete(hdxy);
282 SafeDelete(hdz12);
283
284 pol.WritePosTree();
285
286 return 1;
287}
Long_t Integral()
Definition: EdbCell2.cpp:216
Int_t eNpeak
number of found combinations
Definition: EdbPositionAlignment.h:42

◆ get_run_patterns()

void EdbAlignmentMap::get_run_patterns ( const char *  runfile,
TEnv *  cenv,
EdbPattern p1,
EdbPattern p2 
)

read patterns from the raw files in the root format

665{
667
668 Log(3," EdbAlignmentMap::get_run_patterns","from file: %s",runfile);
669
670 float x = cenv->GetValue("comptonmap.x0" , 40000.);
671 float y = cenv->GetValue("comptonmap.y0" , 40000.);
672 float size = cenv->GetValue("comptonmap.size" , 500000.);
673 float pulsmin = cenv->GetValue("comptonmap.pulsmin" , 8);
674 int side1 = cenv->GetValue("comptonmap.side1" , 1);
675 int side2 = cenv->GetValue("comptonmap.side2" , 2);
676 gEDBDEBUGLEVEL = cenv->GetValue("comptonmap.EdbDebugLevel", 2);
677
678 //cenv->Print();
679
680 EdbRunAccess ra;
681 ra.InitRun(runfile);
682
683 ra.eAFID = cenv->GetValue("fedra.link.AFID" , 1);
684
685 float min[5] = {-500,-500,-1.,-1., pulsmin};
686 float max[5] = { 500, 500, 1., 1., 50};
687 ra.AddSegmentCut(side1, 1, min, max);
688 ra.AddSegmentCut(side2, 1, min, max);
689
690 EdbSegP s(0, x, y, 0,0);
691
692 int n1= ra.GetPatternXY(s, side1, p1, size);
693 int n2= ra.GetPatternXY(s, side2, p2, size);
694
695 p1.SetID(1);
696 p2.SetID(2);
697
698 Log(2," EdbAlignmentMap::get_run_patterns","read segments: n1=%d at z=%f and n2=%d at z=%f",
699 n1,ra.GetLayer(1)->Z(),n2,ra.GetLayer(2)->Z());
700}
float Z() const
Definition: EdbLayer.h:77
void SetID(int id)
Definition: EdbPattern.h:309
helper class for access to the run data
Definition: EdbRunAccess.h:23
bool InitRun(const char *runfile=0, bool do_update=false)
Definition: EdbRunAccess.cxx:112
EdbLayer * GetLayer(int id)
Definition: EdbRunAccess.h:113
void AddSegmentCut(int xi, const char *cutline)
Definition: EdbRunAccess.cxx:1187
Int_t eAFID
if =1 - use affine transformations of the fiducial marks
Definition: EdbRunAccess.h:26
int GetPatternXY(EdbSegP &s, int side, EdbPattern &pat, float rmin=200)
Definition: EdbRunAccess.cxx:377

◆ GetPostreeAsPat()

void EdbAlignmentMap::GetPostreeAsPat ( EdbPattern pat,
const char *  filename 
)
static
725{
726 TFile f(filename);
727 TTree *tree = (TTree*)(f.Get("couples"));
728 int n = tree->GetEntries();
729 if(!n) return;
730 EdbSegP *s = new EdbSegP();
731 tree->SetBranchAddress("s.",&s);
732 for(int i=0; i<n; i++) {
733 tree->GetEntry(i);
734 pat.AddSegment(*s);
735 }
736 Log(2,"GetPostreeAsPat","%d segments read from file %s",n, filename);
737 f.Close();
738}
FILE * f
Definition: RecDispMC.C:150
const char filename[256]
Definition: RecDispNU.C:83
EdbSegP * AddSegment(int i, EdbSegP &s)
Definition: EdbPattern.cxx:72

◆ InitFile()

int EdbAlignmentMap::InitFile ( const char *  file = 0,
const char *  mode = 0 
)
48{
49 Log(2,"EdbAlignmentMap::InitFile","%s for %s",file,mode);
50 if(file) eOutputFile = new TFile(file,mode);
51
52 if( eOutputFile && (strstr("READ",mode) || strstr("UPDATE",mode)) ) {
53 eMapTree = (TTree*)eOutputFile->Get("maptree");
54 TVector3 *v = (TVector3*)eOutputFile->Get("GV1");
55 eGV1 = *v;
56 v = (TVector3*)eOutputFile->Get("GV2");
57 eGV2 = *v;
58 EdbCell2 *c = (EdbCell2*)eOutputFile->Get("Map");
59 eMap.Copy(*c);
61 }
62 else {
63 eMapTree = new TTree("maptree","maptree (izone:x0:y0:dx:dy:dz1:dz2:npeak:Hdxy)");
64 Int_t izone, peak, peakbin;
65 Float_t mean3,mean;
66 Float_t x0,y0,dx,dy,dz1,dz2;
67 TH2F *hdxy=0, *hdz12=0;
68 eMapTree->Branch("izone", &izone,"izone/I"); // zone id
69 eMapTree->Branch("x0", &x0,"x0/F"); // zone position
70 eMapTree->Branch("y0", &y0,"y0/F");
71 eMapTree->Branch("peak", &peak,"peak/I"); // total number of the found coinsidence
72 eMapTree->Branch("peakbin",&peakbin,"peakbin/F"); // peak bin
73 eMapTree->Branch("mean3", &mean3,"mean3/F"); // mean in 3x3 neigbour around the peak bin
74 eMapTree->Branch("mean", &mean,"mean/F"); // mean out of the peak zone
75 eMapTree->Branch("dx", &dx,"dx/F"); // found offsets
76 eMapTree->Branch("dy", &dy,"dy/F");
77 eMapTree->Branch("dz1", &dz1,"dz1/F");
78 eMapTree->Branch("dz2", &dz2,"dz2/F");
79 eMapTree->Branch("Hdxy", "TH2F", &hdxy, 128000,0);
80 eMapTree->Branch("Hdz12","TH2F", &hdz12,128000,0);
81 }
82 return 1;
83}
void Copy(const EdbCell2 &cell)
Definition: EdbCell2.cpp:628

◆ Link() [1/3]

int EdbAlignmentMap::Link ( )
449{
450 return Link( *eEnv, *ePat1, *ePat2, *ePlate);
451}
EdbPlateP * ePlate
Definition: EdbAlignmentMap.h:23
int Link()
Definition: EdbAlignmentMap.cxx:448

◆ Link() [2/3]

void EdbAlignmentMap::Link ( const char *  file,
EdbPlateP plate 
)
704{
705 Log(3," EdbAlignmentMap::Link","file %s",file);
706 if(gEDBDEBUGLEVEL>2) plate.PrintPlateLayout();
707
708 EdbPattern p1, p2;
709 get_run_patterns( file , eEnv, p1, p2);
710 p1.SetZ( plate.GetLayer(1)->Z() ); p1.SetSegmentsZ();
711 p2.SetZ( plate.GetLayer(2)->Z() ); p2.SetSegmentsZ();
712
713 Log(2," EdbAlignmentMap::Link","%d at %f and %d at %f",p1.N(),p1.Z(),p2.N(),p2.Z() );
714
715 eOutputFile->cd();
716 ePat1 = &p1;
717 ePat2 = &p2;
718 ePlate = &plate;
719 Link();
720 SaveAll();
721 }
void get_run_patterns(const char *runfile, TEnv *cenv, EdbPattern &p1, EdbPattern &p2)
Definition: EdbAlignmentMap.cxx:664
Float_t Z() const
Definition: EdbPattern.h:84
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1

◆ Link() [3/3]

int EdbAlignmentMap::Link ( TEnv &  cenv,
EdbPattern p1,
EdbPattern p2,
EdbPlateP plate 
)

-----------------------—
to be implemented:

  • select the subpatterns for shrinkage correction
  • do shrinkage correction
  • apply the results for subpatterns
  • check dz
  • apply the results for full patterns
  • link at the fixed dz
    -----------------------—
455{
465
466 eOutputFile->cd();
467
468 int npeak=0;
469
470 //CheckPattern(p1all, "1");
471 //CheckPattern(p2all, "2");
472
473 float DZ = plate.GetLayer(2)->Z() - plate.GetLayer(1)->Z();
474
475 bool doCorr = cenv.GetValue("fedra.link.DoCorr", false);
476 bool doShr = cenv.GetValue("fedra.link.DoShr", false);
477
478 //atest.FillArrays(p1sel, p2sel); // test linking with the selected segments
479 //atest.DoubletsFilterOut(true);
480 if(doCorr) {
481 float wmin = cenv.GetValue("fedra.link.shr.Wmin", 9. );
482 float wmax = cenv.GetValue("fedra.link.shr.Wmax", 100. );
483 float tmin = cenv.GetValue("fedra.link.shr.Tmin", 0.01 );
484 float tmax = cenv.GetValue("fedra.link.shr.Tmax", 0.5 );
485
486 EdbPattern p1sel,p2sel;
487 SelectSampleForShrinkageCorr(p1all,p1sel, wmin, wmax, tmin, tmax);
488 SelectSampleForShrinkageCorr(p2all,p2sel, wmin, wmax, tmin, tmax);
489
490 float dx1 = p1sel.Xmax() - p1sel.Xmin();
491 float dy1 = p1sel.Ymax() - p1sel.Ymin();
492 float cell1 = Sqrt( dx1*dy1/p1sel.N() );
493 float dx2 = p2sel.Xmax() - p2sel.Xmin();
494 float dy2 = p2sel.Ymax() - p2sel.Ymin();
495 float cell2 = Sqrt( dx2*dy2/p2sel.N() );
496 float cell = Max( (float)Min(cell1,cell2), (float)(cenv.GetValue("fedra.link.CellXY", 25.)) );
497
499 atest.eXcell = atest.eYcell = cell;
500 atest.eBinX = atest.eBinY = cenv.GetValue("fedra.link.BinXY", 4. );
501 atest.eDTXmax = atest.eDTYmax = cenv.GetValue("fedra.link.DTMax", 0.02 );
502
503 atest.FindCorrections( p1sel, p2sel, DZ, doShr);
504 TH2F *hshr1 = atest.eHShr1.DrawH2("hshr1");
505 TH2F *hshr2 = atest.eHShr2.DrawH2("hshr2");
506 hshr1->Write("hshr1"); SafeDelete(hshr1);
507 hshr2->Write("hshr2"); SafeDelete(hshr2);
508
509
510 atest.ePC1.ApplyCorrections(p1all);
511 atest.ePC2.ApplyCorrections(p2all);
512 atest.ePC1.ApplyCorrections( *(plate.GetLayer(1)) );
513 atest.ePC2.ApplyCorrections( *(plate.GetLayer(2)) );
514 plate.Write("plate");
515
516
517// Log(1,"Link","check corrections:");
518// EdbPositionAlignment atest2;
519// atest2.eXcell = atest2.eYcell = cenv.GetValue("fedra.link.CellXY", 10 );
520// atest2.eBinX = atest2.eBinY = cenv.GetValue("fedra.link.BinXY", 4 );
521// atest2.eDTXmax = atest2.eDTYmax = cenv.GetValue("fedra.link.DTMax", 0.02 );
522// atest.ePC1.ApplyCorrections(p1sel);
523// atest.ePC2.ApplyCorrections(p2sel);
524// atest2.FindCorrections(p1sel, p2sel, DZ, false);
525
526 }
527
528 //return 1;
529
530 // --------- end of test linking -----------
531
532
534 aall.eXcell = aall.eYcell = cenv.GetValue("fedra.link.CellXY", 10. );
535 aall.eBinX = aall.eBinY = cenv.GetValue("comptonmap.link.BinXY", 4. );
536 aall.eDTXmax = aall.eDTYmax = cenv.GetValue("comptonmap.link.DTMax", 0.02 );
537 aall.eChi2Max = cenv.GetValue("fedra.link.Chi2Max", 3. );
538
539 aall.ePC1.eDZ = DZ/2;
540 aall.ePC2.eDZ = -DZ/2;
541 aall.ePC1.eApplyCorr = aall.ePC2.eApplyCorr = true;
542
543 aall.FillArrays(p1all, p2all);
544
545 TH2F *hpat1 = aall.ePC1.DrawH2("hpat1");
546 TH2F *hpat2 = aall.ePC2.DrawH2("hpat2");
547 hpat1->Write("hpat1"); SafeDelete(hpat1);
548 hpat2->Write("hpat2"); SafeDelete(hpat2);
549 TH1I *spectr1 = aall.ePC1.DrawSpectrum("spectrum1");
550 spectr1->Write("spectrum1"); SafeDelete(spectr1);
551 TH1I *spectr2 = aall.ePC2.DrawSpectrum("spectrum2");
552 spectr2->Write("spectrum2"); SafeDelete(spectr2);
553
554 aall.DoubletsFilterOut(true);
555 aall.ePC1.eDoubletsRate = new TH1I("DoubletsRate","DoubletsRate",100,-0.5,99.5);
556
557 float dpos = cenv.GetValue("fedra.link.AccPos", 15. );
558 float dang = cenv.GetValue("fedra.link.AccAng", 0.1 );
559 aall.FillCombinations(aall.ePC1, aall.ePC2, dpos, dpos, dang, dang );
560 aall.RankCouples();
561 aall.ePC1.eDoubletsRate->Write(); SafeDelete(aall.ePC1.eDoubletsRate);
562
563 int nfill=0;
564 if(aall.ActivatePosTree("couples")) nfill = aall.FillPosTree(DZ/2, -DZ/2, 10);
565 FillMapTree(aall,10);
566
567 return npeak;
568}
float eDZ
corrections to be applied if eApplyCorr==true
Definition: EdbPatCell2.h:18
TH1I * eDoubletsRate
to be filled in FillCombinations()
Definition: EdbPatCell2.h:23
bool eApplyCorr
Definition: EdbPatCell2.h:17
void ApplyCorrections(EdbPattern &pat)
Definition: EdbPatCell2.cxx:61
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
EdbH2 eHShr2
histogram used for the shrinkage-selection
Definition: EdbPositionAlignment.h:62
EdbH2 eHShr1
histogram used for the shrinkage-selection
Definition: EdbPositionAlignment.h:61
void FindCorrections(EdbPattern &pat1, EdbPattern &pat2, float DZ, bool doShrink)
Definition: EdbPositionAlignment.cxx:1008
Float_t eChi2Max
cut for the chi2 of the basetrack
Definition: EdbPositionAlignment.h:54
int FillPosTree(float dz1, float dz2, int flag)
Definition: EdbPositionAlignment.cxx:189
void RankCouples()
Definition: EdbPositionAlignment.cxx:287

◆ SaveAll()

void EdbAlignmentMap::SaveAll ( )
87{
88 if(!eOutputFile) return;
89 if(!eOutputFile->IsWritable()) return;
90 if(eMapTree) eMapTree->Write();
91 if(eEnv) eEnv->Write();
92 eGV1.Write("GV1");
93 eGV2.Write("GV2");
94 eMap.Write("Map");
95 eOutputFile->Save();
96 eOutputFile->Purge();
97}

◆ SaveMap()

void EdbAlignmentMap::SaveMap ( const char *  file)
291{
292 Log(2,"EdbAlignmentMap::SaveMap","to file %s",file);
293 if(!eMapTree) return;
294 FILE *f = fopen(file,"w");
295 if(!f) return;
296 int N,IZONE;
297 float X, Y, DX, DY, DZ1, DZ2;
298 eMapTree->ResetBranchAddresses();
299 eMapTree->SetBranchAddress("izone",&IZONE);
300 eMapTree->SetBranchAddress("peak" ,&N);
301 eMapTree->SetBranchAddress("x0" ,&X);
302 eMapTree->SetBranchAddress("y0" ,&Y);
303 eMapTree->SetBranchAddress("dx" ,&DX);
304 eMapTree->SetBranchAddress("dy" ,&DY);
305 eMapTree->SetBranchAddress("dz1" ,&DZ1);
306 eMapTree->SetBranchAddress("dz2" ,&DZ2);
307 fprintf(f,"%d %d %f %f %f %f\n",eNx,eNy,eMap.Xbin(),eMap.Ybin(), eGV1.X(), eGV1.Y() );
308 int n = eMapTree->GetEntries();
309 for(int i=0; i<n; i++) {
310 eMapTree->GetEntry(i);
311 fprintf(f,"%5d %13.1f %13.1f %7.2f %7.2f %7.2f %7.2f %5d\n",
312 IZONE,X,Y,DX+eGV1.X(),DY+eGV1.Y(),DZ1,DZ2,N);
313 }
314 fclose(f);
315}
fclose(pFile)

◆ SelectSampleForShrinkageCorr()

int EdbAlignmentMap::SelectSampleForShrinkageCorr ( EdbPattern p,
EdbPattern psel,
float  wmin,
float  wmax,
float  tmin,
float  tmax 
)

to do: analyse patterns before and select the best sample

619{
621 int nsel=0;
622
623 int n = p.N();
624 EdbSegP *s;
625 for(int i=0; i<n; i++) {
626 s = p.GetSegment(i);
627 if( s->W()<wmin ) continue;
628 if( s->Theta() < tmin ) continue;
629 if( s->Theta() > tmax ) continue;
630 if( s->Flag() == -10 ) continue;
631 if( s->W() > wmax ) continue;
632 plnk.AddSegment(*s);
633 }
634 nsel = plnk.N();
635 Log(2,"SelectSampleForShrinkageCorr","%d of %d segments selected",nsel,n);
636 return nsel;
637}

Member Data Documentation

◆ eDensityMax

float EdbAlignmentMap::eDensityMax

the max segments density for patterns selection (in N/100/100 microns)

◆ eEnv

TEnv* EdbAlignmentMap::eEnv

environment used to pass the parameters

◆ eGlobal

EdbPositionAlignment EdbAlignmentMap::eGlobal

service object for the zones selection

◆ eGraphDX

TGraph2D* EdbAlignmentMap::eGraphDX

◆ eGraphDY

TGraph2D* EdbAlignmentMap::eGraphDY

◆ eGraphDZ1

TGraph2D* EdbAlignmentMap::eGraphDZ1

◆ eGraphDZ2

TGraph2D* EdbAlignmentMap::eGraphDZ2

◆ eGV1

TVector3 EdbAlignmentMap::eGV1

◆ eGV2

TVector3 EdbAlignmentMap::eGV2

global offsets found by the test alignment (dx:dy:dz)

◆ eMap

EdbCell2 EdbAlignmentMap::eMap

in each cell TArratF with (x,y, dx,dy,dz1,dz2, n)

◆ eMapTree

TTree* EdbAlignmentMap::eMapTree

◆ eNx

int EdbAlignmentMap::eNx

◆ eNy

int EdbAlignmentMap::eNy

number of divisions calculated using eXcell, eYcell

◆ eOutputFile

TFile* EdbAlignmentMap::eOutputFile

◆ ePat1

EdbPattern* EdbAlignmentMap::ePat1

big patterns to be splitted and aligned piece by piece

◆ ePat2

EdbPattern* EdbAlignmentMap::ePat2

◆ ePlate

EdbPlateP* EdbAlignmentMap::ePlate

◆ eXcell

float EdbAlignmentMap::eXcell

◆ eYcell

float EdbAlignmentMap::eYcell

approximate zones size (for example 10000 microns)


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