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