FEDRA emulsion software from the OPERA Collaboration
EdbTrackAssembler Class Reference

generic class for the tracks assembling from segments More...

#include <EdbScanTracking.h>

Inheritance diagram for EdbTrackAssembler:
Collaboration diagram for EdbTrackAssembler:

Public Member Functions

bool AcceptDZGap (EdbTrackP &t, float z)
 
void AddPattern (EdbPattern &p)
 
EdbTrackPAddSegment (EdbSegP &s)
 owner of the segments!!! More...
 
EdbTrackPAddSegmentAsTrack (EdbSegP &s)
 
void CheckPatternAlignment (EdbPattern &p, EdbPlateP &plate, int nsegmin)
 
void CombTracks (TObjArray &selected)
 
void DoubletsFilterOut (EdbPattern &p)
 
 EdbTrackAssembler ()
 
void ExtrapolateTracksToZ (float z, int nsegmin=0)
 
void FillTrZMap ()
 
void FitTracks ()
 
void InitTrZMap ()
 
void InitTrZMap (const char *str)
 
void InitTrZMap (int nx, float xmi, float xma, int ny, float ymi, float yma, int ncell)
 
float ProbSeg (EdbSegP &s1, EdbSegP &s2)
 
void RecalculateSegmentsProb (EdbTrackP &t)
 
bool SameSegment (EdbSegP &s1, EdbSegP &s2)
 
void SetMomentum (float p)
 
void SetRadLength (float x0)
 
void SetSegmentsErrors ()
 
TObjArray & Tracks ()
 
virtual ~EdbTrackAssembler ()
 

Public Attributes

int eCollisionsRate
 
EdbScanCond eCond
 
bool eDoUseMCS
 flag to use MultipleScattering addition for chi2 More...
 
float eDRmax
 position acceptance for the fast preselection More...
 
float eDTmax
 angular acceptance for the fast preselection More...
 
float eDZGapMax
 maxgap acceptance for the fast preselection More...
 
TH1F * eHistNcnd
 number of candidates after preliminary selection More...
 
TH1F * eHistProbAll
 prob of all candidate More...
 
TH1F * eHistProbBest
 prob of the best candidate More...
 
TH1F * eHistThetaAll
 theta of all candidate More...
 
TH1F * eHistThetaBest
 theta of the best candidate More...
 
float eProbMin
 min acceptable probability for segments preselection More...
 

Private Attributes

Int_t eCellN
 mean cell occupancy More...
 
EdbTrackFitter eFitter
 
Float_t eMapMarg
 margin for the map creation More...
 
EdbPattern eSegments
 all segments of tracks More...
 
TObjArray eTracks
 array of tracks (EdbTrackP) (owner of tracks) More...
 
TObjArray eTrZ
 "predictions" - tracks extrapolated to the given z (not owner) More...
 
EdbCell2 eTrZMap
 map of predictions at given eZ More...
 
Float_t eZ
 the z-position More...
 

Detailed Description

generic class for the tracks assembling from segments

//////////////////////////////////////////////////////////////////////// // EdbScanTracking // // To handle tracking in the scanset (IO) // // ////////////////////////////////////////////////////////////////////////

Constructor & Destructor Documentation

◆ EdbTrackAssembler()

