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