]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterDrawAZ.cxx
Introduced new DE names unique to each det element;
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterDrawAZ.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 // -------------------------------------
19 // Class AliMUONClusterDrawAZ
20 // -------------------------------------
21 // Cluster drawing for AZ cluster finder 
22 // Author: Alexander Zinchenko, JINR Dubna
23
24 #include <stdlib.h>
25 #include <Riostream.h>
26 //#include <TROOT.h>
27 #include <TCanvas.h>
28 #include <TLine.h>
29 //#include <TTree.h>
30 #include <TH2.h>
31 #include <TView.h>
32 #include <TStyle.h>
33
34 #include "AliMUONClusterDrawAZ.h"
35 #include "AliMUONClusterFinderAZ.h"
36 #include "AliMUONGeometryModuleTransformer.h"
37 #include "AliHeader.h"
38 #include "AliRun.h"
39 #include "AliMUON.h"
40 //#include "AliMUONChamber.h"
41 #include "AliMUONDigit.h"
42 #include "AliMUONHit.h"
43 #include "AliMUONRawCluster.h"
44 #include "AliMUONClusterInput.h"
45 #include "AliMUONPixel.h"
46 //#include "AliMC.h"
47 #include "AliMUONLoader.h"
48 #include "AliLog.h"
49
50 ClassImp(AliMUONClusterDrawAZ)
51  
52 //_____________________________________________________________________________
53 AliMUONClusterDrawAZ::AliMUONClusterDrawAZ()
54   : TObject(),
55     fData(0x0),
56     fFind(0x0),
57     fnMu(0),
58     fEvent(0),
59     fChamber(0),
60     fidDE(0),
61     fDebug(0),
62     fModif(0)
63 {
64 /// Default constructor
65   for (Int_t i=0; i<4; i++) fHist[i] = NULL;
66 }
67
68 //_____________________________________________________________________________
69 AliMUONClusterDrawAZ::AliMUONClusterDrawAZ(AliMUONClusterFinderAZ *clusFinder)
70   : TObject(),
71     fData(0x0),
72     fFind(clusFinder),
73     fnMu(0),
74     fEvent(0),
75     fChamber(0),
76     fidDE(0),
77     fDebug(1),
78     fModif(0)
79 {
80 /// Constructor
81   for (Int_t i=0; i<4; i++) fHist[i] = NULL;
82   Init();
83 }
84
85 //_____________________________________________________________________________
86 AliMUONClusterDrawAZ::~AliMUONClusterDrawAZ()
87 {
88 /// Destructor
89 }
90
91 //_____________________________________________________________________________
92 void AliMUONClusterDrawAZ::Init()
93 {
94 /// Initialization
95
96   TCanvas *c1 = new TCanvas("c1","Clusters",0,0,600,700);
97   //c1->SetFillColor(10);
98   c1->Divide(1,2);
99   new TCanvas("c2","Mlem",700,0,600,350);
100
101   // Get pointer to Alice detectors
102   //AliMUON *muon  = (AliMUON*) gAlice->GetModule("MUON");
103   //if (!muon) return;
104   //Loaders
105   AliRunLoader *rl = AliRunLoader::GetRunLoader();
106   AliLoader *gime = rl->GetLoader("MUONLoader");
107   fData = ((AliMUONLoader*)gime)->GetMUONData();
108
109   gime->LoadHits("READ"); 
110   gime->LoadRecPoints("READ"); 
111 }
112
113 //_____________________________________________________________________________
114 Bool_t AliMUONClusterDrawAZ::FindEvCh(Int_t nev, Int_t ch)
115 {
116 /// Find requested event and chamber (skip the ones before the selected)
117
118   if (nev < fEvent) return kFALSE;
119   else if (nev == fEvent) {
120     if (ch < fChamber) return kFALSE;
121     if (AliMUONClusterInput::Instance()->DetElemId() < fidDE) return kFALSE;
122   }
123   fEvent = nev;
124   fChamber = ch;
125   fidDE = AliMUONClusterInput::Instance()->DetElemId();
126   return kTRUE;
127 }
128
129 //_____________________________________________________________________________
130 void AliMUONClusterDrawAZ::DrawCluster()
131 {
132 /// Draw preclusters
133
134   TCanvas *c1 = (TCanvas*) gROOT->GetListOfCanvases()->FindObject("c1");
135
136   cout << " nev         " << fEvent << endl;
137     
138   char hName[4];
139   for (Int_t cath = 0; cath < 2; cath++) {
140     // Build histograms
141     //if (fHist[cath*2]) {fHist[cath*2]->Delete(); fHist[cath*2] = 0;}
142     //if (fHist[cath*2+1]) {fHist[cath*2+1]->Delete(); fHist[cath*2+1] = 0;}
143     if (fHist[cath*2]) fHist[cath*2] = 0;
144     if (fHist[cath*2+1]) fHist[cath*2+1] = 0;
145     if (fFind->GetNPads(cath) == 0) continue; // cluster on one cathode only
146     Float_t wxMin = 999, wxMax = 0, wyMin = 999, wyMax = 0; 
147     Int_t minDx = 0, maxDx = 0, minDy = 0, maxDy = 0;
148     for (Int_t i = 0; i < fFind->GetNPads(0)+fFind->GetNPads(1); i++) {
149       if (fFind->GetIJ(0,i) != cath) continue;
150       if (fFind->GetXyq(3,i) < wxMin) { wxMin = fFind->GetXyq(3,i); minDx = i; }
151       if (fFind->GetXyq(3,i) > wxMax) { wxMax = fFind->GetXyq(3,i); maxDx = i; }
152       if (fFind->GetXyq(4,i) < wyMin) { wyMin = fFind->GetXyq(4,i); minDy = i; }
153       if (fFind->GetXyq(4,i) > wyMax) { wyMax = fFind->GetXyq(4,i); maxDy = i; }
154     }
155     cout << minDx << " " << maxDx << " " << minDy << " " << maxDy << endl;
156     Int_t nx, ny, padSize;
157     Float_t xmin = 9999, xmax = -9999, ymin = 9999, ymax = -9999;
158     if (TMath::Nint(fFind->GetXyq(3,minDx)*1000) == TMath::Nint(fFind->GetXyq(3,maxDx)*1000) &&
159         TMath::Nint(fFind->GetXyq(4,minDy)*1000) == TMath::Nint(fFind->GetXyq(4,maxDy)*1000)) {
160       // the same segmentation
161       cout << " Same" << endl;
162       cout << fFind->GetXyq(3,minDx) << " " << fFind->GetXyq(3,maxDx) << " " 
163            << fFind->GetXyq(4,minDy) << " " << fFind->GetXyq(4,maxDy) << endl;
164       for (Int_t i = 0; i < fFind->GetNPads(0)+fFind->GetNPads(1); i++) {
165         if (fFind->GetIJ(0,i) != cath) continue;
166         if (fFind->GetXyq(0,i) < xmin) xmin = fFind->GetXyq(0,i);
167         if (fFind->GetXyq(0,i) > xmax) xmax = fFind->GetXyq(0,i);
168         if (fFind->GetXyq(1,i) < ymin) ymin = fFind->GetXyq(1,i);
169         if (fFind->GetXyq(1,i) > ymax) ymax = fFind->GetXyq(1,i);
170       }
171       xmin -= fFind->GetXyq(3,minDx); xmax += fFind->GetXyq(3,minDx);
172       ymin -= fFind->GetXyq(4,minDy); ymax += fFind->GetXyq(4,minDy);
173       nx = TMath::Nint ((xmax-xmin)/wxMin/2);
174       ny = TMath::Nint ((ymax-ymin)/wyMin/2);
175       cout << xmin << " " << xmax << " " << nx << " " << ymin << " " << ymax << " " << ny << endl;
176       sprintf(hName,"h%d",cath*2);
177       fHist[cath*2] = new TH2D(hName,"cluster",nx,xmin,xmax,ny,ymin,ymax);
178       for (Int_t i = 0; i < fFind->GetNPads(0)+fFind->GetNPads(1); i++) {
179         if (fFind->GetIJ(0,i) != cath) continue;
180         fHist[cath*2]->Fill(fFind->GetXyq(0,i),fFind->GetXyq(1,i),fFind->GetXyq(2,i));
181       }
182     } else {
183       // different segmentation in the cluster
184       cout << " Different" << endl;
185       cout << fFind->GetXyq(3,minDx) << " " << fFind->GetXyq(3,maxDx) << " " 
186            << fFind->GetXyq(4,minDy) << " " << fFind->GetXyq(4,maxDy) << endl;
187       Int_t nOK = 0;
188       Int_t indx, locMin, locMax;
189       if (TMath::Nint(fFind->GetXyq(3,minDx)*1000) != TMath::Nint(fFind->GetXyq(3,maxDx)*1000)) {
190         // different segmentation along x
191         indx = 0;
192         locMin = minDx;
193         locMax = maxDx;
194       } else {
195         // different segmentation along y
196         indx = 1;
197         locMin = minDy;
198         locMax = maxDy;
199       }
200       Int_t loc = locMin;
201       for (Int_t i = 0; i < 2; i++) {
202         // loop over different pad sizes
203         if (i > 0) loc = locMax;
204         padSize = TMath::Nint(fFind->GetXyq(indx+3,loc)*1000);
205         xmin = 9999; xmax = -9999; ymin = 9999; ymax = -9999;
206         for (Int_t j = 0; j < fFind->GetNPads(0)+fFind->GetNPads(1); j++) {
207           if (fFind->GetIJ(0,j) != cath) continue;
208           if (TMath::Nint(fFind->GetXyq(indx+3,j)*1000) != padSize) continue;
209           nOK++;
210           xmin = TMath::Min (xmin,fFind->GetXyq(0,j));
211           xmax = TMath::Max (xmax,fFind->GetXyq(0,j));
212           ymin = TMath::Min (ymin,fFind->GetXyq(1,j));
213           ymax = TMath::Max (ymax,fFind->GetXyq(1,j));
214         }
215         xmin -= fFind->GetXyq(3,loc); xmax += fFind->GetXyq(3,loc);
216         ymin -= fFind->GetXyq(4,loc); ymax += fFind->GetXyq(4,loc);
217         nx = TMath::Nint ((xmax-xmin)/fFind->GetXyq(3,loc)/2);
218         ny = TMath::Nint ((ymax-ymin)/fFind->GetXyq(4,loc)/2);
219         sprintf(hName,"h%d",cath*2+i);
220         fHist[cath*2+i] = new TH2D(hName,"cluster",nx,xmin,xmax,ny,ymin,ymax);
221         for (Int_t j = 0; j < fFind->GetNPads(0)+fFind->GetNPads(1); j++) {
222           if (fFind->GetIJ(0,j) != cath) continue;
223           if (TMath::Nint(fFind->GetXyq(indx+3,j)*1000) != padSize) continue;
224           fHist[cath*2+i]->Fill(fFind->GetXyq(0,j),fFind->GetXyq(1,j),fFind->GetXyq(2,j));
225         }
226       } // for (Int_t i=0;
227       if (nOK != fFind->GetNPads(cath)) cout << " *** Too many segmentations: nPads, nOK " 
228                                              << fFind->GetNPads(cath) << " " << nOK << endl;
229     } // if (TMath::Nint(fFind->GetXyq(3,minDx)*1000)
230   } // for (Int_t cath = 0;
231         
232   // Draw histograms and coordinates
233   //TH2D *hFake = 0x0;
234   TH2D *hFake = (TH2D*) gROOT->FindObject("hFake");
235   if (hFake) hFake->Delete();
236   if (fModif) {
237     // This part works only after a modification of THistPainter::PaintLego(Option_t *) 
238     // in ROOT
239     Double_t xmin = 9999, ymin = 9999, xmax = -9999, ymax = -9999, aMax = -1;
240     for (Int_t i = 0; i < 4; i++) {
241       if (!fHist[i]) continue;
242       xmin = TMath::Min (xmin, fHist[i]->GetXaxis()->GetXmin());
243       xmax = TMath::Max (xmax, fHist[i]->GetXaxis()->GetXmax());
244       ymin = TMath::Min (ymin, fHist[i]->GetYaxis()->GetXmin());
245       ymax = TMath::Max (ymax, fHist[i]->GetYaxis()->GetXmax());
246       aMax = TMath:: Max (aMax, fHist[i]->GetMaximum());
247     }
248     //if (c1->FindObject("hFake")) delete c1->FindObject("hFake");
249     hFake = new TH2D ("hFake", "hFake", 1, xmin, xmax, 1, ymin, ymax);
250     hFake->SetMaximum(aMax);
251     hFake->SetNdivisions(505,"Z");
252     hFake->SetStats(kFALSE);
253   }
254
255   TObject *stats = 0x0;
256   for (Int_t cath = 0; cath < 2; cath++) {
257     if (cath == 0) ModifyHistos();
258     if (fFind->GetNPads(cath) == 0) continue; // cluster on one cathode only
259     c1->cd(cath+1);
260     gPad->SetTheta(55);
261     gPad->SetPhi(30);
262     if (fModif) {
263       if (fHist[cath*2]) { 
264         fHist[cath*2]->Draw();
265         gPad->Update();
266         stats = fHist[cath*2]->GetListOfFunctions()->FindObject("stats");
267       }
268       hFake->Draw("legoFb");
269     }
270
271     Double_t x, y, x0, y0, r1 = 999, r2 = 0;
272     if (fHist[cath*2+1]) {
273       // 
274       x0 = fHist[cath*2]->GetXaxis()->GetXmin() - 1000*TMath::Cos(30*TMath::Pi()/180);
275       y0 = fHist[cath*2]->GetYaxis()->GetXmin() - 1000*TMath::Sin(30*TMath::Pi()/180);
276       r1 = 0;
277       Int_t ihist=cath*2;
278       for (Int_t iy = 1; iy <= fHist[ihist]->GetNbinsY(); iy++) {
279         y = fHist[ihist]->GetYaxis()->GetBinCenter(iy) 
280           + fHist[ihist]->GetYaxis()->GetBinWidth(iy);
281         for (Int_t ix = 1; ix <= fHist[ihist]->GetNbinsX(); ix++) {
282           if (fHist[ihist]->GetCellContent(ix,iy) > 0.1) {
283             x = fHist[ihist]->GetXaxis()->GetBinCenter(ix)
284               + fHist[ihist]->GetXaxis()->GetBinWidth(ix);
285             r1 = TMath::Max (r1,TMath::Sqrt((x-x0)*(x-x0)+(y-y0)*(y-y0)));
286           }
287         }
288       }
289       ihist = cath*2 + 1 ;
290       for (Int_t iy = 1; iy <= fHist[ihist]->GetNbinsY(); iy++) {
291         y = fHist[ihist]->GetYaxis()->GetBinCenter(iy)
292           + fHist[ihist]->GetYaxis()->GetBinWidth(iy);
293         for (Int_t ix = 1; ix <= fHist[ihist]->GetNbinsX(); ix++) {
294           if (fHist[ihist]->GetCellContent(ix,iy) > 0.1) {
295             x = fHist[ihist]->GetXaxis()->GetBinCenter(ix)
296               + fHist[ihist]->GetXaxis()->GetBinWidth(ix);
297             r2 = TMath::Max (r2,TMath::Sqrt((x-x0)*(x-x0)+(y-y0)*(y-y0)));
298           }
299         }
300       }
301       cout << r1 << " " << r2 << endl;
302     } // if (fHist[cath*2+1])
303     if (r1 > r2) {
304       if (fModif) fHist[cath*2]->Draw("lego1FbSame");
305       else fHist[cath*2]->Draw("lego1Fb");
306       // Draw background contaminated charges
307       //TH2D *hBkg = GetBackground(cath*2);
308       //if (hBkg) hBkg->Draw("lego1FbBbSameAxis");
309       if (fHist[cath*2+1]) fHist[cath*2+1]->Draw("lego1SameAxisBbFb");
310     } else {
311       if (fModif) fHist[cath*2+1]->Draw("lego1FbSame");
312       else fHist[cath*2+1]->Draw("lego1Fb");
313       fHist[cath*2]->Draw("lego1SameAxisFbBb");
314     }
315     c1->Update();
316     if (fModif) stats->Draw();
317   } // for (Int_t cath = 0;
318
319   // Draw simulated and reconstructed hits 
320   DrawHits();
321 }
322
323 //_____________________________________________________________________________
324 void AliMUONClusterDrawAZ::DrawHits()
325 {
326 /// Draw simulated and reconstructed hits 
327
328   TView *view[2] = { 0x0, 0x0 };
329   Double_t p1[3]={0}, p2[3], xNDC[6], xl, yl, zl;
330   TLine *line[199] = {0};
331   TCanvas *c1 = (TCanvas*) gROOT->GetListOfCanvases()->FindObject("c1");
332   if (c1) {
333     c1->cd(1);
334     view[0] = c1->Pad()->GetView();
335     c1->cd(2);
336     view[1] = c1->Pad()->GetView();
337   }
338   fData->SetTreeAddress("H RC");
339
340   TH2D *hist = fHist[0] ? fHist[0] : fHist[2];
341   p2[2] = hist->GetMaximum();
342
343   // Draw simulated hits
344   cout << " *** Simulated hits *** " << endl;
345   Int_t ntracks = 0;
346   if (fData->TreeH()) ntracks = (Int_t) fData->GetNtracks();
347   fnMu = 0;
348   Int_t ix, iy, iok, nLine = 0;
349   TClonesArray *hits = NULL;
350   for (Int_t i = 0; i < ntracks; i++) {
351     fData->GetTrack(i);
352     hits = fData->Hits();
353     Int_t nhits = (Int_t) hits->GetEntriesFast();
354     AliMUONHit* mHit;
355     for (Int_t ihit = 0; ihit < nhits; ihit++) {
356       mHit = (AliMUONHit*) hits->UncheckedAt(ihit);
357       if (mHit->Chamber() != fChamber+1) continue;  // chamber number
358       if (mHit->DetElemId() != fidDE) continue;  // det. elem. Id
359       AliMUONClusterInput::Instance()->Segmentation2(0)->GetTransformer()->
360                          Global2Local(fidDE, mHit->X(), mHit->Y(), mHit->Z(), xl, yl, zl);
361       //if (TMath::Abs(zl-fFind->GetZpad()) > 1) continue; // different slat
362       p2[0] = p1[0] = xl;        // x-pos of hit
363       p2[1] = p1[1] = yl;        // y-pos
364       if (p1[0] < hist->GetXaxis()->GetXmin() || 
365           p1[0] > hist->GetXaxis()->GetXmax()) continue;
366       if (p1[1] < hist->GetYaxis()->GetXmin() || 
367           p1[1] > hist->GetYaxis()->GetXmax()) continue;
368       // Check if track comes thru pads with signal
369       iok = 0;
370       for (Int_t ihist = 0; ihist < 4; ihist++) {
371         if (!fHist[ihist]) continue;
372         ix = fHist[ihist]->GetXaxis()->FindBin(p1[0]);
373         iy = fHist[ihist]->GetYaxis()->FindBin(p1[1]);
374         if (fHist[ihist]->GetCellContent(ix,iy) > 0.5) {iok = 1; break;}
375       }
376       if (!iok) continue;
377       gStyle->SetLineColor(1);
378       if (TMath::Abs((Int_t)mHit->Particle()) == 13) {
379         gStyle->SetLineColor(4);
380         if (fnMu < 2) {
381           fxyMu[fnMu][0] = p1[0];
382           fxyMu[fnMu++][1] = p1[1];
383         }
384       }     
385       printf(" Local coord.:  X=%10.4f, Y=%10.4f\n",p1[0],p1[1]);
386       printf(" Global coord.: X=%10.4f, Y=%10.4f, Z=%10.4f\n",mHit->X(),mHit->Y(),mHit->Z());
387       if (view[0] || view[1]) {
388         // Take into account track angles
389         p2[0] += mHit->Tlength() * TMath::Sin(mHit->Theta()/180*TMath::Pi()) 
390                                  * TMath::Cos(mHit->Phi()/180*TMath::Pi()) / 2;
391         p2[1] += mHit->Tlength() * TMath::Sin(mHit->Theta()/180*TMath::Pi()) 
392                                  * TMath::Sin(mHit->Phi()/180*TMath::Pi()) / 2;
393         for (Int_t ipad = 1; ipad < 3; ipad++) {
394           c1->cd(ipad);
395           if (!view[ipad-1]) continue;
396           view[ipad-1]->WCtoNDC(p1, &xNDC[0]);
397           view[ipad-1]->WCtoNDC(p2, &xNDC[3]);
398           //c1->DrawLine(xpad[0],xpad[1],xpad[3],xpad[4]);
399           line[nLine] = new TLine(xNDC[0],xNDC[1],xNDC[3],xNDC[4]);
400           line[nLine++]->Draw();
401         }
402       }
403     } // for (Int_t ihit = 0; ihit < nhits;
404   } // for (Int_t i = 0; i < ntracks;
405   fData->ResetHits();
406
407   // Draw reconstructed coordinates
408   if (fData->TreeR()) fData->GetRawClusters();
409   TClonesArray *rawclust = fData->RawClusters(fChamber);
410   AliMUONRawCluster *mRaw;
411   gStyle->SetLineColor(3);
412   cout << " *** Reconstructed hits *** " << endl;
413   if (rawclust) {
414     for (Int_t i = 0; i < rawclust ->GetEntriesFast(); i++) {
415       mRaw = (AliMUONRawCluster*)rawclust ->UncheckedAt(i);
416       AliMUONClusterInput::Instance()->Segmentation2(0)->GetTransformer()->
417                     Global2Local(fidDE, mRaw->GetX(0), mRaw->GetY(0), mRaw->GetZ(0), xl, yl, zl);
418       if (TMath::Abs(zl-fFind->GetZpad()) > 1) continue; // different slat
419       p2[0] = p1[0] = xl;        // x-pos of hit
420       p2[1] = p1[1] = yl;        // y-pos
421       if (p1[0] < hist->GetXaxis()->GetXmin() || 
422           p1[0] > hist->GetXaxis()->GetXmax()) continue;
423       if (p1[1] < hist->GetYaxis()->GetXmin() || 
424           p1[1] > hist->GetYaxis()->GetXmax()) continue;
425       /*
426         treeD->GetEvent(cath);
427         cout << mRaw->fMultiplicity[0] << mRaw->fMultiplicity[1] << endl;
428         for (Int_t j=0; j<mRaw->fMultiplicity[cath]; j++) {
429         Int_t digit = mRaw->fIndexMap[j][cath];
430         cout << ((AliMUONDigit*)fMuonDigits->UncheckedAt(digit))->Signal() << endl;
431         }
432       */
433       // Check if track comes thru pads with signal
434       iok = 0;
435       for (Int_t ihist = 0; ihist < 4; ihist++) {
436         if (!fHist[ihist]) continue;
437         ix = fHist[ihist]->GetXaxis()->FindBin(p1[0]);
438         iy = fHist[ihist]->GetYaxis()->FindBin(p1[1]);
439         if (fHist[ihist]->GetCellContent(ix,iy) > 0.5) {iok = 1; break;}
440       }
441       if (!iok) continue;
442       printf(" Local coord.:  X=%10.4f, Y=%10.4f\n",p1[0],p1[1]);
443       printf(" Global coord.: X=%10.4f, Y=%10.4f, Z=%10.4f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
444       if (view[0] || view[1]) {
445         for (Int_t ipad = 1; ipad < 3; ipad++) {
446           c1->cd(ipad);
447           if (!view[ipad-1]) continue;
448           view[ipad-1]->WCtoNDC(p1, &xNDC[0]);
449           view[ipad-1]->WCtoNDC(p2, &xNDC[3]);
450           line[nLine] = new TLine(xNDC[0],xNDC[1],xNDC[3],xNDC[4]);
451           line[nLine++]->Draw();
452         }
453       }
454     } // for (Int_t i = 0; i < rawclust ->GetEntries();
455   } // if (rawclust)
456   if (fData->TreeR()) fData->ResetRawClusters();
457   c1->Update();
458 }
459
460 //_____________________________________________________________________________
461 Int_t AliMUONClusterDrawAZ::Next()
462 {
463 /// What to do next?
464
465   // File
466   FILE *lun = 0;
467   //lun = fopen("pull.dat","w");
468
469   for (Int_t i = 0; i < fnMu; i++) {
470     // Check again if muon comes thru the used pads (due to extra splitting)
471     for (Int_t j = 0; j < fFind->GetNPads(0)+fFind->GetNPads(1); j++) {
472       if (TMath::Abs(fxyMu[i][0]-fFind->GetXyq(0,j))<fFind->GetXyq(3,j) && 
473           TMath::Abs(fxyMu[i][1]-fFind->GetXyq(1,j))<fFind->GetXyq(4,j)) {
474         if (fDebug) printf("%12.3e %12.3e %12.3e %12.3e\n",fxyMu[i][2],fxyMu[i][3],fxyMu[i][4],fxyMu[i][5]);
475         if (lun) fprintf(lun,"%4d %2d %12.3e %12.3e %12.3e %12.3e\n",fEvent,fChamber,fxyMu[i][2],fxyMu[i][3],fxyMu[i][4],fxyMu[i][5]);
476         break;
477       }
478     }
479   } // for (Int_t i=0; i<fnMu;
480
481   // What's next?
482   char command[8];
483   cout << " What is next? " << endl;
484   command[0] = ' '; 
485   gets(command);
486   if (command[0] == 'n' || command[0] == 'N') { fEvent++; fChamber = fidDE = 0; } // next event
487   else if (command[0] == 'q' || command[0] == 'Q') { if (lun) fclose(lun); } // exit display 
488   else if (command[0] == 'c' || command[0] == 'C') sscanf(command+1,"%d",&fChamber); // new chamber
489   else if (command[0] == 'e' || command[0] == 'E') { sscanf(command+1,"%d",&fEvent); fChamber = fidDE = 0; } // new event
490   else if (command[0] == 'd' || command[0] == 'D') { sscanf(command+1,"%d",&fidDE); fChamber = fidDE / 100 - 1; } // new DetElem.
491   else return 1; // Next precluster
492   return 0;
493 }
494
495
496 //_____________________________________________________________________________
497 void AliMUONClusterDrawAZ::ModifyHistos(void)
498 {
499 /// Modify histograms to bring them to (approximately) the same size
500
501   Int_t nhist = 0;
502   Float_t hlim[4][4], hbin[4][4]; // first index - xmin, xmax, ymin, ymax
503
504   Float_t binMin[4] = {999,999,999,999};
505
506   for (Int_t i = 0; i < 4; i++) {
507     if (!fHist[i]) {
508       hlim[0][i] = hlim[2][i] = 999;
509       hlim[1][i] = hlim[3][i] = -999;
510       continue;
511     }
512     hlim[0][i] = fHist[i]->GetXaxis()->GetXmin(); // xmin
513     hlim[1][i] = fHist[i]->GetXaxis()->GetXmax(); // xmax
514     hlim[2][i] = fHist[i]->GetYaxis()->GetXmin(); // ymin
515     hlim[3][i] = fHist[i]->GetYaxis()->GetXmax(); // ymax
516     hbin[0][i] = hbin[1][i] = fHist[i]->GetXaxis()->GetBinWidth(1);
517     hbin[2][i] = hbin[3][i] = fHist[i]->GetYaxis()->GetBinWidth(1);
518     binMin[0] = TMath::Min(binMin[0],hbin[0][i]);
519     binMin[2] = TMath::Min(binMin[2],hbin[2][i]);
520     nhist++;
521   }
522   binMin[1] = binMin[0];
523   binMin[3] = binMin[2];
524   cout << " Nhist: " << nhist << endl;
525
526   // Adjust histo limits for cathode with different segmentation
527   for (Int_t i = 0; i < 4; i+=2) {
528     if (!fHist[i+1]) continue;
529     Int_t imin, imax, i1 = i + 1;
530     for (Int_t lim = 0; lim < 4; lim++) {
531       while (1) {
532         if (hlim[lim][i] < hlim[lim][i1]) {
533           imin = i;
534           imax = i1;
535         } else {
536           imin = i1;
537           imax = i;
538         }
539         if (TMath::Abs(hlim[lim][imin]-hlim[lim][imax])<0.01*binMin[lim]) break;
540         if (lim == 0 || lim == 2) {
541           // find lower limit
542           hlim[lim][imax] -= hbin[lim][imax];
543         } else {
544           // find upper limit
545           hlim[lim][imin] += hbin[lim][imin];
546         }
547       } // while (1)
548     }
549   }
550     
551
552   Int_t imnmx = 0, nExtra = 0;
553   for (Int_t lim = 0; lim < 4; lim++) {
554     if (lim == 0 || lim == 2) imnmx = TMath::LocMin(4,hlim[lim]); // find lower limit
555     else imnmx = TMath::LocMax(4,hlim[lim]); // find upper limit
556
557     // Adjust histogram limit
558     for (Int_t i = 0; i < 4; i++) {
559       if (!fHist[i]) continue;
560       nExtra = TMath::Nint ((hlim[lim][imnmx]-hlim[lim][i]) / hbin[lim][i]);
561       hlim[lim][i] += nExtra * hbin[lim][i];
562     }
563   }
564     
565   // Rebuild histograms 
566   TH2D *hist = 0;
567   Int_t nx, ny;
568   Double_t x, y, cont, cmax=0;
569   char hName[4];
570   for (Int_t ihist = 0; ihist < 4; ihist++) {
571     if (!fHist[ihist]) continue;
572     nx = TMath::Nint((hlim[1][ihist]-hlim[0][ihist])/hbin[0][ihist]);
573     ny = TMath::Nint((hlim[3][ihist]-hlim[2][ihist])/hbin[2][ihist]);
574     cout << ihist << " " << hlim[0][ihist] << " " << hlim[1][ihist] << " " << nx;
575     cout << " " << hlim[2][ihist] << " " << hlim[3][ihist] << " " << ny << endl;
576     sprintf(hName,"hh%d",ihist);
577     hist =  new TH2D(hName,"hist",nx,hlim[0][ihist],hlim[1][ihist],ny,hlim[2][ihist],hlim[3][ihist]);
578     for (Int_t i=1; i<=fHist[ihist]->GetNbinsX(); i++) {
579       x = fHist[ihist]->GetXaxis()->GetBinCenter(i);
580       for (Int_t j=1; j<=fHist[ihist]->GetNbinsY(); j++) {
581         y = fHist[ihist]->GetYaxis()->GetBinCenter(j);
582         cont = fHist[ihist]->GetCellContent(i,j);
583         hist->Fill(x,y,cont);
584       }
585     }
586     cmax = TMath::Max (cmax,hist->GetMaximum());
587     sprintf(hName,"%s%d",fHist[ihist]->GetName(),ihist);
588     fHist[ihist]->Delete();
589     fHist[ihist] = new TH2D(*hist);
590     fHist[ihist]->SetName(hName);
591     fHist[ihist]->SetNdivisions(505,"Z");
592     hist->Delete(); 
593   }
594   if (fDebug) printf("%f \n",cmax);
595
596   for (Int_t ihist = 0; ihist < 4; ihist++) {
597     if (!fHist[ihist]) continue;
598     fHist[ihist]->SetMaximum(cmax);
599     fHist[ihist]->SetMinimum(0);
600   }
601 }
602
603 //_____________________________________________________________________________
604 void AliMUONClusterDrawAZ::AdjustHist(Double_t *xylim, const AliMUONPixel *pixPtr)
605 {
606 /// Adjust histogram limits for pixel drawing
607
608   Float_t xypads[4];
609   if (fHist[0]) {
610     xypads[0] = fHist[0]->GetXaxis()->GetXmin();
611     xypads[1] = -fHist[0]->GetXaxis()->GetXmax();
612     xypads[2] = fHist[0]->GetYaxis()->GetXmin();
613     xypads[3] = -fHist[0]->GetYaxis()->GetXmax();
614     for (Int_t i = 0; i < 4; i++) {
615       while(1) {
616         if (xylim[i] < xypads[i]) break;
617         xylim[i] -= 2*pixPtr->Size(i/2);
618       }
619     }
620   } 
621 }
622
623 //_____________________________________________________________________________
624 void AliMUONClusterDrawAZ::DrawHist(const char* canvas, TH2D *hist)
625 {
626 /// Draw histogram in given canvas 
627
628   Int_t ix = 0;
629   //((TCanvas*)gROOT->FindObject("c2"))->cd();
630   ((TCanvas*)gROOT->FindObject(canvas))->cd();
631   gPad->SetTheta(55);
632   gPad->SetPhi(30);
633   hist->Draw("lego1Fb");
634   gPad->Update();
635   //gets((char*)&ix);
636   if (fnMu) gets((char*)&ix);
637 }
638
639 //_____________________________________________________________________________
640 void AliMUONClusterDrawAZ::FillMuon(Int_t nfit, const Double_t *parOk, const Double_t *errOk)
641 {
642 /// Fill muon information
643
644   Int_t indx, imax;
645   Double_t cmax, rad;
646   for (Int_t i = 0; i < fnMu; i++) {
647     cmax = fxyMu[i][6];
648     for (Int_t j = 0; j < nfit; j++) {
649       indx = j<2 ? j*2 : j*2+1;  
650       rad = (fxyMu[i][0]-parOk[indx])*(fxyMu[i][0]-parOk[indx]) +
651             (fxyMu[i][1]-parOk[indx+1])*(fxyMu[i][1]-parOk[indx+1]);
652       if (rad < cmax) {
653         cmax = rad; 
654         imax = indx;
655         fxyMu[i][6] = cmax;
656         fxyMu[i][2] = parOk[imax] - fxyMu[i][0];
657         fxyMu[i][4] = parOk[imax+1] - fxyMu[i][1];
658         fxyMu[i][3] = errOk[imax];
659         fxyMu[i][5] = errOk[imax+1];
660       }
661     }      
662   }
663 }
664
665 //_____________________________________________________________________________
666 void AliMUONClusterDrawAZ::UpdateCluster(Int_t npad)
667 {
668 /// Update cluster after removing non-overlapped pads
669
670   Int_t cath = 0, ix = 0, iy = 0;
671   cout << " Update cluster " << endl;
672   gets((char*)&ix);
673   for (Int_t i = 0; i < npad; i++) {
674     if (TMath::Nint (fFind->GetXyq(2,i)) != -2) continue;
675     cath = fFind->GetIJ(0,i);
676     for (Int_t j = 0; j < 2; j++) {
677       Int_t ihist = cath * 2 + j;
678       if (!fHist[ihist]) continue;
679       ix = fHist[ihist]->GetXaxis()->FindBin(fFind->GetXyq(0,i));
680       iy = fHist[ihist]->GetYaxis()->FindBin(fFind->GetXyq(1,i));
681       Double_t cont = fHist[ihist]->GetCellContent(ix, iy);
682       if (cont < 0.1) continue;
683       fHist[ihist]->Fill(fFind->GetXyq(0,i), fFind->GetXyq(1,i), -cont);
684     }
685   }
686   TCanvas *c1 = (TCanvas*) gROOT->GetListOfCanvases()->FindObject("c1");
687   if (c1) {
688     c1->cd(1);
689     gPad->Modified();
690     gPad->Update();
691     c1->cd(2);
692     gPad->Modified();
693     gPad->Update();
694   }  
695 }
696
697 //_____________________________________________________________________________
698 TH2D* AliMUONClusterDrawAZ::GetBackground(Int_t iHist)
699 {
700 /// Build histogram with pads from the cluster contaminated by the background
701
702   //return 0x0;
703   Int_t cath = iHist / 2;
704   Double_t xmin = fHist[iHist]->GetXaxis()->GetXmin();
705   Double_t xmax = fHist[iHist]->GetXaxis()->GetXmax();
706   Double_t ymin = fHist[iHist]->GetYaxis()->GetXmin();
707   Double_t ymax = fHist[iHist]->GetYaxis()->GetXmax();
708   
709   // Create histogram
710   char hName[4];
711   sprintf(hName,"Bkg%1d",iHist);
712   TH2D *hist = (TH2D*) gROOT->FindObject(hName);
713   if (hist) hist->Delete();
714   hist = (TH2D*) fHist[iHist]->Clone(hName);
715   hist->Reset();
716
717   // Loop over pads
718   Int_t digit = 0, iok = 0, ix = 0, iy = 0; 
719   AliMUONDigit *mdig = 0x0;
720   Double_t cont = 0, x = 0, y = 0; 
721   for (Int_t i = 0; i < fFind->GetNPads(0)+fFind->GetNPads(1); i++) {
722     if (fFind->GetIJ(0,i) != cath) continue;
723     x = fFind->GetXyq(0,i);
724     y = fFind->GetXyq(1,i);
725     if (x < xmin || x > xmax) continue;
726     if (y < ymin || y > ymax) continue;
727     digit = fFind->GetIJ(4,i);
728     if (digit >= 0) mdig = AliMUONClusterInput::Instance()->Digit(cath, digit);
729     else mdig = AliMUONClusterInput::Instance()->Digit(TMath::Even(cath), -digit-1);
730     if (mdig->Track(1) >= 0 || mdig->Track(0) >= 10000000) {
731       ix = fHist[iHist]->GetXaxis()->FindBin(x);
732       iy = fHist[iHist]->GetYaxis()->FindBin(y);
733       cont = fHist[iHist]->GetCellContent(ix, iy);
734       hist->Fill(x, y, cont);
735       iok = 1;
736     }
737   }
738   if (iok) { 
739     hist->SetFillColor(5); 
740     hist->SetMaximum(fHist[iHist]->GetMaximum());
741     return hist; 
742   }
743   return 0x0;
744 }