EdbTrackAssembler::EdbTrackAssembler ( )
39{
40 eMapMarg = 50.; // [microns]
41 eZ = 0;
42 eCellN=10; //mean n/cell
43 eDTmax=0.07;
44 eDRmax=45.;
45 eDZGapMax = 5000;
46 eProbMin = 0.001;
48
49 eHistNcnd = new TH1F("Ncnd","number of candidates after preliminary selection", 20,0.5,20.5);
50 eHistProbBest = new TH1F("ProbBest","prob for best selected candidate", 250,0,1);
51 eHistProbAll = new TH1F("ProbAll","prob for all candidates", 250,0,1);
52 eHistThetaBest = new TH1F("ThetaBest","angle theta for best selected candidate", 180,0,TMath::PiOver2());
53 eHistThetaAll = new TH1F("ThetaAll","angle theta for all candidates", 180,0,TMath::PiOver2());
54
55 // for basetracks:
57 eCond.SetSigma0( 4, 4, 0.005, 0.005 );
58 eCond.SetPulsRamp0(14., 21.);
59 eCond.SetPulsRamp04(14., 21.);
60}
void SetDefault()
Definition: EdbScanCond.cxx:16
void SetPulsRamp0(float p1, float p2)
Definition: EdbScanCond.h:74
void SetSigma0(float x, float y, float tx, float ty)
Definition: EdbScanCond.h:62
void SetPulsRamp04(float p1, float p2)
Definition: EdbScanCond.h:75
TH1F * eHistNcnd
number of candidates after preliminary selection
Definition: EdbScanTracking.h:46
EdbScanCond eCond
Definition: EdbScanTracking.h:40
TH1F * eHistProbBest
prob of the best candidate
Definition: EdbScanTracking.h:42
Float_t eZ
the z-position
Definition: EdbScanTracking.h:24
Int_t eCellN
mean cell occupancy
Definition: EdbScanTracking.h:28
Float_t eMapMarg
margin for the map creation
Definition: EdbScanTracking.h:27
TH1F * eHistThetaAll
theta of all candidate
Definition: EdbScanTracking.h:45
float eDZGapMax
maxgap acceptance for the fast preselection
Definition: EdbScanTracking.h:35
TH1F * eHistThetaBest
theta of the best candidate
Definition: EdbScanTracking.h:44
TH1F * eHistProbAll
prob of all candidate
Definition: EdbScanTracking.h:43
float eProbMin
min acceptable probability for segments preselection
Definition: EdbScanTracking.h:36
float eDRmax
position acceptance for the fast preselection
Definition: EdbScanTracking.h:34
float eDTmax
angular acceptance for the fast preselection
Definition: EdbScanTracking.h:33
int eCollisionsRate
Definition: EdbScanTracking.h:39

◆ ~EdbTrackAssembler()

virtual EdbTrackAssembler::~EdbTrackAssembler ( )
virtual

Member Function Documentation

◆ AcceptDZGap()

bool EdbTrackAssembler::AcceptDZGap ( EdbTrackP t,
float  z 
)
282{
283 float z1 = t.GetSegmentFirst()->Z();
284 float z2 = t.GetSegmentLast()->Z();
285 if(Min( Abs(z1-z), Abs(z2-z)) > eDZGapMax ) return false;
286 return true;
287}
TTree * t
Definition: check_shower.C:4

◆ AddPattern()

void EdbTrackAssembler::AddPattern ( EdbPattern p)
106{
107 //int ntrBefore=eTracks.GetEntries();
108 //if(ntrBefore>0) ExtrapolateTracksToZ(p.Z());
109
110 //DoubletsFilterOut();
111 int nseg = p.N();
112 Log(3,"EdbTrackAssembler::AddPattern","try to add %d segments",p.N());
113 int attached=0;
114 for(int j=0; j<nseg; j++) {
115 EdbSegP *s = p.GetSegment(j);
116 if(s->Flag()==-10) continue;
117 s->SetErrors();
118 eCond.FillErrorsCov(s->TX(),s->TY(),s->COV());
119 if( !AddSegment( *(p.GetSegment(j)) ) )
120 AddSegmentAsTrack( *(p.GetSegment(j)) );
121 else attached++;
122 }
123 int ntrAfter = eTracks.GetEntries();
124 //int totSegTr=0;
125 //for(int i=0; i<ntrAfter; i++) totSegTr += ((EdbTrackP*)(eTracks.At(i)))->N();
126 Log(2,"EdbTrackAssembler::AddPattern","with z=%10.2f %d/%d attached/tried; total collisions: %d; tracks: %d",
127 p.Z(), attached, nseg, eCollisionsRate, ntrAfter );
128}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
void FillErrorsCov(float tx, float ty, TMatrixD &cov)
Definition: EdbScanCond.cxx:161
Definition: EdbSegP.h:21
TObjArray eTracks
array of tracks (EdbTrackP) (owner of tracks)
Definition: EdbScanTracking.h:21
EdbTrackP * AddSegmentAsTrack(EdbSegP &s)
Definition: EdbScanTracking.cxx:246
EdbTrackP * AddSegment(EdbSegP &s)
owner of the segments!!!
Definition: EdbScanTracking.cxx:131
s
Definition: check_shower.C:55
p
Definition: testBGReduction_AllMethods.C:8

