]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterFinderAZ.cxx
308150d9f83e76d7fa0dc94ea473741e52ec1d16
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderAZ.cxx
1 #include "AliMUONClusterFinderAZ.h"
2
3 #include <stdlib.h>
4 #include <fcntl.h>
5 #include <Riostream.h>
6 #include <TROOT.h>
7 #include <TCanvas.h>
8 #include <TLine.h>
9 #include <TTree.h>
10 #include <TH2.h>
11 #include <TView.h>
12 #include <TStyle.h>
13 #include <TMinuit.h>
14 #include <TMatrixD.h>
15
16 #include "AliHeader.h"
17 #include "AliRun.h"
18 #include "AliMUON.h"
19 #include "AliMUONChamber.h"
20 #include "AliMUONDigit.h"
21 #include "AliMUONHit.h"
22 #include "AliMUONChamber.h"
23 #include "AliMUONRawCluster.h"
24 #include "AliMUONClusterInput.h"
25 #include "AliMUONPixel.h"
26 #include "AliMC.h"
27
28 // This function is used for fitting
29 void fcn1(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag);
30
31 ClassImp(AliMUONClusterFinderAZ)
32
33 AliMUONClusterFinderAZ* AliMUONClusterFinderAZ::fgClusterFinder = NULL;
34 TMinuit* AliMUONClusterFinderAZ::fgMinuit = NULL; 
35
36 //_____________________________________________________________________________
37 AliMUONClusterFinderAZ::AliMUONClusterFinderAZ(Bool_t draw=0, Int_t iReco=0)
38 {
39 // Constructor
40   for (Int_t i=0; i<4; i++) {fHist[i] = 0;}
41   fMuonDigits = 0;
42   fSegmentation[1] = fSegmentation[0] = 0; 
43   if (!fgClusterFinder) fgClusterFinder = this;
44   if (!fgMinuit) fgMinuit = new TMinuit(8);
45   fDraw = draw;
46   fReco = iReco;
47   fPixArray = new TObjArray(20); 
48   /*
49   fPoints = 0;
50   fPhits = 0;
51   fRpoints = 0;
52   fCanvas = 0;
53   fNextCathode = kFALSE; 
54   fColPad = 0;
55   */
56 }
57
58 //_____________________________________________________________________________
59 AliMUONClusterFinderAZ::~AliMUONClusterFinderAZ()
60 {
61   // Destructor
62   delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
63   /*
64   // Delete space point structure
65   if (fPoints) fPoints->Delete();
66   delete fPoints;
67   fPoints     = 0;
68   //
69   if (fPhits) fPhits->Delete();
70   delete fPhits;
71   fPhits     = 0;
72   //
73   if (fRpoints) fRpoints->Delete();
74   delete fRpoints;
75   fRpoints     = 0;
76   */
77 }
78
79 //_____________________________________________________________________________
80 void AliMUONClusterFinderAZ::FindRawClusters()
81 {
82 // To provide the same interface as in AliMUONClusterFinderVS
83
84   EventLoop (gAlice->GetHeader()->GetEvent(), AliMUONClusterInput::Instance()->Chamber());
85 }
86
87 //_____________________________________________________________________________
88 void AliMUONClusterFinderAZ::EventLoop(Int_t nev=0, Int_t ch=0)
89 {
90 // Loop over events
91   
92   FILE *lun = 0;
93   TCanvas *c1 = 0;
94   TView *view = 0;
95   TH2F *hist = 0;
96   Double_t p1[3]={0}, p2[3];
97   TTree *TR = 0;
98   if (fDraw) {
99     // File
100     lun = fopen("pool.dat","w");
101     c1 = new TCanvas("c1","Clusters",0,0,600,700);
102     c1->Divide(1,2);
103     new TCanvas("c2","Mlem",700,0,600,350);
104   }
105
106 newev:
107   Int_t nparticles = 0, nent;
108
109   //Loaders
110   AliRunLoader * rl = AliRunLoader::GetRunLoader();
111   AliLoader * gime  = rl->GetLoader("MUONLoader");
112
113   if (!fReco) nparticles = rl->GetEvent(nev);
114   else nparticles = gAlice->GetMCApp()->GetNtrack();
115   cout << "nev         " << nev <<endl;
116   cout << "nparticles  " << nparticles <<endl;
117   if (nparticles <= 0) return;
118   
119   TTree *TH = gime->TreeH();
120   Int_t ntracks = (Int_t) TH->GetEntries();
121   cout<<"ntracks "<<ntracks<<endl;
122     
123   // Get pointers to Alice detectors and Digits containers
124   AliMUON *MUON  = (AliMUON*) gAlice->GetModule("MUON");
125   if (!MUON) return;
126   //       TClonesArray *Particles = gAlice->Particles();
127   if (!fReco) {
128     TR = gime->TreeR();
129     if (TR) {
130       MUON->ResetRawClusters();
131       nent = (Int_t) TR->GetEntries();
132       if (nent != 1) {
133         cout << "Error in MUONdrawClust" << endl;
134         cout << " nent = " <<  nent << " not equal to 1" << endl;
135         //exit(0);
136       }
137     } // if (TR)
138   } // if (!fReco)
139
140   TTree *TD = gime->TreeD();
141   //MUON->ResetDigits();
142
143   TClonesArray *MUONrawclust;
144   AliMUONChamber*       iChamber = 0;
145   
146   // As default draw the first cluster of the chamber #0
147       
148 newchamber:
149   if (ch > 9) {if (fReco) return; nev++; ch = 0; goto newev;}
150   //gAlice->ResetDigits();
151   fMuonDigits  = MUON->GetMUONData()->Digits(ch);
152   if (fMuonDigits == 0) return;
153   iChamber = &(MUON->Chamber(ch));
154   fSegmentation[0] = iChamber->SegmentationModel(1);
155   fSegmentation[1] = iChamber->SegmentationModel(2);
156   fResponse = iChamber->ResponseModel();
157     
158   nent = 0;
159  
160   if (TD) {
161     nent = (Int_t) TD->GetEntries();
162     //printf(" entries %d \n", nent);
163   }
164
165   Int_t ndigits[2]={9,9}, nShown[2]={0};
166   for (Int_t i=0; i<2; i++) {
167     for (Int_t j=0; j<kDim; j++) {fUsed[i][j]=kFALSE;}
168   }
169
170 next:
171   if (ndigits[0] == nShown[0] && ndigits[1] == nShown[1]) {
172     // No more clusters
173     if (fReco) return;
174     ch++;
175     goto newchamber; // next chamber
176   }
177   Float_t xpad, ypad, zpad, zpad0;
178   TLine *line[99]={0};
179   Int_t nLine = 0;
180   Bool_t first = kTRUE;
181   cout << " *** Event # " << nev << " chamber: " << ch << endl;
182   fnPads[0] = fnPads[1] = 0;
183   for (Int_t i=0; i<kDim; i++) {fPadIJ[1][i] = 0;}
184   //for (Int_t iii = 0; iii<999; iii++) { 
185   for (Int_t iii = 0; iii<2; iii++) { 
186     Int_t cath = TMath::Odd(iii);
187     gAlice->ResetDigits();
188     TD->GetEvent(cath);
189     fMuonDigits  = MUON->GetMUONData()->Digits(ch);
190
191     ndigits[cath] = fMuonDigits->GetEntriesFast();
192     if (!ndigits[0] && !ndigits[1]) {if (fReco) return; ch++; goto newchamber;}
193     if (ndigits[cath] == 0) continue;
194     cout << " ndigits: " << ndigits[cath] << " " << cath << endl;
195
196     AliMUONDigit  *mdig;
197     Int_t digit;
198
199     Bool_t EOC = kTRUE; // end-of-cluster
200     for (digit = 0; digit < ndigits[cath]; digit++) {
201       mdig    = (AliMUONDigit*)fMuonDigits->UncheckedAt(digit);
202       if (mdig->Cathode() != cath) continue;
203       if (first) {
204         // Find first unused pad
205         if (fUsed[cath][digit]) continue;
206         fSegmentation[cath]->GetPadC(mdig->PadX(),mdig->PadY(),xpad,ypad,zpad0);
207       } else {
208         if (fUsed[cath][digit]) continue;
209         fSegmentation[cath]->GetPadC(mdig->PadX(),mdig->PadY(),xpad,ypad,zpad);
210         if (TMath::Abs(zpad-zpad0)>0.1) continue; // different slats
211         // Find a pad overlapping with the cluster
212         if (!Overlap(cath,mdig)) continue;
213       }
214       // Add pad - recursive call
215       AddPad(cath,digit);
216       EOC = kFALSE;
217       if (digit >= 0) break;
218     }
219     if (first && EOC) {
220       // No more unused pads 
221       if (cath == 0) continue; // on cathode #0 - check #1
222       else {
223         // No more clusters
224         if (fReco) return;
225         ch++;
226         goto newchamber; // next chamber
227       }
228     }
229     if (EOC) break; // cluster found
230     first = kFALSE;
231     cout << " nPads: " << fnPads[cath] << " " << nShown[cath]+fnPads[cath] << " " << cath << endl;
232   } // for (Int_t iii = 0;
233
234   
235   if (fReco) goto skip;
236   char hName[4];
237   for (Int_t cath = 0; cath<2; cath++) {
238     // Build histograms
239     if (fHist[cath*2]) {fHist[cath*2]->Delete(); fHist[cath*2] = 0;}
240     if (fHist[cath*2+1]) {fHist[cath*2+1]->Delete(); fHist[cath*2+1] = 0;}
241     if (fnPads[cath] == 0) continue; // cluster on one cathode only
242     Float_t wxMin=999, wxMax=0, wyMin=999, wyMax=0; 
243     Int_t minDx=0, maxDx=0, minDy=0, maxDy=0;
244     for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
245       if (fPadIJ[0][i] != cath) continue;
246       if (fXyq[3][i] < wxMin) {wxMin = fXyq[3][i]; minDx = i;}
247       if (fXyq[3][i] > wxMax) {wxMax = fXyq[3][i]; maxDx = i;}
248       if (fXyq[4][i] < wyMin) {wyMin = fXyq[4][i]; minDy = i;}
249       if (fXyq[4][i] > wyMax) {wyMax = fXyq[4][i]; maxDy = i;}
250     }
251     cout << minDx << maxDx << minDy << maxDy << endl;
252     Int_t nx, ny, padSize;
253     Float_t xmin=9999, xmax=-9999, ymin=9999, ymax=-9999;
254     if (TMath::Nint(fXyq[3][minDx]*1000) == TMath::Nint(fXyq[3][maxDx]*1000) &&
255         TMath::Nint(fXyq[4][minDy]*1000) == TMath::Nint(fXyq[4][maxDy]*1000)) {
256       // the same segmentation
257       cout << " Same" << endl;
258       cout << fXyq[3][minDx] << " " << fXyq[3][maxDx] << " " << fXyq[4][minDy] << " " << fXyq[4][maxDy] << endl;
259       for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
260         if (fPadIJ[0][i] != cath) continue;
261         if (fXyq[0][i] < xmin) xmin = fXyq[0][i];
262         if (fXyq[0][i] > xmax) xmax = fXyq[0][i];
263         if (fXyq[1][i] < ymin) ymin = fXyq[1][i];
264         if (fXyq[1][i] > ymax) ymax = fXyq[1][i];
265       }
266       xmin -= fXyq[3][minDx]; xmax += fXyq[3][minDx];
267       ymin -= fXyq[4][minDy]; ymax += fXyq[4][minDy];
268       nx = TMath::Nint ((xmax-xmin)/wxMin/2);
269       ny = TMath::Nint ((ymax-ymin)/wyMin/2);
270       sprintf(hName,"h%d",cath*2);
271       fHist[cath*2] = new TH2F(hName,"cluster",nx,xmin,xmax,ny,ymin,ymax);
272       cout << fHist[cath*2] << " " << fnPads[cath] << endl;
273       for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
274         if (fPadIJ[0][i] != cath) continue;
275         fHist[cath*2]->Fill(fXyq[0][i],fXyq[1][i],fXyq[2][i]);
276         //cout << fXyq[0][i] << fXyq[1][i] << fXyq[2][i] << endl;
277       }
278     } else {
279       // different segmentation in the cluster
280       cout << " Different" << endl;
281       cout << fXyq[3][minDx] << " " << fXyq[3][maxDx] << " " << fXyq[4][minDy] << " " << fXyq[4][maxDy] << endl;
282       Int_t nOK = 0;
283       Int_t indx, locMin, locMax;
284       if (TMath::Nint(fXyq[3][minDx]*1000) != TMath::Nint(fXyq[3][maxDx]*1000)) {
285         // different segmentation along x
286         indx = 0;
287         locMin = minDx;
288         locMax = maxDx;
289       } else {
290         // different segmentation along y
291         indx = 1;
292         locMin = minDy;
293         locMax = maxDy;
294       }
295       Int_t loc = locMin;
296       for (Int_t i=0; i<2; i++) {
297         // loop over different pad sizes
298         if (i>0) loc = locMax;
299         padSize = TMath::Nint(fXyq[indx+3][loc]*1000);
300         xmin = 9999; xmax = -9999; ymin = 9999; ymax = -9999;
301         for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
302           if (fPadIJ[0][j] != cath) continue;
303           if (TMath::Nint(fXyq[indx+3][j]*1000) != padSize) continue;
304           nOK++;
305           xmin = TMath::Min (xmin,fXyq[0][j]);
306           xmax = TMath::Max (xmax,fXyq[0][j]);
307           ymin = TMath::Min (ymin,fXyq[1][j]);
308           ymax = TMath::Max (ymax,fXyq[1][j]);
309         }
310         xmin -= fXyq[3][loc]; xmax += fXyq[3][loc];
311         ymin -= fXyq[4][loc]; ymax += fXyq[4][loc];
312         nx = TMath::Nint ((xmax-xmin)/fXyq[3][loc]/2);
313         ny = TMath::Nint ((ymax-ymin)/fXyq[4][loc]/2);
314         sprintf(hName,"h%d",cath*2+i);
315         fHist[cath*2+i] = new TH2F(hName,"cluster",nx,xmin,xmax,ny,ymin,ymax);
316         for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
317           if (fPadIJ[0][j] != cath) continue;
318           if (TMath::Nint(fXyq[indx+3][j]*1000) != padSize) continue;
319           fHist[cath*2+i]->Fill(fXyq[0][j],fXyq[1][j],fXyq[2][j]);
320         }
321       } // for (Int_t i=0;
322       if (nOK != fnPads[cath]) cout << " *** Too many segmentations: nPads, nOK " << fnPads[cath] << " " << nOK << endl;
323     } // if (TMath::Nint(fXyq[3][minDx]*1000)
324   } // for (Int_t cath = 0;
325         
326   // Draw histograms and coordinates
327   for (Int_t cath=0; cath<2; cath++) {
328     if (cath == 0) ModifyHistos();
329     if (fnPads[cath] == 0) continue; // cluster on one cathode only
330     if (fDraw) {
331       c1->cd(cath+1);
332       gPad->SetTheta(55);
333       gPad->SetPhi(30);
334       Double_t x, y, x0, y0, r1=999, r2=0;
335       if (fHist[cath*2+1]) {
336         // 
337         x0 = fHist[cath*2]->GetXaxis()->GetXmin() - 1000*TMath::Cos(30*TMath::Pi()/180);
338         y0 = fHist[cath*2]->GetYaxis()->GetXmin() - 1000*TMath::Sin(30*TMath::Pi()/180);
339         r1 = 0;
340         Int_t ihist=cath*2;
341         for (Int_t iy=1; iy<=fHist[ihist]->GetNbinsY(); iy++) {
342           y = fHist[ihist]->GetYaxis()->GetBinCenter(iy) 
343             + fHist[ihist]->GetYaxis()->GetBinWidth(iy);
344           for (Int_t ix=1; ix<=fHist[ihist]->GetNbinsX(); ix++) {
345             if (fHist[ihist]->GetCellContent(ix,iy) > 0.1) {
346               x = fHist[ihist]->GetXaxis()->GetBinCenter(ix)
347                 + fHist[ihist]->GetXaxis()->GetBinWidth(ix);
348               r1 = TMath::Max (r1,TMath::Sqrt((x-x0)*(x-x0)+(y-y0)*(y-y0)));
349             }
350           }
351         }
352         ihist = cath*2 + 1 ;
353         for (Int_t iy=1; iy<=fHist[ihist]->GetNbinsY(); iy++) {
354           y = fHist[ihist]->GetYaxis()->GetBinCenter(iy)
355             + fHist[ihist]->GetYaxis()->GetBinWidth(iy);
356           for (Int_t ix=1; ix<=fHist[ihist]->GetNbinsX(); ix++) {
357             if (fHist[ihist]->GetCellContent(ix,iy) > 0.1) {
358               x = fHist[ihist]->GetXaxis()->GetBinCenter(ix)
359                 + fHist[ihist]->GetXaxis()->GetBinWidth(ix);
360               r2 = TMath::Max (r2,TMath::Sqrt((x-x0)*(x-x0)+(y-y0)*(y-y0)));
361             }
362           }
363         }
364         cout << r1 << " " << r2 << endl;
365       } // if (fHist[cath*2+1])
366       if (r1 > r2) {
367         //fHist[cath*2]->Draw("lego1");
368         fHist[cath*2]->Draw("lego1Fb");
369         //if (fHist[cath*2+1]) fHist[cath*2+1]->Draw("lego1SameAxisBb");
370         if (fHist[cath*2+1]) fHist[cath*2+1]->Draw("lego1SameAxisBbFb");
371       } else {
372         //fHist[cath*2+1]->Draw("lego1");
373         fHist[cath*2+1]->Draw("lego1Fb");
374         //fHist[cath*2]->Draw("lego1SameAxisBb");
375         fHist[cath*2]->Draw("lego1SameAxisFbBb");
376       }
377       c1->Update();
378     } // if (fDraw)
379   } // for (Int_t cath = 0;
380
381   // Draw generated hits
382   Double_t xNDC[6];
383   hist = fHist[0] ? fHist[0] : fHist[2];
384   p2[2] = hist->GetMaximum();
385   view = 0;
386   if (c1) view = c1->Pad()->GetView();
387   cout << " *** GEANT hits *** " << endl;
388   fnMu = 0;
389   Int_t ix, iy, iok;
390   for (Int_t i=0; i<ntracks; i++) {
391     TH->GetEvent(i);
392     for (AliMUONHit* mHit=(AliMUONHit*)MUON->FirstHit(-1); 
393          mHit;
394          mHit=(AliMUONHit*)MUON->NextHit()) {
395       if (mHit->Chamber() != ch+1) continue;  // chamber number
396       if (TMath::Abs(mHit->Z()-zpad0) > 1) continue; // different slat
397       p2[0] = p1[0] = mHit->X();        // x-pos of hit
398       p2[1] = p1[1] = mHit->Y();        // y-pos
399       if (p1[0] < hist->GetXaxis()->GetXmin() || 
400           p1[0] > hist->GetXaxis()->GetXmax()) continue;
401       if (p1[1] < hist->GetYaxis()->GetXmin() || 
402           p1[1] > hist->GetYaxis()->GetXmax()) continue;
403       // Check if track comes thru pads with signal
404       iok = 0;
405       for (Int_t ihist=0; ihist<4; ihist++) {
406         if (!fHist[ihist]) continue;
407         ix = fHist[ihist]->GetXaxis()->FindBin(p1[0]);
408         iy = fHist[ihist]->GetYaxis()->FindBin(p1[1]);
409         if (fHist[ihist]->GetCellContent(ix,iy) > 0.5) {iok = 1; break;}
410       }
411       if (!iok) continue;
412       gStyle->SetLineColor(1);
413       if (TMath::Abs((Int_t)mHit->Particle()) == 13) {
414         gStyle->SetLineColor(4);
415         fnMu++;
416         if (fnMu <= 2) {
417           fxyMu[fnMu-1][0] = p1[0];
418           fxyMu[fnMu-1][1] = p1[1];
419         }
420       }     
421       printf(" X=%10.4f, Y=%10.4f, Z=%10.4f\n",p1[0],p1[1],mHit->Z());
422       if (view) {
423         view->WCtoNDC(p1, &xNDC[0]);
424         view->WCtoNDC(p2, &xNDC[3]);
425         for (Int_t ipad=1; ipad<3; ipad++) {
426           c1->cd(ipad);
427           //c1->DrawLine(xpad[0],xpad[1],xpad[3],xpad[4]);
428           line[nLine] = new TLine(xNDC[0],xNDC[1],xNDC[3],xNDC[4]);
429           line[nLine++]->Draw();
430         }
431       }
432     } // for (AliMUONHit* mHit=
433   } // for (Int_t i=0; i<ntracks;
434
435   // Draw reconstructed coordinates
436   MUONrawclust = MUON->GetMUONData()->RawClusters(ch);
437   TR->GetEvent(ch);
438   //cout << MUONrawclust << " " << MUONrawclust->GetEntries() << endl;
439   AliMUONRawCluster *mRaw;
440   gStyle->SetLineColor(3);
441   cout << " *** Reconstructed hits *** " << endl;
442   for (Int_t i=0; i<MUONrawclust->GetEntries(); i++) {
443     mRaw = (AliMUONRawCluster*)MUONrawclust->UncheckedAt(i);
444     if (TMath::Abs(mRaw->fZ[0]-zpad0) > 1) continue; // different slat
445     p2[0] = p1[0] = mRaw->fX[0];        // x-pos of hit
446     p2[1] = p1[1] = mRaw->fY[0];        // y-pos
447     if (p1[0] < hist->GetXaxis()->GetXmin() || 
448         p1[0] > hist->GetXaxis()->GetXmax()) continue;
449     if (p1[1] < hist->GetYaxis()->GetXmin() || 
450         p1[1] > hist->GetYaxis()->GetXmax()) continue;
451     /*
452       TD->GetEvent(cath);
453       cout << mRaw->fMultiplicity[0] << mRaw->fMultiplicity[1] << endl;
454       for (Int_t j=0; j<mRaw->fMultiplicity[cath]; j++) {
455       Int_t digit = mRaw->fIndexMap[j][cath];
456       cout << ((AliMUONDigit*)fMuonDigits->UncheckedAt(digit))->Signal() << endl;
457       }
458     */
459     // Check if track comes thru pads with signal
460     iok = 0;
461     for (Int_t ihist=0; ihist<4; ihist++) {
462       if (!fHist[ihist]) continue;
463       ix = fHist[ihist]->GetXaxis()->FindBin(p1[0]);
464       iy = fHist[ihist]->GetYaxis()->FindBin(p1[1]);
465       if (fHist[ihist]->GetCellContent(ix,iy) > 0.5) {iok = 1; break;}
466     }
467     if (!iok) continue;
468     printf(" X=%10.4f, Y=%10.4f, Z=%10.4f\n",p1[0],p1[1],mRaw->fZ[0]);
469     if (view) {
470       view->WCtoNDC(p1, &xNDC[0]);
471       view->WCtoNDC(p2, &xNDC[3]);
472       for (Int_t ipad=1; ipad<3; ipad++) {
473         c1->cd(ipad);
474         line[nLine] = new TLine(xNDC[0],xNDC[1],xNDC[3],xNDC[4]);
475         line[nLine++]->Draw();
476       }
477     }
478   } // for (Int_t i=0; i<MUONrawclust->GetEntries();
479   if (fDraw) c1->Update();
480
481 skip:
482   // Use MLEM for cluster finder
483   fZpad = zpad0;
484   Int_t nMax = 1, localMax[100], maxPos[100];
485   Double_t maxVal[100];
486   
487   if (CheckPrecluster(nShown)) {
488     BuildPixArray();
489     if (fnPads[0]+fnPads[1] > 50) nMax = FindLocalMaxima(localMax, maxVal);
490     if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
491     for (Int_t i=0; i<nMax; i++) {
492       if (nMax > 1) FindCluster(localMax, maxPos[i]);
493       if (!MainLoop()) cout << " MainLoop failed " << endl;
494       if (i < nMax-1) {
495         for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
496           if (fPadIJ[1][j] == 0) continue; // pad charge was not modified
497           fPadIJ[1][j] = 0;
498           fXyq[2][j] = fXyq[5][j]; // use backup charge value
499         }
500       }
501     }
502   }
503   if (fReco) goto next;
504
505   for (Int_t i=0; i<fnMu; i++) {
506     // Check again if muon come thru the used pads (due to extra splitting)
507     for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
508       if (TMath::Abs(fxyMu[i][0]-fXyq[0][j])<fXyq[3][j] && 
509           TMath::Abs(fxyMu[i][1]-fXyq[1][j])<fXyq[4][j]) {
510         printf("%12.3e %12.3e %12.3e %12.3e\n",fxyMu[i][2],fxyMu[i][3],fxyMu[i][4],fxyMu[i][5]);
511         if (lun) fprintf(lun,"%4d %2d %12.3e %12.3e %12.3e %12.3e\n",nev,ch,fxyMu[i][2],fxyMu[i][3],fxyMu[i][4],fxyMu[i][5]);
512         break;
513       }
514     }
515   } // for (Int_t i=0; i<fnMu;
516
517   // What's next?
518   char command[8];
519   cout << " What is next? " << endl;
520   command[0] = ' '; 
521   if (fDraw) gets(command);
522   if (command[0] == 'n' || command[0] == 'N') {nev++; goto newev;} // next event 
523   else if (command[0] == 'q' || command[0] == 'Q') {fclose(lun); return;} // exit display 
524   //else if (command[0] == 'r' || command[0] == 'R') goto redraw; // redraw points
525   else if (command[0] == 'c' || command[0] == 'C') {
526     // new chamber
527     sscanf(command+1,"%d",&ch);
528     goto newchamber;
529   } 
530   else if (command[0] == 'e' || command[0] == 'E') {
531     // new event
532     sscanf(command+1,"%d",&nev);
533     goto newev;
534   } 
535   else goto next; // Next cluster
536 }
537
538 //_____________________________________________________________________________
539 void AliMUONClusterFinderAZ::ModifyHistos(void)
540 {
541   // Modify histograms to bring them to the same size
542   Int_t nhist = 0;
543   Float_t hlim[4][4], hbin[4][4]; // first index - xmin, xmax, ymin, ymax
544   Float_t binMin[4] = {999,999,999,999};
545
546   for (Int_t i=0; i<4; i++) {
547     if (!fHist[i]) continue;
548     hlim[0][nhist] = fHist[i]->GetXaxis()->GetXmin(); // xmin
549     hlim[1][nhist] = fHist[i]->GetXaxis()->GetXmax(); // xmax
550     hlim[2][nhist] = fHist[i]->GetYaxis()->GetXmin(); // ymin
551     hlim[3][nhist] = fHist[i]->GetYaxis()->GetXmax(); // ymax
552     hbin[0][nhist] = hbin[1][nhist] = fHist[i]->GetXaxis()->GetBinWidth(1);
553     hbin[2][nhist] = hbin[3][nhist] = fHist[i]->GetYaxis()->GetBinWidth(1);
554     binMin[0] = TMath::Min(binMin[0],hbin[0][nhist]);
555     binMin[2] = TMath::Min(binMin[2],hbin[2][nhist]);
556     nhist++;
557   }
558   binMin[1] = binMin[0];
559   binMin[3] = binMin[2];
560   cout << " Nhist: " << nhist << endl;
561
562   Int_t imin, imax;
563   for (Int_t lim=0; lim<4; lim++) {
564     while (1) {
565       imin = TMath::LocMin(nhist,hlim[lim]);
566       imax = TMath::LocMax(nhist,hlim[lim]);
567       if (TMath::Abs(hlim[lim][imin]-hlim[lim][imax])<0.01*binMin[lim]) break;
568       if (lim == 0 || lim == 2) {
569         // find lower limit
570         hlim[lim][imax] -= hbin[lim][imax];
571       } else {
572         // find upper limit
573         hlim[lim][imin] += hbin[lim][imin];
574       }
575     } // while (1)
576   }
577     
578   // Rebuild histograms 
579   nhist = 0;
580   TH2F *hist = 0;
581   Int_t nx, ny;
582   Double_t x, y, cont, cmax=0;
583   char hName[4];
584   for (Int_t ihist=0; ihist<4; ihist++) {
585     if (!fHist[ihist]) continue;
586     nx = TMath::Nint((hlim[1][nhist]-hlim[0][nhist])/hbin[0][nhist]);
587     ny = TMath::Nint((hlim[3][nhist]-hlim[2][nhist])/hbin[2][nhist]);
588     //hist =  new TH2F("h","hist",nx,hlim[0][nhist],hlim[1][nhist],ny,hlim[2][nhist],hlim[3][nhist]);
589     sprintf(hName,"hh%d",ihist);
590     hist =  new TH2F(hName,"hist",nx,hlim[0][nhist],hlim[1][nhist],ny,hlim[2][nhist],hlim[3][nhist]);
591     for (Int_t i=1; i<=fHist[ihist]->GetNbinsX(); i++) {
592       x = fHist[ihist]->GetXaxis()->GetBinCenter(i);
593       for (Int_t j=1; j<=fHist[ihist]->GetNbinsY(); j++) {
594         y = fHist[ihist]->GetYaxis()->GetBinCenter(j);
595         cont = fHist[ihist]->GetCellContent(i,j);
596         hist->Fill(x,y,cont);
597       }
598     }
599     cmax = TMath::Max (cmax,hist->GetMaximum());
600     fHist[ihist]->Delete();
601     fHist[ihist] = new TH2F(*hist);
602     hist->Delete(); 
603     nhist++;
604   }
605   printf("%f \n",cmax);
606
607   for (Int_t ihist=0; ihist<4; ihist++) {
608     if (!fHist[ihist]) continue;
609     fHist[ihist]->SetMaximum(cmax);
610   }
611 }
612
613 //_____________________________________________________________________________
614 void AliMUONClusterFinderAZ::AddPad(Int_t cath, Int_t digit)
615 {
616   // Add pad to the cluster
617   AliMUONDigit *mdig = (AliMUONDigit*)fMuonDigits->UncheckedAt(digit);
618
619   Int_t charge = mdig->Signal();
620   // get the center of the pad
621   Float_t xpad, ypad, zpad;
622   fSegmentation[cath]->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
623   
624   Int_t   isec = fSegmentation[cath]->Sector(mdig->PadX(), mdig->PadY());
625   Int_t nPads = fnPads[0] + fnPads[1];
626   fXyq[0][nPads] = xpad;
627   fXyq[1][nPads] = ypad;
628   fXyq[2][nPads] = charge;
629   fXyq[3][nPads] = fSegmentation[cath]->Dpx(isec)/2;
630   fXyq[4][nPads] = fSegmentation[cath]->Dpy(isec)/2;
631   fXyq[5][nPads] = digit;
632   fPadIJ[0][nPads] = cath;
633   fPadIJ[1][nPads] = 0;
634   fUsed[cath][digit] = kTRUE;
635   //cout << " bbb " << fXyq[cath][2][nPads] << " " << fXyq[cath][0][nPads] << " " << fXyq[cath][1][nPads] << " " << fXyq[cath][3][nPads] << " " << fXyq[cath][4][nPads] << " " << zpad << " " << nPads << endl;
636   fnPads[cath]++;
637
638   // Check neighbours
639   Int_t nn, ix, iy, xList[10], yList[10];
640   AliMUONDigit  *mdig1;
641
642   Int_t ndigits = fMuonDigits->GetEntriesFast();
643   fSegmentation[cath]->Neighbours(mdig->PadX(),mdig->PadY(),&nn,xList,yList); 
644   for (Int_t in=0; in<nn; in++) {
645     ix=xList[in];
646     iy=yList[in];
647     for (Int_t digit1 = 0; digit1 < ndigits; digit1++) {
648       if (digit1 == digit) continue;
649       mdig1 = (AliMUONDigit*)fMuonDigits->UncheckedAt(digit1);
650       if (mdig1->Cathode() != cath) continue;
651       if (!fUsed[cath][digit1] && mdig1->PadX() == ix && mdig1->PadY() == iy) {
652         fUsed[cath][digit1] = kTRUE;
653         // Add pad - recursive call
654         AddPad(cath,digit1);
655       }
656     } //for (Int_t digit1 = 0;
657   } // for (Int_t in=0;
658 }
659
660 //_____________________________________________________________________________
661 Bool_t AliMUONClusterFinderAZ::Overlap(Int_t cath, TObject *dig)
662 {
663   // Check if the pad from one cathode overlaps with a pad 
664   // in the precluster on the other cathode
665
666   AliMUONDigit *mdig = (AliMUONDigit*) dig;
667
668   Float_t xpad, ypad, zpad;
669   fSegmentation[cath]->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
670   Int_t   isec = fSegmentation[cath]->Sector(mdig->PadX(), mdig->PadY());
671
672   Float_t xy1[4], xy12[4];
673   xy1[0] = xpad - fSegmentation[cath]->Dpx(isec)/2;
674   xy1[1] = xy1[0] + fSegmentation[cath]->Dpx(isec);
675   xy1[2] = ypad - fSegmentation[cath]->Dpy(isec)/2;
676   xy1[3] = xy1[2] + fSegmentation[cath]->Dpy(isec);
677   //cout << " ok " << fnPads[0]+fnPads[1] << xy1[0] << xy1[1] << xy1[2] << xy1[3] << endl;
678
679   Int_t cath1 = TMath::Even(cath);
680   for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
681     if (fPadIJ[0][i] != cath1) continue;
682     if (Overlap(xy1, i, xy12, 0)) return kTRUE;
683   }
684   return kFALSE;
685 }
686
687 //_____________________________________________________________________________
688 Bool_t AliMUONClusterFinderAZ::Overlap(Float_t *xy1, Int_t iPad, Float_t *xy12, Int_t iSkip)
689 {
690   // Check if the pads xy1 and iPad overlap and return overlap area
691
692   Float_t xy2[4];
693   xy2[0] = fXyq[0][iPad] - fXyq[3][iPad];
694   xy2[1] = fXyq[0][iPad] + fXyq[3][iPad];
695   if (xy1[0] > xy2[1]-1.e-4 || xy1[1] < xy2[0]+1.e-4) return kFALSE;
696   xy2[2] = fXyq[1][iPad] - fXyq[4][iPad];
697   xy2[3] = fXyq[1][iPad] + fXyq[4][iPad];
698   if (xy1[2] > xy2[3]-1.e-4 || xy1[3] < xy2[2]+1.e-4) return kFALSE;
699   if (!iSkip) return kTRUE; // just check overlap (w/out computing the area)
700   xy12[0] = TMath::Max (xy1[0],xy2[0]);
701   xy12[1] = TMath::Min (xy1[1],xy2[1]);
702   xy12[2] = TMath::Max (xy1[2],xy2[2]);
703   xy12[3] = TMath::Min (xy1[3],xy2[3]);
704   return kTRUE;
705 }
706
707 //_____________________________________________________________________________
708 /*
709 Bool_t AliMUONClusterFinderAZ::Overlap(Int_t i, Int_t j, Float_t *xy12, Int_t iSkip)
710 {
711   // Check if the pads i and j overlap and return overlap area
712
713   Float_t xy1[4], xy2[4];
714   return Overlap(xy1, xy2, xy12, iSkip);
715 }
716 */
717 //_____________________________________________________________________________
718 Bool_t AliMUONClusterFinderAZ::CheckPrecluster(Int_t *nShown)
719 {
720   // Check precluster in order to attempt to simplify it (mostly for
721   // two-cathode preclusters)
722
723   Int_t i1, i2;
724   Float_t xy1[4], xy12[4];
725   
726   Int_t npad = fnPads[0] + fnPads[1];
727
728   // If pads have the same size take average of pads on both cathodes 
729   Int_t sameSize = (fnPads[0] && fnPads[1]) ? 1 : 0;
730   if (sameSize) {
731     Double_t xSize = -1, ySize = 0;
732     for (Int_t i=0; i<npad; i++) {
733       if (fXyq[2][i] < 0) continue;
734       if (xSize < 0) { xSize = fXyq[3][i]; ySize = fXyq[4][i]; }
735       if (TMath::Abs(xSize-fXyq[3][i]) > 1.e-4 ||  TMath::Abs(ySize-fXyq[4][i]) > 1.e-4) { sameSize = 0; break; }
736     }
737   } // if (sameSize)
738   if (sameSize && (fnPads[0] > 2 || fnPads[1] > 2)) {
739     nShown[0] += fnPads[0];
740     nShown[1] += fnPads[1];
741     fnPads[0] = fnPads[1] = 0;
742     Int_t div;
743     for (Int_t i=0; i<npad; i++) {
744       if (fXyq[2][i] < 0) continue; // used pad
745       fXyq[2][fnPads[0]] = fXyq[2][i];
746       div = 1;
747       for (Int_t j=i+1; j<npad; j++) {
748         if (fPadIJ[0][j] == fPadIJ[0][i]) continue; // same cathode
749         if (TMath::Abs(fXyq[0][j]-fXyq[0][i]) > 1.e-4) continue;
750         if (TMath::Abs(fXyq[1][j]-fXyq[1][i]) > 1.e-4) continue;
751         fXyq[2][fnPads[0]] += fXyq[2][j];
752         div = 2;
753         fXyq[2][j] = -2;
754         break;
755       }
756       fXyq[2][fnPads[0]] /= div;
757       fXyq[0][fnPads[0]] = fXyq[0][i];
758       fXyq[1][fnPads[0]] = fXyq[1][i];
759       fPadIJ[0][fnPads[0]++] = 0;
760     }
761   } // if (sameSize)
762
763   // Check if one-cathode precluster
764   i1 = fnPads[0]!=0 ? 0 : 1;
765   i2 = fnPads[1]!=0 ? 1 : 0;
766
767   if (i1 != i2) { // two-cathode 
768
769     Int_t *flags = new Int_t[npad];
770     for (Int_t i=0; i<npad; i++) { flags[i] = 0; }
771
772     // Check pad overlaps
773     for (Int_t i=0; i<npad; i++) {
774       if (fPadIJ[0][i] != i1) continue;
775       xy1[0] = fXyq[0][i] - fXyq[3][i];
776       xy1[1] = fXyq[0][i] + fXyq[3][i];
777       xy1[2] = fXyq[1][i] - fXyq[4][i];
778       xy1[3] = fXyq[1][i] + fXyq[4][i];
779       for (Int_t j=0; j<npad; j++) {
780         if (fPadIJ[0][j] != i2) continue;
781         if (!Overlap(xy1, j, xy12, 0)) continue;
782         flags[i] = flags[j] = 1; // mark overlapped pads
783       } // for (Int_t j=0;
784     } // for (Int_t i=0;
785
786     // Check if all pads overlap
787     Int_t digit=0, cath, nFlags=0;
788     for (Int_t i=0; i<npad; i++) {nFlags += !flags[i];}
789     if (nFlags) cout << " nFlags = " << nFlags << endl;
790     //if (nFlags > 2 || (Float_t)nFlags / npad > 0.2) { // why 2 ??? - empirical choice
791     if (nFlags > 0) {
792       for (Int_t i=0; i<npad; i++) {
793         if (flags[i]) continue;
794         digit = TMath::Nint (fXyq[5][i]);
795         cath = fPadIJ[0][i];
796         fUsed[cath][digit] = kFALSE; // release pad
797         fXyq[2][i] = -2;
798         fnPads[cath]--;
799       }
800     } // if (nFlags > 2)
801
802     // Check correlations of cathode charges
803     if (fnPads[0] && fnPads[1]) { // two-cathode
804       Double_t sum[2]={0};
805       Int_t over[2] = {1, 1};
806       for (Int_t i=0; i<npad; i++) {
807         cath = fPadIJ[0][i];
808         if (fXyq[2][i] > 0) sum[cath] += fXyq[2][i];
809         if (fXyq[2][i] > fResponse->MaxAdc()-1) over[cath] = 0;
810       }
811       cout << " Total charge: " << sum[0] << " " << sum[1] << endl;
812       if ((over[0] || over[1]) && TMath::Abs(sum[0]-sum[1])/(sum[0]+sum[1])*2 > 1) { // 3 times difference
813         cout << " Release " << endl;
814         // Big difference
815         cath = sum[0]>sum[1] ? 0 : 1;
816         Int_t imax = 0;
817         Double_t cmax=-1;
818         Double_t *dist = new Double_t[npad];
819         for (Int_t i=0; i<npad; i++) {
820           if (fPadIJ[0][i] != cath) continue;
821           if (fXyq[2][i] < cmax) continue;
822           cmax = fXyq[2][i];
823           imax = i;
824         }
825         // Arrange pads according to their distance to the max, 
826         // normalized to the pad size
827         for (Int_t i=0; i<npad; i++) {
828           dist[i] = 0;
829           if (fPadIJ[0][i] != cath) continue;
830           if (i == imax) continue; 
831           if (fXyq[2][i] < 0) continue;
832           dist[i] = (fXyq[0][i]-fXyq[0][imax])*(fXyq[0][i]-fXyq[0][imax])/
833                      fXyq[3][imax]/fXyq[3][imax]/4;
834           dist[i] += (fXyq[1][i]-fXyq[1][imax])*(fXyq[1][i]-fXyq[1][imax])/
835                       fXyq[4][imax]/fXyq[4][imax]/4;
836           dist[i] = TMath::Sqrt (dist[i]);
837         }
838         TMath::Sort(npad, dist, flags, kFALSE); // in increasing order
839         Int_t indx;
840         Double_t xmax = -1;
841         for (Int_t i=0; i<npad; i++) {
842           indx = flags[i];
843           if (fPadIJ[0][indx] != cath) continue;
844           if (fXyq[2][indx] < 0) continue;
845           if (fXyq[2][indx] <= cmax || TMath::Abs(dist[indx]-xmax)<1.e-3) {
846             // Release pads
847             if (TMath::Abs(dist[indx]-xmax)<1.e-3) 
848                 cmax = TMath::Max((Double_t)(fXyq[2][indx]),cmax);
849             else cmax = fXyq[2][indx];
850             xmax = dist[indx];
851             digit = TMath::Nint (fXyq[5][indx]);
852             fUsed[cath][digit] = kFALSE; 
853             fXyq[2][indx] = -2;
854             fnPads[cath]--;
855             // xmax = dist[i]; // Bug?
856           }
857           else break;
858         } 
859         delete [] dist; dist = 0;
860       } // TMath::Abs(sum[0]-sum[1])...
861     } // if (fnPads[0] && fnPads[1])
862     delete [] flags; flags = 0;
863   } // if (i1 != i2) 
864
865   if (!sameSize) { nShown[0] += fnPads[0]; nShown[1] += fnPads[1]; }
866
867   // Move released pads to the right
868   Int_t beg = 0, end = npad-1, padij;
869   Double_t xyq;
870   while (beg < end) {
871     if (fXyq[2][beg] > 0) { beg++; continue; }
872     for (Int_t j=end; j>beg; j--) {
873       if (fXyq[2][j] < 0) continue;
874       end = j - 1;
875       for (Int_t j1=0; j1<2; j1++) {
876         padij = fPadIJ[j1][beg]; 
877         fPadIJ[j1][beg] = fPadIJ[j1][j];
878         fPadIJ[j1][j] = padij;
879       }
880       for (Int_t j1=0; j1<6; j1++) {
881         xyq = fXyq[j1][beg]; 
882         fXyq[j1][beg] = fXyq[j1][j];
883         fXyq[j1][j] = xyq;
884       }
885       break;
886     } // for (Int_t j=end;
887     beg++;
888   } // while
889   npad = fnPads[0] + fnPads[1];
890   if (npad > 500) { cout << " ***** Too large cluster. Give up. " << npad << endl; return kFALSE; }
891   // Back up charge value
892   for (Int_t j=0; j<npad; j++) fXyq[5][j] = fXyq[2][j];
893
894   return kTRUE;
895 }
896
897 //_____________________________________________________________________________
898 void AliMUONClusterFinderAZ::BuildPixArray()
899 {
900   // Build pixel array for MLEM method
901   
902   Int_t nPix=0, i1, i2;
903   Float_t xy1[4], xy12[4];
904   AliMUONPixel *pixPtr=0;
905
906   Int_t npad = fnPads[0] + fnPads[1];
907
908   // One cathode is empty
909   i1 = fnPads[0]!=0 ? 0 : 1;
910   i2 = fnPads[1]!=0 ? 1 : 0;
911
912   // Build array of pixels on anode plane
913   if (i1 == i2) { // one-cathode precluster
914     for (Int_t j=0; j<npad; j++) {
915       pixPtr = new AliMUONPixel();
916       for (Int_t i=0; i<2; i++) {
917         pixPtr->SetCoord(i, fXyq[i][j]); // pixel coordinates
918         pixPtr->SetSize(i, fXyq[i+3][j]); // pixel size
919       }
920       pixPtr->SetCharge(fXyq[2][j]); // charge
921       fPixArray->Add((TObject*)pixPtr);
922       nPix++;
923     }
924   } else { // two-cathode precluster    
925     for (Int_t i=0; i<npad; i++) {
926       if (fPadIJ[0][i] != i1) continue;
927       xy1[0] = fXyq[0][i] - fXyq[3][i];
928       xy1[1] = fXyq[0][i] + fXyq[3][i];
929       xy1[2] = fXyq[1][i] - fXyq[4][i];
930       xy1[3] = fXyq[1][i] + fXyq[4][i];
931       for (Int_t j=0; j<npad; j++) {
932         if (fPadIJ[0][j] != i2) continue;
933         if (!Overlap(xy1, j, xy12, 1)) continue;
934         pixPtr = new AliMUONPixel();
935         for (Int_t k=0; k<2; k++) {
936           pixPtr->SetCoord(k, (xy12[2*k]+xy12[2*k+1])/2); // pixel coordinates
937           pixPtr->SetSize(k, xy12[2*k+1]-pixPtr->Coord(k)); // size
938         }
939         pixPtr->SetCharge(TMath::Min (fXyq[2][i],fXyq[2][j])); //charge
940         fPixArray->Add((TObject*)pixPtr);
941         nPix++;
942       } // for (Int_t j=0;
943     } // for (Int_t i=0;
944   } // else
945
946   Float_t wxmin=999, wymin=999;
947   for (Int_t i=0; i<npad; i++) {
948     if (fPadIJ[0][i] == i1) wymin = TMath::Min (wymin,fXyq[4][i]);
949     if (fPadIJ[0][i] == i2) wxmin = TMath::Min (wxmin,fXyq[3][i]);
950   }
951   cout << wxmin << " " << wymin << endl;
952
953   // Check if small pixel X-size
954   AjustPixel(wxmin, 0);
955   // Check if small pixel Y-size
956   AjustPixel(wymin, 1);
957   // Check if large pixel size
958   AjustPixel(wxmin, wymin);
959
960   // Remove discarded pixels
961   for (Int_t i=0; i<nPix; i++) {
962     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
963     //pixPtr->Print();
964     if (pixPtr->Charge() < 1) { fPixArray->RemoveAt(i); delete pixPtr; }// discarded pixel
965   }
966   fPixArray->Compress();
967   nPix = fPixArray->GetEntriesFast();
968
969   if (nPix > npad) {
970     cout << nPix << endl;
971     // Too many pixels - sort and remove pixels with the lowest signal
972     fPixArray->Sort();
973     for (Int_t i=npad; i<nPix; i++) {
974       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
975       //pixPtr->Print();
976       fPixArray->RemoveAt(i);
977       delete pixPtr;
978     }
979     nPix = npad;
980   } // if (nPix > npad)
981
982   // Set pixel charges to the same value (for MLEM)
983   for (Int_t i=0; i<nPix; i++) {
984     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
985     //pixPtr->SetCharge(10);
986     cout << i+1 << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << " " << pixPtr->Size(0) << " " << pixPtr->Size(1) << endl;
987   }
988 }
989
990 //_____________________________________________________________________________
991 void AliMUONClusterFinderAZ::AjustPixel(Float_t width, Int_t ixy)
992 {
993   // Check if some pixels have small size (ajust if necessary)
994
995   AliMUONPixel *pixPtr, *pixPtr1 = 0;
996   Int_t ixy1 = TMath::Even(ixy);
997   Int_t nPix = fPixArray->GetEntriesFast();
998
999   for (Int_t i=0; i<nPix; i++) {
1000     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1001     if (pixPtr->Charge() < 1) continue; // discarded pixel
1002     if (pixPtr->Size(ixy)-width < -1.e-4) {
1003       // try to merge 
1004       cout << " Small X or Y: " << ixy << " " << pixPtr->Size(ixy) << " " << width << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << endl;
1005       for (Int_t j=i+1; j<nPix; j++) {
1006         pixPtr1 = (AliMUONPixel*) fPixArray->UncheckedAt(j);
1007         if (pixPtr1->Charge() < 1) continue; // discarded pixel
1008         if (TMath::Abs(pixPtr1->Size(ixy)-width) < 1.e-4) continue; // right size 
1009         if (TMath::Abs(pixPtr1->Coord(ixy1)-pixPtr->Coord(ixy1)) > 1.e-4) continue; // different rows/columns
1010         if (TMath::Abs(pixPtr1->Coord(ixy)-pixPtr->Coord(ixy)) < 2*width) {
1011           // merge
1012           pixPtr->SetSize(ixy, width);
1013           pixPtr->SetCoord(ixy, (pixPtr->Coord(ixy)+pixPtr1->Coord(ixy))/2);
1014           pixPtr->SetCharge(TMath::Min (pixPtr->Charge(),pixPtr1->Charge()));
1015           pixPtr1->SetCharge(0);
1016           pixPtr1 = 0;
1017           break;
1018         }
1019       } // for (Int_t j=i+1;
1020       //if (!pixPtr1) { cout << " I am here!" << endl; pixPtr->SetSize(ixy, width); } // ???
1021       //else if (pixPtr1->Charge() > 0.5 || i == nPix-1) {
1022       if (pixPtr1 || i == nPix-1) {
1023         // edge pixel - just increase its size
1024         cout << " Edge ..." << endl; 
1025         for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
1026           // ???if (fPadIJ[0][j] != i1) continue;
1027           if (TMath::Abs(pixPtr->Coord(ixy1)-fXyq[ixy1][j]) > 1.e-4) continue;
1028           if (pixPtr->Coord(ixy) < fXyq[ixy][j]) 
1029             pixPtr->Shift(ixy, -pixPtr->Size(ixy));
1030           else pixPtr->Shift(ixy, pixPtr->Size(ixy));
1031           pixPtr->SetSize(ixy, width);
1032           break;
1033         }
1034       }
1035     } // if (pixPtr->Size(ixy)-width < -1.e-4)
1036   } // for (Int_t i=0; i<nPix;
1037   return;
1038 }
1039   
1040 //_____________________________________________________________________________
1041 void AliMUONClusterFinderAZ::AjustPixel(Float_t wxmin, Float_t wymin)
1042 {
1043   // Check if some pixels have large size (ajust if necessary)
1044
1045   Int_t nx, ny;
1046   Int_t nPix = fPixArray->GetEntriesFast();
1047   AliMUONPixel *pixPtr, *pixPtr1, pix;
1048
1049   // Check if large pixel size
1050   for (Int_t i=0; i<nPix; i++) {
1051     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1052     if (pixPtr->Charge() < 1) continue; // discarded pixel
1053     if (pixPtr->Size(0)-wxmin > 1.e-4 || pixPtr->Size(1)-wymin > 1.e-4) {
1054       cout << " Different " << pixPtr->Size(0) << " " << wxmin << " " << pixPtr->Size(1) << " " << wymin << endl;
1055       pix = *pixPtr;
1056       nx = TMath::Nint (pix.Size(0)/wxmin);
1057       ny = TMath::Nint (pix.Size(1)/wymin);
1058       pix.Shift(0, -pix.Size(0)-wxmin);
1059       pix.Shift(1, -pix.Size(1)-wymin);
1060       pix.SetSize(0, wxmin);
1061       pix.SetSize(1, wymin);
1062       for (Int_t ii=0; ii<nx; ii++) {
1063         pix.Shift(0, wxmin*2);
1064         for (Int_t jj=0; jj<ny; jj++) {
1065           pix.Shift(1, wymin*2);
1066           pixPtr1 = new AliMUONPixel(pix);
1067           fPixArray->Add((TObject*)pixPtr1);
1068         }
1069       }
1070       pixPtr->SetCharge(0);
1071     }
1072   } // for (Int_t i=0; i<nPix;
1073   return;
1074 }
1075
1076 //_____________________________________________________________________________
1077 Bool_t AliMUONClusterFinderAZ::MainLoop()
1078 {
1079   // Repeat MLEM algorithm until pixel size becomes sufficiently small
1080   
1081   TH2D *mlem;
1082
1083   Int_t ix, iy;
1084   //Int_t nn, xList[10], yList[10];
1085   Int_t nPix = fPixArray->GetEntriesFast();
1086   Int_t npadTot = fnPads[0] + fnPads[1], npadOK = 0;
1087   AliMUONPixel *pixPtr = 0;
1088   Double_t *coef = 0, *probi = 0; 
1089   for (Int_t i=0; i<npadTot; i++) if (fPadIJ[1][i] == 0) npadOK++;
1090
1091   while (1) {
1092
1093     mlem = (TH2D*) gROOT->FindObject("mlem");
1094     if (mlem) mlem->Delete();
1095     // Calculate coefficients
1096     cout << " nPix, npadTot, npadOK " << nPix << " " << npadTot << " " << npadOK << endl;
1097
1098     // Calculate coefficients and pixel visibilities
1099     coef = new Double_t [npadTot*nPix];
1100     probi = new Double_t [nPix];
1101     Int_t indx = 0, cath;
1102     for (Int_t ipix=0; ipix<nPix; ipix++) {
1103       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1104       probi[ipix] = 0;
1105       for (Int_t j=0; j<npadTot; j++) {
1106         if (fPadIJ[1][j] < 0) { coef[j*nPix+ipix] = 0; continue; }
1107         cath = fPadIJ[0][j];
1108         fSegmentation[cath]->GetPadI(fXyq[0][j],fXyq[1][j],fZpad,ix,iy);
1109         fSegmentation[cath]->SetPad(ix,iy);
1110         /*
1111           fSegmentation[cath]->Neighbours(ix,iy,&nn,xList,yList); 
1112           if (nn != 4) {
1113           cout << nn << ": ";
1114           for (Int_t i=0; i<nn; i++) {cout << xList[i] << " " << yList[i] << ", ";}
1115           cout << endl;
1116           }
1117         */
1118         Double_t sum = 0;
1119         fSegmentation[cath]->SetHit(pixPtr->Coord(0),pixPtr->Coord(1),fZpad);
1120         sum += fResponse->IntXY(fSegmentation[cath]);
1121         indx = j*nPix + ipix;
1122         coef[indx] = sum; 
1123         probi[ipix] += coef[indx];
1124         //cout << j << " " << ipix << " " << coef[indx] << endl;
1125       } // for (Int_t j=0;
1126       //cout << " prob: " << probi[ipix] << endl;
1127       if (probi[ipix] < 0.01) pixPtr->SetCharge(0); // "invisible" pixel
1128     } // for (Int_t ipix=0;
1129
1130     // MLEM algorithm
1131     Mlem(coef, probi);
1132
1133     Double_t xylim[4] = {999, 999, 999, 999};
1134     for (Int_t ipix=0; ipix<nPix; ipix++) {
1135       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1136       for (Int_t i=0; i<4; i++) 
1137         xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1138       //cout << ipix+1; pixPtr->Print();
1139     }
1140     for (Int_t i=0; i<4; i++) {
1141       xylim[i] -= pixPtr->Size(i/2); cout << (i%2 ? -1 : 1)*xylim[i] << " "; }
1142     cout << endl;
1143
1144     // Ajust histogram to approximately the same limits as for the pads
1145     // (for good presentation)
1146     //*
1147     Float_t xypads[4];
1148     if (fHist[0]) {
1149       xypads[0] = fHist[0]->GetXaxis()->GetXmin();
1150       xypads[1] = -fHist[0]->GetXaxis()->GetXmax();
1151       xypads[2] = fHist[0]->GetYaxis()->GetXmin();
1152       xypads[3] = -fHist[0]->GetYaxis()->GetXmax();
1153       for (Int_t i=0; i<4; i++) {
1154         while(1) {
1155           if (xylim[i] < xypads[i]) break;
1156           xylim[i] -= 2*pixPtr->Size(i/2);
1157         }
1158       }
1159     } // if (fHist[0])
1160     //*/
1161
1162     Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1163     Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1164     mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1165     for (Int_t ipix=0; ipix<nPix; ipix++) {
1166       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1167       mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
1168     }
1169     //gPad->GetCanvas()->cd(3);
1170     if (fDraw) {
1171       ((TCanvas*)gROOT->FindObject("c2"))->cd();
1172       gPad->SetTheta(55);
1173       gPad->SetPhi(30);
1174       mlem->Draw("lego1Fb");
1175       gPad->Update();
1176       gets((char*)&ix);
1177     }
1178
1179     // Check if the total charge of pixels is too low
1180     Double_t qTot = 0;
1181     for (Int_t i=0; i<nPix; i++) {
1182       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1183       qTot += pixPtr->Charge();
1184     }
1185     if (qTot < 1.e-4 || npadOK < 3 && qTot < 50) {
1186       delete [] coef; delete [] probi; coef = 0; probi = 0;
1187       fPixArray->Delete(); 
1188       return kFALSE; 
1189     }
1190
1191     // Plot data - expectation
1192     /*
1193     Double_t x, y, cont;
1194     for (Int_t j=0; j<npadTot; j++) {
1195       Double_t sum1 = 0;
1196       for (Int_t i=0; i<nPix; i++) {
1197         // Caculate expectation
1198         pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1199         sum1 += pixPtr->Charge()*coef[j*nPix+i];
1200       }
1201       sum1 = TMath::Min (sum1,(Double_t)fResponse->MaxAdc());
1202       x = fXyq[0][j];
1203       y = fXyq[1][j];
1204       cath = fPadIJ[0][j];
1205       Int_t ihist = cath*2;
1206       ix = fHist[ihist]->GetXaxis()->FindBin(x);
1207       iy = fHist[ihist]->GetYaxis()->FindBin(y);
1208       cont = fHist[ihist]->GetCellContent(ix,iy);
1209       if (cont == 0 && fHist[ihist+1]) {
1210         ihist += 1;
1211         ix = fHist[ihist]->GetXaxis()->FindBin(x);
1212         iy = fHist[ihist]->GetYaxis()->FindBin(y);
1213       }
1214       fHist[ihist]->SetBinContent(ix,iy,fXyq[2][j]-sum1);
1215     }
1216     ((TCanvas*)gROOT->FindObject("c1"))->cd(1);
1217     //gPad->SetTheta(55);
1218     //gPad->SetPhi(30);
1219     //mlem->Draw("lego1");
1220     gPad->Modified();
1221     ((TCanvas*)gROOT->FindObject("c1"))->cd(2);
1222     gPad->Modified();
1223     */
1224
1225     // Calculate position of the center-of-gravity around the maximum pixel
1226     Double_t xyCOG[2];
1227     FindCOG(mlem, xyCOG);
1228
1229     if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 && pixPtr->Size(0) > pixPtr->Size(1)) break;
1230     //if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) >= 0.07 || pixPtr->Size(0) < pixPtr->Size(1)) {
1231     // Sort pixels according to the charge
1232     fPixArray->Sort();
1233     /*
1234     for (Int_t i=0; i<nPix; i++) {
1235       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1236       cout << i+1; pixPtr->Print();
1237     }
1238     */
1239     Double_t pixMin = 0.01*((AliMUONPixel*)fPixArray->UncheckedAt(0))->Charge();
1240     pixMin = TMath::Min (pixMin,50.);
1241
1242     // Decrease pixel size and shift pixels to make them centered at 
1243     // the maximum one
1244     indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
1245     Double_t width = 0, shift[2]={0};
1246     ix = 1;
1247     for (Int_t i=0; i<4; i++) xylim[i] = 999;
1248     Int_t nPix1 = nPix; nPix = 0;
1249     for (Int_t ipix=0; ipix<nPix1; ipix++) {
1250       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1251       if (nPix >= npadOK) { // too many pixels already
1252         fPixArray->RemoveAt(ipix); 
1253         delete pixPtr; 
1254         continue;
1255       }
1256       if (pixPtr->Charge() < pixMin) { // low charge
1257         fPixArray->RemoveAt(ipix); 
1258         delete pixPtr; 
1259         continue;
1260       }
1261       for (Int_t i=0; i<2; i++) {
1262         if (!i) {
1263           pixPtr->SetCharge(10);
1264           pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
1265           width = -pixPtr->Size(indx);
1266           pixPtr->Shift(indx, width);
1267           // Shift pixel position
1268           if (ix) {
1269             ix = 0;
1270             for (Int_t j=0; j<2; j++) {
1271               shift[j] = pixPtr->Coord(j) - xyCOG[j];
1272               shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
1273             }
1274             //cout << ipix << " " << i << " " << shift[0] << " " << shift[1] << endl;
1275           } // if (ix)
1276           pixPtr->Shift(0, -shift[0]);
1277           pixPtr->Shift(1, -shift[1]);
1278         } else {
1279           pixPtr = new AliMUONPixel(*pixPtr);
1280           pixPtr->Shift(indx, -2*width);
1281           fPixArray->Add((TObject*)pixPtr);
1282         } // else
1283         //pixPtr->Print();
1284         for (Int_t i=0; i<4; i++) 
1285           xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1286       } // for (Int_t i=0; i<2;
1287       nPix += 2;
1288     } // for (Int_t ipix=0;
1289
1290     fPixArray->Compress();
1291     nPix = fPixArray->GetEntriesFast();
1292
1293     // Remove excessive pixels
1294     if (nPix > npadOK) {
1295       for (Int_t ipix=npadOK; ipix<nPix; ipix++) { 
1296         pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1297         fPixArray->RemoveAt(ipix); 
1298         delete pixPtr;
1299       }
1300     } else {
1301       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(0);
1302       // add pixels if the maximum is at the limit of pixel area
1303       // start from Y-direction
1304       Int_t j = 0;
1305       for (Int_t i=3; i>-1; i--) {
1306         if (nPix < npadOK && 
1307             TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2)) {
1308           pixPtr = new AliMUONPixel(*pixPtr);
1309           pixPtr->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1310           j = TMath::Even (i/2);
1311           pixPtr->SetCoord(j, xyCOG[j]);
1312           fPixArray->Add((TObject*)pixPtr);
1313           nPix++;
1314         }
1315       }
1316     } // else    
1317
1318     fPixArray->Compress();
1319     nPix = fPixArray->GetEntriesFast();
1320     delete [] coef; delete [] probi; coef = 0; probi = 0;
1321   } // while (1)
1322
1323   // remove pixels with low signal or low visibility
1324   // Cuts are empirical !!!
1325   Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1326   thresh = TMath::Min (thresh,50.);
1327   Double_t cmax = -1, charge = 0;
1328   for (Int_t i=0; i<nPix; i++) cmax = TMath::Max (cmax,probi[i]); 
1329   // Mark pixels which should be removed
1330   for (Int_t i=0; i<nPix; i++) {
1331     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1332     charge = pixPtr->Charge();
1333     if (charge < thresh) pixPtr->SetCharge(-charge);
1334     else if (cmax > 1.91) {
1335       if (probi[i] < 1.9) pixPtr->SetCharge(-charge);
1336     }
1337     else if (probi[i] < cmax*0.9) pixPtr->SetCharge(-charge);
1338   }
1339   // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1340   Int_t near = 0;
1341   for (Int_t i=0; i<nPix; i++) {
1342     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1343     charge = pixPtr->Charge();
1344     if (charge > 0) continue;
1345     near = FindNearest(pixPtr);
1346     pixPtr->SetCharge(0);
1347     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(near);
1348     pixPtr->SetCharge(pixPtr->Charge() - charge);
1349   }
1350   // Update histogram
1351   for (Int_t i=0; i<nPix; i++) {
1352     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1353     ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1354     iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1355     mlem->SetBinContent(ix, iy, pixPtr->Charge());
1356   }
1357   if (fDraw) {
1358     ((TCanvas*)gROOT->FindObject("c2"))->cd();
1359     gPad->SetTheta(55);
1360     gPad->SetPhi(30);
1361     mlem->Draw("lego1Fb");
1362     gPad->Update();
1363   }
1364
1365   fxyMu[0][6] = fxyMu[1][6] = 9999;
1366   // Try to split into clusters
1367   Bool_t ok = kTRUE;
1368   if (mlem->GetSum() < 1) ok = kFALSE;
1369   else Split(mlem, coef);
1370   delete [] coef; delete [] probi; coef = 0; probi = 0;
1371   fPixArray->Delete(); 
1372   return ok;
1373 }
1374
1375 //_____________________________________________________________________________
1376 void AliMUONClusterFinderAZ::Mlem(Double_t *coef, Double_t *probi)
1377 {
1378   // Use MLEM to find pixel charges
1379   
1380   Int_t nPix = fPixArray->GetEntriesFast();
1381   Int_t npad = fnPads[0] + fnPads[1];
1382   Double_t *probi1 = new Double_t [nPix];
1383   Int_t indx, indx1;
1384   AliMUONPixel *pixPtr;
1385
1386   for (Int_t iter=0; iter<15; iter++) {
1387     // Do iterations
1388     for (Int_t ipix=0; ipix<nPix; ipix++) {
1389       // Correct each pixel
1390       if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1391       Double_t sum = 0;
1392       probi1[ipix] = probi[ipix];
1393       for (Int_t j=0; j<npad; j++) {
1394         if (fPadIJ[1][j] < 0) continue; 
1395         Double_t sum1 = 0;
1396         indx1 = j*nPix;
1397         indx = indx1 + ipix;
1398         for (Int_t i=0; i<nPix; i++) {
1399           // Caculate expectation
1400           pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1401           sum1 += pixPtr->Charge()*coef[indx1+i];
1402         } // for (Int_t i=0;
1403         if (fXyq[2][j] > fResponse->MaxAdc()-1 && sum1 > fResponse->MaxAdc()) { probi1[ipix] -= coef[indx]; continue; } // correct for pad charge overflows
1404         //cout << sum1 << " " << fXyq[2][j] << " " << coef[j*nPix+ipix] << endl;
1405         if (coef[indx] > 1.e-6) sum += fXyq[2][j]*coef[indx]/sum1;
1406       } // for (Int_t j=0;
1407       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1408       if (probi1[ipix] > 1.e-6) pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1409     } // for (Int_t ipix=0;
1410   } // for (Int_t iter=0;
1411   delete [] probi1;
1412   return;
1413 }
1414
1415 //_____________________________________________________________________________
1416 void AliMUONClusterFinderAZ::FindCOG(TH2D *mlem, Double_t *xyc)
1417 {
1418   // Calculate position of the center-of-gravity around the maximum pixel
1419
1420   Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1421   Int_t i1 = -9, j1 = -9;
1422   mlem->GetMaximumBin(ixmax,iymax,ix);
1423   Int_t nx = mlem->GetNbinsX();
1424   Int_t ny = mlem->GetNbinsY();
1425   Double_t thresh = mlem->GetMaximum()/10;
1426   Double_t x, y, cont, xq=0, yq=0, qq=0;
1427
1428   for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1429     y = mlem->GetYaxis()->GetBinCenter(i);
1430     for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1431       cont = mlem->GetCellContent(j,i);
1432       if (cont < thresh) continue;
1433       if (i != i1) {i1 = i; nsumy++;}
1434       if (j != j1) {j1 = j; nsumx++;}
1435       x = mlem->GetXaxis()->GetBinCenter(j);
1436       xq += x*cont;
1437       yq += y*cont;
1438       qq += cont;
1439       nsum++;
1440     }
1441   }
1442
1443   Double_t cmax = 0;
1444   Int_t i2 = 0, j2 = 0;
1445   x = y = 0;
1446   if (nsumy == 1) {
1447     // one bin in Y - add one more (with the largest signal)
1448     for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1449       if (i == iymax) continue;
1450       for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1451         cont = mlem->GetCellContent(j,i);
1452         if (cont > cmax) {
1453           cmax = cont;
1454           x = mlem->GetXaxis()->GetBinCenter(j);
1455           y = mlem->GetYaxis()->GetBinCenter(i);
1456           i2 = i;
1457           j2 = j;
1458         }
1459       }
1460     }
1461     xq += x*cmax;
1462     yq += y*cmax;
1463     qq += cmax;
1464     if (i2 != i1) nsumy++;
1465     if (j2 != j1) nsumx++;
1466     nsum++;
1467   } // if (nsumy == 1)
1468
1469   if (nsumx == 1) {
1470     // one bin in X - add one more (with the largest signal)
1471     cmax = x = y = 0;
1472     for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1473       if (j == ixmax) continue;
1474       for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1475         cont = mlem->GetCellContent(j,i);
1476         if (cont > cmax) {
1477           cmax = cont;
1478           x = mlem->GetXaxis()->GetBinCenter(j);
1479           y = mlem->GetYaxis()->GetBinCenter(i);
1480           i2 = i;
1481           j2 = j;
1482         }
1483       }
1484     }
1485     xq += x*cmax;
1486     yq += y*cmax;
1487     qq += cmax;
1488     if (i2 != i1) nsumy++;
1489     if (j2 != j1) nsumx++;
1490     nsum++;
1491   } // if (nsumx == 1)
1492
1493   xyc[0] = xq/qq; xyc[1] = yq/qq;
1494   cout << xyc[0] << " " << xyc[1] << " " << qq << " " << nsum << " " << nsumx << " " << nsumy << endl;
1495   return;
1496 }
1497
1498 //_____________________________________________________________________________
1499 Int_t AliMUONClusterFinderAZ::FindNearest(AliMUONPixel *pixPtr0)
1500 {
1501   // Find the pixel nearest to the given one
1502   // (algorithm may be not very efficient)
1503
1504   Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1505   Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1506   Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1507   AliMUONPixel *pixPtr;
1508
1509   for (Int_t i=0; i<nPix; i++) {
1510     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1511     if (pixPtr->Charge() < 0.5) continue;
1512     dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1513     dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1514     r = dx *dx + dy * dy;
1515     if (r < rmin) { rmin = r; imin = i; }
1516   }
1517   return imin;
1518 }
1519
1520 //_____________________________________________________________________________
1521 void AliMUONClusterFinderAZ::Split(TH2D *mlem, Double_t *coef)
1522 {
1523   // The main steering function to work with clusters of pixels in anode
1524   // plane (find clusters, decouple them from each other, merge them (if
1525   // necessary), pick up coupled pads, call the fitting function)
1526   
1527   Int_t nx = mlem->GetNbinsX();
1528   Int_t ny = mlem->GetNbinsY();
1529   Int_t nPix = fPixArray->GetEntriesFast();
1530
1531   Bool_t *used = new Bool_t[ny*nx];
1532   Double_t cont;
1533   Int_t nclust = 0, indx, indx1;
1534
1535   for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE; 
1536
1537   TObjArray *clusters[200]={0};
1538   TObjArray *pix;
1539
1540   // Find clusters of histogram bins (easier to work in 2-D space)
1541   for (Int_t i=1; i<=ny; i++) {
1542     for (Int_t j=1; j<=nx; j++) {
1543       indx = (i-1)*nx + j - 1;
1544       if (used[indx]) continue;
1545       cont = mlem->GetCellContent(j,i);
1546       if (cont < 0.5) continue;
1547       pix = new TObjArray(20);
1548       used[indx] = 1;
1549       pix->Add(BinToPix(mlem,j,i));
1550       AddBin(mlem, i, j, 0, used, pix); // recursive call
1551       clusters[nclust++] = pix;
1552       if (nclust > 200) { cout << " Too many clusters " << endl; ::exit(0); }
1553     } // for (Int_t j=1; j<=nx; j++) {
1554   } // for (Int_t i=1; i<=ny;
1555   cout << nclust << endl;
1556   delete [] used; used = 0;
1557   
1558   // Compute couplings between clusters and clusters to pads
1559   Int_t npad = fnPads[0] + fnPads[1];
1560
1561   // Exclude pads with overflows
1562   for (Int_t j=0; j<npad; j++) {
1563     if (fXyq[2][j] > fResponse->MaxAdc()-1) fPadIJ[1][j] = -9;
1564     else fPadIJ[1][j] = 0;
1565   }
1566
1567   // Compute couplings of clusters to pads
1568   TMatrixD *aij_clu_pad = new TMatrixD(nclust,npad);
1569   *aij_clu_pad = 0;
1570   Int_t npxclu;
1571   for (Int_t iclust=0; iclust<nclust; iclust++) {
1572     pix = clusters[iclust];
1573     npxclu = pix->GetEntriesFast();
1574     for (Int_t i=0; i<npxclu; i++) {
1575       indx = fPixArray->IndexOf(pix->UncheckedAt(i));
1576       for (Int_t j=0; j<npad; j++) {
1577         // Exclude overflows
1578         if (fPadIJ[1][j] < 0) continue;
1579         if (coef[j*nPix+indx] < kCouplMin) continue;
1580         (*aij_clu_pad)(iclust,j) += coef[j*nPix+indx];
1581       }
1582     }
1583   }
1584   // Compute couplings between clusters
1585   TMatrixD *aij_clu_clu = new TMatrixD(nclust,nclust);
1586   *aij_clu_clu = 0;
1587   for (Int_t iclust=0; iclust<nclust; iclust++) {
1588     for (Int_t j=0; j<npad; j++) {
1589       // Exclude overflows
1590       if (fPadIJ[1][j] < 0) continue;
1591       if ((*aij_clu_pad)(iclust,j) < kCouplMin) continue;
1592       for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) {
1593         if ((*aij_clu_pad)(iclust1,j) < kCouplMin) continue;
1594         (*aij_clu_clu)(iclust,iclust1) += 
1595           TMath::Sqrt ((*aij_clu_pad)(iclust,j)*(*aij_clu_pad)(iclust1,j));
1596       }
1597     }
1598   }
1599   for (Int_t iclust=0; iclust<nclust; iclust++) {
1600     for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) {
1601       (*aij_clu_clu)(iclust1,iclust) = (*aij_clu_clu)(iclust,iclust1);
1602     }
1603   }
1604
1605   if (nclust > 1) aij_clu_clu->Print();
1606
1607   // Find groups of coupled clusters
1608   used = new Bool_t[nclust];
1609   for (Int_t i=0; i<nclust; i++) used[i] = kFALSE;
1610   Int_t *clustNumb = new Int_t[nclust];
1611   Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
1612   Double_t parOk[8];
1613
1614   for (Int_t igroup=0; igroup<nclust; igroup++) {
1615     if (used[igroup]) continue;
1616     used[igroup] = kTRUE;
1617     clustNumb[0] = igroup;
1618     nCoupled = 1;
1619     // Find group of coupled clusters
1620     AddCluster(igroup, nclust, aij_clu_clu, used, clustNumb, nCoupled); // recursive
1621     cout << " nCoupled: " << nCoupled << endl;
1622     for (Int_t i=0; i<nCoupled; i++) cout << clustNumb[i] << " "; cout << endl; 
1623
1624     while (nCoupled > 0) {
1625
1626       if (nCoupled < 4) {
1627         nForFit = nCoupled;
1628         for (Int_t i=0; i<nCoupled; i++) clustFit[i] = clustNumb[i];
1629       } else {
1630         // Too many coupled clusters to fit - try to decouple them
1631         // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with 
1632         // all the others in the group 
1633         for (Int_t j=0; j<3; j++) minGroup[j] = -1;
1634         Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aij_clu_clu, minGroup);
1635
1636         // Flag clusters for fit
1637         nForFit = 0;
1638         while (minGroup[nForFit] >= 0 && nForFit < 3) {
1639           cout << clustNumb[minGroup[nForFit]] << " ";
1640           clustFit[nForFit] = clustNumb[minGroup[nForFit]];
1641           clustNumb[minGroup[nForFit]] -= 999;
1642           nForFit++;
1643         }
1644         cout << nForFit << " " << coupl << endl;
1645       } // else
1646
1647       // Select pads for fit. 
1648       if (SelectPad(nCoupled, nForFit, clustNumb, clustFit, aij_clu_pad) < 3 && nCoupled > 1) {
1649         // Deselect pads
1650         for (Int_t j=0; j<npad; j++) if (TMath::Abs(fPadIJ[1][j]) == 1) fPadIJ[1][j] = 0;
1651         // Merge the failed cluster candidates (with too few pads to fit) with 
1652         // the one with the strongest coupling
1653         Merge(nForFit, nCoupled, clustNumb, clustFit, clusters, aij_clu_clu, aij_clu_pad);
1654       } else {
1655         // Do the fit
1656         nfit = Fit(nForFit, clustFit, clusters, parOk);
1657       }
1658
1659       // Subtract the fitted charges from pads with strong coupling and/or
1660       // return pads for further use
1661       UpdatePads(nfit, parOk);
1662
1663       // Mark used pads
1664       for (Int_t j=0; j<npad; j++) {if (fPadIJ[1][j] == 1) fPadIJ[1][j] = -1;}
1665
1666       // Sort the clusters (move to the right the used ones)
1667       Int_t beg = 0, end = nCoupled - 1;
1668       while (beg < end) {
1669         if (clustNumb[beg] >= 0) { beg++; continue; }
1670         for (Int_t j=end; j>beg; j--) {
1671           if (clustNumb[j] < 0) continue;
1672           end = j - 1;
1673           indx = clustNumb[beg];
1674           clustNumb[beg] = clustNumb[j];
1675           clustNumb[j] = indx;
1676           break;
1677         }
1678         beg++;
1679       }
1680
1681       nCoupled -= nForFit;
1682       if (nCoupled > 3) {
1683         // Remove couplings of used clusters
1684         for (Int_t iclust=nCoupled; iclust<nCoupled+nForFit; iclust++) {
1685           indx = clustNumb[iclust] + 999;
1686           for (Int_t iclust1=0; iclust1<nCoupled; iclust1++) {
1687             indx1 = clustNumb[iclust1];
1688             (*aij_clu_clu)(indx,indx1) = (*aij_clu_clu)(indx1,indx) = 0;
1689           }
1690         }
1691
1692         // Update the remaining clusters couplings (exclude couplings from 
1693         // the used pads)
1694         for (Int_t j=0; j<npad; j++) {
1695           if (fPadIJ[1][j] != -1) continue;
1696           for (Int_t iclust=0; iclust<nCoupled; iclust++) {
1697             indx = clustNumb[iclust];
1698             if ((*aij_clu_pad)(indx,j) < kCouplMin) continue;
1699             for (Int_t iclust1=iclust+1; iclust1<nCoupled; iclust1++) {
1700               indx1 = clustNumb[iclust1];
1701               if ((*aij_clu_pad)(indx1,j) < kCouplMin) continue;
1702               // Check this
1703               (*aij_clu_clu)(indx,indx1) -= 
1704                 TMath::Sqrt ((*aij_clu_pad)(indx,j)*(*aij_clu_pad)(indx1,j));
1705               (*aij_clu_clu)(indx1,indx) = (*aij_clu_clu)(indx,indx1);
1706             }
1707           }
1708           fPadIJ[1][j] = -9;
1709         } // for (Int_t j=0; j<npad;
1710       } // if (nCoupled > 3)
1711     } // while (nCoupled > 0)
1712   } // for (Int_t igroup=0; igroup<nclust;
1713
1714   //delete aij_clu; aij_clu = 0; delete aij_clu_pad; aij_clu_pad = 0;
1715   aij_clu_clu->Delete(); aij_clu_pad->Delete();
1716   for (Int_t iclust=0; iclust<nclust; iclust++) {
1717     pix = clusters[iclust]; 
1718     pix->Clear();
1719     delete pix; pix = 0;
1720   }
1721   delete [] clustNumb; clustNumb = 0; delete [] used; used = 0;
1722 }
1723
1724 //_____________________________________________________________________________
1725 void AliMUONClusterFinderAZ::AddBin(TH2D *mlem, Int_t ic, Int_t jc, Int_t mode, Bool_t *used, TObjArray *pix)
1726 {
1727   // Add a bin to the cluster
1728
1729   Int_t nx = mlem->GetNbinsX();
1730   Int_t ny = mlem->GetNbinsY();
1731   Double_t cont1, cont = mlem->GetCellContent(jc,ic);
1732   AliMUONPixel *pixPtr = 0;
1733
1734   for (Int_t i=TMath::Max(ic-1,1); i<=TMath::Min(ic+1,ny); i++) {
1735     for (Int_t j=TMath::Max(jc-1,1); j<=TMath::Min(jc+1,nx); j++) {
1736       if (i != ic && j != jc) continue;
1737       if (used[(i-1)*nx+j-1]) continue;
1738       cont1 = mlem->GetCellContent(j,i);
1739       if (mode && cont1 > cont) continue;
1740       used[(i-1)*nx+j-1] = kTRUE;
1741       if (cont1 < 0.5) continue;
1742       if (pix) pix->Add(BinToPix(mlem,j,i)); 
1743       else {
1744         pixPtr = new AliMUONPixel (mlem->GetXaxis()->GetBinCenter(j), 
1745                                    mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1746         fPixArray->Add((TObject*)pixPtr);
1747       }
1748       AddBin(mlem, i, j, mode, used, pix); // recursive call
1749     }
1750   }
1751 }
1752
1753 //_____________________________________________________________________________
1754 TObject* AliMUONClusterFinderAZ::BinToPix(TH2D *mlem, Int_t jc, Int_t ic)
1755 {
1756   // Translate histogram bin to pixel 
1757   
1758   Double_t yc = mlem->GetYaxis()->GetBinCenter(ic);
1759   Double_t xc = mlem->GetXaxis()->GetBinCenter(jc);
1760   
1761   Int_t nPix = fPixArray->GetEntriesFast();
1762   AliMUONPixel *pixPtr;
1763
1764   // Compare pixel and bin positions
1765   for (Int_t i=0; i<nPix; i++) {
1766     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1767     if (pixPtr->Charge() < 0.5) continue;
1768     if (TMath::Abs(pixPtr->Coord(0)-xc)<1.e-4 && TMath::Abs(pixPtr->Coord(1)-yc)<1.e-4) return (TObject*) pixPtr;
1769   }
1770   cout << " Something wrong ??? " << endl;
1771   return NULL;
1772 }
1773
1774 //_____________________________________________________________________________
1775 void AliMUONClusterFinderAZ::AddCluster(Int_t ic, Int_t nclust, TMatrixD *aij_clu_clu, Bool_t *used, Int_t *clustNumb, Int_t &nCoupled)
1776 {
1777   // Add a cluster to the group of coupled clusters
1778
1779   for (Int_t i=0; i<nclust; i++) {
1780     if (used[i]) continue;
1781     if ((*aij_clu_clu)(i,ic) < kCouplMin) continue;
1782     used[i] = kTRUE;
1783     clustNumb[nCoupled++] = i;
1784     AddCluster(i, nclust, aij_clu_clu, used, clustNumb, nCoupled);
1785   }
1786 }
1787
1788 //_____________________________________________________________________________
1789 Double_t AliMUONClusterFinderAZ::MinGroupCoupl(Int_t nCoupled, Int_t *clustNumb, TMatrixD *aij_clu_clu, Int_t *minGroup)
1790 {
1791   // Find group of clusters with minimum coupling to all the others
1792
1793   Int_t i123max = TMath::Min(3,nCoupled/2); 
1794   Int_t indx, indx1, indx2, indx3, nTot = 0;
1795   Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1796
1797   for (Int_t i123=1; i123<=i123max; i123++) {
1798
1799     if (i123 == 1) {
1800       coupl1 = new Double_t [nCoupled];
1801       for (Int_t i=0; i<nCoupled; i++) coupl1[i] = 0;
1802     }
1803     else if (i123 == 2) {
1804       nTot = nCoupled*nCoupled;
1805       coupl2 = new Double_t [nTot];
1806       for (Int_t i=0; i<nTot; i++) coupl2[i] = 9999;
1807     } else {
1808       nTot = nTot*nCoupled;
1809       coupl3 = new Double_t [nTot];
1810       for (Int_t i=0; i<nTot; i++) coupl3[i] = 9999;
1811     } // else
1812
1813     for (Int_t i=0; i<nCoupled; i++) {
1814       indx1 = clustNumb[i];
1815       for (Int_t j=i+1; j<nCoupled; j++) {
1816         indx2 = clustNumb[j];
1817         if (i123 == 1) {
1818           coupl1[i] += (*aij_clu_clu)(indx1,indx2);
1819           coupl1[j] += (*aij_clu_clu)(indx1,indx2);
1820         } 
1821         else if (i123 == 2) {
1822           indx = i*nCoupled + j;
1823           coupl2[indx] = coupl1[i] + coupl1[j];
1824           coupl2[indx] -= 2 * ((*aij_clu_clu)(indx1,indx2));
1825         } else {
1826           for (Int_t k=j+1; k<nCoupled; k++) {
1827             indx3 = clustNumb[k];
1828             indx = i*nCoupled*nCoupled + j*nCoupled + k;
1829             coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1830             coupl3[indx] -= 2 * ((*aij_clu_clu)(indx1,indx3)+(*aij_clu_clu)(indx2,indx3));
1831           }
1832         } // else
1833       } // for (Int_t j=i+1;
1834     } // for (Int_t i=0;
1835   } // for (Int_t i123=1;
1836
1837   // Find minimum coupling
1838   Double_t couplMin = 9999;
1839   Int_t locMin = 0;
1840
1841   for (Int_t i123=1; i123<=i123max; i123++) {
1842     if (i123 == 1) {
1843       locMin = TMath::LocMin(nCoupled, coupl1);
1844       couplMin = coupl1[locMin];
1845       minGroup[0] = locMin;
1846       delete [] coupl1; coupl1 = 0;
1847     } 
1848     else if (i123 == 2) {
1849       locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
1850       if (coupl2[locMin] < couplMin) {
1851         couplMin = coupl2[locMin];
1852         minGroup[0] = locMin/nCoupled;
1853         minGroup[1] = locMin%nCoupled;
1854       }
1855       delete [] coupl2; coupl2 = 0;
1856     } else {
1857       locMin = TMath::LocMin(nTot, coupl3);
1858       if (coupl3[locMin] < couplMin) {
1859         couplMin = coupl3[locMin];
1860         minGroup[0] = locMin/nCoupled/nCoupled;
1861         minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
1862         minGroup[2] = locMin%nCoupled;
1863       }
1864       delete [] coupl3; coupl3 = 0;
1865     } // else
1866   } // for (Int_t i123=1;
1867   return couplMin;
1868 }
1869
1870 //_____________________________________________________________________________
1871 Int_t AliMUONClusterFinderAZ::SelectPad(Int_t nCoupled, Int_t nForFit, Int_t *clustNumb, Int_t *clustFit, TMatrixD *aij_clu_pad)
1872 {
1873   // Select pads for fit. If too many coupled clusters, find pads giving 
1874   // the strongest coupling with the rest of clusters and exclude them from the fit.
1875
1876   Int_t npad = fnPads[0] + fnPads[1];
1877   Double_t *pad_pix = 0;
1878
1879   if (nCoupled > 3) {
1880     pad_pix = new Double_t[npad];
1881     for (Int_t i=0; i<npad; i++) pad_pix[i] = 0; 
1882   }
1883
1884   Int_t nOK = 0, indx, indx1;
1885   for (Int_t iclust=0; iclust<nForFit; iclust++) {
1886     indx = clustFit[iclust];
1887     for (Int_t j=0; j<npad; j++) {
1888       if (fPadIJ[1][j] < 0) continue; // exclude overflows and used pads
1889       if ((*aij_clu_pad)(indx,j) < kCouplMin) continue;
1890       fPadIJ[1][j] = 1; // pad to be used in fit
1891       nOK++;
1892       if (nCoupled > 3) {
1893         // Check other clusters
1894         for (Int_t iclust1=0; iclust1<nCoupled; iclust1++) {
1895           indx1 = clustNumb[iclust1];
1896           if (indx1 < 0) continue;
1897           if ((*aij_clu_pad)(indx1,j) < kCouplMin) continue;
1898           pad_pix[j] += (*aij_clu_pad)(indx1,j);
1899         }
1900       } // if (nCoupled > 3)
1901     } // for (Int_t j=0; j<npad;
1902   } // for (Int_t iclust=0; iclust<nForFit
1903   if (nCoupled < 4) return nOK;
1904
1905   Double_t aaa = 0;
1906   for (Int_t j=0; j<npad; j++) {
1907     if (pad_pix[j] < kCouplMin) continue;
1908     cout << j << " " << pad_pix[j] << " "; 
1909     cout << fXyq[0][j] << " " << fXyq[1][j] << endl;
1910     aaa += pad_pix[j];
1911     fPadIJ[1][j] = -1; // exclude pads with strong coupling to the other clusters
1912     nOK--;
1913   }
1914   delete [] pad_pix; pad_pix = 0;
1915   return nOK;
1916 }
1917   
1918 //_____________________________________________________________________________
1919 void AliMUONClusterFinderAZ::Merge(Int_t nForFit, Int_t nCoupled, Int_t *clustNumb, Int_t *clustFit, TObjArray **clusters, TMatrixD *aij_clu_clu, TMatrixD *aij_clu_pad)
1920 {
1921   // Merge the group of clusters with the one having the strongest coupling with them
1922
1923   Int_t indx, indx1, npxclu, npxclu1, imax=0;
1924   TObjArray *pix, *pix1;
1925   Double_t couplMax;
1926
1927   for (Int_t icl=0; icl<nForFit; icl++) {
1928     indx = clustFit[icl];
1929     pix = clusters[indx];
1930     npxclu = pix->GetEntriesFast();
1931     couplMax = -1;
1932     for (Int_t icl1=0; icl1<nCoupled; icl1++) {
1933       indx1 = clustNumb[icl1];
1934       if (indx1 < 0) continue;
1935       if ((*aij_clu_clu)(indx,indx1) > couplMax) {
1936         couplMax = (*aij_clu_clu)(indx,indx1);
1937         imax = indx1;
1938       }
1939     } // for (Int_t icl1=0;
1940     /*if (couplMax < kCouplMin) {
1941       cout << " Oops " << couplMax << endl;
1942       aij_clu_clu->Print();
1943       cout << icl << " " << indx << " " << npxclu << " " << nLinks << endl;
1944       ::exit(0);
1945       }*/
1946     // Add to it
1947     pix1 = clusters[imax];
1948     npxclu1 = pix1->GetEntriesFast();
1949     // Add pixels 
1950     for (Int_t i=0; i<npxclu; i++) { pix1->Add(pix->UncheckedAt(i)); pix->RemoveAt(i); }
1951     cout << " New number of pixels: " << npxclu1 << " " << pix1->GetEntriesFast() << endl;
1952     //Add cluster-to-cluster couplings
1953     //aij_clu_clu->Print();
1954     for (Int_t icl1=0; icl1<nCoupled; icl1++) {
1955       indx1 = clustNumb[icl1];
1956       if (indx1 < 0 || indx1 == imax) continue;
1957       (*aij_clu_clu)(indx1,imax) += (*aij_clu_clu)(indx,indx1);
1958       (*aij_clu_clu)(imax,indx1) = (*aij_clu_clu)(indx1,imax);
1959     }
1960     (*aij_clu_clu)(indx,imax) = (*aij_clu_clu)(imax,indx) = 0;
1961     //aij_clu_clu->Print();
1962     //Add cluster-to-pad couplings
1963     for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
1964       if (fPadIJ[1][j] < 0) continue; // exclude overflows and used pads
1965       (*aij_clu_pad)(imax,j) += (*aij_clu_pad)(indx,j);
1966       (*aij_clu_pad)(indx,j) = 0;
1967     }
1968   } // for (Int_t icl=0; icl<nForFit;
1969 }
1970
1971 //_____________________________________________________________________________
1972 Int_t AliMUONClusterFinderAZ::Fit(Int_t nfit, Int_t *clustFit, TObjArray **clusters, Double_t *parOk)
1973 {
1974   // Find selected clusters to selected pad charges
1975   
1976   TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
1977   //Int_t nx = mlem->GetNbinsX();
1978   //Int_t ny = mlem->GetNbinsY();
1979   Double_t xmin = mlem->GetXaxis()->GetXmin() - mlem->GetXaxis()->GetBinWidth(1);
1980   Double_t xmax = mlem->GetXaxis()->GetXmax() + mlem->GetXaxis()->GetBinWidth(1);
1981   Double_t ymin = mlem->GetYaxis()->GetXmin() - mlem->GetYaxis()->GetBinWidth(1);
1982   Double_t ymax = mlem->GetYaxis()->GetXmax() + mlem->GetYaxis()->GetBinWidth(1);
1983   //Double_t qmin = 0, qmax = 1;
1984   Double_t step[3]={0.01,0.002,0.02};
1985
1986   Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8];
1987   TObjArray *pix;
1988   Int_t npxclu;
1989
1990   // Number of pads to use
1991   Int_t npads = 0;
1992   for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {if (fPadIJ[1][i] == 1) npads++;}
1993   for (Int_t i=0; i<nfit; i++) {cout << i+1 << " " << clustFit[i] << " ";}
1994   cout << nfit << endl;
1995   cout << " Number of pads to fit: " << npads << endl;
1996   fNpar = 0;
1997   fQtot = 0;
1998   if (npads < 2) return 0; 
1999
2000   // Take cluster maxima as fitting seeds
2001   AliMUONPixel *pixPtr;
2002   Double_t xyseed[3][2], qseed[3];
2003   for (Int_t ifit=1; ifit<=nfit; ifit++) {
2004     cmax = 0;
2005     pix = clusters[clustFit[ifit-1]];
2006     npxclu = pix->GetEntriesFast();
2007     for (Int_t clu=0; clu<npxclu; clu++) {
2008       pixPtr = (AliMUONPixel*) pix->UncheckedAt(clu);
2009       cont = pixPtr->Charge();
2010       fQtot += cont;
2011       if (cont > cmax) { 
2012         cmax = cont; 
2013         xseed = pixPtr->Coord(0);
2014         yseed = pixPtr->Coord(1);
2015       }
2016     }
2017     xyseed[ifit-1][0] = xseed;
2018     xyseed[ifit-1][1] = yseed;
2019     qseed[ifit-1] = cmax;
2020   } // for (Int_t ifit=1;
2021
2022   Int_t nDof, maxSeed[3];
2023   Double_t fmin, chi2o = 9999, chi2n;
2024
2025   // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is 
2026   // lower, try 3-track (if number of pads is sufficient).
2027   
2028   TMath::Sort(nfit, qseed, maxSeed, kTRUE); // in decreasing order
2029   nfit = TMath::Min (nfit, (npads + 1) / 3);
2030
2031   Double_t *gin = 0, func0, func1, param[8], param0[2][8], deriv[2][8], step0[8];
2032   Double_t shift[8], stepMax, derMax, parmin[8], parmax[8], func2[2], shift0;
2033   Double_t delta[8], scMax, dder[8], estim, shiftSave = 0;
2034   Int_t min, max, nCall = 0, memory[8] = {0}, nLoop, idMax = 0, iestMax = 0, nFail;
2035
2036   for (Int_t iseed=0; iseed<nfit; iseed++) {
2037
2038     for (Int_t j=0; j<3; j++) step0[fNpar+j] = shift[fNpar+j] = step[j];
2039     param[fNpar] = xyseed[maxSeed[iseed]][0];
2040     parmin[fNpar] = xmin; 
2041     parmax[fNpar++] = xmax; 
2042     param[fNpar] = xyseed[maxSeed[iseed]][1];
2043     parmin[fNpar] = ymin; 
2044     parmax[fNpar++] = ymax; 
2045     if (fNpar > 2) {
2046       param[fNpar] = fNpar == 4 ? 0.5 : 0.3;
2047       parmin[fNpar] = 0; 
2048       parmax[fNpar++] = 1; 
2049     }
2050
2051     // Try new algorithm
2052     min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
2053
2054     while (1) {
2055       max = !min;
2056       fcn1(fNpar, gin, func0, param, 1); nCall++;
2057       //cout << " Func: " << func0 << endl;
2058
2059       func2[max] = func0;
2060       for (Int_t j=0; j<fNpar; j++) {
2061         param0[max][j] = param[j];
2062         delta[j] = step0[j];
2063         param[j] += delta[j] / 10;
2064         if (j > 0) param[j-1] -= delta[j-1] / 10;
2065         fcn1(fNpar, gin, func1, param, 1); nCall++;
2066         deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
2067         //cout << j << " " << deriv[max][j] << endl;
2068         dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) / 
2069                                                  (param0[0][j] - param0[1][j]) : 0; // second derivative
2070       }
2071       param[fNpar-1] -= delta[fNpar-1] / 10;
2072       if (nCall > 2000) ::exit(0);
2073
2074       min = func2[0] < func2[1] ? 0 : 1;
2075       nFail = min == max ? 0 : nFail + 1;
2076
2077       stepMax = derMax = estim = 0;
2078       for (Int_t j=0; j<fNpar; j++) { 
2079         // Estimated distance to minimum
2080         shift0 = shift[j];
2081         if (nLoop == 1) shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
2082         else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3) shift[j] = 0;
2083         else if (deriv[min][j]*deriv[!min][j] > 0 && TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])
2084               || TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3) {
2085           shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
2086           if (min == max) { 
2087             if (memory[j] > 1) { shift[j] *= 2; } //cout << " Memory " << memory[j] << " " << shift[j] << endl; }
2088             memory[j]++;
2089           }
2090         } else {
2091           shift[j] = -deriv[min][j] / dder[j];
2092           memory[j] = 0;
2093         }
2094         if (TMath::Abs(shift[j])/step0[j] > estim) { 
2095           estim = TMath::Abs(shift[j])/step0[j];
2096           iestMax = j;
2097         }
2098
2099         // Too big step
2100         if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; // 
2101
2102         // Failed to improve minimum
2103         if (min != max) {
2104           memory[j] = 0;
2105           param[j] = param0[min][j];
2106           if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j]) shift[j] = (shift[j] + shift0) / 2;
2107           else shift[j] /= -2;
2108         } 
2109
2110         // Too big step
2111         if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min]) 
2112           shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
2113
2114         // Introduce step relaxation factor
2115         if (memory[j] < 3) {
2116           scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
2117           if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax) 
2118             shift[j] = TMath::Sign (shift0*scMax, shift[j]);
2119         }
2120         param[j] += shift[j]; 
2121           
2122         //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
2123         stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
2124         if (TMath::Abs(deriv[min][j]) > derMax) {
2125           idMax = j;
2126           derMax = TMath::Abs (deriv[min][j]);
2127         }
2128       } // for (Int_t j=0; j<fNpar;
2129       //cout << max << " " << func2[min] << " " << derMax << " " << stepMax << " " << estim << " " << iestMax << " " << nCall << endl;
2130       if (estim < 1 && derMax < 2 || nLoop > 100) break; // minimum was found
2131
2132       nLoop++;
2133       // Check for small step
2134       if (shift[idMax] == 0) { shift[idMax] = step0[idMax]/10; param[idMax] += shift[idMax]; continue; }
2135       if (!memory[idMax] && derMax > 0.5 && nLoop > 10) {
2136         //cout << " ok " << deriv[min][idMax] << " " << deriv[!min][idMax] << " " << dder[idMax]*shift[idMax] << " " << shift[idMax] << endl;
2137         if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10) {
2138           if (min == max) dder[idMax] = -dder[idMax];
2139           shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10; 
2140           param[idMax] += shift[idMax];
2141           stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
2142           //cout << shift[idMax] << " " << param[idMax] << endl;
2143           if (min == max) shiftSave = shift[idMax];
2144         }
2145         if (nFail > 10) {
2146           param[idMax] -= shift[idMax];
2147           shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
2148           param[idMax] += shift[idMax];
2149           //cout << shift[idMax] << endl;
2150         }
2151       }      
2152     } // while (1)
2153     fmin = func2[min];
2154
2155     nDof = npads - fNpar;
2156     chi2n = nDof ? fmin/nDof : 0;
2157
2158     if (chi2n*1.2+1.e-6 > chi2o ) { fNpar -= 3; break; }
2159     // Save parameters and errors
2160     for (Int_t i=0; i<fNpar; i++) {
2161       parOk[i] = param0[min][i];
2162       errOk[i] = fmin;
2163     }
2164
2165     cout << chi2o << " " << chi2n << endl;
2166     chi2o = chi2n;
2167     if (fmin < 0.1) break; // !!!???
2168   } // for (Int_t iseed=0; 
2169
2170   for (Int_t i=0; i<fNpar; i++) {
2171     if (i == 4 || i == 7) continue;
2172     cout << parOk[i] << " " << errOk[i] << endl;
2173   }
2174   nfit = (fNpar + 1) / 3;
2175   Double_t rad;
2176   Int_t indx, imax;
2177   if (fReco) {
2178     for (Int_t j=0; j<nfit; j++) {
2179       indx = j<2 ? j*2 : j*2+1;  
2180       AddRawCluster (parOk[indx], parOk[indx+1], errOk[indx]);
2181     }
2182     return nfit;
2183   } 
2184   for (Int_t i=0; i<fnMu; i++) {
2185     cmax = fxyMu[i][6];
2186     for (Int_t j=0; j<nfit; j++) {
2187       indx = j<2 ? j*2 : j*2+1;  
2188       rad = (fxyMu[i][0]-parOk[indx])*(fxyMu[i][0]-parOk[indx]) +
2189             (fxyMu[i][1]-parOk[indx+1])*(fxyMu[i][1]-parOk[indx+1]);
2190       if (rad < cmax) {
2191         cmax = rad; 
2192         imax = indx;
2193         fxyMu[i][6] = cmax;
2194         fxyMu[i][2] = parOk[imax] - fxyMu[i][0];
2195         fxyMu[i][4] = parOk[imax+1] - fxyMu[i][1];
2196         fxyMu[i][3] = errOk[imax];
2197         fxyMu[i][5] = errOk[imax+1];
2198       }
2199     }      
2200   }
2201   return nfit;
2202 }  
2203
2204 //_____________________________________________________________________________
2205 void fcn1(Int_t & /*npar*/, Double_t * /*gin*/, Double_t &f, Double_t *par, Int_t /*iflag*/)
2206 {
2207   // Fit for one track
2208   AliMUONClusterFinderAZ& c = *(AliMUONClusterFinderAZ::fgClusterFinder);    
2209   
2210   Int_t cath, ix, iy, indx, npads=0;
2211   Double_t charge, delta, coef=0, chi2=0;
2212   for (Int_t j=0; j<c.fnPads[0]+c.fnPads[1]; j++) {
2213     if (c.fPadIJ[1][j] != 1) continue;
2214     cath = c.fPadIJ[0][j];
2215     npads++;
2216     c.fSegmentation[cath]->GetPadI(c.fXyq[0][j],c.fXyq[1][j],c.fZpad,ix,iy);
2217     c.fSegmentation[cath]->SetPad(ix,iy);
2218     charge = 0;
2219     for (Int_t i=c.fNpar/3; i>=0; i--) { // sum over tracks
2220       indx = i<2 ? 2*i : 2*i+1;
2221       c.fSegmentation[cath]->SetHit(par[indx],par[indx+1],c.fZpad);
2222       //charge += c.fResponse->IntXY(c.fSegmentation[cath])*par[icl*3+2];
2223       if (c.fNpar == 2) coef = 1;
2224       else coef = i==c.fNpar/3 ? par[indx+2] : 1-coef;
2225       //coef = TMath::Max (coef, 0.);
2226       if (c.fNpar == 8 && i < 2) coef = i==1 ? coef*par[indx+2] : coef - par[7];
2227       //coef = TMath::Max (coef, 0.);
2228       charge += c.fResponse->IntXY(c.fSegmentation[cath])*coef;
2229     }
2230     charge *= c.fQtot;
2231     //if (c.fXyq[2][j] > c.fResponse->MaxAdc()-1 && charge > 
2232     //  c.fResponse->MaxAdc()) charge = c.fResponse->MaxAdc(); 
2233     delta = charge - c.fXyq[2][j];
2234     delta /= TMath::Sqrt ((Double_t)c.fXyq[2][j]);
2235     //chi2 += TMath::Abs(delta);
2236     chi2 += delta*delta;
2237   } // for (Int_t j=0;
2238   f = chi2; 
2239   Double_t qAver = c.fQtot/npads; //(c.fnPads[0]+c.fnPads[1]);
2240   f = chi2/qAver;
2241 }
2242
2243 //_____________________________________________________________________________
2244 void AliMUONClusterFinderAZ::UpdatePads(Int_t /*nfit*/, Double_t *par)
2245 {
2246   // Subtract the fitted charges from pads with strong coupling
2247
2248   Int_t cath, ix, iy, indx;
2249   Double_t charge, coef=0;
2250   for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2251     if (fPadIJ[1][j] != -1) continue;
2252     if (fNpar != 0) {
2253       cath = fPadIJ[0][j];
2254       fSegmentation[cath]->GetPadI(fXyq[0][j],fXyq[1][j],fZpad,ix,iy);
2255       fSegmentation[cath]->SetPad(ix,iy);
2256       charge = 0;
2257       for (Int_t i=fNpar/3; i>=0; i--) { // sum over tracks
2258         indx = i<2 ? 2*i : 2*i+1;
2259         fSegmentation[cath]->SetHit(par[indx],par[indx+1],fZpad);
2260         if (fNpar == 2) coef = 1;
2261         else coef = i==fNpar/3 ? par[indx+2] : 1-coef;
2262         if (fNpar == 8 && i < 2) coef = i==1 ? coef*par[indx+2] : coef - par[7];
2263         charge += fResponse->IntXY(fSegmentation[cath])*coef;
2264       }
2265       charge *= fQtot;
2266       fXyq[2][j] -= charge;
2267     } // if (fNpar != 0)
2268     if (fXyq[2][j] > fResponse->ZeroSuppression()) fPadIJ[1][j] = 0; // return pad for further using
2269   } // for (Int_t j=0;
2270 }  
2271
2272 //_____________________________________________________________________________
2273 Bool_t AliMUONClusterFinderAZ::TestTrack(Int_t /*t*/) {
2274 // Test if track was user selected
2275   return kTRUE;
2276   /*
2277     if (fTrack[0]==-1 || fTrack[1]==-1) {
2278         return kTRUE;
2279     } else if (t==fTrack[0] || t==fTrack[1]) {
2280         return kTRUE;
2281     } else {
2282         return kFALSE;
2283     }
2284   */
2285 }
2286
2287 //_____________________________________________________________________________
2288 void AliMUONClusterFinderAZ::AddRawCluster(Double_t x, Double_t y, Double_t fmin)
2289 {
2290   //
2291   // Add a raw cluster copy to the list
2292   //
2293   AliMUONRawCluster cnew;
2294   AliMUON *pMUON=(AliMUON*)gAlice->GetModule("MUON");
2295   //pMUON->AddRawCluster(fInput->Chamber(),c); 
2296
2297   Int_t cath;    
2298   for (cath=0; cath<2; cath++) {
2299     cnew.fX[cath] = x;
2300     cnew.fY[cath] = y;
2301     cnew.fZ[cath] = fZpad;
2302     cnew.fQ[cath] = 100;
2303     cnew.fPeakSignal[cath] = 20;
2304     cnew.fMultiplicity[cath] = 5;
2305     cnew.fNcluster[cath] = 1;
2306     cnew.fChi2[cath] = fmin; //0.1;
2307     /*
2308     cnew.fMultiplicity[cath]=c->fMultiplicity[cath];
2309     for (i=0; i<fMul[cath]; i++) {
2310       cnew.fIndexMap[i][cath]=c->fIndexMap[i][cath];
2311       fSeg[cath]->SetPad(fIx[i][cath], fIy[i][cath]);
2312     }
2313     fprintf(stderr,"\nRawCluster %d cath %d\n",ico,cath);
2314     fprintf(stderr,"mult_av %d\n",c->fMultiplicity[cath]);
2315     FillCluster(&cnew,cath);
2316     */
2317   } 
2318   //cnew.fClusterType=cnew.PhysicsContribution();
2319   pMUON->GetMUONData()->AddRawCluster(AliMUONClusterInput::Instance()->Chamber(),cnew); 
2320   //fNPeaks++;
2321 }
2322
2323 //_____________________________________________________________________________
2324 Int_t AliMUONClusterFinderAZ::FindLocalMaxima(Int_t *localMax, Double_t *maxVal)
2325 {
2326   // Find local maxima in pixel space for large preclusters in order to
2327   // try to split them into smaller pieces (to speed up the MLEM procedure)
2328
2329   TH2D *hist = (TH2D*) gROOT->FindObject("anode");
2330   if (hist) hist->Delete();
2331
2332   Double_t xylim[4] = {999, 999, 999, 999};
2333   Int_t nPix = fPixArray->GetEntriesFast();
2334   AliMUONPixel *pixPtr = 0;
2335   for (Int_t ipix=0; ipix<nPix; ipix++) {
2336     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
2337     for (Int_t i=0; i<4; i++) 
2338          xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
2339   }
2340   for (Int_t i=0; i<4; i++) xylim[i] -= pixPtr->Size(i/2); 
2341
2342   Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
2343   Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
2344   hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
2345   for (Int_t ipix=0; ipix<nPix; ipix++) {
2346     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
2347     hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
2348   }
2349   if (fDraw) {
2350     ((TCanvas*)gROOT->FindObject("c2"))->cd();
2351     gPad->SetTheta(55);
2352     gPad->SetPhi(30);
2353     hist->Draw("lego1Fb");
2354     gPad->Update();
2355     int ia;
2356     cin >> ia;
2357   }
2358
2359   Int_t nMax = 0, indx;
2360   Int_t *isLocalMax = new Int_t[ny*nx];
2361   for (Int_t i=0; i<ny*nx; i++) isLocalMax[i] = 0;
2362
2363   for (Int_t i=1; i<=ny; i++) {
2364     indx = (i-1) * nx;
2365     for (Int_t j=1; j<=nx; j++) {
2366       if (hist->GetCellContent(j,i) < 0.5) continue;
2367       //if (isLocalMax[indx+j-1] < 0) continue;
2368       if (isLocalMax[indx+j-1] != 0) continue;
2369       FlagLocalMax(hist, i, j, isLocalMax);
2370     }
2371   }
2372
2373   for (Int_t i=1; i<=ny; i++) {
2374     indx = (i-1) * nx;
2375     for (Int_t j=1; j<=nx; j++) {
2376       if (isLocalMax[indx+j-1] > 0) { 
2377         localMax[nMax] = indx + j - 1; 
2378         maxVal[nMax++] = hist->GetCellContent(j,i);
2379       }
2380       if (nMax > 99) { cout << " Too many local maxima !!!" << endl; ::exit(0); }
2381     }
2382   }
2383   cout << " Local max: " << nMax << endl;
2384   delete [] isLocalMax; isLocalMax = 0;
2385   return nMax;
2386 }
2387
2388 //_____________________________________________________________________________
2389 void AliMUONClusterFinderAZ::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
2390 {
2391   // Flag pixels (whether or not local maxima)
2392
2393   Int_t nx = hist->GetNbinsX();
2394   Int_t ny = hist->GetNbinsY();
2395   Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
2396   Int_t cont1 = 0;
2397
2398   for (Int_t i1=i-1; i1<i+2; i1++) {
2399     if (i1 < 1 || i1 > ny) continue;
2400     for (Int_t j1=j-1; j1<j+2; j1++) {
2401       if (j1 < 1 || j1 > nx) continue;
2402       if (i == i1 && j == j1) continue;
2403       cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
2404       if (cont < cont1) { isLocalMax[(i-1)*nx+j-1] = -1; return; }
2405       else if (cont > cont1) isLocalMax[(i1-1)*nx+j1-1] = -1;
2406       else { // the same charge
2407         isLocalMax[(i-1)*nx+j-1] = 1; 
2408         if (isLocalMax[(i1-1)*nx+j1-1] == 0) {
2409           FlagLocalMax(hist, i1, j1, isLocalMax);
2410           if (isLocalMax[(i1-1)*nx+j1-1] < 0) { isLocalMax[(i-1)*nx+j-1] = -1; return; }
2411           else isLocalMax[(i1-1)*nx+j1-1] = -1;
2412         }
2413       } 
2414     }
2415   }
2416   isLocalMax[(i-1)*nx+j-1] = 1; // local maximum
2417 }
2418
2419 //_____________________________________________________________________________
2420 void AliMUONClusterFinderAZ::FindCluster(Int_t *localMax, Int_t iMax)
2421 {
2422   // Find pixel cluster around local maximum #iMax and pick up pads
2423   // overlapping with it
2424
2425   TH2D *hist = (TH2D*) gROOT->FindObject("anode");
2426   Int_t nx = hist->GetNbinsX();
2427   Int_t ny = hist->GetNbinsY();
2428   Int_t ic = localMax[iMax] / nx + 1;
2429   Int_t jc = localMax[iMax] % nx + 1;
2430   Bool_t *used = new Bool_t[ny*nx];
2431   for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE;
2432
2433   // Drop all pixels from the array - pick up only the ones from the cluster
2434   fPixArray->Delete();
2435
2436   Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2; 
2437   Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;  
2438   Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
2439   Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
2440   Double_t cont = hist->GetCellContent(jc,ic);
2441   AliMUONPixel *pixPtr = new AliMUONPixel (xc, yc, wx, wy, cont);
2442   fPixArray->Add((TObject*)pixPtr);
2443   used[(ic-1)*nx+jc-1] = kTRUE;
2444   AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
2445
2446   Int_t nPix = fPixArray->GetEntriesFast(), npad = fnPads[0] + fnPads[1];
2447   for (Int_t i=0; i<nPix; i++) {
2448     ((AliMUONPixel*)fPixArray->UncheckedAt(i))->SetSize(0,wx); 
2449     ((AliMUONPixel*)fPixArray->UncheckedAt(i))->SetSize(1,wy); 
2450   }
2451   cout << iMax << " " << nPix << endl;
2452
2453   Float_t xy[4], xy12[4];
2454   // Pick up pads which overlap with found pixels
2455   for (Int_t i=0; i<npad; i++) fPadIJ[1][i] = -1;
2456   for (Int_t i=0; i<nPix; i++) {
2457     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
2458     for (Int_t j=0; j<4; j++) 
2459       xy[j] = pixPtr->Coord(j/2) + (j%2 ? 1 : -1)*pixPtr->Size(j/2);
2460     for (Int_t j=0; j<npad; j++) 
2461       if (Overlap(xy, j, xy12, 0)) fPadIJ[1][j] = 0; // flag for use
2462   }
2463
2464   delete [] used; used = 0;
2465 }