1 |
/** |
2 |
* \file TrkLevel2.cpp |
3 |
* \author Elena Vannuccini |
4 |
*/ |
5 |
#include <TrkLevel2.h> |
6 |
#include <iostream> |
7 |
#include <math.h> |
8 |
using namespace std; |
9 |
//...................................... |
10 |
// F77 routines |
11 |
//...................................... |
12 |
extern "C" { |
13 |
void dotrack_(int*, double*, double*, double*, double*, int*); |
14 |
void dotrack2_(int*, double*, double*, double*, double*,double*, double*, double*,int*); |
15 |
void dotrack3_(int*, double*, double*, double*, double*,double*, double*, double*,double*,int*); |
16 |
void mini2_(int*,int*,int*); |
17 |
void guess_(); |
18 |
void gufld_(float*, float*); |
19 |
float risxeta2_(float *); |
20 |
float risxeta3_(float *); |
21 |
float risxeta4_(float *); |
22 |
float risyeta2_(float *); |
23 |
} |
24 |
|
25 |
//-------------------------------------- |
26 |
// |
27 |
// |
28 |
//-------------------------------------- |
29 |
TrkTrack::TrkTrack(){ |
30 |
// cout << "TrkTrack::TrkTrack()" << endl; |
31 |
seqno = -1; |
32 |
image = -1; |
33 |
chi2 = 0; |
34 |
nstep = 0; |
35 |
for(int it1=0;it1<5;it1++){ |
36 |
al[it1] = 0; |
37 |
for(int it2=0;it2<5;it2++)coval[it1][it2] = 0; |
38 |
}; |
39 |
for(int ip=0;ip<6;ip++){ |
40 |
xgood[ip] = 0; |
41 |
ygood[ip] = 0; |
42 |
xm[ip] = 0; |
43 |
ym[ip] = 0; |
44 |
zm[ip] = 0; |
45 |
resx[ip] = 0; |
46 |
resy[ip] = 0; |
47 |
tailx[ip] = 0; |
48 |
taily[ip] = 0; |
49 |
xv[ip] = 0; |
50 |
yv[ip] = 0; |
51 |
zv[ip] = 0; |
52 |
axv[ip] = 0; |
53 |
ayv[ip] = 0; |
54 |
dedx_x[ip] = 0; |
55 |
dedx_y[ip] = 0; |
56 |
multmaxx[ip] = 0; |
57 |
multmaxy[ip] = 0; |
58 |
seedx[ip] = 0; |
59 |
seedy[ip] = 0; |
60 |
xpu[ip] = 0; |
61 |
ypu[ip] = 0; |
62 |
|
63 |
}; |
64 |
|
65 |
// TrkParams::SetTrackingMode(); |
66 |
// TrkParams::SetPrecisionFactor(); |
67 |
// TrkParams::SetStepMin(); |
68 |
TrkParams::SetMiniDefault(); |
69 |
TrkParams::SetPFA(); |
70 |
|
71 |
int ngf = TrkParams::nGF; |
72 |
for(int i=0; i<ngf; i++){ |
73 |
xGF[i] = 0.; |
74 |
yGF[i] = 0.; |
75 |
} |
76 |
|
77 |
|
78 |
}; |
79 |
//-------------------------------------- |
80 |
// |
81 |
// |
82 |
//-------------------------------------- |
83 |
TrkTrack::TrkTrack(const TrkTrack& t){ |
84 |
seqno = t.seqno; |
85 |
image = t.image; |
86 |
chi2 = t.chi2; |
87 |
nstep = t.nstep; |
88 |
for(int it1=0;it1<5;it1++){ |
89 |
al[it1] = t.al[it1]; |
90 |
for(int it2=0;it2<5;it2++)coval[it1][it2] = t.coval[it1][it2]; |
91 |
}; |
92 |
for(int ip=0;ip<6;ip++){ |
93 |
xgood[ip] = t.xgood[ip]; |
94 |
ygood[ip] = t.ygood[ip]; |
95 |
xm[ip] = t.xm[ip]; |
96 |
ym[ip] = t.ym[ip]; |
97 |
zm[ip] = t.zm[ip]; |
98 |
resx[ip] = t.resx[ip]; |
99 |
resy[ip] = t.resy[ip]; |
100 |
tailx[ip] = t.tailx[ip]; |
101 |
taily[ip] = t.taily[ip]; |
102 |
xv[ip] = t.xv[ip]; |
103 |
yv[ip] = t.yv[ip]; |
104 |
zv[ip] = t.zv[ip]; |
105 |
axv[ip] = t.axv[ip]; |
106 |
ayv[ip] = t.ayv[ip]; |
107 |
dedx_x[ip] = t.dedx_x[ip]; |
108 |
dedx_y[ip] = t.dedx_y[ip]; |
109 |
multmaxx[ip] = t.multmaxx[ip]; |
110 |
multmaxy[ip] = t.multmaxy[ip]; |
111 |
seedx[ip] = t.seedx[ip]; |
112 |
seedy[ip] = t.seedy[ip]; |
113 |
xpu[ip] = t.xpu[ip]; |
114 |
ypu[ip] = t.ypu[ip]; |
115 |
}; |
116 |
|
117 |
// TrkParams::SetTrackingMode(); |
118 |
// TrkParams::SetPrecisionFactor(); |
119 |
// TrkParams::SetStepMin(); |
120 |
TrkParams::SetMiniDefault(); |
121 |
TrkParams::SetPFA(); |
122 |
|
123 |
int ngf = TrkParams::nGF; |
124 |
for(int i=0; i<ngf; i++){ |
125 |
xGF[i] = t.xGF[i]; |
126 |
yGF[i] = t.yGF[i]; |
127 |
} |
128 |
}; |
129 |
//-------------------------------------- |
130 |
// |
131 |
// |
132 |
//-------------------------------------- |
133 |
void TrkTrack::Copy(TrkTrack& t){ |
134 |
|
135 |
t.seqno = seqno; |
136 |
t.image = image; |
137 |
t.chi2 = chi2; |
138 |
t.nstep = nstep; |
139 |
for(int it1=0;it1<5;it1++){ |
140 |
t.al[it1] = al[it1]; |
141 |
for(int it2=0;it2<5;it2++)t.coval[it1][it2] = coval[it1][it2]; |
142 |
}; |
143 |
for(int ip=0;ip<6;ip++){ |
144 |
t.xgood[ip] = xgood[ip]; |
145 |
t.ygood[ip] = ygood[ip]; |
146 |
t.xm[ip] = xm[ip]; |
147 |
t.ym[ip] = ym[ip]; |
148 |
t.zm[ip] = zm[ip]; |
149 |
t.resx[ip] = resx[ip]; |
150 |
t.resy[ip] = resy[ip]; |
151 |
t.tailx[ip] = tailx[ip]; |
152 |
t.taily[ip] = taily[ip]; |
153 |
t.xv[ip] = xv[ip]; |
154 |
t.yv[ip] = yv[ip]; |
155 |
t.zv[ip] = zv[ip]; |
156 |
t.axv[ip] = axv[ip]; |
157 |
t.ayv[ip] = ayv[ip]; |
158 |
t.dedx_x[ip] = dedx_x[ip]; |
159 |
t.dedx_y[ip] = dedx_y[ip]; |
160 |
t.multmaxx[ip] = multmaxx[ip]; |
161 |
t.multmaxy[ip] = multmaxy[ip]; |
162 |
t.seedx[ip] = seedx[ip]; |
163 |
t.seedy[ip] = seedy[ip]; |
164 |
t.xpu[ip] = xpu[ip]; |
165 |
t.ypu[ip] = ypu[ip]; |
166 |
|
167 |
}; |
168 |
int ngf = TrkParams::nGF; |
169 |
for(int i=0; i<ngf; i++){ |
170 |
t.xGF[i] = xGF[i]; |
171 |
t.yGF[i] = yGF[i]; |
172 |
} |
173 |
|
174 |
|
175 |
}; |
176 |
//-------------------------------------- |
177 |
// |
178 |
// |
179 |
//-------------------------------------- |
180 |
/** |
181 |
* |
182 |
* >>> OBSOLETE !!! use TrkTrack::DoTrack(Trajectory* t) instead |
183 |
* |
184 |
*/ |
185 |
int TrkTrack::DoTrack2(Trajectory* t){ |
186 |
|
187 |
cout << endl; |
188 |
cout << " int TrkTrack::DoTrack2(Trajectory* t) --->> NB NB !! this method is going to be eliminated !!! "<<endl; |
189 |
cout << " >>>> replace it with TrkTrack::DoTrack(Trajectory* t) <<<<"<<endl; |
190 |
cout << " (Sorry Wolfgang!! Don't be totally confused!! By Elena)"<<endl; |
191 |
cout << endl; |
192 |
|
193 |
return DoTrack(t); |
194 |
|
195 |
}; |
196 |
//-------------------------------------- |
197 |
// |
198 |
// |
199 |
//-------------------------------------- |
200 |
/** |
201 |
* Evaluates the trajectory in the apparatus associated to the track state-vector. |
202 |
* It integrates the equations of motion in the magnetic field. |
203 |
* @param t pointer to an object of the class Trajectory, |
204 |
* which z coordinates should be previously assigned. |
205 |
* @return error flag. |
206 |
*/ |
207 |
int TrkTrack::DoTrack(Trajectory* t){ |
208 |
|
209 |
double *dxout = new double[t->npoint]; |
210 |
double *dyout = new double[t->npoint]; |
211 |
double *dthxout = new double[t->npoint]; |
212 |
double *dthyout = new double[t->npoint]; |
213 |
double *dtlout = new double[t->npoint]; |
214 |
double *dzin = new double[t->npoint]; |
215 |
double dal[5]; |
216 |
|
217 |
int ifail = 0; |
218 |
|
219 |
for (int i=0; i<5; i++) dal[i] = (double)al[i]; |
220 |
for (int i=0; i<t->npoint; i++) dzin[i] = (double)t->z[i]; |
221 |
|
222 |
TrkParams::Load(1); |
223 |
if( !TrkParams::IsLoaded(1) ){ |
224 |
cout << "int TrkTrack::DoTrack(Trajectory* t) --- ERROR --- m.field not loaded"<<endl; |
225 |
return 0; |
226 |
} |
227 |
dotrack2_(&(t->npoint),dzin,dxout,dyout,dthxout,dthyout,dtlout,dal,&ifail); |
228 |
|
229 |
for (int i=0; i<t->npoint; i++){ |
230 |
t->x[i] = (float)*(dxout+i); |
231 |
t->y[i] = (float)*(dyout+i); |
232 |
t->thx[i] = (float)*(dthxout+i); |
233 |
t->thy[i] = (float)*(dthyout+i); |
234 |
t->tl[i] = (float)*(dtlout+i); |
235 |
} |
236 |
|
237 |
delete [] dxout; |
238 |
delete [] dyout; |
239 |
delete [] dzin; |
240 |
delete [] dthxout; |
241 |
delete [] dthyout; |
242 |
delete [] dtlout; |
243 |
|
244 |
return ifail; |
245 |
}; |
246 |
//-------------------------------------- |
247 |
// |
248 |
// |
249 |
//-------------------------------------- |
250 |
//float TrkTrack::BdL(){ |
251 |
//}; |
252 |
//-------------------------------------- |
253 |
// |
254 |
// |
255 |
//-------------------------------------- |
256 |
Float_t TrkTrack::GetRigidity(){ |
257 |
Float_t rig=0; |
258 |
if(chi2>0)rig=1./al[4]; |
259 |
if(rig<0) rig=-rig; |
260 |
return rig; |
261 |
}; |
262 |
// |
263 |
Float_t TrkTrack::GetDeflection(){ |
264 |
Float_t def=0; |
265 |
if(chi2>0)def=al[4]; |
266 |
return def; |
267 |
}; |
268 |
// |
269 |
/** |
270 |
* Method to retrieve the dE/dx measured on a tracker view. |
271 |
* @param ip plane (0-5) |
272 |
* @param iv view (0=x 1=y) |
273 |
*/ |
274 |
Float_t TrkTrack::GetDEDX(int ip, int iv){ |
275 |
if(iv==0 && ip>=0 && ip<6)return fabs(dedx_x[ip]); |
276 |
else if(iv==1 && ip>=0 && ip<6)return fabs(dedx_y[ip]); |
277 |
else { |
278 |
cout << "TrkTrack::GetDEDX(int ip, int iv) -- wrong input parameters "<<ip<<iv<<endl; |
279 |
return 0.; |
280 |
} |
281 |
} |
282 |
/** |
283 |
* Method to evaluate the dE/dx measured on a tracker plane. |
284 |
* The two measurements on x- and y-view are averaged. |
285 |
* @param ip plane (0-5) |
286 |
*/ |
287 |
Float_t TrkTrack::GetDEDX(int ip){ |
288 |
if( (Int_t)XGood(ip)+(Int_t)YGood(ip) == 0 ) return 0; |
289 |
return (GetDEDX(ip,0)+GetDEDX(ip,1))/((Int_t)XGood(ip)+(Int_t)YGood(ip)); |
290 |
}; |
291 |
|
292 |
/** |
293 |
* Method to evaluate the dE/dx averaged over all planes. |
294 |
*/ |
295 |
Float_t TrkTrack::GetDEDX(){ |
296 |
Float_t dedx=0; |
297 |
for(Int_t ip=0; ip<6; ip++)dedx+=GetDEDX(ip,0)*XGood(ip)+GetDEDX(ip,1)*YGood(ip); |
298 |
dedx = dedx/(GetNX()+GetNY()); |
299 |
return dedx; |
300 |
}; |
301 |
/** |
302 |
* Returns 1 if the cluster on a tracker view includes bad strips |
303 |
* (at least one bad strip among the four strip used by p.f.a.) |
304 |
* @param ip plane (0-5) |
305 |
* @param iv view (0=x 1=y) |
306 |
*/ |
307 |
Bool_t TrkTrack::IsBad(int ip,int iv){ |
308 |
if(iv==0 && ip>=0 && ip<6)return (xgood[ip]<0) ; |
309 |
else if(iv==1 && ip>=0 && ip<6)return (ygood[ip]<0) ; |
310 |
else { |
311 |
cout << "TrkTrack::IsBad(int ip, int iv) -- wrong input parameters "<<ip<<iv<<endl; |
312 |
return 0.; |
313 |
} |
314 |
}; |
315 |
/** |
316 |
* Returns 1 if the signal on a tracker view is saturated. |
317 |
* @param ip plane (0-5) |
318 |
* @param iv view (0=x 1=y) |
319 |
*/ |
320 |
Bool_t TrkTrack::IsSaturated(int ip,int iv){ |
321 |
if(iv==0 && ip>=0 && ip<6)return (dedx_x[ip]<0) ; |
322 |
else if(iv==1 && ip>=0 && ip<6)return (dedx_y[ip]<0) ; |
323 |
else { |
324 |
cout << "TrkTrack::IsSaturated(int ip, int iv) -- wrong input parameters "<<ip<<iv<<endl; |
325 |
return 0.; |
326 |
} |
327 |
}; |
328 |
/** |
329 |
* Returns 1 if either the x or the y signal on a tracker plane is saturated. |
330 |
* @param ip plane (0-5) |
331 |
*/ |
332 |
Bool_t TrkTrack::IsSaturated(int ip){ |
333 |
return (IsSaturated(ip,0)||IsSaturated(ip,1)); |
334 |
}; |
335 |
/** |
336 |
* Returns 1 if there is at least a saturated signal along the track. |
337 |
*/ |
338 |
Bool_t TrkTrack::IsSaturated(){ |
339 |
for(int ip=0; ip<6; ip++)for(int iv=0; iv<2; iv++)if(IsSaturated(ip,iv))return true; |
340 |
return false; |
341 |
} |
342 |
/** |
343 |
* Returns the track "lever-arm" on the x view, defined as the distance (in planes) between |
344 |
* the upper and lower x measurements (the maximum value of lever-arm is 6). |
345 |
*/ |
346 |
Int_t TrkTrack::GetLeverArmX(){ |
347 |
int first_plane = -1; |
348 |
int last_plane = -1; |
349 |
for(Int_t ip=0; ip<6; ip++){ |
350 |
if( XGood(ip) && first_plane == -1 )first_plane = ip; |
351 |
if( XGood(ip) && first_plane != -1 )last_plane = ip; |
352 |
} |
353 |
if( first_plane == -1 || last_plane == -1){ |
354 |
cout<< "Int_t TrkTrack::GetLeverArmX() -- XGood(ip) always false ??? "<<endl; |
355 |
return 0; |
356 |
} |
357 |
return (last_plane-first_plane+1); |
358 |
} |
359 |
/** |
360 |
* Returns the track "lever-arm" on the y view, defined as the distance (in planes) between |
361 |
* the upper and lower y measurements (the maximum value of lever-arm is 6). |
362 |
*/ |
363 |
Int_t TrkTrack::GetLeverArmY(){ |
364 |
int first_plane = -1; |
365 |
int last_plane = -1; |
366 |
for(Int_t ip=0; ip<6; ip++){ |
367 |
if( YGood(ip) && first_plane == -1 )first_plane = ip; |
368 |
if( YGood(ip) && first_plane != -1 )last_plane = ip; |
369 |
} |
370 |
if( first_plane == -1 || last_plane == -1){ |
371 |
cout<< "Int_t TrkTrack::GetLeverArmY() -- YGood(ip) always false ??? "<<endl; |
372 |
return 0; |
373 |
} |
374 |
return (last_plane-first_plane+1); |
375 |
} |
376 |
/** |
377 |
* Returns the track "lever-arm" on the x+y view, defined as the distance (in planes) between |
378 |
* the upper and lower x,y (couple) measurements (the maximum value of lever-arm is 6). |
379 |
*/ |
380 |
Int_t TrkTrack::GetLeverArmXY(){ |
381 |
int first_plane = -1; |
382 |
int last_plane = -1; |
383 |
for(Int_t ip=0; ip<6; ip++){ |
384 |
if( XGood(ip)*YGood(ip) && first_plane == -1 )first_plane = ip; |
385 |
if( XGood(ip)*YGood(ip) && first_plane != -1 )last_plane = ip; |
386 |
} |
387 |
if( first_plane == -1 || last_plane == -1){ |
388 |
cout<< "Int_t TrkTrack::GetLeverArmXY() -- XGood(ip)*YGood(ip) always false ??? "<<endl; |
389 |
return 0; |
390 |
} |
391 |
return (last_plane-first_plane+1); |
392 |
} |
393 |
/** |
394 |
* Returns the number of hit planes |
395 |
*/ |
396 |
Int_t TrkTrack::GetNhit() { |
397 |
int np=0; |
398 |
for(Int_t ip=0; ip<6; ip++) np += (XGood(ip)||YGood(ip)) ; |
399 |
return np; |
400 |
}; |
401 |
/** |
402 |
* Returns the reduced chi-square of track x-projection |
403 |
*/ |
404 |
Float_t TrkTrack::GetChi2X(){ |
405 |
float chiq=0; |
406 |
for(int ip=0; ip<6; ip++)if(XGood(ip))chiq+= pow((xv[ip]-xm[ip])/resx[ip],2.); |
407 |
if(GetNX()>3)chiq=chiq/(GetNX()-3); |
408 |
else chiq=0; |
409 |
if(chiq==0)cout << " Float_t TrkTrack::GetChi2X() -- WARNING -- value not defined "<<chiq<<endl; |
410 |
return chiq; |
411 |
} |
412 |
/** |
413 |
* Returns the reduced chi-square of track y-projection |
414 |
*/ |
415 |
Float_t TrkTrack::GetChi2Y(){ |
416 |
float chiq=0; |
417 |
for(int ip=0; ip<6; ip++)if(YGood(ip))chiq+= pow((yv[ip]-ym[ip])/resy[ip],2.); |
418 |
if(GetNY()>2)chiq=chiq/(GetNY()-2); |
419 |
else chiq=0; |
420 |
if(chiq==0)cout << " Float_t TrkTrack::GetChi2Y() -- WARNING -- value not defined "<<chiq<<endl; |
421 |
return chiq; |
422 |
} |
423 |
/** |
424 |
* Returns the logarythm of the likeliwood-function of track x-projection |
425 |
*/ |
426 |
Float_t TrkTrack::GetLnLX(){ |
427 |
float lnl=0; |
428 |
for(int ip=0; ip<6; ip++) |
429 |
if( XGood(ip) && tailx[ip]!=0 ) |
430 |
lnl += (tailx[ip]+1.) * log( (tailx[ip]*pow(resx[ip],2.) + pow(xv[ip]-xm[ip],2.)) / (tailx[ip]*pow(resx[ip],2)) ); |
431 |
if(GetNX()>3)lnl=lnl/(GetNX()-3); |
432 |
else lnl=0; |
433 |
if(lnl==0){ |
434 |
cout << " Float_t TrkTrack::GetLnLX() -- WARNING -- value not defined "<<lnl<<endl; |
435 |
Dump(); |
436 |
} |
437 |
return lnl; |
438 |
|
439 |
} |
440 |
/** |
441 |
* Returns the logarythm of the likeliwood-function of track y-projection |
442 |
*/ |
443 |
Float_t TrkTrack::GetLnLY(){ |
444 |
float lnl=0; |
445 |
for(int ip=0; ip<6; ip++) |
446 |
if( YGood(ip) && taily[ip]!=0 ) |
447 |
lnl += (taily[ip]+1.) * log( (taily[ip]*pow(resy[ip],2.) + pow(yv[ip]-ym[ip],2.)) / (taily[ip]*pow(resy[ip],2)) ); |
448 |
if(GetNY()>2)lnl=lnl/(GetNY()-2); |
449 |
else lnl=0; |
450 |
if(lnl==0){ |
451 |
cout << " Float_t TrkTrack::GetLnLY() -- WARNING -- value not defined "<<lnl<<endl; |
452 |
Dump(); |
453 |
} |
454 |
return lnl; |
455 |
|
456 |
} |
457 |
/** |
458 |
* Returns the effective angle, relative to the sensor, on each plane. |
459 |
* @param ip plane (0-5) |
460 |
* @param iv view (0=x 1=y) |
461 |
*/ |
462 |
Float_t TrkTrack::GetEffectiveAngle(int ip, int iv){ |
463 |
|
464 |
if(ip<0 || ip>5){ |
465 |
cout << "Float_t TrkTrack::GetEffectiveAngle(int "<<ip<<", int "<<iv<<") ==> wrong input"<<endl; |
466 |
return 0.; |
467 |
} |
468 |
|
469 |
float v[3]={xv[ip],yv[ip],zv[ip]}; |
470 |
//----------------------------------------- |
471 |
// effective angle (relative to the sensor) |
472 |
//----------------------------------------- |
473 |
float axv_geo = axv[ip]; |
474 |
float muhall_h = 297.61; //cm**2/Vs |
475 |
float BY = TrkParams::GetBY(v); |
476 |
float axv_eff = 0; |
477 |
if(ip==5) axv_geo = -1*axv_geo; |
478 |
if(ip==5) BY = -1*BY; |
479 |
axv_eff = 180.*atan( tan(axv_geo*acos(-1.)/180.) + muhall_h * BY * 0.0001)/acos(-1.); |
480 |
//----------------------------------------- |
481 |
// effective angle (relative to the sensor) |
482 |
//----------------------------------------- |
483 |
float ayv_geo = ayv[ip]; |
484 |
float muhall_e = 1258.18; //cm**2/Vs |
485 |
float BX = TrkParams::GetBX(v); |
486 |
float ayv_eff = 0; |
487 |
ayv_eff = 180.*atan( tan(ayv_geo*acos(-1.)/180.) + muhall_e * BX * 0.0001)/acos(-1.); |
488 |
|
489 |
if (iv==0)return axv_eff; |
490 |
else if(iv==1)return ayv_eff; |
491 |
else{ |
492 |
cout << "Float_t TrkTrack::GetEffectiveAngle(int "<<ip<<", int "<<iv<<") ==> wrong input"<<endl; |
493 |
return 0.; |
494 |
} |
495 |
|
496 |
}; |
497 |
|
498 |
//-------------------------------------- |
499 |
// |
500 |
// |
501 |
//-------------------------------------- |
502 |
void TrkTrack::Dump(){ |
503 |
cout << endl << "========== Track " ; |
504 |
cout << endl << "seq. n. : "<< seqno; |
505 |
cout << endl << "image n. : "<< image; |
506 |
cout << endl << "al : "; for(int i=0; i<5; i++)cout << al[i] << " "; |
507 |
cout << endl << "chi^2 : "<< chi2; |
508 |
cout << endl << "n.step : "<< nstep; |
509 |
cout << endl << "xgood : "; for(int i=0; i<6; i++)cout << XGood(i) ; |
510 |
cout << endl << "ygood : "; for(int i=0; i<6; i++)cout << YGood(i) ; |
511 |
cout << endl << "xm : "; for(int i=0; i<6; i++)cout << xm[i] << " "; |
512 |
cout << endl << "ym : "; for(int i=0; i<6; i++)cout << ym[i] << " "; |
513 |
cout << endl << "zm : "; for(int i=0; i<6; i++)cout << zm[i] << " "; |
514 |
cout << endl << "xv : "; for(int i=0; i<6; i++)cout << xv[i] << " "; |
515 |
cout << endl << "yv : "; for(int i=0; i<6; i++)cout << yv[i] << " "; |
516 |
cout << endl << "zv : "; for(int i=0; i<6; i++)cout << zv[i] << " "; |
517 |
cout << endl << "resx : "; for(int i=0; i<6; i++)cout << resx[i] << " "; |
518 |
cout << endl << "resy : "; for(int i=0; i<6; i++)cout << resy[i] << " "; |
519 |
cout << endl << "tailx : "; for(int i=0; i<6; i++)cout << tailx[i] << " "; |
520 |
cout << endl << "taily : "; for(int i=0; i<6; i++)cout << taily[i] << " "; |
521 |
cout << endl << "coval : "; for(int i=0; i<5; i++)cout << coval[0][i]<<" "; |
522 |
cout << endl << " "; for(int i=0; i<5; i++)cout << coval[1][i]<<" "; |
523 |
cout << endl << " "; for(int i=0; i<5; i++)cout << coval[2][i]<<" "; |
524 |
cout << endl << " "; for(int i=0; i<5; i++)cout << coval[3][i]<<" "; |
525 |
cout << endl << " "; for(int i=0; i<5; i++)cout << coval[4][i]<<" "; |
526 |
cout << endl << "dedx_x : "; for(int i=0; i<6; i++)cout << dedx_x[i] << " "; |
527 |
cout << endl << "dedx_y : "; for(int i=0; i<6; i++)cout << dedx_y[i] << " "; |
528 |
cout << endl << "maxs x : "; for(int i=0; i<6; i++)cout << GetClusterX_MaxStrip(i) << " "; |
529 |
cout << endl << "maxs y : "; for(int i=0; i<6; i++)cout << GetClusterY_MaxStrip(i) << " "; |
530 |
cout << endl << "mult x : "; for(int i=0; i<6; i++)cout << GetClusterX_Multiplicity(i) << " "; |
531 |
cout << endl << "mult y : "; for(int i=0; i<6; i++)cout << GetClusterY_Multiplicity(i) << " "; |
532 |
cout << endl << "seed x : "; for(int i=0; i<6; i++)cout << GetClusterX_Seed(i) << " "; |
533 |
cout << endl << "seed y : "; for(int i=0; i<6; i++)cout << GetClusterY_Seed(i) << " "; |
534 |
cout << endl << "xpu : "; for(int i=0; i<6; i++)cout << xpu[i] << " "; |
535 |
cout << endl << "ypu : "; for(int i=0; i<6; i++)cout << ypu[i] << " "; |
536 |
|
537 |
cout << endl; |
538 |
} |
539 |
/** |
540 |
* Set the TrkTrack position measurements |
541 |
*/ |
542 |
void TrkTrack::SetMeasure(double *xmeas, double *ymeas, double *zmeas){ |
543 |
for(int i=0; i<6; i++) xm[i]=*xmeas++; |
544 |
for(int i=0; i<6; i++) ym[i]=*ymeas++; |
545 |
for(int i=0; i<6; i++) zm[i]=*zmeas++; |
546 |
} |
547 |
/** |
548 |
* Set the TrkTrack position resolution |
549 |
*/ |
550 |
void TrkTrack::SetResolution(double *rx, double *ry){ |
551 |
for(int i=0; i<6; i++) resx[i]=*rx++; |
552 |
for(int i=0; i<6; i++) resy[i]=*ry++; |
553 |
} |
554 |
/** |
555 |
* Set the TrkTrack tails position resolution |
556 |
*/ |
557 |
void TrkTrack::SetTail(double *tx, double *ty, double factor){ |
558 |
for(int i=0; i<6; i++) tailx[i]=factor*(*tx++); |
559 |
for(int i=0; i<6; i++) taily[i]=factor*(*ty++); |
560 |
} |
561 |
/** |
562 |
* Set the TrkTrack Student parameter (resx,resy,tailx,taily) |
563 |
* from previous gausian fit |
564 |
*@param flag =0 standard, =1 with noise correction |
565 |
*/ |
566 |
void TrkTrack::SetStudentParam(int flag){ |
567 |
float sx[11]={0.000128242, |
568 |
0.000136942, |
569 |
0.000162718, |
570 |
0.000202644, |
571 |
0.00025597, |
572 |
0.000317456, |
573 |
0.000349048, |
574 |
0.000384638, |
575 |
0.000457295, |
576 |
0.000512319, |
577 |
0.000538573}; |
578 |
float tx[11]={1.79402, |
579 |
2.04876, |
580 |
2.88376, |
581 |
3.3, |
582 |
3.14084, |
583 |
4.07686, |
584 |
4.44736, |
585 |
3.5179, |
586 |
3.38697, |
587 |
3.45739, |
588 |
3.18627}; |
589 |
float sy[11]={0.000483075, |
590 |
0.000466925, |
591 |
0.000431658, |
592 |
0.000428317, |
593 |
0.000433854, |
594 |
0.000444044, |
595 |
0.000482098, |
596 |
0.000537579, |
597 |
0.000636279, |
598 |
0.000741998, |
599 |
0.000864261}; |
600 |
float ty[11]={0.997032, |
601 |
1.11147, |
602 |
1.18526, |
603 |
1.61404, |
604 |
2.21908, |
605 |
3.08959, |
606 |
4.48833, |
607 |
4.42687, |
608 |
4.65253, |
609 |
4.52043, |
610 |
4.29926}; |
611 |
int index; |
612 |
float fact=0.; |
613 |
for(int i=0; i<6; i++) { |
614 |
index = int((fabs(axv[i])+1.)/2.); |
615 |
if(index>10) index=10; |
616 |
tailx[i]=tx[index]; |
617 |
if(flag==1) { |
618 |
if(fabs(axv[i])<=10.) fact = resx[i]/risxeta2_(&(axv[i])); |
619 |
if(fabs(axv[i])>10.&&fabs(axv[i])<=15.) fact = resx[i]/risxeta3_(&(axv[i])); |
620 |
if(fabs(axv[i])>15.) fact = resx[i]/risxeta4_(&(axv[i])); |
621 |
} else fact = 1.; |
622 |
resx[i] = sx[index]*fact; |
623 |
} |
624 |
for(int i=0; i<6; i++) { |
625 |
index = int((fabs(ayv[i])+1.)/2.); |
626 |
if(index>10) index=10; |
627 |
taily[i]=ty[index]; |
628 |
if(flag==1) fact = resy[i]/risyeta2_(&(ayv[i])); |
629 |
else fact = 1.; |
630 |
resy[i] = sy[index]*fact; |
631 |
} |
632 |
} |
633 |
/** |
634 |
* Set the TrkTrack good measurement |
635 |
*/ |
636 |
void TrkTrack::SetGood(int *xg, int *yg){ |
637 |
|
638 |
for(int i=0; i<6; i++) xgood[i]=*xg++; |
639 |
for(int i=0; i<6; i++) ygood[i]=*yg++; |
640 |
} |
641 |
|
642 |
/** |
643 |
* Load the magnetic field |
644 |
*/ |
645 |
void TrkTrack::LoadField(TString path){ |
646 |
|
647 |
// strcpy(path_.path,path.Data()); |
648 |
// path_.pathlen = path.Length(); |
649 |
// path_.error = 0; |
650 |
// readb_(); |
651 |
|
652 |
// TrkParams::SetTrackingMode(); |
653 |
// TrkParams::SetPrecisionFactor(); |
654 |
// TrkParams::SetStepMin(); |
655 |
TrkParams::SetMiniDefault(); |
656 |
|
657 |
TrkParams::Set(path,1); |
658 |
TrkParams::Load(1); |
659 |
if( !TrkParams::IsLoaded(1) ){ |
660 |
cout << "void TrkTrack::LoadField(TString path) --- ERROR --- m.field not loaded"<<endl; |
661 |
} |
662 |
|
663 |
}; |
664 |
|
665 |
|
666 |
/** |
667 |
* Method to fill minimization-routine common |
668 |
*/ |
669 |
void TrkTrack::FillMiniStruct(cMini2track& track){ |
670 |
|
671 |
for(int i=0; i<6; i++){ |
672 |
|
673 |
// cout << i<<" - "<<xgood[i]<<" "<<XGood(i)<<endl; |
674 |
// cout << i<<" - "<<ygood[i]<<" "<<YGood(i)<<endl; |
675 |
track.xgood[i] = (double) XGood(i); |
676 |
track.ygood[i] = (double) YGood(i); |
677 |
|
678 |
track.xm[i]= (double) xm[i]; |
679 |
track.ym[i]= (double) ym[i]; |
680 |
track.zm[i]= (double) zm[i]; |
681 |
|
682 |
// --- temporaneo ---------------------------- |
683 |
// float segment = 100.; |
684 |
// track.xm_a[i]=xm[i]; |
685 |
// track.xm_b[i]=xm[i]; |
686 |
// track.ym_a[i]=ym[i]; |
687 |
// track.ym_b[i]=ym[i]; |
688 |
// if( XGood(i) && !YGood(i) ){ |
689 |
// track.ym_a[i] = track.ym_a[i]+segment; |
690 |
// track.ym_b[i] = track.ym_b[i]-segment; |
691 |
// }else if( !XGood(i) && YGood(i)){ |
692 |
// track.xm_a[i] = track.xm_a[i]+segment; |
693 |
// track.xm_b[i] = track.xm_b[i]-segment; |
694 |
// } |
695 |
// --- temporaneo ---------------------------- |
696 |
|
697 |
if( XGood(i) || YGood(i) ){ |
698 |
//NB!! the length of the sensor is not exactely taken into account |
699 |
double segment = 7.;// 2.;//cm //Elena 10th |
700 |
// NB: i parametri di allineamento hanno una notazione particolare!!! |
701 |
// sensor = 0 (hybrid side), 1 |
702 |
// ladder = 0-2 (increasing x) |
703 |
// plane = 0-5 (from bottom to top!!!) |
704 |
int is = (int)GetSensor(i); if(i==5)is=1-is; |
705 |
int ip = 5-i; |
706 |
int il = (int)GetLadder(i); |
707 |
|
708 |
double omega = 0.; |
709 |
// double beta = 0.;// EM GCC 4.7 |
710 |
// double gamma = 0.; |
711 |
if( |
712 |
(is < 0 || is > 1 || ip < 0 || ip > 5 || il < 0 || il > 2) && |
713 |
true){ |
714 |
// se il piano risulta colpito, ladder e sensore devono essere |
715 |
// assegnati correttamente |
716 |
cout << " void TrkTrack::FillMiniStruct(cMini2track&) --- WARNING --- sensor not defined, cannot read alignment parameters "<<endl; |
717 |
cout << " is ip il = "<<is<<" "<<ip<<" "<<il<<endl; |
718 |
}else{ |
719 |
omega = alignparameters_.omega[is][il][ip]; |
720 |
// beta = alignparameters_.beta[is][il][ip];// EM GCC 4.7 unused |
721 |
// gamma = alignparameters_.gamma[is][il][ip];// EM GCC 4.7 unused |
722 |
} |
723 |
|
724 |
if( XGood(i) && !YGood(i) ){ |
725 |
track.xm_a[i] = (double) xm[i] - omega * segment; |
726 |
track.ym_a[i] = (double) ym[i] + segment; |
727 |
// track.zm_a[i] = zm[i] + beta * segment;//not used yet |
728 |
track.xm_b[i] = (double) xm[i] + omega * segment; |
729 |
track.ym_b[i] = (double) ym[i] - segment; |
730 |
// track.zm_b[i] = zm[i] - beta * segment;//not used yet |
731 |
}else if( !XGood(i) && YGood(i) ){ |
732 |
track.xm_a[i] = (double) xm[i] + segment; |
733 |
track.ym_a[i] = (double) ym[i] + omega * segment; |
734 |
// track.zm_a[i] = zm[i] - gamma * segment;//not used yet |
735 |
track.xm_b[i] = (double) xm[i] - segment; |
736 |
track.ym_b[i] = (double) ym[i] - omega * segment; |
737 |
// track.zm_b[i] = zm[i] + gamma * segment;//not used yet |
738 |
} |
739 |
} |
740 |
|
741 |
track.resx[i]= (double) resx[i]; |
742 |
track.resy[i]= (double) resy[i]; |
743 |
track.tailx[i]= (double) tailx[i]; |
744 |
track.taily[i]= (double) taily[i]; |
745 |
} |
746 |
|
747 |
for(int i=0; i<5; i++) track.al[i]= (double) al[i]; |
748 |
track.zini = 23.5; |
749 |
// ZINI = 23.5 !!! it should be the same parameter in all codes |
750 |
|
751 |
} |
752 |
/** |
753 |
* Method to set values from minimization-routine common |
754 |
*/ |
755 |
void TrkTrack::SetFromMiniStruct(cMini2track *track){ |
756 |
|
757 |
for(int i=0; i<5; i++) { |
758 |
al[i]= (float) (track->al[i]); |
759 |
for(int j=0; j<5; j++) coval[i][j]= (float) (track->cov[i][j]); |
760 |
} |
761 |
chi2 = (float) (track->chi2); |
762 |
nstep = (float) (track->nstep); |
763 |
for(int i=0; i<6; i++){ |
764 |
xv[i] = (float) (track->xv[i]); |
765 |
yv[i] = (float) (track->yv[i]); |
766 |
zv[i] = (float) (track->zv[i]); |
767 |
xm[i] = (float) (track->xm[i]); |
768 |
ym[i] = (float) (track->ym[i]); |
769 |
zm[i] = (float) (track->zm[i]); |
770 |
axv[i] = (float) (track->axv[i]); |
771 |
ayv[i] = (float) (track->ayv[i]); |
772 |
resx[i] = (float) (track->resx[i]); //Elena 10th |
773 |
resy[i] = (float) (track->resy[i]); |
774 |
} |
775 |
|
776 |
} |
777 |
/** |
778 |
* \brief Method to re-evaluate coordinates of clusters associated with a track. |
779 |
* |
780 |
* The method can be applied only after recovering level1 information |
781 |
* (either by reprocessing single events from level0 or from |
782 |
* the TrkLevel1 branch, if present); it calls F77 subroutines that |
783 |
* read the level1 common and fill the minimization-routine common. |
784 |
* Some clusters can be excluded or added by means of the methods: |
785 |
* |
786 |
* TrkTrack::ResetXGood(int ip) |
787 |
* TrkTrack::ResetYGood(int ip) |
788 |
* TrkTrack::SetXGood(int ip, int cid, int is) |
789 |
* TrkTrack::SetYGood(int ip, int cid, int is) |
790 |
* |
791 |
* NB! The method TrkTrack::SetGood(int *xg, int *yg) set the plane-mask (0-1) |
792 |
* for the minimization-routine common. It deletes the cluster information |
793 |
* (at least for the moment...) thus cannot be applied before |
794 |
* TrkTrack::EvaluateClusterPositions(). |
795 |
* |
796 |
* Different p.f.a. can be applied by calling (once) the method: |
797 |
* |
798 |
* TrkParams::SetPFA(0); //Set ETA p.f.a. |
799 |
* |
800 |
* @see TrkParams::SetPFA(int) |
801 |
*/ |
802 |
Bool_t TrkTrack::EvaluateClusterPositions(){ |
803 |
|
804 |
// cout << "void TrkTrack::GetClusterositions() "<<endl; |
805 |
|
806 |
TrkParams::Load(1); |
807 |
if( !TrkParams::IsLoaded(1) ){ |
808 |
cout << "Bool_t TrkTrack::EvaluateClusterPositions() ---ERROR--- m.field not loaded "<<endl; |
809 |
return false; |
810 |
} |
811 |
TrkParams::Load(4); |
812 |
if( !TrkParams::IsLoaded(4) ){ |
813 |
cout << "Bool_t TrkTrack::EvaluateClusterPositions() ---ERROR--- p.f.a. par. not loaded "<<endl; |
814 |
return false; |
815 |
} |
816 |
TrkParams::Load(5); |
817 |
if( !TrkParams::IsLoaded(5) ){ |
818 |
cout << "Bool_t TrkTrack::EvaluateClusterPositions() ---ERROR--- alignment par. not loaded "<<endl; |
819 |
return false; |
820 |
} |
821 |
|
822 |
for(int ip=0; ip<6; ip++){ |
823 |
// cout << ip<<" ** "<<xm[ip]<<" / "<<ym[ip]<<endl;; |
824 |
int icx = GetClusterX_ID(ip)+1;//0=no-cluster,1-N |
825 |
int icy = GetClusterY_ID(ip)+1;//0=no-cluster,1-N |
826 |
int sensor = GetSensor(ip)+1;//<< convenzione "Paolo" |
827 |
if(ip==5 && sensor!=0)sensor=3-sensor;//<< convenzione "Elena" |
828 |
int ladder = GetLadder(ip)+1; |
829 |
float ax = axv[ip]; |
830 |
float ay = ayv[ip]; |
831 |
float v[3]; |
832 |
v[0]=xv[ip]; |
833 |
v[1]=yv[ip]; |
834 |
v[2]=zv[ip]; |
835 |
float bfx = 10*TrkParams::GetBX(v);//Tesla |
836 |
float bfy = 10*TrkParams::GetBY(v);//Tesla |
837 |
int ipp=ip+1; |
838 |
xyzpam_(&ipp,&icx,&icy,&ladder,&sensor,&ax,&ay,&bfx,&bfy); |
839 |
// if(icx<0 || icy<0)return false; |
840 |
} |
841 |
return true; |
842 |
} |
843 |
/** |
844 |
* \brief Tracking method. It calls F77 mini routine. |
845 |
* |
846 |
* @param pfixed Particle momentum. If pfixed=0 the momentum |
847 |
* is left as a free parameter, otherwise it is fixed to the input value. |
848 |
* @param fail Output flag (!=0 if the fit failed). |
849 |
* @param iprint Flag to set debug mode ( 0 = no output; 1 = verbose; 2 = debug). |
850 |
* @param froml1 Flag to re-evaluate positions (see TrkTrack::GetClusterPositions()). |
851 |
* |
852 |
* The option to re-evaluate positions can be used only after recovering |
853 |
* level1 information, eg. by reprocessing the single event. |
854 |
* |
855 |
* Example: |
856 |
* |
857 |
* if( !event->GetTrkLevel0() )return false; |
858 |
* event->GetTrkLevel0()->ProcessEvent(); // re-processing level0->level1 |
859 |
* int fail=0; |
860 |
* event->GetTrkLevel2()->GetTrack(0)->Fit(0.,fail,0,1); |
861 |
* |
862 |
* @see EvaluateClusterPositions() |
863 |
* |
864 |
* The fitting procedure can be varied by changing the tracking mode, |
865 |
* the fit-precision factor, the minimum number of step, etc. |
866 |
* @see SetTrackingMode(int) |
867 |
* @see SetPrecisionFactor(double) |
868 |
* @see SetStepMin(int) |
869 |
* @see SetDeltaB(int,double) |
870 |
*/ |
871 |
void TrkTrack::Fit(double pfixed, int& fail, int iprint, int froml1){ |
872 |
|
873 |
TrkParams::Load(1); |
874 |
if( !TrkParams::IsLoaded(1) ){ |
875 |
cout << "void TrkTrack::Fit(double,int&,int,int) ---ERROR--- m.field not loaded "<<endl; |
876 |
return; |
877 |
} |
878 |
TrkParams::Load(5); |
879 |
if( !TrkParams::IsLoaded(5) ){ |
880 |
cout << "void TrkTrack::Fit(double,int&,int,int) ---ERROR--- align.param. not loaded "<<endl; |
881 |
return; |
882 |
} |
883 |
|
884 |
float al_ini[] = {0.,0.,0.,0.,0.}; |
885 |
|
886 |
extern cMini2track track_; |
887 |
fail = 0; |
888 |
|
889 |
// FillMiniStruct(track_); |
890 |
|
891 |
if(froml1!=0){ |
892 |
if( !EvaluateClusterPositions() ){ |
893 |
cout << "void TrkTrack::Fit("<<pfixed<<","<<fail<<","<<iprint<<","<<froml1<<") --- ERROR evaluating cluster positions "<<endl; |
894 |
FillMiniStruct(track_) ; |
895 |
fail = 1; |
896 |
return; |
897 |
} |
898 |
}else{ |
899 |
FillMiniStruct(track_); |
900 |
} |
901 |
|
902 |
// if fit variables have been reset, evaluate the initial guess |
903 |
if(al[0]==-9999.&&al[1]==-9999.&&al[2]==-9999.&&al[3]==-9999.&&al[4]==-9999.)guess_(); |
904 |
|
905 |
// --------------------- free momentum |
906 |
if(pfixed==0.) { |
907 |
track_.pfixed=0.; |
908 |
} |
909 |
// --------------------- fixed momentum |
910 |
if(pfixed!=0.) { |
911 |
al[4]=1./pfixed; |
912 |
track_.pfixed=pfixed; |
913 |
} |
914 |
|
915 |
// store temporarily the initial guess |
916 |
for(int i=0; i<5; i++) al_ini[i]=track_.al[i]; |
917 |
|
918 |
// ------------------------------------------ |
919 |
// call mini routine |
920 |
// ------------------------------------------ |
921 |
int istep=0; |
922 |
int ifail=0; |
923 |
mini2_(&istep,&ifail, &iprint); |
924 |
if(ifail!=0) { |
925 |
if(iprint)cout << "ERROR: ifail= " << ifail << endl; |
926 |
fail = 1; |
927 |
} |
928 |
if(chi2!=chi2){ |
929 |
if(iprint)cout << "ERROR: chi2= " << chi2 << endl; |
930 |
FitReset(); |
931 |
fail = 1; |
932 |
} |
933 |
// ------------------------------------------ |
934 |
|
935 |
SetFromMiniStruct(&track_); |
936 |
|
937 |
if(fail){ |
938 |
if(iprint)cout << " >>>> fit failed "<<endl; |
939 |
for(int i=0; i<5; i++) al[i]=al_ini[i]; |
940 |
} |
941 |
|
942 |
}; |
943 |
/** |
944 |
* Reset the fit parameters |
945 |
*/ |
946 |
void TrkTrack::FitReset(){ |
947 |
for(int i=0; i<5; i++) al[i]=-9999.; |
948 |
chi2=0.; |
949 |
nstep=0; |
950 |
// for(int i=0; i<6; i++) xv[i]=0.; |
951 |
// for(int i=0; i<6; i++) yv[i]=0.; |
952 |
// for(int i=0; i<6; i++) zv[i]=0.; |
953 |
// for(int i=0; i<6; i++) axv[i]=0.; |
954 |
// for(int i=0; i<6; i++) ayv[i]=0.; |
955 |
for(int i=0; i<5; i++) { |
956 |
for(int j=0; j<5; j++) coval[i][j]=0.; |
957 |
} |
958 |
} |
959 |
/** |
960 |
* Set the tracking mode |
961 |
*/ |
962 |
void TrkTrack::SetTrackingMode(int trackmode){ |
963 |
extern cMini2track track_; |
964 |
track_.trackmode = trackmode; |
965 |
} |
966 |
/** |
967 |
* Set the factor scale for tracking precision |
968 |
*/ |
969 |
void TrkTrack::SetPrecisionFactor(double fact){ |
970 |
extern cMini2track track_; |
971 |
track_.fact = fact; |
972 |
} |
973 |
/** |
974 |
* Set the minimum number of steps for tracking precision |
975 |
*/ |
976 |
void TrkTrack::SetStepMin(int istepmin){ |
977 |
extern cMini2track track_; |
978 |
track_.istepmin = istepmin; |
979 |
} |
980 |
/** |
981 |
* Set deltaB parameters (id=0,1). By default they are set to zero. |
982 |
*/ |
983 |
void TrkTrack::SetDeltaB(int id, double db){ |
984 |
if(id!=0 && id!=1)cout << "void TrkTrack::SetDeltaB(int id,double db) -- wrong input parameters: "<<id<<" "<<db<<endl; |
985 |
TrkParams::SetDeltaB(id,db); |
986 |
} |
987 |
|
988 |
/** |
989 |
* Returns true if the track is inside the magnet cavity. |
990 |
* @param toll Tolerance around the nominal volume (toll>0 define an inner fiducial volume) |
991 |
*/ |
992 |
Bool_t TrkTrack::IsInsideCavity(float toll){ |
993 |
|
994 |
// float xmagntop, ymagntop, xmagnbottom, ymagnbottom; |
995 |
// xmagntop = xv[0] + (ZMAGNHIGH-zv[0])*tan(acos(-1.0)*axv[0]/180.); |
996 |
// ymagntop = yv[0] + (ZMAGNHIGH-zv[0])*tan(acos(-1.0)*ayv[0]/180.); |
997 |
// xmagnbottom = xv[5] + (ZMAGNLOW-zv[5])*tan(acos(-1.0)*axv[5]/180.); |
998 |
// ymagnbottom = yv[5] + (ZMAGNLOW-zv[5])*tan(acos(-1.0)*ayv[5]/180.); |
999 |
// if( xmagntop>XMAGNLOW && xmagntop<XMAGNHIGH && |
1000 |
// ymagntop>YMAGNLOW && ymagntop<YMAGNHIGH && |
1001 |
// xmagnbottom>XMAGNLOW && xmagnbottom<XMAGNHIGH && |
1002 |
// ymagnbottom>YMAGNLOW && ymagnbottom<YMAGNHIGH ) return(true); |
1003 |
// else return(false); |
1004 |
|
1005 |
int ngf = TrkParams::nGF; |
1006 |
for(int i=0; i<ngf; i++){ |
1007 |
// |
1008 |
// cout << endl << TrkParams::GF_element[i]; |
1009 |
if( |
1010 |
TrkParams::GF_element[i].CompareTo("CUF") && |
1011 |
TrkParams::GF_element[i].CompareTo("T2") && |
1012 |
TrkParams::GF_element[i].CompareTo("T3") && |
1013 |
TrkParams::GF_element[i].CompareTo("T4") && |
1014 |
TrkParams::GF_element[i].CompareTo("T5") && |
1015 |
TrkParams::GF_element[i].CompareTo("CLF") && |
1016 |
true)continue; |
1017 |
// apply condition only within the cavity |
1018 |
// cout << " -- "<<xGF[i]<<" "<<yGF[i]; |
1019 |
if( |
1020 |
xGF[i] <= TrkParams::xGF_min[i] + toll || |
1021 |
xGF[i] >= TrkParams::xGF_max[i] - toll || |
1022 |
yGF[i] <= TrkParams::yGF_min[i] + toll || |
1023 |
yGF[i] >= TrkParams::yGF_max[i] - toll || |
1024 |
false){ |
1025 |
|
1026 |
return false; |
1027 |
} |
1028 |
} |
1029 |
return true; |
1030 |
|
1031 |
|
1032 |
} |
1033 |
/** |
1034 |
* Returns true if the track is inside the nominal acceptance, which is defined |
1035 |
* by the intersection among magnet cavity, silicon-plane sensitive area and |
1036 |
* ToF sensitive area (nominal values from the official document used to |
1037 |
* calculate the geometrical factor) |
1038 |
* @param toll Tolerance around the nominal volume (toll>0 define an inner fiducial volume) |
1039 |
*/ |
1040 |
// Bool_t TrkTrack::IsInsideAcceptance(){ |
1041 |
|
1042 |
// int ngf = TrkParams::nGF; |
1043 |
// for(int i=0; i<ngf; i++){ |
1044 |
// if( |
1045 |
// xGF[i] <= TrkParams::xGF_min[i] || |
1046 |
// xGF[i] >= TrkParams::xGF_max[i] || |
1047 |
// yGF[i] <= TrkParams::yGF_min[i] || |
1048 |
// yGF[i] >= TrkParams::yGF_max[i] || |
1049 |
// false)return false; |
1050 |
// } |
1051 |
// return true; |
1052 |
|
1053 |
// } |
1054 |
Bool_t TrkTrack::IsInsideAcceptance(float toll){ |
1055 |
|
1056 |
|
1057 |
int ngf = TrkParams::nGF; |
1058 |
for(int i=0; i<ngf; i++){ |
1059 |
// |
1060 |
// cout << endl << TrkParams::GF_element[i]; |
1061 |
if( |
1062 |
TrkParams::GF_element[i].CompareTo("S11") && |
1063 |
TrkParams::GF_element[i].CompareTo("S12") && |
1064 |
TrkParams::GF_element[i].CompareTo("S21") && |
1065 |
TrkParams::GF_element[i].CompareTo("S22") && |
1066 |
TrkParams::GF_element[i].CompareTo("T1") && |
1067 |
TrkParams::GF_element[i].CompareTo("CUF") && |
1068 |
TrkParams::GF_element[i].CompareTo("T2") && |
1069 |
TrkParams::GF_element[i].CompareTo("T3") && |
1070 |
TrkParams::GF_element[i].CompareTo("T4") && |
1071 |
TrkParams::GF_element[i].CompareTo("T5") && |
1072 |
TrkParams::GF_element[i].CompareTo("CLF") && |
1073 |
TrkParams::GF_element[i].CompareTo("T6") && |
1074 |
TrkParams::GF_element[i].CompareTo("S31") && |
1075 |
TrkParams::GF_element[i].CompareTo("S32") && |
1076 |
true)continue; |
1077 |
// apply condition only within the cavity |
1078 |
// cout << " -- "<<xGF[i]<<" "<<yGF[i]; |
1079 |
if( |
1080 |
xGF[i] <= TrkParams::xGF_min[i] + toll || |
1081 |
xGF[i] >= TrkParams::xGF_max[i] - toll || |
1082 |
yGF[i] <= TrkParams::yGF_min[i] + toll || |
1083 |
yGF[i] >= TrkParams::yGF_max[i] - toll || |
1084 |
false){ |
1085 |
|
1086 |
return false; |
1087 |
} |
1088 |
} |
1089 |
return true; |
1090 |
} |
1091 |
|
1092 |
/** |
1093 |
* Returns true if the track is inside one of the surfaces which define the |
1094 |
* geometrical acceptance. |
1095 |
* @param surf tag of the surface (possible values are: S11 S12 S21 S22 T1 |
1096 |
* CUF T2 T3 T4 T5 CLF T6 S31 S32). |
1097 |
* @param toll Tolerance around the nominal surface (toll>0 define an inner |
1098 |
* fiducial surface) |
1099 |
*/ |
1100 |
Bool_t TrkTrack::IsInsideGFSurface(const char* surf, float toll){ |
1101 |
|
1102 |
|
1103 |
int ngf = TrkParams::nGF; |
1104 |
bool SURFOK = false; |
1105 |
for(int i=0; i<ngf; i++){ |
1106 |
if( !TrkParams::GF_element[i].CompareTo(surf) ){ |
1107 |
SURFOK=true; |
1108 |
if( |
1109 |
xGF[i] > TrkParams::xGF_min[i] + toll && |
1110 |
xGF[i] < TrkParams::xGF_max[i] - toll && |
1111 |
yGF[i] > TrkParams::yGF_min[i] + toll && |
1112 |
yGF[i] < TrkParams::yGF_max[i] - toll && |
1113 |
true)return true; |
1114 |
} |
1115 |
} |
1116 |
if( !SURFOK )cout << " Bool_t TrkTrack::IsInsideGFSurface(char* surf, float toll) --> suface "<<surf<<" not defined "<<endl; |
1117 |
return false; |
1118 |
|
1119 |
} |
1120 |
|
1121 |
/** |
1122 |
* Method to retrieve ID (0,1,...) of x-cluster (if any) associated to this track. |
1123 |
* If no cluster is associated, ID=-1. |
1124 |
* @param ip Tracker plane (0-5) |
1125 |
*/ |
1126 |
Int_t TrkTrack::GetClusterX_ID(int ip){ |
1127 |
return ((Int_t)fabs(xgood[ip]))%10000000-1; |
1128 |
}; |
1129 |
/** |
1130 |
* Method to retrieve ID (0-xxx) of y-cluster (if any) associated to this track. |
1131 |
* If no cluster is associated, ID=-1. |
1132 |
* @param ip Tracker plane (0-5) |
1133 |
*/ |
1134 |
Int_t TrkTrack::GetClusterY_ID(int ip){ |
1135 |
return ((Int_t)fabs(ygood[ip]))%10000000-1; |
1136 |
}; |
1137 |
|
1138 |
/** |
1139 |
* Method to retrieve the ladder (0-2, increasing x) traversed by the track on this plane. |
1140 |
* If no ladder is traversed (dead area) the metod retuns -1. |
1141 |
* @param ip Tracker plane (0-5) |
1142 |
*/ |
1143 |
Int_t TrkTrack::GetLadder(int ip){ |
1144 |
if(XGood(ip))return (Int_t)fabs(xgood[ip]/100000000)-1; |
1145 |
if(YGood(ip))return (Int_t)fabs(ygood[ip]/100000000)-1; |
1146 |
return -1; |
1147 |
}; |
1148 |
/** |
1149 |
* Method to retrieve the sensor (0-1, increasing y) traversed by the track on this plane. |
1150 |
* If no sensor is traversed (dead area) the metod retuns -1. |
1151 |
* @param ip Tracker plane (0-5) |
1152 |
*/ |
1153 |
Int_t TrkTrack::GetSensor(int ip){ |
1154 |
if(XGood(ip))return (Int_t)((Int_t)fabs(xgood[ip]/10000000)%10)-1; |
1155 |
if(YGood(ip))return (Int_t)((Int_t)fabs(ygood[ip]/10000000)%10)-1; |
1156 |
return -1; |
1157 |
}; |
1158 |
|
1159 |
/** |
1160 |
* \brief Method to include a x-cluster to the track. |
1161 |
* @param ip Tracker plane (0-5) |
1162 |
* @param clid Cluster ID (0 = no-cluster, 1,2,... otherwise ) |
1163 |
* @param il Ladder (0-2, increasing x, -1 if no sensitive area is hit) |
1164 |
* @param is Sensor (0-1, increasing y, -1 if no sensitive area is hit) |
1165 |
* @param bad True if the cluster contains bad strips |
1166 |
* @see Fit(double pfixed, int& fail, int iprint, int froml1) |
1167 |
*/ |
1168 |
void TrkTrack::SetXGood(int ip, int clid, int il, int is, bool bad){ |
1169 |
// int il=0; //ladder (temporary) |
1170 |
// bool bad=false; //ladder (temporary) |
1171 |
if(ip<0||ip>5||clid<1||il<-1||il>2||is<-1||is>1) |
1172 |
cout << " void TrkTrack::SetXGood(int,int,int,int,bool) --> MA SEI DI COCCIO?!?!"<<endl; |
1173 |
xgood[ip]=(il+1)*100000000+(is+1)*10000000+clid; |
1174 |
if(bad)xgood[ip]=-xgood[ip]; |
1175 |
}; |
1176 |
/** |
1177 |
* \brief Method to include a y-cluster to the track. |
1178 |
* @param ip Tracker plane (0-5) |
1179 |
* @param clid Cluster ID (0 = no-cluster, 1,2,... otherwise ) |
1180 |
* @param il Ladder (0-2, increasing x, -1 if no sensitive area is hit) |
1181 |
* @param is Sensor (0-1, increasing y, -1 if no sensitive area is hit) |
1182 |
* @param bad True if the cluster contains bad strips |
1183 |
* @see Fit(double pfixed, int& fail, int iprint, int froml1) |
1184 |
*/ |
1185 |
void TrkTrack::SetYGood(int ip, int clid, int il, int is, bool bad){ |
1186 |
// int il=0; //ladder (temporary) |
1187 |
// bool bad=false; //ladder (temporary) |
1188 |
if(ip<0||ip>5||clid<1||il<-1||il>2||is<-1||is>1) |
1189 |
cout << " void TrkTrack::SetYGood(int,int,int,int,bool) --> MA SEI DI COCCIO?!?!"<<endl; |
1190 |
ygood[ip]=(il+1)*100000000+(is+1)*10000000+clid; |
1191 |
if(bad)ygood[ip]=-ygood[ip]; |
1192 |
}; |
1193 |
|
1194 |
/** |
1195 |
* \brief Average X |
1196 |
* Average value of <xv>, evaluated from the first to the last hit x view. |
1197 |
*/ |
1198 |
Float_t TrkTrack::GetXav(){ |
1199 |
|
1200 |
int first_plane = -1; |
1201 |
int last_plane = -1; |
1202 |
for(Int_t ip=0; ip<6; ip++){ |
1203 |
if( XGood(ip) && first_plane == -1 )first_plane = ip; |
1204 |
if( XGood(ip) && first_plane != -1 )last_plane = ip; |
1205 |
} |
1206 |
if( first_plane == -1 || last_plane == -1){ |
1207 |
return -100; |
1208 |
} |
1209 |
if( last_plane-first_plane+1 ==0 )return -100; |
1210 |
|
1211 |
Float_t av = 0; |
1212 |
for(int ip=first_plane; ip<=last_plane; ip++)av+=xv[ip]; |
1213 |
|
1214 |
return (av/(last_plane-first_plane+1)); |
1215 |
} |
1216 |
/** |
1217 |
* \brief Average Y |
1218 |
* Average value of <yv>, evaluated from the first to the last hit x view. |
1219 |
*/ |
1220 |
Float_t TrkTrack::GetYav(){ |
1221 |
|
1222 |
int first_plane = -1; |
1223 |
int last_plane = -1; |
1224 |
for(Int_t ip=0; ip<6; ip++){ |
1225 |
if( XGood(ip) && first_plane == -1 )first_plane = ip; |
1226 |
if( XGood(ip) && first_plane != -1 )last_plane = ip; |
1227 |
} |
1228 |
if( first_plane == -1 || last_plane == -1){ |
1229 |
return -100; |
1230 |
} |
1231 |
if( last_plane-first_plane+1 ==0 )return -100; |
1232 |
|
1233 |
Float_t av = 0; |
1234 |
for(int ip=first_plane; ip<=last_plane; ip++)av+=yv[ip]; |
1235 |
|
1236 |
return (av/(last_plane-first_plane+1)); |
1237 |
} |
1238 |
/** |
1239 |
* \brief Average Z |
1240 |
* Average value of <zv>, evaluated from the first to the last hit x view. |
1241 |
*/ |
1242 |
Float_t TrkTrack::GetZav(){ |
1243 |
|
1244 |
int first_plane = -1; |
1245 |
int last_plane = -1; |
1246 |
for(Int_t ip=0; ip<6; ip++){ |
1247 |
if( XGood(ip) && first_plane == -1 )first_plane = ip; |
1248 |
if( XGood(ip) && first_plane != -1 )last_plane = ip; |
1249 |
} |
1250 |
if( first_plane == -1 || last_plane == -1){ |
1251 |
return -100; |
1252 |
} |
1253 |
if( last_plane-first_plane+1 ==0 )return -100; |
1254 |
|
1255 |
Float_t av = 0; |
1256 |
for(int ip=first_plane; ip<=last_plane; ip++)av+=zv[ip]; |
1257 |
|
1258 |
return (av/(last_plane-first_plane+1)); |
1259 |
} |
1260 |
|
1261 |
/** |
1262 |
* \brief Number of column traversed |
1263 |
*/ |
1264 |
Int_t TrkTrack::GetNColumns(){ |
1265 |
int sensors[] = {0,0,0,0,0,0}; |
1266 |
for(int ip=0; ip<6; ip++){ |
1267 |
int sensorid = GetLadder(ip)+3*GetSensor(ip); |
1268 |
if(XGood(ip)||YGood(ip)) |
1269 |
if(sensorid>=0 && sensorid<6)sensors[sensorid]=1; |
1270 |
} |
1271 |
int nsensors=0; |
1272 |
for(int is=0; is<6; is++)nsensors += sensors[is]; |
1273 |
return nsensors; |
1274 |
}; |
1275 |
/** |
1276 |
* \brief Give the maximum energy release |
1277 |
*/ |
1278 |
Float_t TrkTrack::GetDEDX_max(int ip, int iv){ |
1279 |
Float_t max=0; |
1280 |
int pfrom = 0; |
1281 |
int pto = 6; |
1282 |
int vfrom = 0; |
1283 |
int vto = 2; |
1284 |
if(ip>=0&&ip<6){ |
1285 |
pfrom = ip; |
1286 |
pto = ip+1; |
1287 |
} |
1288 |
if(iv>=0&&iv<2){ |
1289 |
vfrom = iv; |
1290 |
vto = iv+1; |
1291 |
} |
1292 |
for(int i=pfrom; i<pto; i++) |
1293 |
for(int j=vfrom; j<vto; j++){ |
1294 |
if(j==0 && XGood(i) && GetDEDX(i,j)>max)max=GetDEDX(i,j); |
1295 |
if(j==1 && YGood(i) && GetDEDX(i,j)>max)max=GetDEDX(i,j); |
1296 |
} |
1297 |
return max; |
1298 |
|
1299 |
}; |
1300 |
|
1301 |
/** |
1302 |
* \brief Give the minimum energy release |
1303 |
*/ |
1304 |
Float_t TrkTrack::GetDEDX_min(int ip, int iv){ |
1305 |
Float_t min=100000000; |
1306 |
int pfrom = 0; |
1307 |
int pto = 6; |
1308 |
int vfrom = 0; |
1309 |
int vto = 2; |
1310 |
if(ip>=0&&ip<6){ |
1311 |
pfrom = ip; |
1312 |
pto = ip+1; |
1313 |
} |
1314 |
if(iv>=0&&iv<2){ |
1315 |
vfrom = iv; |
1316 |
vto = iv+1; |
1317 |
} |
1318 |
for(int i=pfrom; i<pto; i++) |
1319 |
for(int j=vfrom; j<vto; j++){ |
1320 |
if(j==0 && XGood(i) && GetDEDX(i,j)<min)min=GetDEDX(i,j); |
1321 |
if(j==1 && YGood(i) && GetDEDX(i,j)<min)min=GetDEDX(i,j); |
1322 |
} |
1323 |
return min; |
1324 |
|
1325 |
}; |
1326 |
|
1327 |
/** |
1328 |
* \brief Give the maximum spatial residual |
1329 |
*/ |
1330 |
Float_t TrkTrack::GetResidual_max(int ip, int iv){ |
1331 |
Float_t max=0; |
1332 |
int pfrom = 0; |
1333 |
int pto = 6; |
1334 |
int vfrom = 0; |
1335 |
int vto = 2; |
1336 |
if(ip>=0&&ip<6){ |
1337 |
pfrom = ip; |
1338 |
pto = ip+1; |
1339 |
} |
1340 |
if(iv>=0&&iv<2){ |
1341 |
vfrom = iv; |
1342 |
vto = iv+1; |
1343 |
} |
1344 |
for(int i=pfrom; i<pto; i++){ |
1345 |
for(int j=vfrom; j<vto; j++){ |
1346 |
if(j==0 && XGood(i) && fabs(xm[i]-xv[i])>fabs(max))max=xm[i]-xv[i]; |
1347 |
if(j==1 && YGood(i) && fabs(ym[i]-yv[i])>fabs(max))max=ym[i]-yv[i]; |
1348 |
} |
1349 |
} |
1350 |
return max; |
1351 |
|
1352 |
}; |
1353 |
/** |
1354 |
* \brief Give the anerage spatial residual |
1355 |
*/ |
1356 |
Float_t TrkTrack::GetResidual_av(int ip, int iv){ |
1357 |
// |
1358 |
//Sum$((xm>-50)*(xm-xv)/resx)/sqrt(TrkTrack.GetNX()*TrkTrack.GetChi2X())<0.3 |
1359 |
|
1360 |
Float_t av = 0.; |
1361 |
int nav = 0; |
1362 |
// |
1363 |
int pfrom = 0; |
1364 |
int pto = 6; |
1365 |
int vfrom = 0; |
1366 |
int vto = 2; |
1367 |
if(ip>=0&&ip<6){ |
1368 |
pfrom = ip; |
1369 |
pto = ip+1; |
1370 |
} |
1371 |
if(iv>=0&&iv<2){ |
1372 |
vfrom = iv; |
1373 |
vto = iv+1; |
1374 |
} |
1375 |
for(int i=pfrom; i<pto; i++){ |
1376 |
for(int j=vfrom; j<vto; j++){ |
1377 |
nav++; |
1378 |
if(j==0 && XGood(i)) av += (xm[i]-xv[i])/resx[i]; |
1379 |
if(j==1 && YGood(i)) av += (ym[i]-yv[i])/resy[i]; |
1380 |
} |
1381 |
} |
1382 |
if(nav==0)return -100.; |
1383 |
return av/nav; |
1384 |
|
1385 |
}; |
1386 |
|
1387 |
|
1388 |
/** |
1389 |
* \brief Give the maximum multiplicity on the x view |
1390 |
*/ |
1391 |
Int_t TrkTrack::GetClusterX_Multiplicity_max(){ |
1392 |
int max=0; |
1393 |
for(int ip=0; ip<6; ip++) |
1394 |
if(GetClusterX_Multiplicity(ip)>max)max=GetClusterX_Multiplicity(ip); |
1395 |
return max; |
1396 |
}; |
1397 |
/** |
1398 |
* \brief Give the minimum multiplicity on the x view |
1399 |
*/ |
1400 |
Int_t TrkTrack::GetClusterX_Multiplicity_min(){ |
1401 |
int min=50; |
1402 |
for(int ip=0; ip<6; ip++) |
1403 |
if(GetClusterX_Multiplicity(ip)<min)min=GetClusterX_Multiplicity(ip); |
1404 |
return min; |
1405 |
}; |
1406 |
/** |
1407 |
* \brief Give the maximum multiplicity on the x view |
1408 |
*/ |
1409 |
Int_t TrkTrack::GetClusterY_Multiplicity_max(){ |
1410 |
int max=0; |
1411 |
for(int ip=0; ip<6; ip++) |
1412 |
if(GetClusterY_Multiplicity(ip)>max)max=GetClusterY_Multiplicity(ip); |
1413 |
return max; |
1414 |
}; |
1415 |
/** |
1416 |
* \brief Give the minimum multiplicity on the x view |
1417 |
*/ |
1418 |
Int_t TrkTrack::GetClusterY_Multiplicity_min(){ |
1419 |
int min=50; |
1420 |
for(int ip=0; ip<6; ip++) |
1421 |
if(GetClusterY_Multiplicity(ip)<min)min=GetClusterY_Multiplicity(ip); |
1422 |
return min; |
1423 |
}; |
1424 |
|
1425 |
/** |
1426 |
* \brief Give the minimum seed on the x view |
1427 |
*/ |
1428 |
Float_t TrkTrack::GetClusterX_Seed_min(){ |
1429 |
Float_t min=100000; |
1430 |
for(int ip=0; ip<6; ip++) |
1431 |
if(XGood(ip) && GetClusterX_Seed(ip)<min)min=GetClusterX_Seed(ip); |
1432 |
return min; |
1433 |
}; |
1434 |
/** |
1435 |
* \brief Give the minimum seed on the x view |
1436 |
*/ |
1437 |
Float_t TrkTrack::GetClusterY_Seed_min(){ |
1438 |
Float_t min=100000; |
1439 |
for(int ip=0; ip<6; ip++) |
1440 |
if(YGood(ip) && GetClusterY_Seed(ip)<min)min=GetClusterY_Seed(ip); |
1441 |
return min; |
1442 |
}; |
1443 |
|
1444 |
|
1445 |
//-------------------------------------- |
1446 |
// |
1447 |
// |
1448 |
//-------------------------------------- |
1449 |
void TrkTrack::Clear(){ |
1450 |
// cout << "TrkTrack::Clear()"<<endl; |
1451 |
seqno = -1; |
1452 |
image = -1; |
1453 |
chi2 = 0; |
1454 |
nstep = 0; |
1455 |
for(int it1=0;it1<5;it1++){ |
1456 |
al[it1] = 0; |
1457 |
for(int it2=0;it2<5;it2++)coval[it1][it2] = 0; |
1458 |
}; |
1459 |
for(int ip=0;ip<6;ip++){ |
1460 |
xgood[ip] = 0; |
1461 |
ygood[ip] = 0; |
1462 |
xm[ip] = 0; |
1463 |
ym[ip] = 0; |
1464 |
zm[ip] = 0; |
1465 |
resx[ip] = 0; |
1466 |
resy[ip] = 0; |
1467 |
tailx[ip] = 0; |
1468 |
taily[ip] = 0; |
1469 |
xv[ip] = 0; |
1470 |
yv[ip] = 0; |
1471 |
zv[ip] = 0; |
1472 |
axv[ip] = 0; |
1473 |
ayv[ip] = 0; |
1474 |
dedx_x[ip] = 0; |
1475 |
dedx_y[ip] = 0; |
1476 |
|
1477 |
}; |
1478 |
int ngf = TrkParams::nGF; |
1479 |
for(int i=0; i<ngf; i++){ |
1480 |
xGF[i] = 0.; |
1481 |
yGF[i] = 0.; |
1482 |
} |
1483 |
// if(clx)clx->Clear(); |
1484 |
// if(cly)cly->Clear(); |
1485 |
// clx.Clear(); |
1486 |
// cly.Clear(); |
1487 |
}; |
1488 |
//-------------------------------------- |
1489 |
// |
1490 |
// |
1491 |
//-------------------------------------- |
1492 |
void TrkTrack::Delete(){ |
1493 |
// cout << "TrkTrack::Delete()"<<endl; |
1494 |
Clear(); |
1495 |
// if(clx)delete clx; |
1496 |
// if(cly)delete cly; |
1497 |
}; |
1498 |
//-------------------------------------- |
1499 |
// |
1500 |
// |
1501 |
//-------------------------------------- |
1502 |
|
1503 |
//-------------------------------------- |
1504 |
// |
1505 |
// |
1506 |
//-------------------------------------- |
1507 |
TrkSinglet::TrkSinglet(){ |
1508 |
// cout << "TrkSinglet::TrkSinglet() " << GetUniqueID()<<endl; |
1509 |
// plane = 0; |
1510 |
// coord[0] = 0; |
1511 |
// coord[1] = 0; |
1512 |
// sgnl = 0; |
1513 |
// multmax = 0; |
1514 |
// cls = 0; |
1515 |
Clear(); |
1516 |
}; |
1517 |
//-------------------------------------- |
1518 |
// |
1519 |
// |
1520 |
//-------------------------------------- |
1521 |
TrkSinglet::TrkSinglet(const TrkSinglet& s){ |
1522 |
// cout << "TrkSinglet::TrkSinglet(const TrkSinglet& s) " << GetUniqueID()<<endl; |
1523 |
plane = s.plane; |
1524 |
coord[0] = s.coord[0]; |
1525 |
coord[1] = s.coord[1]; |
1526 |
sgnl = s.sgnl; |
1527 |
multmax = s.multmax; |
1528 |
// cls = 0;//<<<<pointer |
1529 |
// cls = TRef(s.cls); |
1530 |
}; |
1531 |
//-------------------------------------- |
1532 |
// |
1533 |
// |
1534 |
//-------------------------------------- |
1535 |
void TrkSinglet::Dump(){ |
1536 |
int i=0; |
1537 |
cout << endl << "========== Singlet " ; |
1538 |
cout << endl << "plane : " << plane; |
1539 |
cout << endl << "coord[2] : "; while( i<2 && cout << coord[i] << " ") i++; |
1540 |
cout << endl << "sgnl : " << sgnl; |
1541 |
cout << endl << "max.strip : "; |
1542 |
cout << endl << "multiplicity : "; |
1543 |
} |
1544 |
//-------------------------------------- |
1545 |
// |
1546 |
// |
1547 |
//-------------------------------------- |
1548 |
void TrkSinglet::Clear(){ |
1549 |
// cout << "TrkSinglet::Clear() " << GetUniqueID()<<endl; |
1550 |
// cls=0; |
1551 |
plane=-1; |
1552 |
coord[0]=-999; |
1553 |
coord[1]=-999; |
1554 |
sgnl=0; |
1555 |
multmax = 0; |
1556 |
|
1557 |
} |
1558 |
//-------------------------------------- |
1559 |
// |
1560 |
// |
1561 |
//-------------------------------------- |
1562 |
TrkLevel2::TrkLevel2(){ |
1563 |
// cout <<"TrkLevel2::TrkLevel2()"<<endl; |
1564 |
for(Int_t i=0; i<12 ; i++){ |
1565 |
good[i] = -1; |
1566 |
VKmask[i] = 0; |
1567 |
VKflag[i] = 0; |
1568 |
}; |
1569 |
Track = 0; |
1570 |
SingletX = 0; |
1571 |
SingletY = 0; |
1572 |
|
1573 |
} |
1574 |
//-------------------------------------- |
1575 |
// |
1576 |
// |
1577 |
//-------------------------------------- |
1578 |
void TrkLevel2::Set(){ |
1579 |
if(!Track)Track = new TClonesArray("TrkTrack"); |
1580 |
if(!SingletX)SingletX = new TClonesArray("TrkSinglet"); |
1581 |
if(!SingletY)SingletY = new TClonesArray("TrkSinglet"); |
1582 |
} |
1583 |
//-------------------------------------- |
1584 |
// |
1585 |
// |
1586 |
//-------------------------------------- |
1587 |
void TrkLevel2::Dump(){ |
1588 |
|
1589 |
// |
1590 |
cout << endl << endl << "=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-"; |
1591 |
cout << endl << "good : "; for(int i=0; i<12; i++) cout << hex <<" 0x"<< good[i]<<dec; |
1592 |
cout << endl << "ntrk() : " << ntrk() ; |
1593 |
cout << endl << "nclsx() : " << nclsx(); |
1594 |
cout << endl << "nclsy() : " << nclsy(); |
1595 |
if(Track){ |
1596 |
TClonesArray &t = *Track; |
1597 |
for(int i=0; i<ntrk(); i++) ((TrkTrack *)t[i])->Dump(); |
1598 |
} |
1599 |
// if(SingletX){ |
1600 |
// TClonesArray &sx = *SingletX; |
1601 |
// for(int i=0; i<nclsx(); i++) ((TrkSinglet *)sx[i])->Dump(); |
1602 |
// } |
1603 |
// if(SingletY){ |
1604 |
// TClonesArray &sy = *SingletY; |
1605 |
// for(int i=0; i<nclsy(); i++) ((TrkSinglet *)sy[i])->Dump(); |
1606 |
// } |
1607 |
cout << endl; |
1608 |
} |
1609 |
/** |
1610 |
* \brief Dump processing status |
1611 |
*/ |
1612 |
void TrkLevel2::StatusDump(int view){ |
1613 |
cout << "DSP n. "<<view+1<<" status: "<<hex<<good[view]<<endl; |
1614 |
}; |
1615 |
/** |
1616 |
* \brief Check event status |
1617 |
* |
1618 |
* Check the event status, according to a flag-mask given as input. |
1619 |
* Return true if the view passes the check. |
1620 |
* |
1621 |
* @param view View number (0-11) |
1622 |
* @param flagmask Mask of flags to check (eg. flagmask=0x111 no missing packet, |
1623 |
* no crc error, no software alarm) |
1624 |
* |
1625 |
* @see TrkLevel2 class definition to know how the status flag is defined |
1626 |
* |
1627 |
*/ |
1628 |
Bool_t TrkLevel2::StatusCheck(int view, int flagmask){ |
1629 |
|
1630 |
if( view<0 || view >= 12)return false; |
1631 |
return !(good[view]&flagmask); |
1632 |
|
1633 |
}; |
1634 |
|
1635 |
|
1636 |
//-------------------------------------- |
1637 |
// |
1638 |
// |
1639 |
//-------------------------------------- |
1640 |
/** |
1641 |
* The method returns false if the viking-chip was masked |
1642 |
* either apriori ,on the basis of the mask read from the DB, |
1643 |
* or run-by-run, on the basis of the calibration parameters) |
1644 |
* @param iv Tracker view (0-11) |
1645 |
* @param ivk Viking-chip number (0-23) |
1646 |
*/ |
1647 |
Bool_t TrkLevel2::GetVKMask(int iv, int ivk){ |
1648 |
Int_t whichbit = (Int_t)pow(2,ivk); |
1649 |
return (whichbit&VKmask[iv])!=0; |
1650 |
} |
1651 |
/** |
1652 |
* The method returns false if the viking-chip was masked |
1653 |
* for this event due to common-noise computation failure. |
1654 |
* @param iv Tracker view (0-11) |
1655 |
* @param ivk Viking-chip number (0-23) |
1656 |
*/ |
1657 |
Bool_t TrkLevel2::GetVKFlag(int iv, int ivk){ |
1658 |
Int_t whichbit = (Int_t)pow(2,ivk); |
1659 |
return (whichbit&VKflag[iv])!=0; |
1660 |
} |
1661 |
/** |
1662 |
* The method returns true if the viking-chip was masked, either |
1663 |
* forced (see TrkLevel2::GetVKMask(int,int)) or |
1664 |
* for this event only (TrkLevel2::GetVKFlag(int,int)). |
1665 |
* @param iv Tracker view (0-11) |
1666 |
* @param ivk Viking-chip number (0-23) |
1667 |
*/ |
1668 |
Bool_t TrkLevel2::IsMaskedVK(int iv, int ivk){ |
1669 |
return !(GetVKMask(iv,ivk)&&GetVKFlag(iv,ivk) ); |
1670 |
}; |
1671 |
|
1672 |
//-------------------------------------- |
1673 |
// |
1674 |
// |
1675 |
//-------------------------------------- |
1676 |
/** |
1677 |
* Fills a TrkLevel2 object with values from a struct cTrkLevel2 (to get data from F77 common). |
1678 |
* Ref to Level1 data (clusters) is also set. If l1==NULL no references are set. |
1679 |
* (NB It make sense to set references only if events are stored in a tree that contains also the Level1 branch) |
1680 |
*/ |
1681 |
void TrkLevel2::SetFromLevel2Struct(cTrkLevel2 *l2, TrkLevel1 *l1){ |
1682 |
|
1683 |
// cout << "void TrkLevel2::SetFromLevel2Struct(cTrkLevel2 *l2, TrkLevel1 *l1)"<<endl; |
1684 |
Clear(); |
1685 |
|
1686 |
// temporary objects: |
1687 |
TrkSinglet* t_singlet = new TrkSinglet(); |
1688 |
TrkTrack* t_track = new TrkTrack(); |
1689 |
|
1690 |
// ----------------- |
1691 |
// general variables |
1692 |
// ----------------- |
1693 |
for(Int_t i=0; i<12 ; i++){ |
1694 |
good[i] = l2->good[i]; |
1695 |
VKmask[i]=0; |
1696 |
VKflag[i]=0; |
1697 |
for(Int_t ii=0; ii<24 ; ii++){ |
1698 |
Int_t setbit = (Int_t)pow(2,ii); |
1699 |
if( l2->vkflag[ii][i]!=-1 )VKmask[i]=VKmask[i]|setbit; |
1700 |
if( l2->vkflag[ii][i]!=0 )VKflag[i]=VKflag[i]|setbit; |
1701 |
}; |
1702 |
}; |
1703 |
// -------------- |
1704 |
// *** TRACKS *** |
1705 |
// -------------- |
1706 |
if(!Track) Track = new TClonesArray("TrkTrack"); |
1707 |
TClonesArray &t = *Track; |
1708 |
|
1709 |
for(int i=0; i<l2->ntrk; i++){ |
1710 |
t_track->seqno = i;// NBNBNBNB deve sempre essere = i |
1711 |
t_track->image = l2->image[i]-1; |
1712 |
t_track->chi2 = l2->chi2_nt[i]; |
1713 |
t_track->nstep = l2->nstep_nt[i]; |
1714 |
for(int it1=0;it1<5;it1++){ |
1715 |
t_track->al[it1] = l2->al_nt[i][it1]; |
1716 |
for(int it2=0;it2<5;it2++) |
1717 |
t_track->coval[it1][it2] = l2->coval[i][it2][it1]; |
1718 |
}; |
1719 |
for(int ip=0;ip<6;ip++){ |
1720 |
// --------------------------------- |
1721 |
// new implementation of xgood/ygood |
1722 |
// --------------------------------- |
1723 |
t_track->xgood[ip] = l2->cltrx[i][ip]; //cluster ID |
1724 |
t_track->ygood[ip] = l2->cltry[i][ip]; //cluster ID |
1725 |
t_track->xgood[ip] += 10000000*l2->ls[i][ip]; // ladder+sensor |
1726 |
t_track->ygood[ip] += 10000000*l2->ls[i][ip]; // ladder+sensor |
1727 |
if(l2->xbad[i][ip]>0)t_track->xgood[ip]=-t_track->xgood[ip]; |
1728 |
if(l2->ybad[i][ip]>0)t_track->ygood[ip]=-t_track->ygood[ip]; |
1729 |
// if(l2->xbad[i][ip]>0 || l2->ybad[i][ip]>0){ |
1730 |
// if(l2->dedx_x[i][ip]<0 || l2->dedx_y[i][ip]<0){ |
1731 |
// cout << ip << " - "<< l2->cltrx[i][ip] << " "<<l2->cltry[i][ip]<<" "<<l2->ls[i][ip]<<endl; |
1732 |
// cout << ip << " - "<<t_track->xgood[ip]<<" "<<t_track->ygood[ip]<<endl; |
1733 |
// cout << ip << " - "<<t_track->GetClusterX_ID(ip)<<" "<<t_track->GetClusterY_ID(ip)<<" "<<t_track->GetLadder(ip)<<" "<<t_track->GetSensor(ip)<<endl; |
1734 |
// cout << ip << " - "<<t_track->BadClusterX(ip)<<" "<<t_track->BadClusterY(ip)<<endl; |
1735 |
// cout << ip << " - "<<t_track->SaturatedClusterX(ip)<<" "<<t_track->SaturatedClusterY(ip)<<endl; |
1736 |
// } |
1737 |
t_track->xm[ip] = l2->xm_nt[i][ip]; |
1738 |
t_track->ym[ip] = l2->ym_nt[i][ip]; |
1739 |
t_track->zm[ip] = l2->zm_nt[i][ip]; |
1740 |
t_track->resx[ip] = l2->resx_nt[i][ip]; |
1741 |
t_track->resy[ip] = l2->resy_nt[i][ip]; |
1742 |
t_track->tailx[ip] = l2->tailx[i][ip]; |
1743 |
t_track->taily[ip] = l2->taily[i][ip]; |
1744 |
t_track->xv[ip] = l2->xv_nt[i][ip]; |
1745 |
t_track->yv[ip] = l2->yv_nt[i][ip]; |
1746 |
t_track->zv[ip] = l2->zv_nt[i][ip]; |
1747 |
t_track->axv[ip] = l2->axv_nt[i][ip]; |
1748 |
t_track->ayv[ip] = l2->ayv_nt[i][ip]; |
1749 |
t_track->dedx_x[ip] = l2->dedx_x[i][ip]; |
1750 |
t_track->dedx_y[ip] = l2->dedx_y[i][ip]; |
1751 |
t_track->multmaxx[ip] = l2->multmaxx[i][ip]; |
1752 |
t_track->multmaxy[ip] = l2->multmaxy[i][ip]; |
1753 |
t_track->seedx[ip] = l2->seedx[i][ip]; |
1754 |
t_track->seedy[ip] = l2->seedy[i][ip]; |
1755 |
t_track->xpu[ip] = l2->xpu[i][ip]; |
1756 |
t_track->ypu[ip] = l2->ypu[i][ip]; |
1757 |
//----------------------------------------------------- |
1758 |
//----------------------------------------------------- |
1759 |
//----------------------------------------------------- |
1760 |
//----------------------------------------------------- |
1761 |
}; |
1762 |
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
1763 |
// evaluated coordinates (to define GF) |
1764 |
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
1765 |
int ngf = TrkParams::nGF; |
1766 |
float *zgf = TrkParams::zGF; |
1767 |
Trajectory tgf = Trajectory(ngf,zgf); |
1768 |
tgf.DoTrack(t_track->al);//<<<< integrate the trajectory |
1769 |
for(int ip=0; ip<ngf; ip++){ |
1770 |
t_track->xGF[ip] = tgf.x[ip]; |
1771 |
t_track->yGF[ip] = tgf.y[ip]; |
1772 |
} |
1773 |
|
1774 |
// if(t_track->IsSaturated())t_track->Dump(); |
1775 |
new(t[i]) TrkTrack(*t_track); |
1776 |
t_track->Clear(); |
1777 |
};//end loop over track |
1778 |
|
1779 |
// ---------------- |
1780 |
// *** SINGLETS *** |
1781 |
// ---------------- |
1782 |
if(!SingletX)SingletX = new TClonesArray("TrkSinglet"); |
1783 |
TClonesArray &sx = *SingletX; |
1784 |
for(int i=0; i<l2->nclsx; i++){ |
1785 |
t_singlet->plane = l2->planex[i]; |
1786 |
t_singlet->coord[0] = l2->xs[i][0]; |
1787 |
t_singlet->coord[1] = l2->xs[i][1]; |
1788 |
t_singlet->sgnl = l2->signlxs[i]; |
1789 |
t_singlet->multmax = l2->multmaxsx[i]; |
1790 |
if(l2->sxbad[i]>0) t_singlet->multmax = -1*t_singlet->multmax; |
1791 |
//----------------------------------------------------- |
1792 |
// if(l1) t_singlet->cls = l1->GetCluster(l2->clsx[i]-1); |
1793 |
//----------------------------------------------------- |
1794 |
new(sx[i]) TrkSinglet(*t_singlet); |
1795 |
t_singlet->Clear(); |
1796 |
} |
1797 |
if(!SingletY)SingletY = new TClonesArray("TrkSinglet"); |
1798 |
TClonesArray &sy = *SingletY; |
1799 |
for(int i=0; i<l2->nclsy; i++){ |
1800 |
t_singlet->plane = l2->planey[i]; |
1801 |
t_singlet->coord[0] = l2->ys[i][0]; |
1802 |
t_singlet->coord[1] = l2->ys[i][1]; |
1803 |
t_singlet->sgnl = l2->signlys[i]; |
1804 |
t_singlet->multmax = l2->multmaxsy[i]; |
1805 |
if(l2->sybad[i]>0) t_singlet->multmax = -1*t_singlet->multmax; |
1806 |
//----------------------------------------------------- |
1807 |
// if(l1) t_singlet->cls = l1->GetCluster(l2->clsy[i]-1); |
1808 |
//----------------------------------------------------- |
1809 |
new(sy[i]) TrkSinglet(*t_singlet); |
1810 |
t_singlet->Clear(); |
1811 |
}; |
1812 |
|
1813 |
|
1814 |
|
1815 |
delete t_track; |
1816 |
delete t_singlet; |
1817 |
} |
1818 |
/** |
1819 |
* Fills a struct cTrkLevel2 with values from a TrkLevel2 object (to put data into a F77 common). |
1820 |
*/ |
1821 |
|
1822 |
void TrkLevel2::GetLevel2Struct(cTrkLevel2 *l2) const { |
1823 |
|
1824 |
// general variables |
1825 |
// l2->good2 = good2 ; |
1826 |
for(Int_t i=0; i<12 ; i++){ |
1827 |
// l2->crc[i] = crc[i]; |
1828 |
l2->good[i] = good[i]; |
1829 |
}; |
1830 |
// *** TRACKS *** |
1831 |
|
1832 |
if(Track){ |
1833 |
l2->ntrk = Track->GetEntries(); |
1834 |
for(Int_t i=0;i<l2->ntrk;i++){ |
1835 |
l2->image[i] = 1 + ((TrkTrack *)Track->At(i))->image; |
1836 |
l2->chi2_nt[i] = ((TrkTrack *)Track->At(i))->chi2; |
1837 |
l2->nstep_nt[i] = ((TrkTrack *)Track->At(i))->nstep; |
1838 |
for(int it1=0;it1<5;it1++){ |
1839 |
l2->al_nt[i][it1] = ((TrkTrack *)Track->At(i))->al[it1]; |
1840 |
for(int it2=0;it2<5;it2++) |
1841 |
l2->coval[i][it2][it1] = ((TrkTrack *)Track->At(i))->coval[it1][it2]; |
1842 |
}; |
1843 |
for(int ip=0;ip<6;ip++){ |
1844 |
l2->xgood_nt[i][ip] = ((TrkTrack *)Track->At(i))->XGood(ip); |
1845 |
l2->ygood_nt[i][ip] = ((TrkTrack *)Track->At(i))->YGood(ip); |
1846 |
l2->xm_nt[i][ip] = ((TrkTrack *)Track->At(i))->xm[ip]; |
1847 |
l2->ym_nt[i][ip] = ((TrkTrack *)Track->At(i))->ym[ip]; |
1848 |
l2->zm_nt[i][ip] = ((TrkTrack *)Track->At(i))->zm[ip]; |
1849 |
l2->resx_nt[i][ip] = ((TrkTrack *)Track->At(i))->resx[ip]; |
1850 |
l2->resy_nt[i][ip] = ((TrkTrack *)Track->At(i))->resy[ip]; |
1851 |
l2->tailx[i][ip] = ((TrkTrack *)Track->At(i))->tailx[ip]; |
1852 |
l2->taily[i][ip] = ((TrkTrack *)Track->At(i))->taily[ip]; |
1853 |
l2->xv_nt[i][ip] = ((TrkTrack *)Track->At(i))->xv[ip]; |
1854 |
l2->yv_nt[i][ip] = ((TrkTrack *)Track->At(i))->yv[ip]; |
1855 |
l2->zv_nt[i][ip] = ((TrkTrack *)Track->At(i))->zv[ip]; |
1856 |
l2->axv_nt[i][ip] = ((TrkTrack *)Track->At(i))->axv[ip]; |
1857 |
l2->ayv_nt[i][ip] = ((TrkTrack *)Track->At(i))->ayv[ip]; |
1858 |
l2->dedx_x[i][ip] = ((TrkTrack *)Track->At(i))->dedx_x[ip]; |
1859 |
l2->dedx_y[i][ip] = ((TrkTrack *)Track->At(i))->dedx_y[ip]; |
1860 |
}; |
1861 |
} |
1862 |
} |
1863 |
// *** SINGLETS *** |
1864 |
if(SingletX){ |
1865 |
l2->nclsx = SingletX->GetEntries(); |
1866 |
for(Int_t i=0;i<l2->nclsx;i++){ |
1867 |
l2->planex[i] = ((TrkSinglet *)SingletX->At(i))->plane; |
1868 |
l2->xs[i][0] = ((TrkSinglet *)SingletX->At(i))->coord[0]; |
1869 |
l2->xs[i][1] = ((TrkSinglet *)SingletX->At(i))->coord[1]; |
1870 |
l2->signlxs[i] = ((TrkSinglet *)SingletX->At(i))->sgnl; |
1871 |
} |
1872 |
} |
1873 |
|
1874 |
if(SingletY){ |
1875 |
l2->nclsy = SingletY->GetEntries(); |
1876 |
for(Int_t i=0;i<l2->nclsy;i++){ |
1877 |
l2->planey[i] = ((TrkSinglet *)SingletY->At(i))->plane; |
1878 |
l2->ys[i][0] = ((TrkSinglet *)SingletY->At(i))->coord[0]; |
1879 |
l2->ys[i][1] = ((TrkSinglet *)SingletY->At(i))->coord[1]; |
1880 |
l2->signlys[i] = ((TrkSinglet *)SingletY->At(i))->sgnl; |
1881 |
} |
1882 |
} |
1883 |
} |
1884 |
//-------------------------------------- |
1885 |
// |
1886 |
// |
1887 |
//-------------------------------------- |
1888 |
void TrkLevel2::Clear(){ |
1889 |
for(Int_t i=0; i<12 ; i++){ |
1890 |
good[i] = -1; |
1891 |
VKflag[i] = 0; |
1892 |
VKmask[i] = 0; |
1893 |
}; |
1894 |
// if(Track)Track->Clear("C"); |
1895 |
// if(SingletX)SingletX->Clear("C"); |
1896 |
// if(SingletY)SingletY->Clear("C"); |
1897 |
if(Track)Track->Delete(); |
1898 |
if(SingletX)SingletX->Delete(); |
1899 |
if(SingletY)SingletY->Delete(); |
1900 |
} |
1901 |
// //-------------------------------------- |
1902 |
// // |
1903 |
// // |
1904 |
// //-------------------------------------- |
1905 |
void TrkLevel2::Delete(){ |
1906 |
|
1907 |
// cout << "void TrkLevel2::Delete()"<<endl; |
1908 |
Clear(); |
1909 |
if(Track)delete Track; |
1910 |
if(SingletX)delete SingletX; |
1911 |
if(SingletY)delete SingletY; |
1912 |
|
1913 |
} |
1914 |
//-------------------------------------- |
1915 |
// |
1916 |
// |
1917 |
//-------------------------------------- |
1918 |
/** |
1919 |
* Sort physical tracks and stores them in a TObjectArray, ordering by increasing chi**2 value (in case of track image, it selects the one with lower chi**2). The total number of physical tracks is given by GetNTracks() and the it-th physical track can be retrieved by means of the method GetTrack(int it). |
1920 |
* This method is overridden by PamLevel2::GetTracks(), where calorimeter and TOF information is used. |
1921 |
*/ |
1922 |
TRefArray *TrkLevel2::GetTracks_NFitSorted(){ |
1923 |
|
1924 |
if(!Track)return 0; |
1925 |
|
1926 |
// TRefArray *sorted = new TRefArray(); |
1927 |
TRefArray *sorted = NULL; |
1928 |
|
1929 |
TClonesArray &t = *Track; |
1930 |
// TClonesArray &ts = *PhysicalTrack; |
1931 |
int N = ntrk(); |
1932 |
vector<int> m(N); for(int i=0; i<N; i++)m[i]=1; |
1933 |
// int m[50]; for(int i=0; i<N; i++)m[i]=1; |
1934 |
|
1935 |
int indo=0; |
1936 |
int indi=0; |
1937 |
while(N > 0){ |
1938 |
// while(N != 0){ |
1939 |
int nfit =0; |
1940 |
float chi2ref = numeric_limits<float>::max(); |
1941 |
|
1942 |
// first loop to search maximum num. of fit points |
1943 |
for(int i=0; i < ntrk(); i++){ |
1944 |
if( ((TrkTrack *)t[i])->GetNtot() >= nfit && m[i]==1){ |
1945 |
nfit = ((TrkTrack *)t[i])->GetNtot(); |
1946 |
} |
1947 |
} |
1948 |
//second loop to search minimum chi2 among selected |
1949 |
for(int i=0; i<ntrk(); i++){ |
1950 |
Float_t chi2 = ((TrkTrack *)t[i])->chi2; |
1951 |
if(chi2 < 0) chi2 = -chi2*1000; |
1952 |
if( chi2 < chi2ref |
1953 |
&& ((TrkTrack *)t[i])->GetNtot() == nfit |
1954 |
&& m[i]==1){ |
1955 |
chi2ref = ((TrkTrack *)t[i])->chi2; |
1956 |
indi = i; |
1957 |
}; |
1958 |
}; |
1959 |
if( ((TrkTrack *)t[indi])->HasImage() ){ |
1960 |
m[((TrkTrack *)t[indi])->image] = 0; |
1961 |
N--; |
1962 |
|
1963 |
// cout << "i** "<< ((TrkTrack *)t[indi])->image << " " << nfiti <<" "<<chi2i<<endl; |
1964 |
}; |
1965 |
if(!sorted)sorted = new TRefArray( TProcessID::GetProcessWithUID(t[indi])); |
1966 |
sorted->Add( (TrkTrack*)t[indi] ); |
1967 |
|
1968 |
m[indi] = 0; |
1969 |
// cout << "SORTED "<< indo << " "<< indi << " "<< N << " "<<((TrkTrack *)t[indi])->image<<" "<<chi2ref<<endl; |
1970 |
N--; |
1971 |
indo++; |
1972 |
} |
1973 |
m.clear(); |
1974 |
// cout << "GetTracks_NFitSorted(it): Done"<< endl; |
1975 |
|
1976 |
return sorted; |
1977 |
// return PhysicalTrack; |
1978 |
} |
1979 |
//-------------------------------------- |
1980 |
// |
1981 |
// |
1982 |
//-------------------------------------- |
1983 |
/** |
1984 |
* Retrieves the is-th stored track. |
1985 |
* @param it Track number, ranging from 0 to ntrk(). |
1986 |
* Fitted tracks ( images included ) are stored in a TObjectArray ( TrkLevel2::Track ) in the same order they are returned by the F77 fitting routine. |
1987 |
*/ |
1988 |
TrkTrack *TrkLevel2::GetStoredTrack(int is){ |
1989 |
|
1990 |
if(is >= this->ntrk()){ |
1991 |
cout << "TrkTrack *TrkLevel2::GetStoredTrack(int) >> Track "<< is << "doen not exits! " << endl; |
1992 |
cout << "Stored tracks ntrk() = "<< this->ntrk() << endl; |
1993 |
return 0; |
1994 |
} |
1995 |
if(!Track){ |
1996 |
cout << "TrkTrack *TrkLevel2::GetStoredTrack(int is) >> (TClonesArray*) Track ==0 "<<endl; |
1997 |
}; |
1998 |
TClonesArray &t = *(Track); |
1999 |
TrkTrack *track = (TrkTrack*)t[is]; |
2000 |
return track; |
2001 |
} |
2002 |
//-------------------------------------- |
2003 |
// |
2004 |
// |
2005 |
//-------------------------------------- |
2006 |
/** |
2007 |
* Retrieves the is-th stored X singlet. |
2008 |
* @param it Singlet number, ranging from 0 to nclsx(). |
2009 |
*/ |
2010 |
TrkSinglet *TrkLevel2::GetSingletX(int is){ |
2011 |
|
2012 |
if(is >= this->nclsx()){ |
2013 |
cout << "TrkSinglet *TrkLevel2::GetSingletX(int) >> Singlet "<< is << "doen not exits! " << endl; |
2014 |
cout << "Stored x-singlets nclsx() = "<< this->nclsx() << endl; |
2015 |
return 0; |
2016 |
} |
2017 |
if(!SingletX)return 0; |
2018 |
TClonesArray &t = *(SingletX); |
2019 |
TrkSinglet *singlet = (TrkSinglet*)t[is]; |
2020 |
return singlet; |
2021 |
} |
2022 |
//-------------------------------------- |
2023 |
// |
2024 |
// |
2025 |
//-------------------------------------- |
2026 |
/** |
2027 |
* Retrieves the is-th stored Y singlet. |
2028 |
* @param it Singlet number, ranging from 0 to nclsx(). |
2029 |
*/ |
2030 |
TrkSinglet *TrkLevel2::GetSingletY(int is){ |
2031 |
|
2032 |
if(is >= this->nclsy()){ |
2033 |
cout << "TrkSinglet *TrkLevel2::GetSingletY(int) >> Singlet "<< is << "doen not exits! " << endl; |
2034 |
cout << "Stored y-singlets nclsx() = "<< this->nclsx() << endl; |
2035 |
return 0; |
2036 |
} |
2037 |
if(!SingletY)return 0; |
2038 |
TClonesArray &t = *(SingletY); |
2039 |
TrkSinglet *singlet = (TrkSinglet*)t[is]; |
2040 |
return singlet; |
2041 |
} |
2042 |
//-------------------------------------- |
2043 |
// |
2044 |
// |
2045 |
//-------------------------------------- |
2046 |
/** |
2047 |
* Retrieves the it-th "physical" track, sorted by the method GetNTracks(). |
2048 |
* @param it Track number, ranging from 0 to GetNTracks(). |
2049 |
*/ |
2050 |
|
2051 |
TrkTrack *TrkLevel2::GetTrack(int it){ |
2052 |
|
2053 |
if(it >= this->GetNTracks()){ |
2054 |
cout << "TrkTrack *TrkLevel2::GetTrack(int) >> Track "<< it << "does not exits! " << endl; |
2055 |
cout << "Physical tracks GetNTracks() = "<< this->ntrk() << endl; |
2056 |
return 0; |
2057 |
} |
2058 |
|
2059 |
TRefArray *sorted = GetTracks(); //TEMPORANEO |
2060 |
if(!sorted)return 0; |
2061 |
TrkTrack *track = (TrkTrack*)sorted->At(it); |
2062 |
sorted->Clear(); |
2063 |
delete sorted; |
2064 |
return track; |
2065 |
} |
2066 |
/** |
2067 |
* Give the number of "physical" tracks, sorted by the method GetTracks(). |
2068 |
*/ |
2069 |
Int_t TrkLevel2::GetNTracks(){ |
2070 |
|
2071 |
Float_t ntot=0; |
2072 |
if(!Track)return 0; |
2073 |
TClonesArray &t = *Track; |
2074 |
for(int i=0; i<ntrk(); i++) { |
2075 |
if( ((TrkTrack *)t[i])->GetImageSeqNo() == -1 ) ntot+=1.; |
2076 |
else ntot+=0.5; |
2077 |
} |
2078 |
return (Int_t)ntot; |
2079 |
|
2080 |
}; |
2081 |
//-------------------------------------- |
2082 |
// |
2083 |
// |
2084 |
//-------------------------------------- |
2085 |
/** |
2086 |
* Retrieves (if present) the image of the it-th "physical" track, sorted by the method GetNTracks(). |
2087 |
* @param it Track number, ranging from 0 to GetNTracks(). |
2088 |
*/ |
2089 |
TrkTrack *TrkLevel2::GetTrackImage(int it){ |
2090 |
|
2091 |
if(it >= this->GetNTracks()){ |
2092 |
cout << "TrkTrack *TrkLevel2::GetTrackImage(int) >> Track "<< it << "does not exits! " << endl; |
2093 |
cout << "Physical tracks GetNTracks() = "<< this->ntrk() << endl; |
2094 |
return 0; |
2095 |
} |
2096 |
|
2097 |
TRefArray* sorted = GetTracks(); //TEMPORANEO |
2098 |
if(!sorted)return 0; |
2099 |
TrkTrack *track = (TrkTrack*)sorted->At(it); |
2100 |
|
2101 |
if(!track->HasImage()){ |
2102 |
cout << "TrkTrack *TrkLevel2::GetTrackImage(int) >> Track "<< it << "does not have image! " << endl; |
2103 |
return 0; |
2104 |
} |
2105 |
if(!Track)return 0; |
2106 |
TrkTrack *image = (TrkTrack*)(*Track)[track->image]; |
2107 |
|
2108 |
sorted->Delete(); |
2109 |
delete sorted; |
2110 |
|
2111 |
return image; |
2112 |
|
2113 |
} |
2114 |
//-------------------------------------- |
2115 |
// |
2116 |
// |
2117 |
//-------------------------------------- |
2118 |
/** |
2119 |
* Loads the magnetic field. |
2120 |
* @param s Path of the magnetic-field files. |
2121 |
*/ |
2122 |
void TrkLevel2::LoadField(TString path){ |
2123 |
// |
2124 |
// strcpy(path_.path,path.Data()); |
2125 |
// path_.pathlen = path.Length(); |
2126 |
// path_.error = 0; |
2127 |
// readb_(); |
2128 |
|
2129 |
// TrkParams::SetTrackingMode(); |
2130 |
// TrkParams::SetPrecisionFactor(); |
2131 |
// TrkParams::SetStepMin(); |
2132 |
TrkParams::SetMiniDefault(); |
2133 |
|
2134 |
TrkParams::Set(path,1); |
2135 |
TrkParams::Load(1); |
2136 |
if( !TrkParams::IsLoaded(1) ){ |
2137 |
cout << "void TrkLevel2::LoadField(TString path) --- ERROR --- m.field not loaded"<<endl; |
2138 |
} |
2139 |
|
2140 |
// |
2141 |
}; |
2142 |
// /** |
2143 |
// * Get BY (kGauss) |
2144 |
// * @param v (x,y,z) coordinates in cm |
2145 |
// */ |
2146 |
// float TrkLevel2::GetBX(float* v){ |
2147 |
// float b[3]; |
2148 |
// gufld_(v,b); |
2149 |
// return b[0]/10.; |
2150 |
// } |
2151 |
// /** |
2152 |
// * Get BY (kGauss) |
2153 |
// * @param v (x,y,z) coordinates in cm |
2154 |
// */ |
2155 |
// float TrkLevel2::GetBY(float* v){ |
2156 |
// float b[3]; |
2157 |
// gufld_(v,b); |
2158 |
// return b[1]/10.; |
2159 |
// } |
2160 |
// /** |
2161 |
// * Get BY (kGauss) |
2162 |
// * @param v (x,y,z) coordinates in cm |
2163 |
// */ |
2164 |
// float TrkLevel2::GetBZ(float* v){ |
2165 |
// float b[3]; |
2166 |
// gufld_(v,b); |
2167 |
// return b[2]/10.; |
2168 |
// } |
2169 |
//-------------------------------------- |
2170 |
// |
2171 |
// |
2172 |
//-------------------------------------- |
2173 |
/** |
2174 |
* Get tracker-plane (mechanical) z-coordinate |
2175 |
* @param plane_id plane index (1=TOP,2,3,4,5,6=BOTTOM) |
2176 |
*/ |
2177 |
Float_t TrkLevel2::GetZTrk(Int_t plane_id){ |
2178 |
switch(plane_id){ |
2179 |
case 1: return ZTRK1; |
2180 |
case 2: return ZTRK2; |
2181 |
case 3: return ZTRK3; |
2182 |
case 4: return ZTRK4; |
2183 |
case 5: return ZTRK5; |
2184 |
case 6: return ZTRK6; |
2185 |
default: return 0.; |
2186 |
}; |
2187 |
}; |
2188 |
//-------------------------------------- |
2189 |
// |
2190 |
// |
2191 |
//-------------------------------------- |
2192 |
/** |
2193 |
* Trajectory default constructor. |
2194 |
* (By default is created with z-coordinates inside the tracking volume) |
2195 |
*/ |
2196 |
Trajectory::Trajectory(){ |
2197 |
npoint = 6; |
2198 |
x = new float[npoint]; |
2199 |
y = new float[npoint]; |
2200 |
z = new float[npoint]; |
2201 |
thx = new float[npoint]; |
2202 |
thy = new float[npoint]; |
2203 |
tl = new float[npoint]; |
2204 |
float dz = ((ZTRK1)-(ZTRK6))/(npoint-1); |
2205 |
for(int i=0; i<npoint; i++){ |
2206 |
x[i] = 0; |
2207 |
y[i] = 0; |
2208 |
z[i] = (ZTRK1) - i*dz; |
2209 |
thx[i] = 0; |
2210 |
thy[i] = 0; |
2211 |
tl[i] = 0; |
2212 |
} |
2213 |
} |
2214 |
//-------------------------------------- |
2215 |
// |
2216 |
// |
2217 |
//-------------------------------------- |
2218 |
/** |
2219 |
* Trajectory constructor. |
2220 |
* (By default is created with z-coordinates inside the tracking volume) |
2221 |
* \param n Number of points |
2222 |
*/ |
2223 |
Trajectory::Trajectory(int n){ |
2224 |
if(n<=0){ |
2225 |
cout << "NB! Trajectory must have at least 1 point >>> created with 10 points" << endl; |
2226 |
n=10; |
2227 |
} |
2228 |
npoint = n; |
2229 |
x = new float[npoint]; |
2230 |
y = new float[npoint]; |
2231 |
z = new float[npoint]; |
2232 |
thx = new float[npoint]; |
2233 |
thy = new float[npoint]; |
2234 |
tl = new float[npoint]; |
2235 |
float dz = ((ZTRK1)-(ZTRK6))/(npoint-1); |
2236 |
for(int i=0; i<npoint; i++){ |
2237 |
x[i] = 0; |
2238 |
y[i] = 0; |
2239 |
z[i] = (ZTRK1) - i*dz; |
2240 |
thx[i] = 0; |
2241 |
thy[i] = 0; |
2242 |
tl[i] = 0; |
2243 |
} |
2244 |
} |
2245 |
//-------------------------------------- |
2246 |
// |
2247 |
// |
2248 |
//-------------------------------------- |
2249 |
/** |
2250 |
* Trajectory constructor. |
2251 |
* \param n Number of points |
2252 |
* \param pz Pointer to float array, defining z coordinates |
2253 |
*/ |
2254 |
Trajectory::Trajectory(int n, float* zin){ |
2255 |
npoint = 10; |
2256 |
if(n>0)npoint = n; |
2257 |
x = new float[npoint]; |
2258 |
y = new float[npoint]; |
2259 |
z = new float[npoint]; |
2260 |
thx = new float[npoint]; |
2261 |
thy = new float[npoint]; |
2262 |
tl = new float[npoint]; |
2263 |
int i=0; |
2264 |
do{ |
2265 |
x[i] = 0.; |
2266 |
y[i] = 0.; |
2267 |
z[i] = zin[i]; |
2268 |
thx[i] = 0.; |
2269 |
thy[i] = 0.; |
2270 |
tl[i] = 0.; |
2271 |
i++; |
2272 |
}while(zin[i-1] > zin[i] && i < npoint); |
2273 |
npoint=i; |
2274 |
if(npoint != n)cout << "NB! Trajectory created with "<<npoint<<" points instean of "<<n<<endl; |
2275 |
// Dump(); |
2276 |
} |
2277 |
void Trajectory::Delete(){ |
2278 |
|
2279 |
if(x) delete [] x; |
2280 |
if(y) delete [] y; |
2281 |
if(z) delete [] z; |
2282 |
if(thx) delete [] thx; |
2283 |
if(thy) delete [] thy; |
2284 |
if(tl) delete [] tl; |
2285 |
|
2286 |
} |
2287 |
//-------------------------------------- |
2288 |
// |
2289 |
// |
2290 |
//-------------------------------------- |
2291 |
/** |
2292 |
* Dump the trajectory coordinates. |
2293 |
*/ |
2294 |
void Trajectory::Dump(){ |
2295 |
cout <<endl<< "Trajectory ========== "<<endl; |
2296 |
for (int i=0; i<npoint; i++){ |
2297 |
cout << i <<" >> " << x[i] <<" "<< y[i] <<" "<< z[i] ; |
2298 |
cout <<" -- " << thx[i] <<" "<< thy[i] ; |
2299 |
cout <<" -- " << tl[i] << endl; |
2300 |
}; |
2301 |
} |
2302 |
//-------------------------------------- |
2303 |
// |
2304 |
// |
2305 |
//-------------------------------------- |
2306 |
/** |
2307 |
* Get trajectory length between two points |
2308 |
* @param ifirst first point (default 0) |
2309 |
* @param ilast last point (default npoint) |
2310 |
*/ |
2311 |
float Trajectory::GetLength(int ifirst, int ilast){ |
2312 |
if( ifirst<0 ) ifirst = 0; |
2313 |
if( ilast>=npoint) ilast = npoint-1; |
2314 |
float l=0; |
2315 |
for(int i=ifirst;i<=ilast;i++){ |
2316 |
l=l+tl[i]; |
2317 |
}; |
2318 |
if(z[ilast] > ZINI)l=l-tl[ilast]; |
2319 |
if(z[ifirst] < ZINI) l=l-tl[ifirst]; |
2320 |
|
2321 |
return l; |
2322 |
|
2323 |
} |
2324 |
|
2325 |
/** |
2326 |
* Evaluates the trajectory in the apparatus associated to the track. |
2327 |
* It integrates the equations of motion in the magnetic field. |
2328 |
* @param al Track state-vector (X0,Y0,sin(theta),phi,deflection). |
2329 |
* @param zini z-coordinate of the reference plane (Z0). |
2330 |
* @return error flag. |
2331 |
* |
2332 |
* This method is needed when you want to integrate the particle trajectory |
2333 |
* starting from a track state-vector relative to an arbitrary reference plane. |
2334 |
* The default reference plane, used by the tracker routines, is at zini=23.5. |
2335 |
* If you give as input the track state-vector from a TrkTrack object, |
2336 |
* you can use Trajectory::DoTrack(float* al) instead. |
2337 |
*/ |
2338 |
int Trajectory::DoTrack(float* al, float zini){ |
2339 |
|
2340 |
// double *dxout = new double[npoint]; |
2341 |
// double *dyout = new double[npoint]; |
2342 |
// double *dthxout = new double[npoint]; |
2343 |
// double *dthyout = new double[npoint]; |
2344 |
// double *dtlout = new double[npoint]; |
2345 |
// double *dzin = new double[npoint]; |
2346 |
|
2347 |
double *dxout; |
2348 |
double *dyout; |
2349 |
double *dthxout; |
2350 |
double *dthyout; |
2351 |
double *dtlout; |
2352 |
double *dzin; |
2353 |
|
2354 |
dxout = (double*) malloc(npoint*sizeof(double)); |
2355 |
dyout = (double*) malloc(npoint*sizeof(double)); |
2356 |
dthxout = (double*) malloc(npoint*sizeof(double)); |
2357 |
dthyout = (double*) malloc(npoint*sizeof(double)); |
2358 |
dtlout = (double*) malloc(npoint*sizeof(double)); |
2359 |
dzin = (double*) malloc(npoint*sizeof(double)); |
2360 |
|
2361 |
double dal[5]; |
2362 |
|
2363 |
double dzini = (double)zini; |
2364 |
|
2365 |
int ifail = 0; |
2366 |
|
2367 |
for (int i=0; i<5; i++) dal[i] = (double)al[i]; |
2368 |
for (int i=0; i<npoint; i++) dzin[i] = (double)z[i]; |
2369 |
|
2370 |
TrkParams::Load(1); |
2371 |
if( !TrkParams::IsLoaded(1) ){ |
2372 |
cout << "int Trajectory::DoTrack(float* al) --- ERROR --- m.field not loaded"<<endl; |
2373 |
return 0; |
2374 |
} |
2375 |
// dotrack2_(&(npoint),dzin,dxout,dyout,dthxout,dthyout,dtlout,dal,&ifail); |
2376 |
dotrack3_(&(npoint),dzin,dxout,dyout,dthxout,dthyout,dtlout,dal,&dzini,&ifail); |
2377 |
|
2378 |
for (int i=0; i<npoint; i++){ |
2379 |
x[i] = (float)*(dxout+i); |
2380 |
y[i] = (float)*(dyout+i); |
2381 |
thx[i] = (float)*(dthxout+i); |
2382 |
thy[i] = (float)*(dthyout+i); |
2383 |
tl[i] = (float)*(dtlout+i); |
2384 |
} |
2385 |
|
2386 |
if(dxout) free( dxout ); |
2387 |
if(dyout) free( dyout ); |
2388 |
if(dthxout)free( dthxout ); |
2389 |
if(dthyout)free( dthyout ); |
2390 |
if(dtlout) free( dtlout ); |
2391 |
if(dzin) free( dzin ); |
2392 |
|
2393 |
// delete [] dxout; |
2394 |
// delete [] dyout; |
2395 |
// delete [] dthxout; |
2396 |
// delete [] dthyout; |
2397 |
// delete [] dtlout; |
2398 |
// delete [] dzin; |
2399 |
|
2400 |
|
2401 |
return ifail; |
2402 |
}; |
2403 |
|
2404 |
/** |
2405 |
* |
2406 |
* >>> OBSOLETE !!! use Trajectory::DoTrack(float* al, float zini) instead |
2407 |
* |
2408 |
*/ |
2409 |
int Trajectory::DoTrack2(float* al, float zini){ |
2410 |
|
2411 |
cout << endl; |
2412 |
cout << " int Trajectory::DoTrack2(float* al, float zini) --->> NB NB !! this method is going to be eliminated !!! "<<endl; |
2413 |
cout << " >>>> replace it with TrkTrack::DoTrack(Trajectory* t) <<<<"<<endl; |
2414 |
cout << " (Sorry Wolfgang!! Don't be totally confused!! By Elena)"<<endl; |
2415 |
cout << endl; |
2416 |
|
2417 |
return DoTrack(al,zini); |
2418 |
|
2419 |
}; |
2420 |
|
2421 |
|
2422 |
|
2423 |
ClassImp(TrkLevel2); |
2424 |
ClassImp(TrkSinglet); |
2425 |
ClassImp(TrkTrack); |
2426 |
ClassImp(Trajectory); |