◆ AddSegment()

EdbTrackP * EdbTrackAssembler::AddSegment ( EdbSegP s)

owner of the segments!!!

132{
133 TObjArray trsel;
134 float v[2] = { s.X(), s.Y() };
135 int nsel = eTrZMap.SelectObjectsC( v, eDRmax+50 , trsel );
136 Log(3,"EdbTrackAssembler::AddSegment", "nsel = %d",nsel);
137 if(!nsel) {
138 return 0; }
139 float prob, probmax = eProbMin;
140 EdbSegP *ssbest = 0;
141 int ncnd = 0;
142 for(int i=0; i<nsel; i++) {
143 EdbSegP *ss = (EdbSegP*)(trsel.At(i));
144 prob = ProbSeg( *ss, s );
145 Log(3,"EdbTrackAssembler::AddSegment", "prob(probmin) = %f (%f) ",prob, eProbMin);
146 if(prob<eProbMin) continue;
147 ncnd++;
148 if(eHistProbAll) eHistProbAll->Fill(prob);
149 if(eHistThetaAll) eHistThetaAll->Fill(ATan(ss->Theta()));
150 if( prob > probmax ) { ssbest = ss; probmax=prob; }
151 }
152 if(!ssbest) return 0;
153 s.SetProb(probmax);
154 if(eHistNcnd) eHistNcnd->Fill(ncnd);
155 if(eHistProbBest) eHistProbBest->Fill(probmax);
156 if(eHistThetaBest) eHistThetaBest->Fill(ATan(ssbest->Theta()));
157 EdbTrackP *t = (EdbTrackP*)(ssbest);
158 EdbSegP *sz = t->GetSegmentWithClosestZ( t->Z(), 45. );
159 if(!sz) t->AddSegment( eSegments.AddSegment(s) );
160 else {
161 if( !SameSegment(s,*sz) ) {
162 if( s.Prob() > sz->Prob() ) t->SubstituteSegment( sz , eSegments.AddSegment(s) );
164 }
165 }
166 return t;
167}
int SelectObjectsC(int iv[2], int ir[2], TObjArray &arr)
Definition: EdbCell2.cpp:707
Float_t Prob() const
Definition: EdbSegP.h:156
Float_t Theta() const
Definition: EdbSegP.h:184
EdbSegP * AddSegment(int i, EdbSegP &s)
Definition: EdbPattern.cxx:72
EdbCell2 eTrZMap
map of predictions at given eZ
Definition: EdbScanTracking.h:26
EdbPattern eSegments
all segments of tracks
Definition: EdbScanTracking.h:20
bool SameSegment(EdbSegP &s1, EdbSegP &s2)
Definition: EdbScanTracking.cxx:170
float ProbSeg(EdbSegP &s1, EdbSegP &s2)
Definition: EdbScanTracking.cxx:190
Definition: EdbPattern.h:113
ss
Definition: energy.C:62

◆ AddSegmentAsTrack()

EdbTrackP * EdbTrackAssembler::AddSegmentAsTrack ( EdbSegP s)
247{
248 if(s.W()<16 ) return 0;
249 if(s.Chi2()>2.5) return 0;
250 EdbTrackP *t = new EdbTrackP( eSegments.AddSegment(s), 0.139); // EdbTrackAssembler is owner of segments
251 eTracks.Add(t);
252 //EdbSegP ss;
253 //t->MakePredictionTo(eZ,ss);
254 //eTrZ.AddSegment(ss);
255 return t;
256}

