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