/[PAMELA software]/tracker/flight/PamelaViewer/src/TrkViewer.cpp
ViewVC logotype

Contents of /tracker/flight/PamelaViewer/src/TrkViewer.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1.1.1 - (show annotations) (download) (vendor branch)
Tue Jan 20 10:29:57 2009 UTC (15 years, 10 months ago) by pam-fi
Branch: tracker
CVS Tags: V00
Changes since 1.1: +0 -0 lines
Error occurred while calculating annotation data.
tracker utilities

1 /**
2 * \file TrkViewer.cpp
3 * \author Elena Vannuccini
4 */
5 #include <TrkViewer.h>
6 // (mechanical) z-coordinate of the tracker planes
7 #define ZTRK6 -22.23
8 #define ZTRK5 -13.32
9 #define ZTRK4 -4.42
10 #define ZTRK3 4.48
11 #define ZTRK2 13.38
12 #define ZTRK1 22.28
13 // (mechanical) x/y-coordinates of magnet cavity
14 #define XTRKL -8.1
15 #define XTRKR 8.1
16 #define YTRKL -6.6
17 #define YTRKR 6.6
18 // (approximate) position of dead area among sensors
19 #define DX1 -2.7
20 #define DX2 2.7
21 #define DY -0.2
22
23 #define VCOL 18//29 // volume color
24 #define PCOL 1 // plane color
25 #define CCOL 15 // column color
26
27 #define HTCOL 2 // track-hit color
28 #define HSCOL 11 // singlet-hit color
29
30 #define TCCOL 40 // color of track-candidates
31
32 /**
33 * Tracker volume (x/y/top projections) constructor
34 */
35 GTrkVolume::GTrkVolume(){
36
37 cavity_l = new TBox(XTRKL,ZTRK6,XTRKR,ZTRK1);
38 cavity_l->SetFillColor(VCOL);
39 plane_l[0] = new TLine(XTRKL,ZTRK1,XTRKR,ZTRK1);
40 plane_l[1] = new TLine(XTRKL,ZTRK2,XTRKR,ZTRK2);
41 plane_l[2] = new TLine(XTRKL,ZTRK3,XTRKR,ZTRK3);
42 plane_l[3] = new TLine(XTRKL,ZTRK4,XTRKR,ZTRK4);
43 plane_l[4] = new TLine(XTRKL,ZTRK5,XTRKR,ZTRK5);
44 plane_l[5] = new TLine(XTRKL,ZTRK6,XTRKR,ZTRK6);
45 for(Int_t ip=1; ip<=6;ip++)plane_l[ip-1]->SetLineColor(PCOL);
46
47 cavity_r = new TBox(YTRKL,ZTRK6,YTRKR,ZTRK1);
48 cavity_r->SetFillColor(VCOL);
49 plane_r[0] = new TLine(YTRKL,ZTRK1,YTRKR,ZTRK1);
50 plane_r[1] = new TLine(YTRKL,ZTRK2,YTRKR,ZTRK2);
51 plane_r[2] = new TLine(YTRKL,ZTRK3,YTRKR,ZTRK3);
52 plane_r[3] = new TLine(YTRKL,ZTRK4,YTRKR,ZTRK4);
53 plane_r[4] = new TLine(YTRKL,ZTRK5,YTRKR,ZTRK5);
54 plane_r[5] = new TLine(YTRKL,ZTRK6,YTRKR,ZTRK6);
55 for(Int_t ip=1; ip<=6;ip++)plane_r[ip-1]->SetLineColor(PCOL);
56
57 cavity_c = new TBox(XTRKL,YTRKL,XTRKR,YTRKR);
58 cavity_c->SetFillColor(VCOL);
59
60 column_l[0] = new TLine(DX1,ZTRK6,DX1,ZTRK1);
61 column_l[1] = new TLine(DX2,ZTRK6,DX2,ZTRK1);
62 column_l[0]->SetLineColor(CCOL);
63 column_l[1]->SetLineColor(CCOL);
64 column_l[0]->SetLineStyle(2);
65 column_l[1]->SetLineStyle(2);
66
67 column_r[0] = new TLine(DY,ZTRK6,DY,ZTRK1);
68 column_r[0]->SetLineColor(CCOL);
69 column_r[0]->SetLineStyle(2);
70
71 column_c[0] = new TLine(DX1,YTRKL,DX1,YTRKR);
72 column_c[1] = new TLine(DX2,YTRKL,DX2,YTRKR);
73 column_c[2] = new TLine(XTRKL,DY,XTRKR,DY);
74 column_c[0]->SetLineColor(CCOL);
75 column_c[1]->SetLineColor(CCOL);
76 column_c[2]->SetLineColor(CCOL);
77 column_c[0]->SetLineStyle(2);
78 column_c[1]->SetLineStyle(2);
79 column_c[2]->SetLineStyle(2);
80
81
82 }
83 //
84 /**
85 *Draw x projection
86 */
87 void GTrkVolume::DrawProjectionX(){
88 cavity_l->Draw();
89 for(Int_t ip=1; ip<=6;ip++)plane_l[ip-1]->Draw();
90 column_l[0]->Draw();
91 column_l[1]->Draw();
92 }
93 //
94 /**
95 *Draw y projection
96 */
97 void GTrkVolume::DrawProjectionY(){
98 cavity_r->Draw();
99 for(Int_t ip=1; ip<=6;ip++)plane_r[ip-1]->Draw();
100 column_r[0]->Draw();
101 }
102 //
103 /**
104 *Draw top projection
105 */
106 void GTrkVolume::DrawProjectionT(){
107 cavity_c->Draw();
108 column_c[0]->Draw();
109 column_c[1]->Draw();
110 column_c[2]->Draw();
111 }
112 //
113 void GTrkVolume::Delete(){
114 cavity_l->Delete();
115 cavity_r->Delete();
116 cavity_c->Delete();
117 for(Int_t ip=1; ip<=6;ip++){
118 plane_l[ip-1]->Delete();
119 plane_r[ip-1]->Delete();
120 };
121 column_l[0]->Delete();
122 column_l[1]->Delete();
123 column_r[0]->Delete();
124 column_c[0]->Delete();
125 column_c[1]->Delete();
126 column_c[2]->Delete();
127
128
129 };
130
131 //===============================================================
132 //
133 //
134 //===============================================================
135 GTrkLevel2::GTrkLevel2(){
136 Track = new TClonesArray("GTrkTrack");
137 Image = new TClonesArray("GTrkTrack");
138 SingletX = new TClonesArray("GTrkHit");
139 SingletY = new TClonesArray("GTrkHit");
140
141 DRAWIMAGE=false;
142
143 }
144 //
145 //Int_t GTrkLevel2::Set(TrkLevel2* trk_event){
146 // return 0;
147 //}
148 //
149 /**
150 * Constructor
151 */
152 GTrkLevel2::GTrkLevel2(TrkLevel2* trk_event){
153
154 // cout << "Tracker: view status ";
155 for(Int_t iv=0; iv<12; iv++)cout << trk_event->good[iv];
156 cout << endl;
157
158 // cout << "Tracker: n.tracks "<<trk_event->GetNTracks()<<endl;
159 // cout << "Tracker: n.singlet X "<<trk_event->nclsx()<<endl;
160 // cout << "Tracker: n.singlet Y "<<trk_event->nclsy()<<endl;
161
162 Track = new TClonesArray("GTrkTrack");
163 Image = new TClonesArray("GTrkTrack");
164 SingletX = new TClonesArray("GTrkHit");
165 SingletY = new TClonesArray("GTrkHit");
166
167 DRAWIMAGE=false;
168
169 SetTracks(trk_event);
170 SetSinglets(trk_event);
171
172 }
173
174 /**
175 * Constructor
176 */
177 GTrkLevel2::GTrkLevel2(PamLevel2* pam_event){
178
179 TrkLevel2* trk_event = pam_event->GetTrkLevel2();
180
181 // cout << "Tracker: view status ";
182 // for(Int_t iv=0; iv<12; iv++)cout << trk_event->good[iv];
183 // cout << endl;
184
185 // cout << "Tracker: n.tracks "<<trk_event->GetNTracks()<<endl;
186 // cout << "Tracker: n.singlet X "<<trk_event->nclsx()<<endl;
187 // cout << "Tracker: n.singlet Y "<<trk_event->nclsy()<<endl;
188
189 Track = new TClonesArray("GTrkTrack");
190 Image = new TClonesArray("GTrkTrack");
191 SingletX = new TClonesArray("GTrkHit");
192 SingletY = new TClonesArray("GTrkHit");
193
194 DRAWIMAGE=false;
195
196 SetTracks(pam_event);
197 SetSinglets(trk_event);
198
199 }
200
201 void GTrkLevel2::SetTracks(TrkLevel2* trk_event){
202 TClonesArray &gt = *Track;
203 TClonesArray &gi = *Image;
204 TClonesArray &gsx = *SingletX;
205 TClonesArray &gsy = *SingletY;
206 Int_t isx = SingletX->GetEntries();
207 Int_t isy = SingletY->GetEntries();
208 for(Int_t it=0; it<trk_event->GetNTracks(); it++){
209 TrkTrack *tk = trk_event->GetTrack(it);
210 new(gt[it]) GTrkTrack(tk);
211 for(Int_t i=0; i<6; i++){
212 if(tk->XGood(i)) new(gsx[isx++]) GTrkHit(tk->xm[i],tk->zm[i],tk->dedx_x[i]);
213 if(tk->YGood(i)) new(gsy[isy++]) GTrkHit(tk->ym[i],tk->zm[i],tk->dedx_y[i]);
214 }
215
216 if( tk->HasImage() ){
217 TrkTrack *tk = trk_event->GetTrackImage(it);
218 new(gi[it]) GTrkTrack(tk);
219 }
220
221
222 }
223
224
225 }
226 void GTrkLevel2::SetTracks(PamLevel2* pam_event){
227 TClonesArray &gt = *Track;
228 TClonesArray &gi = *Image;
229 TClonesArray &gsx = *SingletX;
230 TClonesArray &gsy = *SingletY;
231 Int_t isx = SingletX->GetEntries();
232 Int_t isy = SingletY->GetEntries();
233
234 // ==== NON FUNZIONA ===
235 // TClonesArray *sorted = pam_event->GetTracks();
236 // cout << "==== "<<sorted->GetEntries()<<endl;
237 // for(Int_t it=0; it<pam_event->GetTrkLevel2()->GetNTracks(); it++){
238
239 // TrkTrack *tk = (TrkTrack*)(*sorted)[it];
240 // new(gt[it]) GTrkTrack(tk);
241 // for(Int_t i=0; i<6; i++){
242 // if(tk->XGood(i)) new(gsx[isx++]) GTrkHit(tk->xm[i],tk->zm[i],tk->dedx_x[i]);
243 // if(tk->YGood(i)) new(gsy[isy++]) GTrkHit(tk->ym[i],tk->zm[i],tk->dedx_y[i]);
244 // }
245 // }
246
247 for(Int_t it=0; it<pam_event->GetTrkLevel2()->GetNTracks(); it++){
248 TrkTrack *tk = pam_event->GetTrack(it)->GetTrkTrack();
249 new(gt[it]) GTrkTrack(tk);
250 for(Int_t i=0; i<6; i++){
251 if(tk->XGood(i)) new(gsx[isx++]) GTrkHit(tk->xm[i],tk->zm[i],tk->dedx_x[i]);
252 if(tk->YGood(i)) new(gsy[isy++]) GTrkHit(tk->ym[i],tk->zm[i],tk->dedx_y[i]);
253 }
254 if( tk->HasImage() ){
255 // cout << "*** HAS IMAGE ***"<<endl;
256 TrkTrack *tk = pam_event->GetTrackImage(it)->GetTrkTrack();
257 new(gi[it]) GTrkTrack(tk);
258 }
259 }
260
261
262 }
263
264 void GTrkLevel2::SetSinglets(TrkLevel2* trk_event){
265 TClonesArray &gsx = *SingletX;
266 Int_t isx = SingletX->GetEntries();
267 for(Int_t is=0; is<trk_event->nclsx(); is++){
268
269 TrkSinglet *s = trk_event->GetSingletX(is);
270
271 GTrkHit h = GTrkHit(s->coord[0],trk_event->GetZTrk(s->plane),s->sgnl);
272 // h.SetCluster(s->GetCluster());
273 new(gsx[isx++]) GTrkHit(&h);
274 }
275 TClonesArray &gsy = *SingletY;
276 Int_t isy = SingletY->GetEntries();
277 Int_t iss=0;
278 for(Int_t is=0; is<trk_event->nclsy(); is++){
279
280 for(Int_t i=0; i<2; i++){
281 TrkSinglet *s = trk_event->GetSingletY(is);
282
283 GTrkHit h = GTrkHit(s->coord[i],trk_event->GetZTrk(s->plane),s->sgnl);
284 // h.SetCluster(s->GetCluster());
285 new(gsy[isy++]) GTrkHit(&h);
286 iss++;
287 }
288 }
289 }
290
291 //
292 void GTrkLevel2::DrawProjectionX(){
293
294
295 TClonesArray &s = *(SingletX);
296 for(Int_t is=0; is<SingletX->GetEntries(); is++){
297 ((GTrkHit*)s[is])->Draw();
298 };
299
300 TClonesArray &t = *(Track);
301 for(Int_t it=0; it<Track->GetEntries(); it++){
302 ((GTrkTrack*)t[it])->DrawProjectionX();
303 };
304
305 };
306 void GTrkLevel2::DrawProjectionY(){
307
308 TClonesArray &s = *(SingletY);
309 for(Int_t is=0; is<SingletY->GetEntries(); is++){
310 ((GTrkHit*)s[is])->Draw();
311 };
312
313 TClonesArray &t = *(Track);
314 for(Int_t it=0; it<Track->GetEntries(); it++){
315 ((GTrkTrack*)t[it])->DrawProjectionY();
316 };
317
318 if(DRAWIMAGE){
319 TClonesArray &i = *(Image);
320 for(Int_t it=0; it<Image->GetEntries(); it++){
321 if( i[it] ){
322 ((GTrkTrack*)i[it])->GetProjectionY()->SetLineStyle(2);
323 ((GTrkTrack*)i[it])->DrawProjectionY();
324 }
325 };
326 }
327
328 };
329 void GTrkLevel2::DrawProjectionT(){
330
331 TClonesArray &t = *(Track);
332 for(Int_t it=0; it<Track->GetEntries(); it++){
333 ((GTrkTrack*)t[it])->DrawProjectionT();
334 };
335
336 if(DRAWIMAGE){
337 TClonesArray &i = *(Image);
338 for(Int_t it=0; it<Image->GetEntries(); it++){
339 if( i[it] ){
340 ((GTrkTrack*)i[it])->GetProjectionT()->SetLineStyle(2);
341 ((GTrkTrack*)i[it])->DrawProjectionT();
342 }
343 };
344 }
345
346 };
347 //
348 void GTrkLevel2::Delete(){
349 Track->Delete();
350 Image->Delete();
351 SingletX->Delete();
352 SingletY->Delete();
353 }
354 void GTrkLevel2::Clear(){
355 Track->Clear("C");
356 Image->Clear("C");
357 SingletX->Clear("C");
358 SingletY->Clear("C");
359 }
360 //===============================================================
361 //
362 //
363 //===============================================================
364 GTrkTrack::GTrkTrack(){
365
366 tr_l = new TPolyLine();
367 tr_r = new TPolyLine();
368 tr_c = new TPolyLine();
369
370
371 for (Int_t i=0; i<6; i++){
372 // m_l[i] = m_r[i] = m_c[i] = NULL;
373 mt_l[i] = new TMarker();
374 mt_r[i] = new TMarker();
375 mt_c[i] = new TMarker();
376 }
377
378 }
379 GTrkTrack::~GTrkTrack(){
380
381 // delete tr_l;
382 // delete tr_r;
383 // delete tr_c;
384
385 tr_l->Delete();
386 tr_r->Delete();
387 tr_c->Delete();
388
389
390 for (Int_t i=0; i<6; i++){
391 // cout << " GTrkTrack::~GTrkTrack() ==> m_l[i] "<<m_l[i]<<endl;
392 // if( m_l[i] )m_l[i]->Delete();
393 // if( m_r[i] )m_r[i]->Delete();
394 // if( m_c[i] )m_c[i]->Delete();
395
396 if( mt_l[i] )mt_l[i]->Delete();
397 if( mt_r[i] )mt_r[i]->Delete();
398 if( mt_c[i] )mt_c[i]->Delete();
399 }
400
401 };
402 //
403 GTrkTrack::GTrkTrack(TrkTrack* tr){
404
405 // --------------------------------------------
406 // create trajectory in the tracker
407
408 float pz[100];
409 float ztop = 55.;
410 float zbot = -27.;
411 float dz = (ztop-zbot)/100;
412 for(int i=0; i<100; i++)pz[i]=ztop-dz*i;
413 // Trajectory *tj = new Trajectory(100);
414 Trajectory *tj = new Trajectory(100,pz);
415 tr->DoTrack2(tj);
416
417 // cout << tr->GetRigidity()<<" GV"<<endl;
418
419 // for (Int_t ii=0;ii<100;ii+=10)cout<<tj->x[ii]<<" - "<<tj->y[ii]<<" - "<<tj->z[ii]<<endl;
420
421 Int_t iwall=100;
422 // Int_t iwall=0;
423 // for(iwall=0; iwall<tj->npoint; iwall++){
424 // if(tj->x[iwall] <= XTRKL || tj->x[iwall] >= XTRKR || tj->y[iwall]<= YTRKL || tj->y[iwall]>= YTRKR)break;
425 // }
426 // create track projections
427 tr_l = new TPolyLine( iwall, tj->x, tj->z );
428 tr_r = new TPolyLine( iwall, tj->y, tj->z );
429 tr_c = new TPolyLine( iwall, tj->x, tj->y );
430 tr_l->SetLineColor(HTCOL);
431 tr_r->SetLineColor(HTCOL);
432 tr_c->SetLineColor(HTCOL);
433 tj->Delete();
434
435 // --------------------------------------------
436 // create
437 for(Int_t i=0; i<6; i++){
438 if(tr->XGood(i)){
439 mt_l[i] = new TMarker(tr->xm[i],tr->zm[i],20);
440 mt_l[i]->SetMarkerSize(0.6);
441 mt_l[i]->SetMarkerColor(HTCOL);
442 // cout << "X-hit "<<i<<" "<<tr->xm[i]<<endl;
443 }else{
444 mt_l[i] = NULL;
445 };
446 if(tr->YGood(i)){
447 mt_r[i] = new TMarker(tr->ym[i],tr->zm[i],20);
448 mt_r[i]->SetMarkerSize(0.6);
449 mt_r[i]->SetMarkerColor(HTCOL);
450 // cout << "Y-hit "<<i<<" "<<tr->ym[i]<<endl;
451 }else{
452 mt_r[i] = NULL;
453 };
454 if(tr->XGood(i) || tr->YGood(i)){
455 mt_c[i] = new TMarker(tr->xm[i],tr->ym[i],20);
456 mt_c[i]->SetMarkerSize(0.9-i*0.1);
457 mt_c[i]->SetMarkerColor(HTCOL);
458 }else{
459 mt_c[i] = NULL;
460 };
461 }
462
463
464 }
465
466 //
467 GTrkTrack::GTrkTrack(GTrkTrack* gt){
468
469 tr_l = new TPolyLine( *(gt->tr_l) );
470 tr_r = new TPolyLine( *(gt->tr_r) );
471 tr_c = new TPolyLine( *(gt->tr_c) );
472
473
474 for(Int_t i=0; i<6; i++){
475 // m_l[i] = m_r[i] = m_c[i] = NULL;
476 // if( gt->m_l[i] ) m_l[i] = new GTrkHit( gt->m_l[i] );
477 // if( gt->m_r[i] ) m_r[i] = new GTrkHit( gt->m_r[i] );
478 // if( gt->m_c[i] ) m_c[i] = new GTrkHit( gt->m_c[i] );
479
480 mt_l[i] = new TMarker( *(gt->mt_l[i]) );
481 mt_r[i] = new TMarker( *(gt->mt_r[i]) );
482 mt_c[i] = new TMarker( *(gt->mt_c[i]) );
483 }
484
485 }
486 //
487 void GTrkTrack::Delete(){
488
489 tr_l->Delete();
490 tr_r->Delete();
491 tr_c->Delete();
492 for (Int_t i=0; i<6; i++){
493 if( mt_l[i] )mt_l[i]->Delete();
494 if( mt_r[i] )mt_r[i]->Delete();
495 if( mt_c[i] )mt_c[i]->Delete();
496 }
497 };
498
499 void GTrkTrack::Clear(){
500
501
502 };
503
504 void GTrkTrack::DrawProjectionX(){
505 // cout << "draw track x\n";
506 for(Int_t i=0; i<6; i++){
507 if( mt_l[i] )mt_l[i]->Draw();
508 };
509 tr_l->Draw();
510 };
511 void GTrkTrack::DrawProjectionY(){
512 // cout << "draw track y\n";
513 for(Int_t i=0; i<6; i++){
514 if( mt_r[i] )mt_r[i]->Draw();
515 };
516 tr_r->Draw();
517 };
518 void GTrkTrack::DrawProjectionT(){
519 // cout << "draw track t\n";
520 for(Int_t i=0; i<6; i++){
521 if( mt_c[i] )mt_c[i]->Draw();
522 };
523 tr_c->Draw();
524 };
525 //===============================================================
526 //
527 //
528 //===============================================================
529 GTrkHit::GTrkHit(){
530 m1 = new TMarker();
531 m2 = new TEllipse();
532 cl = 0;
533 };
534
535 GTrkHit::GTrkHit(GTrkHit* h){
536 m1 = new TMarker(*(h->m1));
537 m2 = new TEllipse(*(h->m2));
538 cl = h->cl;
539 };
540 //
541 GTrkHit::GTrkHit(Float_t x, Float_t y, Float_t w){
542
543 m1 = new TMarker(x,y,20);
544 // m2 = new TMarker(x,y,20);
545 m2 = new TEllipse(x,y,0.,0.);
546 cl = 0;
547
548
549
550 // create a palette vvvvvvvvvvvvvvvvvvvvvvvvvvv
551 // const Int_t colNum = 30;//12;
552 // const Int_t colBase = 500;
553 // // Int_t palette[colNum];
554 // TColor *color;
555 // // for (Int_t i=0;i<colNum;i++) {
556 // // /* Int_t ii = 255*i/colNum;
557 // // Int_t R=255;
558 // // Int_t G=255;
559 // // Int_t B=255-ii;*/
560 // // Int_t ii = 255*i/colNum;
561 // // Int_t R=255-51*((Int_t)(i/6))%6;
562 // // Int_t G=255;
563 // // Int_t B=255-51*i%6;
564 // // R=R/255;
565 // // G=G/255;
566 // // B=B/255;
567 // // if( !gROOT->GetColor(colBase+i) ){
568 // // color = new TColor(colBase+i,R,G,B,"");
569 // // // , pow(i/((colNum)*1.0),0.3)
570 // // // , pow(i/((colNum)*1.0),0.3)
571 // // // ,0.5*(i/((colNum)*1.0)),"");
572 // // }else{
573 // // color = gROOT->GetColor(colBase+i);
574 // // color->SetRGB(R,G,B);
575 // // };
576 // // // palette[i] = 301+i;
577 // // }
578 // // // gStyle->SetPalette(colNum,palette);
579 // // // create a palette ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
580
581 // // Int_t col = colNum*(Int_t)sqrt(w)/8 + colBase;
582 // m2->SetMarkerSize(w*0.2+1.1);
583 // m2->SetMarkerSize(sqrt(w)*.9);
584 // m2->SetMarkerColor(col);
585 Int_t colo = 10;
586 if ( w > 0.7 ) colo = 38;
587 if ( w > 2. ) colo = 4;
588 if ( w > 10. ) colo = 3;
589 if ( w > 100. ) colo = 2;
590 if ( w > 500. ) colo = 6;
591
592 m1->SetMarkerSize(0.7);
593 // m1->SetMarkerColor(HTCOL);
594 m1->SetMarkerColor(colo);
595
596 m2->SetR1(sqrt(w)*.8);
597 m2->SetR2(m2->GetR1());
598 m2->SetLineColor(colo);
599 m2->SetLineWidth(1);
600 // cout << w << " "<< col << endl;
601
602 // cout << "creato \n";
603 };
604 //
605 void GTrkHit::Draw(){
606
607 m2->Draw();
608 m1->Draw();
609
610 // gPad->GetCanvas()->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", "GTrkHit", this, "ExecuteEvent(Int_t,Int_t,Int_t,TObject*)");
611
612 };
613 void GTrkHit::Delete(){
614 if( m2 )m2->Delete();
615 if( m1 )m1->Delete();
616 // if( cl )cl->Delete();
617 };
618
619 void GTrkHit::ExecuteEvent(Int_t event, Int_t px, Int_t py, TObject* obj)
620 {
621 // Execute action corresponding to one event.
622 // This member function is called when a marker is clicked with the locator
623 // If Left button is clicked on a marker, the marker is moved to
624 // a new position when the mouse button is released.
625
626 if(obj != m1 && obj != m2) return;
627
628 if(!cl) return;
629
630 TCanvas *tempc = 0;
631 TPad *tempp = 0;
632 TCanvas *clc = 0;
633 TPad *clpx = 0;
634 TPad *clpy = 0;
635 TPad *clp = 0;
636 TPad *from = 0;
637 TObject *tempo = 0;
638
639 TIter nextcanv(gROOT->GetListOfCanvases());
640 while ( (tempc = (TCanvas*) nextcanv()) ) {
641 // cout << " canvas - "<<tempc->GetName()<<endl;
642 TIter nextpad( tempc->GetListOfPrimitives() );
643 while( (tempp = (TPad*) nextpad()) ){
644 // cout << " pad - "<<tempp->GetName()<<endl;
645 TString name = tempp->GetName();
646 // if( !name.CompareTo("clx") )clpx = tempp;
647 // if( !name.CompareTo("cly") )clpy = tempp;
648 if( !name.CompareTo("cl") )clp = tempp;
649 TIter nextobj( tempp->GetListOfPrimitives() );
650 while( (tempo = nextobj()) ){
651 // cout << " pad - "<<tempp->GetName()<<endl;
652 if( tempo == obj ) from = tempp;
653 };
654 };
655 }
656 // if(clp) clp->cd();
657 // if( (from->GetName()).Contains("px") && !clpx)
658 if(clp){
659 // clp->cd();
660 clp->Divide(1,2);
661 clp->Update();
662 TString name = from->GetName();
663 // cout << "obj selected from: "<<from->GetName()<<endl;
664 if( !name.CompareTo("px") ){cout<< "clp->cd(1)\n";clp->cd(1);}
665 if( !name.CompareTo("py") ){cout<< "clp->cd(2)\n";clp->cd(2);}
666
667 }else {
668 cout << "cl pad not fownd -- create new canvas \n";
669 clc = new TCanvas("cl","Cluster");
670 clp = (TPad*) clc;
671 clp->cd();
672 };
673
674 cout << "dove siamo: "<<gPad->GetName()<<endl;
675 switch (event) {
676
677
678 case kButton1Down:
679 // gVirtualX->SetTextColor(-1); // invalidate current text color (use xor mode)
680 // TAttMarker::Modify(); //Change marker attributes only if necessary
681 // No break !!!
682
683 cout << "kButton1Down\n";
684
685 break;
686
687 case kMouseMotion:
688 // pxold = px; pyold = py;
689 // gPad->SetCursor(kMove);
690 cout << "kMouseMotion\n";
691
692 cl->Draw();
693 gPad->Update();
694 // clp->Update();
695 // cout << clp->GetName() << endl;
696
697 break;
698
699 case kButton1Motion:
700 // p.fX = pxold; p.fY = pyold;
701 // gVirtualX->DrawPolyMarker(1, &p);
702 // p.fX = px; p.fY = py;
703 // gVirtualX->DrawPolyMarker(1, &p);
704 // pxold = px; pyold = py;
705 cout << "kButton1Motion\n";
706 break;
707
708 case kButton1Up:
709 // Double_t dpx, dpy, xp1,yp1;
710 // if (TestBit(kMarkerNDC)) {
711 // dpx = gPad->GetX2() - gPad->GetX1();
712 // dpy = gPad->GetY2() - gPad->GetY1();
713 // xp1 = gPad->GetX1();
714 // yp1 = gPad->GetY1();
715 // fX = (gPad->AbsPixeltoX(pxold)-xp1)/dpx;
716 // fY = (gPad->AbsPixeltoY(pyold)-yp1)/dpy;
717 // } else {
718 // fX = gPad->PadtoX(gPad->AbsPixeltoX(px));
719 // fY = gPad->PadtoY(gPad->AbsPixeltoY(py));
720 // }
721 // gPad->Modified(kTRUE);
722 // gVirtualX->SetTextColor(-1);
723 cout << "kButton1Up\n";
724 break;
725 }
726 }
727
728 //===============================================================
729 //
730 //
731 //===============================================================
732 GTrkCluster::GTrkCluster(){
733 sgnl = new TH1F();
734 cut = new TH1F();
735 seed = new TH1F();
736 sat = new TH1F();
737 maxs = new TLine();
738 zero = new TLine();
739 }
740 //
741 GTrkCluster::~GTrkCluster(){
742 if(sgnl) sgnl->Delete();
743 if(cut) cut->Delete();
744 if(seed) seed->Delete();
745 if(sat) sat->Delete();
746 if(maxs) maxs->Delete();
747 if(zero) zero->Delete();
748 }
749 //
750 GTrkCluster::GTrkCluster(GTrkCluster* gcl){
751 sgnl = gcl->sgnl;
752 cut = gcl->cut;
753 seed = gcl->seed;
754 sat = gcl->sat;
755 maxs = gcl->maxs;
756 zero = gcl->zero;
757 }
758 //
759 GTrkCluster::GTrkCluster(TrkCluster* cl){
760
761 maxs = new TLine(cl->maxs, 0, cl->maxs, cl->clsignal[cl->indmax]);
762
763 Float_t min = cl->maxs - cl->indmax - 0.5;
764 Float_t max = min + cl->CLlength;
765
766 zero = new TLine(min,0.,max,0.);
767
768 TString title = "";
769 title.Form("Cluster - view %i",cl->view);
770
771 TString h = "";
772 h = ""; h.Form("sgnl-%i",cl->maxs);
773 sgnl = new TH1F(h.Data(),title.Data(),cl->CLlength,min,max);
774 h = ""; h.Form("cut-%i",cl->maxs);
775 cut = new TH1F(h.Data(),title.Data(),cl->CLlength,min,max);
776 h = ""; h.Form("seed-%i",cl->maxs);
777 seed = new TH1F(h.Data(),title.Data(),cl->CLlength,min,max);
778 h = ""; h.Form("sat-%i",cl->maxs);
779 sat = new TH1F(h.Data(),title.Data(),cl->CLlength,min,max);
780 for (Int_t is = 0; is < cl->CLlength; is++){
781 sgnl->Fill( (Float_t)(min+is), cl->clsignal[is] );
782 cut->Fill( (Float_t)(min+is), 4*cl->clsigma[is] );
783 if( (cl->view)%2 )seed->Fill( (Float_t)(min+is), 6*cl->clsigma[is] );
784 else seed->Fill( (Float_t)(min+is), 7*cl->clsigma[is] );
785
786 // if( (cl->view)%2 && cl->cladc[is] < 100 ) sat->Fill( (Float_t)(min+is), cl->clsignal[is] );
787 // if( !((cl->view)%2) && cl->cladc[is] > 3100) sat->Fill( (Float_t)(min+is), cl->clsignal[is] );
788 if( (cl->view)%2 )sat->Fill( (Float_t)(min+is), -80+cl->clsignal[is]+cl->cladc[is] );
789 else sat->Fill( (Float_t)(min+is), 2980+cl->clsignal[is]-cl->cladc[is] );
790
791 // if( (cl->view)%2 )cout << -100+cl->clsignal[is]+cl->cladc[is] <<endl;
792 // else cout << 3100+cl->clsignal[is]-cl->cladc[is]<<endl ;
793
794 };
795
796 Float_t mi = -5;
797 Float_t ma = 1.1*sgnl->GetMaximum(5000);
798
799 sgnl->SetStats(0);
800 sgnl->SetLineWidth(2);
801 sgnl->SetFillColor(19);
802 sgnl->SetMinimum(mi);
803 sgnl->SetMaximum(ma);
804
805 cut->SetStats(0);
806 cut->SetLineStyle(2);
807 cut->SetLineWidth(2);
808 cut->SetLineColor(9);
809 cut->SetMinimum(mi);
810 cut->SetMaximum(ma);
811
812
813 sat->SetStats(0);
814 // sat->SetFillColor(18);
815 // sat->SetLineWidth(0);
816 // sat->SetLineColor(18);
817 sat->SetLineWidth(2);
818 sat->SetLineStyle(2);
819 sat->SetLineColor(2);
820 sat->SetMinimum(mi);
821 sat->SetMaximum(ma);
822
823 seed->SetStats(0);
824 // seed->SetFillColor(19);
825 seed->SetLineWidth(2);
826 seed->SetLineStyle(2);
827 seed->SetLineColor(3);
828 seed->SetMinimum(mi);
829 seed->SetMaximum(ma);
830
831 maxs->SetLineStyle(2);
832 maxs->SetLineWidth(2);
833 maxs->SetLineColor(2);
834 zero->SetLineStyle(0);
835 zero->SetLineWidth(1);
836 zero->SetLineColor(0);
837
838 }
839 //
840 void GTrkCluster::Draw(){
841
842 sgnl->Draw();
843 sat->Draw("same");
844
845 cut->Draw("same");
846 seed->Draw("same");
847
848 maxs->Draw();
849
850
851 }
852 //
853 void GTrkCluster::Delete(){
854 if(sgnl) sgnl->Delete();
855 if(cut) cut->Delete();
856 if(seed) seed->Delete();
857 if(sat) sat->Delete();
858 if(maxs) maxs->Delete();
859 if(zero) zero->Delete();
860 }
861 //===============================================================
862 //
863 //
864 //===============================================================
865 GTrkHough::GTrkHough(){
866
867 TrackCandidates = new TClonesArray("GTrkTrack");
868
869 }
870 /**
871 * Constructor
872 */
873 GTrkHough::GTrkHough(TrkHough *lh){
874
875
876 // cout << "Tracker: n.track-candidates "<<lh->GetNCandidates()<<endl;
877
878 TrackCandidates = new TClonesArray("GTrkTrack");
879 for(int i=0; i<100; i++)gtrip[i]=0;
880 for(int i=0; i<100; i++)gdoub[i]=0;
881
882 if(!lh)return;
883 SetCandidates(lh);
884 SetHough(lh);
885
886 }
887 void GTrkHough::SetCandidates(TrkHough* lh){
888
889 if(!lh)return;
890 TClonesArray &gt = *TrackCandidates;
891 for(Int_t it=0; it<lh->GetNCandidates(); it++){
892 TrkTrack *tk = lh->GetCandidate(it);
893 new(gt[it]) GTrkTrack(tk);
894 }
895
896 }
897
898 void GTrkHough::SetHough(TrkHough* lh){
899
900 if(!lh)return;
901
902 trip = new TMultiGraph();
903 doub = new TMultiGraph();
904
905 // cout << "Tracker: n.doublets "<<lh->ndblt<<endl;
906 // ALL --------------------------------------------------
907 int ndblt=0;
908 gdoub[ndblt] = new TGraph();
909 for(int id=0; id<lh->ndblt; id++){
910 // cout << "**Y** "<<lh->alfayz1[id]<<" "<<lh->alfayz2[id]<<endl;
911 // cout << "**Y** "<<lh->db_cloud[id]<<endl;
912 gdoub[ndblt]->SetPoint(id+1,lh->alfayz1[id],lh->alfayz2[id]);
913 }
914 doub->Add(gdoub[ndblt],"p");
915 ndblt++;
916 // CLOUDS --------------------------------------------------
917
918
919 // cout << "Tracker: n.triplets "<<lh->ntrpt<<endl;
920 int ntrpt=0;
921 gtrip[ntrpt] = new TGraph();
922 for(int id=0; id<lh->ntrpt; id++){
923 // cout << "**X** "<<lh->alfaxz1[id]<<" "<<lh->alfaxz2[id]<<endl;
924 // cout << "**X** "<<lh->tr_cloud[id]<<endl;
925 gtrip[ntrpt]->SetPoint(id+1,lh->alfaxz1[id],lh->alfaxz2[id]);
926 }
927 trip->Add(gtrip[ntrpt],"p");
928 ntrpt++;
929
930 // cout << "done"<<endl;
931
932 }
933
934 void GTrkHough::DrawProjectionX(){
935
936
937 if(!TrackCandidates)return;
938 TClonesArray &t = *(TrackCandidates);
939 for(Int_t it=0; it<TrackCandidates->GetEntries(); it++){
940 ((GTrkTrack*)t[it])->GetProjectionX()->SetLineStyle(2);
941 ((GTrkTrack*)t[it])->GetProjectionX()->SetLineColor(TCCOL);
942 for(int ip=0; ip<6; ip++) if( ((GTrkTrack*)t[it])->GetProjectionX_Points(ip) ) ((GTrkTrack*)t[it])->GetProjectionX_Points(ip)->SetMarkerColor(TCCOL);
943 ((GTrkTrack*)t[it])->DrawProjectionX();
944 };
945
946 // TCanvas *c = new TCanvas();
947 // trip->Draw("ap");
948
949 };
950 void GTrkHough::DrawProjectionY(){
951
952 if(!TrackCandidates)return;
953 TClonesArray &t = *(TrackCandidates);
954 for(Int_t it=0; it<TrackCandidates->GetEntries(); it++){
955 ((GTrkTrack*)t[it])->GetProjectionY()->SetLineStyle(2);
956 ((GTrkTrack*)t[it])->GetProjectionY()->SetLineColor(TCCOL);
957 for(int ip=0; ip<6; ip++)if(((GTrkTrack*)t[it])->GetProjectionY_Points(ip))((GTrkTrack*)t[it])->GetProjectionY_Points(ip)->SetMarkerColor(TCCOL);
958 ((GTrkTrack*)t[it])->DrawProjectionY();
959 };
960
961 // TCanvas *c = new TCanvas();
962 // doub->Draw("ap");
963
964 };
965 void GTrkHough::DrawProjectionT(){
966
967 if(!TrackCandidates)return;
968 TClonesArray &t = *(TrackCandidates);
969 for(Int_t it=0; it<TrackCandidates->GetEntries(); it++){
970 ((GTrkTrack*)t[it])->GetProjectionT()->SetLineStyle(2);
971 ((GTrkTrack*)t[it])->GetProjectionT()->SetLineColor(TCCOL);
972 for(int ip=0; ip<6; ip++)if(((GTrkTrack*)t[it])->GetProjectionT_Points(ip))((GTrkTrack*)t[it])->GetProjectionT_Points(ip)->SetMarkerColor(TCCOL);
973 ((GTrkTrack*)t[it])->DrawProjectionT();
974 };
975
976
977 };
978 void GTrkHough::Delete(){
979
980 // for(int i=0; i<100; i++)cout << gtrip[i] << endl;
981 // for(int i=0; i<100; i++)cout << gdoub[i] << endl;
982
983 // for(int i=0; i<100; i++)if(gtrip[i])trip->Delete();
984 // if(trip)trip->Delete();
985
986 // for(int i=0; i<100; i++)if(gdoub[i])doub->Delete();
987 // if(doub)doub->Delete();
988
989 TrackCandidates->Delete();
990
991 }
992 GTrkHough::~GTrkHough(){
993
994 // cout << "distruttore"<<endl;
995
996 // for(int i=0; i<100; i++)cout << gtrip[i] << endl;
997 // for(int i=0; i<100; i++)cout << gdoub[i] << endl;
998
999 // for(int i=0; i<100; i++)if(gtrip[i])trip->Delete();
1000 // if(trip)trip->Delete();
1001
1002 // for(int i=0; i<100; i++)if(gdoub[i])doub->Delete();
1003 // if(doub)doub->Delete();
1004
1005 TrackCandidates->Delete();
1006
1007 }
1008 void GTrkHough::Clear(){
1009
1010 TrackCandidates->Clear("C");
1011
1012 }
1013
1014
1015
1016 ClassImp(GTrkVolume);
1017 ClassImp(GTrkTrack);
1018 ClassImp(GTrkLevel2);
1019 ClassImp(GTrkHit);
1020 ClassImp(GTrkCluster);
1021 ClassImp(GTrkHough);
1022

  ViewVC Help
Powered by ViewVC 1.1.23