◆ CheckPatternAlignment()

void EdbTrackAssembler::CheckPatternAlignment ( EdbPattern p,
EdbPlateP plate,
int  nsegmin 
)
75{
76 Log(0,"EdbTrackAssembler::CheckPatternAlignment","");
77
79 int ntr = eTrZ.GetEntries();
80 EdbPattern ptr( 0, 0, p.Z(), ntr );
81 for(int i=0; i<ntr; i++) ptr.AddSegment( *((EdbSegP*)(eTrZ.UncheckedAt(i))) );
82
84 al.SetSigma( 25, 0.015 );
85 al.eOffsetMax = 500.;
86 al.eDZ = 0;
87 al.eDPHI = 0.00;
88 //al.eDoCoarse=1;
89 al.Align(ptr,p,0);
90
91 EdbAffine2D *aff = al.eCorrL[0].GetAffineXY();
92 aff->Invert();
93 aff->Print();
94 p.Transform(aff);
95 plate.GetAffineXY()->Transform(aff);
96
97 EdbAffine2D *afftxty = al.eCorrL[0].GetAffineTXTY();
98 afftxty->Invert();
99 afftxty->Print();
100 p.TransformA(afftxty);
101 plate.GetAffineTXTY()->Transform(afftxty);
102}
Definition: EdbAffine.h:17
void Invert()
Definition: EdbAffine.cxx:103
void Print(Option_t *opt="") const
Definition: EdbAffine.cxx:52
EdbLayer eCorrL[2]
corrections in form of affine transformations - the final output
Definition: EdbAlignmentV.h:25
EdbAffine2D * GetAffineTXTY()
Definition: EdbLayer.h:120
EdbAffine2D * GetAffineXY()
Definition: EdbLayer.h:119
Definition: EdbPattern.h:273
plate-to-plate alignment
Definition: EdbPlateAlignment.h:8
Float_t eOffsetMax
the maximal offset to be looked for
Definition: EdbPlateAlignment.h:12
void SetSigma(float spos, float sang)
Definition: EdbPlateAlignment.h:53
void Align(EdbPattern &p1, EdbPattern &p2, float dz)
Definition: EdbPlateAlignment.cxx:61
Float_t eDZ
the range +- dz will be scanned by coarce align
Definition: EdbPlateAlignment.h:14
Float_t eDPHI
the range +- dphi will be scanned by coarce align
Definition: EdbPlateAlignment.h:15
void ExtrapolateTracksToZ(float z, int nsegmin=0)
Definition: EdbScanTracking.cxx:259
TObjArray eTrZ
"predictions" - tracks extrapolated to the given z (not owner)
Definition: EdbScanTracking.h:22
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1
int nsegmin
Definition: check_vertex.C:23

◆ CombTracks()

void EdbTrackAssembler::CombTracks ( TObjArray &  selected)

eliminate crossing&overlapping tracks with multiple segments usage

