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