371{
373
374 int nsegMin=2;
375 int nGapMax=50;
376
377 int ntr = eTracks.GetEntries();
378 Log(3,"EdbTrackAssembler::CombTracks","Comb %d tracks");
379
380 // *** sort tracks by quality
381
382 TIndexCell cn; //"nseg:prob:entry"
383 Long_t v[3];
384
385 int nsegtot=0;
386 EdbTrackP *tr=0;
387 for(int i=0; i<ntr; i++) {
388 tr = (EdbTrackP*)(eTracks.At(i));
389 if( tr->Flag() == -10 ) continue;
390 if( tr->N() < nsegMin ) continue;
391 tr->SetID(i);
392 tr->SetCounters();
393 nsegtot += tr->SetSegmentsTrack(-1);
394 v[0]= -(tr->N());
395 v[1]= (Long_t)((1.-tr->Prob())*100);
396 v[2]= i;
397 cn.Add(3,v);
398 }
399 cn.Sort();
400
401 Log(3,"EdbTrackAssembler::CombTracks","%d tracks with %d segments for processing...",ntr,nsegtot);
402
403 // *** set track ID for segments attached to
404
405 TIndexCell *cp=0, *c=0;
406 int nn=cn.GetEntries();
407 for(int i=nn-1; i>=0; i--) {
408 cp = cn.At(i); // tracks with fixed npl
409 int np = cp->GetEntries();
410 for(int ip=np-1; ip>=0; ip--) {
411 c = cp->At(ip); // tracks with fixed Npl & Prob
412 int nt = c->GetEntries();
413 for(int it=0; it<nt; it++) {
414 tr = (EdbTrackP*)(eTracks.At( c->At(it)->Value() ) );
415 tr->SetSegmentsTrack();
416 }
417 }
418 }
419
420
421 cp=0; c=0;
422 nn=cn.GetEntries();
423 for(int i=0; i<nn; i++) {
424 cp = cn.At(i); // tracks with fixed npl
425
426 int np = cp->GetEntries();
427 for(int ip=0; ip<np; ip++) {
428 c = cp->At(ip); // tracks with fixed Npl & Prob
429
430 int nt = c->GetEntries();
431 for(int it=0; it<nt; it++) {
432
433 tr = (EdbTrackP*)(eTracks.At( c->At(it)->Value() ) );
434
435 if(tr->RemoveAliasSegments()>0){
436 if(tr->N()<nsegMin) tr->SetFlag(-10);
437 if(tr->CheckMaxGap()>nGapMax) tr->SetFlag(-10);
438 }
439
440 if( tr->Flag() != -10 ) selected.Add(tr);
441 }
442 }
443 }
444
445}
int nsegMin
Definition: RecDispEX.C:19
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
sort collection with attributes
Definition: TIndexCell.h:19
Int_t GetEntries() const
Definition: TIndexCell.h:82
TIndexCell const * At(Int_t narg, Int_t vind[]) const
Definition: TIndexCell.cpp:519
void Sort(Int_t upto=kMaxInt)
Definition: TIndexCell.cpp:539
Int_t Add(Int_t narg, Long_t varg[])
Definition: TIndexCell.cpp:602
EdbSegCouple * cp
Definition: tlg2couples.C:28

◆ DoubletsFilterOut()

void EdbTrackAssembler::DoubletsFilterOut ( EdbPattern p)
64{
65 EdbAlignmentV adup;
66 adup.eDVsame[0]=adup.eDVsame[1]=10.;
67 adup.eDVsame[2]=adup.eDVsame[3]=0.08;
68 adup.FillGuessCell(p,p,1.);
69 adup.FillCombinations();
70 adup.DoubletsFilterOut(0); // assign flag -10 to the duplicated segments
71}
universal basic alignment class
Definition: EdbAlignmentV.h:13
Float_t eDVsame[4]
(dx,dy,dtx,dty) condition for the coinsidence
Definition: EdbAlignmentV.h:16
int DoubletsFilterOut(int checkview, TH2F *hxy=0, TH2F *htxty=0)
Definition: EdbAlignmentV.cxx:67
void FillGuessCell(EdbPattern &p1, EdbPattern &p2, float binOK=1., float offsetMax=2000.)
Definition: EdbAlignmentV.cxx:880
int FillCombinations()
Definition: EdbAlignmentV.cxx:208

◆ ExtrapolateTracksToZ()

void EdbTrackAssembler::ExtrapolateTracksToZ ( float  z,
int  nsegmin = 0 
)
260{
261 eTrZ.Clear();
263
264 eZ=z;
265 int n=eTracks.GetEntries();
266 for(int i=0; i<n; i++) {
267 EdbTrackP *t = (EdbTrackP*)(eTracks.At(i));
268
269 if( t->N() < nsegmin ) continue;
270 if( !AcceptDZGap(*t, z) ) continue;
271 t->MakePredictionTo(eZ,*t); // extrapolation of tracks itself
272 //((EdbSegP *)t)->PrintNice();
273 eTrZ.Add(t);
274 }
275
276 //FillTrZMap();
277 //if(gEDBDEBUGLEVEL>2) eTrZMap.PrintStat();
278}
void CleanCells()
Definition: EdbCell2.cpp:114
bool AcceptDZGap(EdbTrackP &t, float z)
Definition: EdbScanTracking.cxx:281

◆ FillTrZMap()

void EdbTrackAssembler::FillTrZMap ( )
291{
292 int n=eTrZ.GetEntries();
293 Log(2,"EdbTrackAssembler::FillTrZMap", "with %d tracks",n);
294 for(int i=0; i<n; i++) {
295 EdbSegP *s = (EdbSegP*)(eTrZ.At(i));
296 eTrZMap.AddObject( s->X(), s->Y(), s );
297 }
298}
bool AddObject(float v[2], TObject *obj)
Definition: EdbCell2.h:173

◆ FitTracks()

void EdbTrackAssembler::FitTracks ( )
355{
356 EdbTrackFitter fit;
357
358 int ntr = eTracks.GetEntries();
359 for( int i=0; i<ntr; i++ ) {
360 EdbTrackP *t = (EdbTrackP*)(eTracks.At(i));
361 if(t->Flag()==-10) continue;
362 int nseg=t->N();
363 t->FitTrackKFS(0,10000);
364 //fit.FitTrackLine(*t);
365 if(nseg>1) RecalculateSegmentsProb(*t);
366 }
367}
void RecalculateSegmentsProb(EdbTrackP &t)
Definition: EdbScanTracking.cxx:181
Definition: EdbTrackFitter.h:17

◆ InitTrZMap() [1/3]

void EdbTrackAssembler::InitTrZMap ( )
321{
322 /* float mi[2] = { eTrZ.Xmin()-eMapMarg, eTrZ.Ymin()-eMapMarg };
323 float ma[2] = { eTrZ.Xmax()+eMapMarg, eTrZ.Ymax()+eMapMarg };
324 float dens = eTrZ.N()/( (ma[0]-mi[0])*(ma[1]-mi[1]));
325 float step=10000;
326 if(dens>0.00000001) step = Sqrt( eCellN/dens );
327 int n[2] = { int((ma[0]-mi[0])/step)+1, int((ma[1]-mi[1])/step)+1 };
328 float stepX = (ma[0]-mi[0])/n[0];
329 float stepY = (ma[1]-mi[1])/n[1];
330 n[0] = int((ma[0]-mi[0]+1.)/stepX);
331 n[1] = int((ma[1]-mi[1]+1.)/stepY);
332 eTrZMap.InitCell(3*eCellN, n, mi, ma);*/
333}

◆ InitTrZMap() [2/3]

void EdbTrackAssembler::InitTrZMap ( const char *  str)
302{
303 int nx=0, ny=0, ncell=0;
304 float xmi,xma, ymi, yma;
305 sscanf(str,"%d %f %f %d %f %f %d",&nx,&xmi,&xma,&ny,&ymi,&yma,&ncell);
306 InitTrZMap( nx,xmi,xma,ny,ymi,yma,ncell );
307}
void InitTrZMap()
Definition: EdbScanTracking.cxx:320

◆ InitTrZMap() [3/3]

void EdbTrackAssembler::InitTrZMap ( int  nx,
float  xmi,
float  xma,
int  ny,
float  ymi,
float  yma,
int  ncell 
)
312{
313 float mi[2] = { xmi, ymi };
314 float ma[2] = { xma, yma };
315 int n[2] = { nx, ny };
316 eTrZMap.InitCell(ncell, n, mi, ma);
317}
int InitCell(EdbCell2 &c)
Definition: EdbCell2.h:167

◆ ProbSeg()

float EdbTrackAssembler::ProbSeg ( EdbSegP s1,
EdbSegP s2 
)

return the probability that the second segment can belong to track defined by s1

workaround to get previous(not propagated) segment

191{
193 float dtx = s1.TX() - s2.TX();
194 if( Abs( dtx ) > eDTmax ) return 0;
195 float dty = s1.TY() - s2.TY();
196 if( Abs( dty ) > eDTmax ) return 0;
197 double dt2 = dtx*dtx + dty*dty;
198 if(dt2>eDTmax*eDTmax) return 0;
199
200 float dz = s2.Z()-s1.Z();
201 float dx = s2.X() - (s1.X() + dz*s1.TX());
202 if( Abs( dx ) > eDRmax ) return 0;
203 float dy = s2.Y() - (s1.Y() + dz*s1.TY());
204 if( Abs( dy ) > eDRmax ) return 0;
205 double dr2 = dx*dx + dy*dy;
206 if(dr2>eDRmax*eDRmax) return 0;
207
208 float prob=0;
209 if(eDoUseMCS==1){
211 EdbTrackP* t = dynamic_cast<EdbTrackP*> (&s1);
212 EdbSegP* seg = t?(const_cast<EdbSegP*>(t->TrackEnd())):0;
213 prob = seg?(eFitter.ProbSegMCS(seg, &s2)):0;
214 }
215 else if(eDoUseMCS==2)
216 {
217 EdbTrackP* t = dynamic_cast<EdbTrackP*> (&s1);
218 EdbSegP* seg = t?(const_cast<EdbSegP*>(t->TrackEnd())):0;
219 if(seg) {
220 float chi = eFitter.Chi2Seg( seg, &s2 );
221 prob = (float)TMath::Prob( chi*chi, 4);
222 }
223 if(prob<0.001&&gEDBDEBUGLEVEL>2) {
224 Log(3,"EdbTrackAssembler::ProbSeg", "dx,dy,dz,dtx,dty: %f %f %f %f %f prob= %f",dx,dy,dz,dtx,dty, prob);
225 seg->Print();
226 s2.Print();
227 eCond.Print();
228 }
229
230 }
231 else
232 {
233 EdbSegP s;
234 float chi = eFitter.Chi2SegM( s1, s2, s, eCond, eCond );
235 prob = (float)TMath::Prob( chi*chi, 4);
236 }
237
238
239 prob *= eCond.ProbSeg( s2.Theta(), s2.W() ); // the proability component depending on the grains number
240 prob *= (float)TMath::Prob( s2.Chi2()*s2.Chi2(), 4 ); // the proability component depending on the segment strength
241
242 return prob;
243}
T Prob(const T &rhs, int n)
Definition: Prob.hh:37
brick dz
Definition: RecDispMC.C:107
void Print() const
Definition: EdbScanCond.cxx:50
float ProbSeg(float tx, float ty, float puls) const
Definition: EdbScanCond.cxx:119
void Print(Option_t *opt="") const
Definition: EdbSegP.cxx:379
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t X() const
Definition: EdbSegP.h:173
Float_t Chi2() const
Definition: EdbSegP.h:157
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174
Float_t W() const
Definition: EdbSegP.h:151
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
bool eDoUseMCS
flag to use MultipleScattering addition for chi2
Definition: EdbScanTracking.h:37
EdbTrackFitter eFitter
Definition: EdbScanTracking.h:30
double ProbSegMCS(EdbSegP *s1, EdbSegP *s2)
Definition: EdbTrackFitter.cxx:288
float Chi2SegM(EdbSegP s1, EdbSegP s2, EdbSegP &s, EdbScanCond &cond1, EdbScanCond &cond2)
Definition: EdbTrackFitter.cxx:315
static float Chi2Seg(EdbSegP *s1, EdbSegP *s2)
Definition: EdbTrackFitter.cxx:62
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ RecalculateSegmentsProb()

void EdbTrackAssembler::RecalculateSegmentsProb ( EdbTrackP t)

assumed that track is fitted: reset the segmets probabilities

182{
184 int n=tr.N();
185 for(int i=0; i<n; i++)
186 tr.GetSegment(i)->SetProb( ProbSeg( *(tr.GetSegmentF(i)), *(tr.GetSegment(i)) ) );
187}

◆ SameSegment()

bool EdbTrackAssembler::SameSegment ( EdbSegP s1,
EdbSegP s2 
)
171{
172 if( Abs( s1.X() - s2.X() ) <0.000001 &&
173 Abs( s1.Y() - s2.Y() ) <0.000001 &&
174 Abs( s1.TX()- s2.TX() ) <0.000001 &&
175 Abs( s1.TY()- s2.TY() ) <0.000001 &&
176 Abs( s1.W() - s2.W() ) <0.000001 ) return true;
177 return false;
178}

◆ SetMomentum()

void EdbTrackAssembler::SetMomentum ( float  p)
inline
float ePdef
default momentum
Definition: EdbTrackFitter.h:25

◆ SetRadLength()

void EdbTrackAssembler::SetRadLength ( float  x0)
inline
53{eFitter.eX0=x0; eCond.SetRadX0(x0);}
void SetRadX0(float x0)
Definition: EdbScanCond.h:57
float eX0
rad length of the media [microns]
Definition: EdbTrackFitter.h:23

◆ SetSegmentsErrors()

void EdbTrackAssembler::SetSegmentsErrors ( )
337{
338 int ntr = eTracks.GetEntries();
339 for( int i=0; i<ntr; i++ ) {
340 EdbTrackP *t = (EdbTrackP*)(eTracks.At(i));
341 if(t->Flag()==-10) continue;
342 int nseg=t->N();
343 if(nseg>0) {
344 for(int j=0; j<nseg; j++) {
345 EdbSegP *s = t->GetSegment(j);
346 s->SetErrors();
347 eCond.FillErrorsCov(s->TX(),s->TY(),s->COV());
348 }
349 }
350 }
351}

◆ Tracks()

TObjArray & EdbTrackAssembler::Tracks ( )
inline
75{return eTracks;}

Member Data Documentation

◆ eCellN

Int_t EdbTrackAssembler::eCellN
private

mean cell occupancy

◆ eCollisionsRate

int EdbTrackAssembler::eCollisionsRate

◆ eCond

EdbScanCond EdbTrackAssembler::eCond

◆ eDoUseMCS

bool EdbTrackAssembler::eDoUseMCS

flag to use MultipleScattering addition for chi2

◆ eDRmax

float EdbTrackAssembler::eDRmax

position acceptance for the fast preselection

◆ eDTmax

float EdbTrackAssembler::eDTmax

angular acceptance for the fast preselection

◆ eDZGapMax

float EdbTrackAssembler::eDZGapMax

maxgap acceptance for the fast preselection

◆ eFitter

EdbTrackFitter EdbTrackAssembler::eFitter
private

◆ eHistNcnd

TH1F* EdbTrackAssembler::eHistNcnd

number of candidates after preliminary selection

◆ eHistProbAll

TH1F* EdbTrackAssembler::eHistProbAll

prob of all candidate

◆ eHistProbBest

TH1F* EdbTrackAssembler::eHistProbBest

prob of the best candidate

◆ eHistThetaAll

TH1F* EdbTrackAssembler::eHistThetaAll

theta of all candidate

◆ eHistThetaBest

TH1F* EdbTrackAssembler::eHistThetaBest

theta of the best candidate

◆ eMapMarg

Float_t EdbTrackAssembler::eMapMarg
private

margin for the map creation

◆ eProbMin

float EdbTrackAssembler::eProbMin

min acceptable probability for segments preselection

◆ eSegments

EdbPattern EdbTrackAssembler::eSegments
private

all segments of tracks

◆ eTracks

TObjArray EdbTrackAssembler::eTracks
private

array of tracks (EdbTrackP) (owner of tracks)

◆ eTrZ

TObjArray EdbTrackAssembler::eTrZ
private

"predictions" - tracks extrapolated to the given z (not owner)

◆ eTrZMap

EdbCell2 EdbTrackAssembler::eTrZMap
private

map of predictions at given eZ

◆ eZ

Float_t EdbTrackAssembler::eZ
private

the z-position


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