]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterFinderAZ.cxx
Raw2SDigits method implemented
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderAZ.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 // -------------------------------
19 // Class AliMUONClusterFinderAZ
20 // -------------------------------
21 // Clusterizer class based on the Expectation-Maximization algorithm
22 // Author: Alexander Zinchenko, JINR Dubna
23
24 #include <stdlib.h>
25 #include <Riostream.h>
26 #include <TH2.h>
27 #include <TMinuit.h>
28 #include <TMatrixD.h>
29 #include <TRandom.h>
30 #include <TROOT.h>
31
32 #include "AliMUONClusterFinderAZ.h"
33 #include "AliMUONClusterDrawAZ.h"
34 #include "AliMUONVGeometryDESegmentation.h"
35 #include "AliMUONGeometryModuleTransformer.h"
36 #include "AliMUON.h"
37 #include "AliMUONDigit.h"
38 #include "AliMUONRawCluster.h"
39 #include "AliMUONClusterInput.h"
40 #include "AliMUONPixel.h"
41 #include "AliMUONMathieson.h"
42 #include "AliLog.h"
43
44 /// \cond CLASSIMP
45 ClassImp(AliMUONClusterFinderAZ)
46 /// \endcond
47  
48  const Double_t AliMUONClusterFinderAZ::fgkCouplMin = 1.e-3; // threshold on coupling 
49  const Double_t AliMUONClusterFinderAZ::fgkZeroSuppression = 6; // average zero suppression value
50  const Double_t AliMUONClusterFinderAZ::fgkSaturation = 3000; // average saturation level
51  AliMUONClusterFinderAZ* AliMUONClusterFinderAZ::fgClusterFinder = 0x0;
52  TMinuit* AliMUONClusterFinderAZ::fgMinuit = 0x0;
53 //FILE *lun1 = fopen("nxny.dat","w");
54
55 //_____________________________________________________________________________
56 AliMUONClusterFinderAZ::AliMUONClusterFinderAZ(Bool_t draw)
57   : AliMUONClusterFinderVS(),
58     fZpad(0),
59     fNpar(0),
60     fQtot(0),
61     fReco(1),
62     fCathBeg(0),
63     fDraw(0x0),
64     fPixArray(0x0),
65     fnCoupled(0),
66     fDebug(0)
67 {
68 /// Constructor
69   fnPads[0]=fnPads[1]=0;
70   
71   for (Int_t i=0; i<7; i++)
72     for (Int_t j=0; j<fgkDim; j++)
73       fXyq[i][j]= 9999.;
74
75   for (Int_t i=0; i<4; i++)
76     for (Int_t j=0; j<fgkDim; j++) 
77       fPadIJ[i][j]=-1;
78
79   for (Int_t i=0; i<2; i++)
80     for (Int_t j=0; j<fgkDim; j++) 
81       fUsed[i][j] = 0;
82
83   fSegmentation[1] = fSegmentation[0] = 0x0; 
84
85   fPadBeg[0] = fPadBeg[1] = 0;
86
87   if (!fgMinuit) fgMinuit = new TMinuit(8);
88   if (!fgClusterFinder) fgClusterFinder = this;
89   fPixArray = new TObjArray(20); 
90
91   if (draw) {
92     fDebug = 1;
93     fReco = 0;
94     fDraw = new AliMUONClusterDrawAZ(this);
95   }
96   cout << " *** Running AZ cluster finder *** " << endl;
97 }
98
99 //_____________________________________________________________________________
100 AliMUONClusterFinderAZ::~AliMUONClusterFinderAZ()
101 {
102 /// Destructor
103   delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
104   delete fDraw;
105 }
106
107 //_____________________________________________________________________________
108 void AliMUONClusterFinderAZ::FindRawClusters()
109 {
110 /// To provide the same interface as in AliMUONClusterFinderVS
111
112   ResetRawClusters(); 
113   EventLoop (fEvtNumber, fInput->Chamber());
114 }
115
116 //_____________________________________________________________________________
117 void AliMUONClusterFinderAZ::EventLoop(Int_t nev, Int_t ch)
118 {
119 /// Loop over digits
120   
121   if (fDraw && !fDraw->FindEvCh(nev, ch)) return;
122
123   fSegmentation[0] = (AliMUONVGeometryDESegmentation*) fInput->
124     Segmentation2(0)->GetDESegmentation(fInput->DetElemId());
125   fSegmentation[1] = (AliMUONVGeometryDESegmentation*) fInput->
126     Segmentation2(1)->GetDESegmentation(fInput->DetElemId());
127     
128   Int_t ndigits[2] = {9,9}, nShown[2] = {0};
129   if (fReco != 2) { // skip initialization for the combined cluster / track
130     fCathBeg = fPadBeg[0] = fPadBeg[1] = 0;
131     for (Int_t i = 0; i < 2; i++) {
132       for (Int_t j = 0; j < fgkDim; j++) { fUsed[i][j] = kFALSE; }
133     }
134   }
135
136 next:
137   if (fReco == 2 && (nShown[0] || nShown[1])) return; // only one precluster for the combined finder
138   if (ndigits[0] == nShown[0] && ndigits[1] == nShown[1]) return;
139
140   Float_t xpad, ypad, zpad, zpad0;
141   Bool_t first = kTRUE;
142   if (fDebug) cout << " *** Event # " << nev << " chamber: " << ch << endl;
143   fnPads[0] = fnPads[1] = 0;
144   for (Int_t i = 0; i < fgkDim; i++) fPadIJ[1][i] = 0; 
145
146   for (Int_t iii = fCathBeg; iii < 2; iii++) { 
147     Int_t cath = TMath::Odd(iii);
148     ndigits[cath] = fInput->NDigits(cath); 
149     if (!ndigits[0] && !ndigits[1]) return;
150     if (ndigits[cath] == 0) continue;
151     if (fDebug) cout << " ndigits: " << ndigits[cath] << " " << cath << endl;
152
153     AliMUONDigit  *mdig;
154     Int_t digit;
155
156     Bool_t eEOC = kTRUE; // end-of-cluster
157     for (digit = fPadBeg[cath]; digit < ndigits[cath]; digit++) {
158       mdig = AliMUONClusterInput::Instance()->Digit(cath,digit); 
159       if (first) {
160         // Find first unused pad
161         if (fUsed[cath][digit]) continue;
162         //if (!fSegmentation[cath]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad0)) { 
163         if (!fSegmentation[cath]->HasPad(mdig->PadX(), mdig->PadY())) {
164           // Handle "non-existing" pads
165           fUsed[cath][digit] = kTRUE; 
166           continue; 
167         } 
168         fSegmentation[cath]->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad0); 
169       } else {
170         if (fUsed[cath][digit]) continue;
171         //if (!fSegmentation[cath]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad)) {
172         if (!fSegmentation[cath]->HasPad(mdig->PadX(), mdig->PadY())) {
173           // Handle "non-existing" pads
174           fUsed[cath][digit] = kTRUE; 
175           continue; 
176         } 
177         fSegmentation[cath]->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad); 
178         //if (TMath::Abs(zpad-zpad0) > 0.1) continue; // different slats
179         // Find a pad overlapping with the cluster
180         if (!Overlap(cath,mdig)) continue;
181       }
182       // Add pad - recursive call
183       AddPad(cath,digit);
184       //AZ !!!!!! Temporary fix of St1 overlap regions !!!!!!!! 
185       /*
186       if (cath && ch < 2) {
187         Int_t npads = fnPads[0] + fnPads[1] - 1;
188         Int_t cath1 = fPadIJ[0][npads];
189         Int_t idig = TMath::Nint (fXyq[5][npads]);
190         mdig = AliMUONClusterInput::Instance()->Digit(cath1,idig);
191         //fSegmentation[cath1]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad);
192         fSegmentation[cath1]->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
193         if (TMath::Abs(zpad-zpad0) > 0.1) zpad0 = zpad;
194       } 
195       */
196       eEOC = kFALSE;
197       if (digit >= 0) break;
198     }
199     if (first && eEOC) {
200       // No more unused pads 
201       if (cath == 0) continue; // on cathode #0 - check #1
202       else return; // No more clusters
203     }
204     if (eEOC) break; // cluster found
205     first = kFALSE;
206     if (fDebug) cout << " nPads: " << fnPads[cath] << " " << nShown[cath]+fnPads[cath] << " " << cath << endl;
207   } // for (Int_t iii = 0;
208
209   fZpad = zpad0;
210   if (fDraw) fDraw->DrawCluster(); 
211
212   // Use MLEM for cluster finder
213   Int_t nMax = 1, localMax[100], maxPos[100];
214   Double_t maxVal[100];
215   
216   if (CheckPrecluster(nShown)) {
217     BuildPixArray();
218     //*
219     if (fnPads[0]+fnPads[1] > 50) nMax = FindLocalMaxima(fPixArray, localMax, maxVal);
220     if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
221     Int_t iSimple = 0, nInX = -1, nInY;
222     PadsInXandY(nInX, nInY);
223     if (fDebug) cout << "Pads in X and Y: " << nInX << " " << nInY << endl;
224     if (nMax == 1 && nInX < 4 && nInY < 4) iSimple = 1; //1; // simple cluster
225     //*/
226     /* For test
227     Int_t iSimple = 0, nInX = -1, nInY;
228     PadsInXandY(nInX, nInY);
229     if (fDebug) cout << "Pads in X and Y: " << nInX << " " << nInY << endl;
230     if (nMax == 1 && nInX < 4 && nInY < 4) iSimple = 1; //1; // simple cluster
231     if (!iSimple) nMax = FindLocalMaxima(fPixArray, localMax, maxVal);
232     nMax = 1;
233     if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
234     */
235     for (Int_t i=0; i<nMax; i++) {
236       if (nMax > 1) FindCluster(localMax, maxPos[i]);
237       MainLoop(iSimple);
238       if (i < nMax-1) {
239         for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
240           if (fPadIJ[1][j] == 0) continue; // pad charge was not modified
241           fPadIJ[1][j] = 0;
242           fXyq[2][j] = fXyq[6][j]; // use backup charge value
243         }
244       }
245     } // for (Int_t i=0; i<nMax;
246     if (nMax > 1) ((TH2D*) gROOT->FindObject("anode"))->Delete();
247     TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
248     if (mlem) mlem->Delete();
249   }
250   if (!fDraw || fDraw->Next()) goto next;
251 }
252
253 //_____________________________________________________________________________
254 void AliMUONClusterFinderAZ::AddPad(Int_t cath, Int_t digit)
255 {
256 /// Add pad to the cluster
257
258   AliMUONDigit *mdig = fInput->Digit(cath,digit); 
259
260   Float_t charge = mdig->Signal();
261   // get the center of the pad
262   Float_t xpad, ypad, zpad0;
263   //if (!fSegmentation[cath]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad0)) {   // Handle "non-existing" pads
264   if (!fSegmentation[cath]->HasPad(mdig->PadX(), mdig->PadY())) {
265     fUsed[cath][digit] = kTRUE; 
266     return; 
267   } 
268   fSegmentation[cath]->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad0);
269   Int_t isec = fSegmentation[cath]->Sector(mdig->PadX(), mdig->PadY());
270   Int_t nPads = fnPads[0] + fnPads[1];
271   fXyq[0][nPads] = xpad;
272   fXyq[1][nPads] = ypad;
273   fXyq[2][nPads] = charge;
274   fXyq[3][nPads] = fSegmentation[cath]->Dpx(isec)/2;
275   fXyq[4][nPads] = fSegmentation[cath]->Dpy(isec)/2;
276   fXyq[5][nPads] = digit;
277   fXyq[6][nPads] = 0;
278   fPadIJ[0][nPads] = cath;
279   fPadIJ[1][nPads] = 0;
280   fPadIJ[2][nPads] = mdig->PadX();
281   fPadIJ[3][nPads] = mdig->PadY();
282   fUsed[cath][digit] = kTRUE;
283   if (fDebug) printf(" bbb %d %d %f %f %f %f %f %f %3d %3d \n", nPads, cath, 
284                      xpad, ypad, zpad0, fXyq[3][nPads]*2, fXyq[4][nPads]*2, 
285                      charge, mdig->PadX(), mdig->PadY());
286   fnPads[cath]++;
287
288   // Check neighbours
289   Int_t nn, ix, iy, xList[10], yList[10];
290   AliMUONDigit  *mdig1;
291
292   Int_t ndigits = fInput->NDigits(cath); 
293   fSegmentation[cath]->Neighbours(mdig->PadX(), mdig->PadY(), &nn, xList, yList); 
294   for (Int_t in = 0; in < nn; in++) {
295     ix = xList[in];
296     iy = yList[in];
297     for (Int_t digit1 = 0; digit1 < ndigits; digit1++) {
298       if (digit1 == digit) continue;
299       mdig1 = fInput->Digit(cath,digit1); 
300       if (!fUsed[cath][digit1] && mdig1->PadX() == ix && mdig1->PadY() == iy) {
301         fUsed[cath][digit1] = kTRUE;
302         // Add pad - recursive call
303         AddPad(cath,digit1);
304       }
305     } //for (Int_t digit1 = 0;
306   } // for (Int_t in = 0;
307 }
308
309 //_____________________________________________________________________________
310 Bool_t AliMUONClusterFinderAZ::Overlap(Int_t cath, AliMUONDigit *mdig)
311 {
312 /// Check if the pad from one cathode overlaps with a pad 
313 /// in the precluster on the other cathode
314
315   Float_t xpad, ypad, zpad;
316   fSegmentation[cath]->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
317   Int_t isec = fSegmentation[cath]->Sector(mdig->PadX(), mdig->PadY());
318
319   Float_t xy1[4], xy12[4];
320   xy1[0] = xpad - fSegmentation[cath]->Dpx(isec)/2;
321   xy1[1] = xy1[0] + fSegmentation[cath]->Dpx(isec);
322   xy1[2] = ypad - fSegmentation[cath]->Dpy(isec)/2;
323   xy1[3] = xy1[2] + fSegmentation[cath]->Dpy(isec);
324   //cout << " ok " << fnPads[0]+fnPads[1] << xy1[0] << xy1[1] << xy1[2] << xy1[3] << endl;
325
326   Int_t cath1 = TMath::Even(cath);
327   for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
328     if (fPadIJ[0][i] != cath1) continue;
329     if (Overlap(xy1, i, xy12, 0)) return kTRUE;
330   }
331   return kFALSE;
332 }
333
334 //_____________________________________________________________________________
335 Bool_t AliMUONClusterFinderAZ::Overlap(Float_t *xy1, Int_t iPad, Float_t *xy12, Int_t iSkip)
336 {
337 /// Check if the pads xy1 and iPad overlap and return overlap area
338
339   Float_t xy2[4];
340   xy2[0] = fXyq[0][iPad] - fXyq[3][iPad];
341   xy2[1] = fXyq[0][iPad] + fXyq[3][iPad];
342   if (xy1[0] > xy2[1]-1.e-4 || xy1[1] < xy2[0]+1.e-4) return kFALSE;
343   xy2[2] = fXyq[1][iPad] - fXyq[4][iPad];
344   xy2[3] = fXyq[1][iPad] + fXyq[4][iPad];
345   if (xy1[2] > xy2[3]-1.e-4 || xy1[3] < xy2[2]+1.e-4) return kFALSE;
346   if (!iSkip) return kTRUE; // just check overlap (w/out computing the area)
347   xy12[0] = TMath::Max (xy1[0],xy2[0]);
348   xy12[1] = TMath::Min (xy1[1],xy2[1]);
349   xy12[2] = TMath::Max (xy1[2],xy2[2]);
350   xy12[3] = TMath::Min (xy1[3],xy2[3]);
351   return kTRUE;
352 }
353
354 //_____________________________________________________________________________
355 Bool_t AliMUONClusterFinderAZ::CheckPrecluster(Int_t *nShown)
356 {
357 /// Check precluster in order to attempt to simplify it (mostly for
358 /// two-cathode preclusters)
359
360   Int_t i1, i2, cath=0, digit=0;
361   Float_t xy1[4], xy12[4];
362   
363   Int_t npad = fnPads[0] + fnPads[1];
364   if (npad == 1) { 
365     // Disregard one-pad clusters (leftovers from splitting)
366     nShown[0] += fnPads[0]; 
367     nShown[1] += fnPads[1]; 
368     return kFALSE;
369   }
370
371   // If pads have the same size take average of pads on both cathodes 
372   //Int_t sameSize = (fnPads[0] && fnPads[1]) ? 1 : 0;
373   Int_t sameSize = 0; //AZ - 17-01-06
374   
375   if (sameSize) {
376     Double_t xSize = -1, ySize = 0;
377     for (Int_t i=0; i<npad; i++) {
378       if (fXyq[2][i] < 0) continue;
379       if (xSize < 0) { xSize = fXyq[3][i]; ySize = fXyq[4][i]; }
380       if (TMath::Abs(xSize-fXyq[3][i]) > 1.e-4 ||  TMath::Abs(ySize-fXyq[4][i]) > 1.e-4) { sameSize = 0; break; }
381     }
382   } // if (sameSize)
383   if (sameSize && fnPads[0] == 1 && fnPads[1] == 1) sameSize = 0; //AZ
384   // Handle shift by half a pad in Station 1
385   if (sameSize) {
386     Int_t cath0 = fPadIJ[0][0];
387     for (Int_t i = 1; i < npad; i++) {
388       if (fPadIJ[0][i] == cath0) continue;
389       Double_t dx = TMath::Abs ((fXyq[0][i] - fXyq[0][0]) / fXyq[3][i] / 2);
390       Int_t idx = (Int_t) TMath::Abs ((fXyq[0][i] - fXyq[0][0]) / fXyq[3][i] / 2);
391       if (TMath::Abs (dx - idx) > 0.001) sameSize = 0;
392       break;
393     }
394   } // if (sameSize)
395    
396   if (sameSize && (fnPads[0] >= 2 || fnPads[1] >= 2)) {
397     nShown[0] += fnPads[0];
398     nShown[1] += fnPads[1];
399     fnPads[0] = fnPads[1] = 0;
400     Int_t div;
401     for (Int_t i=0; i<npad; i++) {
402       if (fXyq[2][i] < 0) continue; // used pad
403       fXyq[2][fnPads[0]] = fXyq[2][i];
404       div = 1;
405       cath = fPadIJ[0][i];
406       for (Int_t j=i+1; j<npad; j++) {
407         if (fPadIJ[0][j] == fPadIJ[0][i]) continue; // same cathode
408         if (TMath::Abs(fXyq[0][j]-fXyq[0][i]) > 1.e-4) continue;
409         if (TMath::Abs(fXyq[1][j]-fXyq[1][i]) > 1.e-4) continue;
410         fXyq[2][fnPads[0]] += fXyq[2][j];
411         div = 2;
412         fXyq[2][j] = -2;
413         if (cath) fXyq[5][fnPads[0]] = fXyq[5][j]; // save digit number for cath 0
414         break;
415       }
416       // Flag that the digit from the other cathode
417       if (cath && div == 1) fXyq[5][fnPads[0]] = -fXyq[5][i] - 1; 
418       // If low pad charge take the other equal to 0 
419       //if (div == 1 && fXyq[2][fnPads[0]] < fgkZeroSuppression + 1.5*3) div = 2;
420       fXyq[2][fnPads[0]] /= div;
421       fXyq[0][fnPads[0]] = fXyq[0][i];
422       fXyq[1][fnPads[0]] = fXyq[1][i];
423       fPadIJ[2][fnPads[0]] = fPadIJ[2][i];
424       fPadIJ[3][fnPads[0]] = fPadIJ[3][i];
425       fPadIJ[0][fnPads[0]++] = 0;
426     }
427   } // if (sameSize)
428
429   // Check if one-cathode precluster
430   i1 = fnPads[0]!=0 ? 0 : 1;
431   i2 = fnPads[1]!=0 ? 1 : 0;
432
433   if (i1 != i2) { // two-cathode 
434
435     Int_t *flags = new Int_t[npad];
436     for (Int_t i=0; i<npad; i++) { flags[i] = 0; }
437
438     // Check pad overlaps
439     for (Int_t i=0; i<npad; i++) {
440       if (fPadIJ[0][i] != i1) continue;
441       xy1[0] = fXyq[0][i] - fXyq[3][i];
442       xy1[1] = fXyq[0][i] + fXyq[3][i];
443       xy1[2] = fXyq[1][i] - fXyq[4][i];
444       xy1[3] = fXyq[1][i] + fXyq[4][i];
445       for (Int_t j=0; j<npad; j++) {
446         if (fPadIJ[0][j] != i2) continue;
447         if (!Overlap(xy1, j, xy12, 0)) continue;
448         flags[i] = flags[j] = 1; // mark overlapped pads
449       } // for (Int_t j=0;
450     } // for (Int_t i=0;
451
452     // Check if all pads overlap
453     Int_t nFlags=0;
454     for (Int_t i=0; i<npad; i++) {
455       if (flags[i]) continue;
456       nFlags ++;
457       if (fDebug) cout << i << " " << fPadIJ[0][i] << " " << fXyq[0][i] << " " << fXyq[1][i] << endl;
458     }
459     if (fDebug && nFlags) cout << " nFlags = " << nFlags << endl;
460     //if (nFlags > 2 || (Float_t)nFlags / npad > 0.2) { // why 2 ??? - empirical choice
461     if (nFlags > 0) {
462       for (Int_t i=0; i<npad; i++) {
463         if (flags[i]) continue;
464         digit = TMath::Nint (fXyq[5][i]);
465         cath = fPadIJ[0][i];
466         // Check for edge effect (missing pads on the other cathode)
467         Int_t cath1 = TMath::Even(cath), ix, iy;
468         ix = iy = 0;
469         //if (!fSegmentation[cath1]->GetPadI(fInput->DetElemId(),fXyq[0][i],fXyq[1][i],fZpad,ix,iy)) continue;
470         if (!fSegmentation[cath1]->HasPad(fXyq[0][i], fXyq[1][i], fZpad)) continue;
471         if (nFlags == 1 && fXyq[2][i] < fgkZeroSuppression * 3) continue;
472         fUsed[cath][digit] = kFALSE; // release pad
473         fXyq[2][i] = -2;
474         fnPads[cath]--;
475       }
476       if (fDraw) fDraw->UpdateCluster(npad);
477     } // if (nFlags > 2)
478
479     // Check correlations of cathode charges
480     if (fnPads[0] && fnPads[1]) { // two-cathode
481       Double_t sum[2]={0};
482       Int_t over[2] = {1, 1};
483       for (Int_t i=0; i<npad; i++) {
484         cath = fPadIJ[0][i];
485         if (fXyq[2][i] > 0) sum[cath] += fXyq[2][i];
486         if (fXyq[2][i] > fgkSaturation-1) over[cath] = 0;
487       }
488       if (fDebug) cout << " Total charge: " << sum[0] << " " << sum[1] << endl;
489       if ((over[0] || over[1]) && TMath::Abs(sum[0]-sum[1])/(sum[0]+sum[1])*2 > 1) { // 3 times difference
490         if (fDebug) cout << " Release " << endl;
491         // Big difference
492         cath = sum[0] > sum[1] ? 0 : 1;
493         Int_t imax = 0, imin = 0;
494         Double_t cmax = -1, cmin = 9999, dxMin = 0, dyMin = 0;
495         Double_t *dist = new Double_t[npad];
496         for (Int_t i = 0; i < npad; i++) {
497           if (fPadIJ[0][i] != cath || fXyq[2][i] < 0) continue;
498           if (fXyq[2][i] < cmin) {
499             cmin = fXyq[2][i];
500             imin = i;
501           }
502           if (fXyq[2][i] < cmax) continue;
503           cmax = fXyq[2][i];
504           imax = i;
505         }
506         // Arrange pads according to their distance to the max, 
507         // normalized to the pad size
508         for (Int_t i = 0; i < npad; i++) {
509           dist[i] = 0;
510           if (fPadIJ[0][i] != cath || fXyq[2][i] < 0) continue;
511           if (i == imax) continue; 
512           Double_t dx = (fXyq[0][i] - fXyq[0][imax]) / fXyq[3][imax] / 2;
513           Double_t dy = (fXyq[1][i] - fXyq[1][imax]) / fXyq[4][imax] / 2;
514           dist[i] = TMath::Sqrt (dx * dx + dy * dy);
515           if (i == imin) {
516             cmin = dist[i] + 0.001; // distance to the pad with minimum charge 
517             dxMin = dx;
518             dyMin = dy;
519           }
520         }
521         TMath::Sort(npad, dist, flags, kFALSE); // in increasing order
522         Int_t indx;
523         Double_t xmax = -1;
524         for (Int_t i = 0; i < npad; i++) {
525           indx = flags[i];
526           if (fPadIJ[0][indx] != cath || fXyq[2][indx] < 0) continue;
527           if (dist[indx] > cmin) {
528             // Farther than the minimum pad
529             Double_t dx = (fXyq[0][indx] - fXyq[0][imax]) / fXyq[3][imax] / 2;
530             Double_t dy = (fXyq[1][indx] - fXyq[1][imax]) / fXyq[4][imax] / 2;
531             dx *= dxMin;
532             dy *= dyMin;
533             if (dx >= 0 && dy >= 0) continue;
534             if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
535             if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
536           }
537           if (fXyq[2][indx] <= cmax || TMath::Abs(dist[indx]-xmax) < 1.e-3) {
538             // Release pads
539             if (TMath::Abs(dist[indx]-xmax) < 1.e-3) 
540                 cmax = TMath::Max((Double_t)(fXyq[2][indx]),cmax);
541             else cmax = fXyq[2][indx];
542             xmax = dist[indx];
543             digit = TMath::Nint (fXyq[5][indx]);
544             fUsed[cath][digit] = kFALSE; 
545             fXyq[2][indx] = -2;
546             fnPads[cath]--;
547           } 
548         } // for (Int_t i = 0; i < npad;
549
550         // Check pad overlaps once more
551         for (Int_t j = 0; j < npad; j++) flags[j] = 0; 
552         for (Int_t k = 0; k < npad; k++) {
553           if (fXyq[2][k] < 0 || fPadIJ[0][k] != i1) continue;
554           xy1[0] = fXyq[0][k] - fXyq[3][k];
555           xy1[1] = fXyq[0][k] + fXyq[3][k];
556           xy1[2] = fXyq[1][k] - fXyq[4][k];
557           xy1[3] = fXyq[1][k] + fXyq[4][k];
558           for (Int_t j = 0; j < npad; j++) {
559             if (fXyq[2][j] < 0) continue;
560             if (fPadIJ[0][j] != i2) continue;
561             if (!Overlap(xy1, j, xy12, 0)) continue;
562             flags[k] = flags[j] = 1; // mark overlapped pads
563           } // for (Int_t j = 0;
564         } // for (Int_t k = 0;
565         nFlags = 0;
566         for (Int_t j = 0; j < npad; j++) {
567           if (fXyq[2][j] < 0 || flags[j]) continue;
568           nFlags ++;
569         }
570         if (nFlags == fnPads[0] + fnPads[1]) {
571           // No overlap
572           for (Int_t j = 0; j < npad; j++) {
573             if (fXyq[2][j] < 0 || fPadIJ[0][j] != cath) continue;
574             fXyq[2][j] = -2;
575             fnPads[cath]--;
576           }
577         }
578         delete [] dist; dist = 0;
579         if (fDraw) fDraw->UpdateCluster(npad);
580       } // TMath::Abs(sum[0]-sum[1])...
581     } // if (fnPads[0] && fnPads[1])
582     delete [] flags; flags = 0;
583   } // if (i1 != i2) 
584
585   if (!sameSize) { nShown[0] += fnPads[0]; nShown[1] += fnPads[1]; }
586
587   // Move released pads to the right
588   Int_t beg = 0, end = npad-1, padij;
589   Double_t xyq;
590   while (beg < end) {
591     if (fXyq[2][beg] > 0) { beg++; continue; }
592     for (Int_t j=end; j>beg; j--) {
593       if (fXyq[2][j] < 0) continue;
594       end = j - 1;
595       for (Int_t j1=0; j1<4; j1++) {
596         padij = fPadIJ[j1][beg]; 
597         fPadIJ[j1][beg] = fPadIJ[j1][j];
598         fPadIJ[j1][j] = padij;
599       }
600       for (Int_t j1=0; j1<6; j1++) {
601         xyq = fXyq[j1][beg]; 
602         fXyq[j1][beg] = fXyq[j1][j];
603         fXyq[j1][j] = xyq;
604       }
605       break;
606     } // for (Int_t j=end;
607     beg++;
608   } // while
609   npad = fnPads[0] + fnPads[1];
610   if (npad > 500) { 
611     AliWarning(Form(" *** Too large cluster. Give up. %d ", npad));
612     return kFALSE; 
613   }
614   // Back up charge value
615   for (Int_t j = 0; j < npad; j++) fXyq[6][j] = fXyq[2][j];
616
617   return kTRUE;
618 }
619
620 //_____________________________________________________________________________
621 void AliMUONClusterFinderAZ::BuildPixArray()
622 {
623 /// Build pixel array for MLEM method
624   
625   Int_t nPix=0, i1, i2;
626   Float_t xy1[4], xy12[4];
627   AliMUONPixel *pixPtr=0;
628
629   Int_t npad = fnPads[0] + fnPads[1];
630
631   // One cathode is empty
632   i1 = fnPads[0]!=0 ? 0 : 1;
633   i2 = fnPads[1]!=0 ? 1 : 0;
634
635   // Build array of pixels on anode plane
636   if (i1 == i2) { // one-cathode precluster
637     for (Int_t j=0; j<npad; j++) {
638       pixPtr = new AliMUONPixel();
639       for (Int_t i=0; i<2; i++) {
640         pixPtr->SetCoord(i, fXyq[i][j]); // pixel coordinates
641         pixPtr->SetSize(i, fXyq[i+3][j]); // pixel size
642       }
643       pixPtr->SetCharge(fXyq[2][j]); // charge
644       fPixArray->Add((TObject*)pixPtr);
645       nPix++;
646     }
647   } else { // two-cathode precluster    
648     i1 = fPadIJ[0][0];
649     i2 = TMath::Even (i1);
650     for (Int_t i = 0; i < npad; i++) {
651       if (fPadIJ[0][i] != i1) continue;
652       xy1[0] = fXyq[0][i] - fXyq[3][i];
653       xy1[1] = fXyq[0][i] + fXyq[3][i];
654       xy1[2] = fXyq[1][i] - fXyq[4][i];
655       xy1[3] = fXyq[1][i] + fXyq[4][i];
656       for (Int_t j = 1; j < npad; j++) {
657         if (fPadIJ[0][j] != i2) continue;
658         if (!Overlap(xy1, j, xy12, 1)) continue;
659         pixPtr = new AliMUONPixel();
660         for (Int_t k=0; k<2; k++) {
661           pixPtr->SetCoord(k, (xy12[2*k]+xy12[2*k+1])/2); // pixel coordinates
662           pixPtr->SetSize(k, xy12[2*k+1]-pixPtr->Coord(k)); // size
663         }
664         pixPtr->SetCharge(TMath::Min (fXyq[2][i],fXyq[2][j])); //charge
665         fPixArray->Add((TObject*)pixPtr);
666         //cout << nPix << " " << pixPtr->Coord(0) << " " << pixPtr->Size(0) << " " << pixPtr->Coord(1) << " " << pixPtr->Size(1) << " " << pixPtr->Charge() << endl;
667         nPix++;
668       } // for (Int_t j=0;
669     } // for (Int_t i=0;
670   } // else
671
672   Float_t xPadMin = 999, yPadMin = 999;
673   for (Int_t i = 0; i < npad; i++) {
674     xPadMin = TMath::Min (xPadMin, fXyq[3][i]);
675     yPadMin = TMath::Min (yPadMin, fXyq[4][i]);
676   }
677   if (fDebug) cout << xPadMin << " " << yPadMin << endl;
678
679   Float_t wxmin = 999, wymin = 999;
680   for (Int_t i = 0; i < nPix; i++) {
681     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
682     wxmin = TMath::Min ((Double_t)wxmin, pixPtr->Size(0));
683     wymin = TMath::Min ((Double_t)wymin, pixPtr->Size(1));
684   }
685   if (fDebug) cout << wxmin << " " << wymin << endl;
686   wxmin = TMath::Abs (wxmin - xPadMin/2) > 0.001 ? xPadMin : xPadMin / 2;
687   wymin = TMath::Abs (wymin - yPadMin/2) > 0.001 ? yPadMin : yPadMin / 2;
688   //wxmin = xPadMin; wymin = yPadMin;
689
690   // Check if small pixel X-size
691   AdjustPixel(wxmin, 0);
692   // Check if small pixel Y-size
693   AdjustPixel(wymin, 1);
694   // Check if large pixel size
695   AdjustPixel(wxmin, wymin);
696
697   // Remove discarded pixels
698   for (Int_t i=0; i<nPix; i++) {
699     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
700     //pixPtr->Print();
701     if (pixPtr->Charge() < 1) { fPixArray->RemoveAt(i); delete pixPtr; }// discarded pixel
702   }
703   fPixArray->Compress();
704   nPix = fPixArray->GetEntriesFast();
705
706   if (nPix > npad) {
707     if (fDebug) cout << nPix << endl;
708     // Too many pixels - sort and remove pixels with the lowest signal
709     fPixArray->Sort();
710     for (Int_t i=npad; i<nPix; i++) {
711       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
712       //pixPtr->Print();
713       fPixArray->RemoveAt(i);
714       delete pixPtr;
715     }
716     nPix = npad;
717   } // if (nPix > npad)
718
719   // Set pixel charges to the same value (for MLEM)
720   for (Int_t i=0; i<nPix; i++) {
721     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
722     //pixPtr->SetCharge(10);
723     if (fDebug) cout << i+1 << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << " " << pixPtr->Size(0) << " " << pixPtr->Size(1) << endl;
724   }
725 }
726
727 //_____________________________________________________________________________
728 void AliMUONClusterFinderAZ::AdjustPixel(Float_t width, Int_t ixy)
729 {
730 /// Check if some pixels have small size (adjust if necessary)
731
732   AliMUONPixel *pixPtr, *pixPtr1 = 0;
733   Int_t ixy1 = TMath::Even(ixy);
734   Int_t nPix = fPixArray->GetEntriesFast();
735
736   for (Int_t i=0; i<nPix; i++) {
737     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
738     if (pixPtr->Charge() < 1) continue; // discarded pixel
739     if (pixPtr->Size(ixy)-width < -1.e-4) {
740       // try to merge 
741       if (fDebug) cout << i << " Small X or Y: " << ixy << " " << pixPtr->Size(ixy) << " " << width << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << endl;
742       for (Int_t j=i+1; j<nPix; j++) {
743         pixPtr1 = (AliMUONPixel*) fPixArray->UncheckedAt(j);
744         if (pixPtr1->Charge() < 1) continue; // discarded pixel
745         if (TMath::Abs(pixPtr1->Size(ixy)-width) < 1.e-4) continue; // right size 
746         if (TMath::Abs(pixPtr1->Coord(ixy1)-pixPtr->Coord(ixy1)) > 1.e-4) continue; // different rows/columns
747         if (TMath::Abs(pixPtr1->Coord(ixy)-pixPtr->Coord(ixy)) < 2*width) {
748           // merge
749           Double_t tmp = pixPtr->Coord(ixy) + pixPtr1->Size(ixy) * 
750             TMath::Sign (1., pixPtr1->Coord(ixy) - pixPtr->Coord(ixy));
751           pixPtr->SetCoord(ixy, tmp);
752           pixPtr->SetSize(ixy, width);
753           pixPtr->SetCharge(TMath::Min (pixPtr->Charge(),pixPtr1->Charge()));
754           pixPtr1->SetCharge(0);
755           pixPtr1 = 0;
756           break;
757         }
758       } // for (Int_t j=i+1;
759       //if (!pixPtr1) { cout << " I am here!" << endl; pixPtr->SetSize(ixy, width); } // ???
760       //else if (pixPtr1->Charge() > 0.5 || i == nPix-1) {
761       if (pixPtr1 || i == nPix-1) {
762         // edge pixel - just increase its size
763         if (fDebug) cout << " Edge ..." << endl; 
764         for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
765           //if (fPadIJ[0][j] != ixy1) continue;
766           //???-check if (TMath::Abs(pixPtr->Coord(ixy1)-fXyq[ixy1][j]) > 1.e-4) continue;
767           if (pixPtr->Coord(ixy) < fXyq[ixy][j]) 
768             //pixPtr->Shift(ixy, -pixPtr->Size(ixy));
769             pixPtr->Shift(ixy, pixPtr->Size(ixy)-width);
770           //else pixPtr->Shift(ixy, pixPtr->Size(ixy));
771           else pixPtr->Shift(ixy, -pixPtr->Size(ixy)+width);
772           pixPtr->SetSize(ixy, width);
773           break;
774         }
775       }
776     } // if (pixPtr->Size(ixy)-width < -1.e-4)
777   } // for (Int_t i=0; i<nPix;
778   return;
779 }
780   
781 //_____________________________________________________________________________
782 void AliMUONClusterFinderAZ::AdjustPixel(Float_t wxmin, Float_t wymin)
783 {
784 /// Check if some pixels have large size (adjust if necessary)
785
786   Int_t n1[2], n2[2], iOK = 1, nPix = fPixArray->GetEntriesFast();
787   AliMUONPixel *pixPtr, pix;
788   Double_t xy0[2] = {9999, 9999}, wxy[2], dist[2] = {0};
789
790   // Check if large pixel size
791   for (Int_t i = 0; i < nPix; i++) {
792     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
793     if (pixPtr->Charge() < 1) continue; // discarded pixel
794     if (pixPtr->Size(0) - wxmin < 1.e-4) {
795       if (xy0[0] > 9998) xy0[0] = pixPtr->Coord(0); // position of a "normal" pixel
796       if (pixPtr->Size(1) - wymin < 1.e-4) { 
797         if (xy0[1] > 9998) xy0[1] = pixPtr->Coord(1); // position of a "normal" pixel
798         continue;
799       } else iOK = 0; // large pixel
800     } else {
801       iOK = 0; // large pixel
802       if (xy0[1] > 9998 && pixPtr->Size(1) - wymin < 1.e-4) xy0[1] = pixPtr->Coord(1); // "normal" pixel
803     }      
804     if (xy0[0] < 9998 && xy0[1] < 9998) break;
805   }
806   if (iOK) return;
807
808   wxy[0] = wxmin;
809   wxy[1] = wymin;
810   //cout << xy0[0] << " " << xy0[1] << endl;
811   for (Int_t i = 0; i < nPix; i++) {
812     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
813     if (pixPtr->Charge() < 1) continue; // discarded pixel
814     n1[0] = n1[1] = 999;
815     n2[0] = n2[1] = 1;
816     for (Int_t j = 0; j < 2; j++) {
817       if (pixPtr->Size(j) - wxy[j] < 1.e-4) continue;
818       dist[j] = (pixPtr->Coord(j) - xy0[j]) / wxy[j] / 2; // normalized distance to "normal" pixel
819       n2[j] = TMath::Nint (pixPtr->Size(j) / wxy[j]);
820       n1[j] = n2[j] == 1 ? TMath::Nint(dist[j]) : (Int_t)dist[j];
821     }
822     if (n1[0] > 998 && n1[1] > 998) continue;
823     if (fDebug) cout << " Different " << pixPtr->Size(0) << " " << wxy[0] << " "
824                      << pixPtr->Size(1) << " " << wxy[1] <<endl;
825     
826     if (n2[0] > 2 || n2[1] > 2) { 
827       //cout << n2[0] << " " << n2[1] << endl; 
828       if (n2[0] > 2 && n1[0] < 999) n1[0]--;
829       if (n2[1] > 2 && n1[1] < 999) n1[1]--;
830     }
831     //cout << n1[0] << " " << n2[0] << " " << n1[1] << " " << n2[1] << endl;
832     pix = *pixPtr;
833     pix.SetSize(0, wxy[0]); pix.SetSize(1, wxy[1]);
834     //pixPtr->Print();
835     for (Int_t ii = 0; ii < n2[0]; ii++) {
836       if (n1[0] < 999) pix.SetCoord(0, xy0[0] + (n1[0] + TMath::Sign(1.,dist[0]) * ii) * 2 * wxy[0]);
837       for (Int_t jj = 0; jj < n2[1]; jj++) {
838         if (n1[1] < 999) pix.SetCoord(1, xy0[1] + (n1[1] + TMath::Sign(1.,dist[1]) * jj) * 2 * wxy[1]);
839         fPixArray->Add(new AliMUONPixel(pix));
840         //pix.Print();
841       }
842     }
843     pixPtr->SetCharge(0);
844   } // for (Int_t i = 0; i < nPix;
845 }
846
847 //_____________________________________________________________________________
848 Bool_t AliMUONClusterFinderAZ::MainLoop(Int_t iSimple)
849 {
850 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
851   
852   TH2D *mlem;
853
854   Int_t ix, iy;
855   //Int_t nn, xList[10], yList[10];
856   Int_t nPix = fPixArray->GetEntriesFast();
857   AliMUONPixel *pixPtr = 0;
858   Double_t *coef = 0, *probi = 0; 
859   AddVirtualPad(); // add virtual pads if necessary
860   Int_t npadTot = fnPads[0] + fnPads[1], npadOK = 0;
861   for (Int_t i = 0; i < npadTot; i++) if (fPadIJ[1][i] == 0) npadOK++;
862   if (fDraw) fDraw->ResetMuon();
863
864   while (1) {
865
866     mlem = (TH2D*) gROOT->FindObject("mlem");
867     if (mlem) mlem->Delete();
868     // Calculate coefficients
869     if (fDebug) cout << " nPix, npadTot, npadOK " << nPix << " " << npadTot << " " << npadOK << endl;
870
871     // Calculate coefficients and pixel visibilities
872     coef = new Double_t [npadTot*nPix];
873     probi = new Double_t [nPix];
874     for (Int_t ipix=0; ipix<nPix; ipix++) probi[ipix] = 0;
875     Int_t indx = 0, indx1 = 0, cath = 0;
876
877     for (Int_t j=0; j<npadTot; j++) {
878       indx = j*nPix;
879       if (fPadIJ[1][j] == 0) {
880         cath = fPadIJ[0][j];
881         ix = fPadIJ[2][j];
882         iy = fPadIJ[3][j];
883         fSegmentation[cath]->SetPad(ix, iy);
884         /*
885           fSegmentation[cath]->Neighbours(fInput->DetElemId(),ix,iy,&nn,xList,yList); 
886           if (nn != 4) {
887           cout << nn << ": ";
888           for (Int_t i=0; i<nn; i++) {cout << xList[i] << " " << yList[i] << ", ";}
889           cout << endl;
890           }
891         */
892       }
893
894       for (Int_t ipix=0; ipix<nPix; ipix++) {
895         indx1 = indx + ipix;
896         if (fPadIJ[1][j] < 0) { coef[indx1] = 0; continue; }
897         pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
898         fSegmentation[cath]->SetHit(pixPtr->Coord(0), pixPtr->Coord(1), fZpad);
899         coef[indx1] = fInput->Mathieson()->IntXY(fInput->DetElemId(),fInput->Segmentation2(cath));
900         probi[ipix] += coef[indx1];
901       } // for (Int_t ipix=0;
902     } // for (Int_t j=0;
903     for (Int_t ipix=0; ipix<nPix; ipix++) if (probi[ipix] < 0.01) pixPtr->SetCharge(0); // "invisible" pixel
904
905     // MLEM algorithm
906     Mlem(coef, probi, 15);
907
908     Double_t xylim[4] = {999, 999, 999, 999};
909     for (Int_t ipix=0; ipix<nPix; ipix++) {
910       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
911       //cout << ipix+1; pixPtr->Print();
912       for (Int_t i=0; i<4; i++) 
913         xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
914     }
915     for (Int_t i=0; i<4; i++) {
916       xylim[i] -= pixPtr->Size(i/2); if (fDebug) cout << (i%2 ? -1 : 1)*xylim[i] << " "; }
917     if (fDebug) cout << endl;
918
919     // Adjust histogram to approximately the same limits as for the pads
920     // (for good presentation)
921     if (fDraw) fDraw->AdjustHist(xylim, pixPtr);
922
923     Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
924     Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
925     
926     mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
927     for (Int_t ipix=0; ipix<nPix; ipix++) {
928       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
929       mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
930     }
931     if (fDraw) fDraw->DrawHist("c2", mlem);
932
933     // Check if the total charge of pixels is too low
934     Double_t qTot = 0;
935     for (Int_t i=0; i<nPix; i++) {
936       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
937       qTot += pixPtr->Charge();
938     }
939     if (qTot < 1.e-4 || npadOK < 3 && qTot < 7) {
940       delete [] coef; delete [] probi; coef = 0; probi = 0;
941       fPixArray->Delete(); 
942       for (Int_t i=0; i<npadTot; i++) if (fPadIJ[1][i] == 0) fPadIJ[1][i] = -1;
943       return kFALSE; 
944     }
945
946     // Plot data - expectation
947     /*
948     Double_t x, y, cont;
949     for (Int_t j=0; j<npadTot; j++) {
950       Double_t sum1 = 0;
951       for (Int_t i=0; i<nPix; i++) {
952         // Caculate expectation
953         pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
954         sum1 += pixPtr->Charge()*coef[j*nPix+i];
955       }
956       sum1 = TMath::Min (sum1,fgkSaturation);
957       x = fXyq[0][j];
958       y = fXyq[1][j];
959       cath = fPadIJ[0][j];
960       Int_t ihist = cath*2;
961       ix = fHist[ihist]->GetXaxis()->FindBin(x);
962       iy = fHist[ihist]->GetYaxis()->FindBin(y);
963       cont = fHist[ihist]->GetCellContent(ix,iy);
964       if (cont == 0 && fHist[ihist+1]) {
965         ihist += 1;
966         ix = fHist[ihist]->GetXaxis()->FindBin(x);
967         iy = fHist[ihist]->GetYaxis()->FindBin(y);
968       }
969       fHist[ihist]->SetBinContent(ix,iy,fXyq[2][j]-sum1);
970     }
971     ((TCanvas*)gROOT->FindObject("c1"))->cd(1);
972     //gPad->SetTheta(55);
973     //gPad->SetPhi(30);
974     //mlem->Draw("lego1");
975     gPad->Modified();
976     ((TCanvas*)gROOT->FindObject("c1"))->cd(2);
977     gPad->Modified();
978     */
979
980     if (iSimple) {
981       // Simple cluster - skip further passes thru EM-procedure
982       Simple();
983       delete [] coef; delete [] probi; coef = 0; probi = 0;
984       fPixArray->Delete(); 
985       return kTRUE;
986     }
987
988     // Calculate position of the center-of-gravity around the maximum pixel
989     Double_t xyCOG[2];
990     FindCOG(mlem, xyCOG);
991
992     if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 && pixPtr->Size(0) > pixPtr->Size(1)) break;
993     //if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.007 && pixPtr->Size(0) > pixPtr->Size(1)) break;
994     //if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) >= 0.07 || pixPtr->Size(0) < pixPtr->Size(1)) {
995     // Sort pixels according to the charge
996     fPixArray->Sort();
997     /*
998     for (Int_t i=0; i<nPix; i++) {
999       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1000       cout << i+1; pixPtr->Print();
1001     }
1002     */
1003     Double_t pixMin = 0.01*((AliMUONPixel*)fPixArray->UncheckedAt(0))->Charge();
1004     pixMin = TMath::Min (pixMin,50.);
1005
1006     // Decrease pixel size and shift pixels to make them centered at 
1007     // the maximum one
1008     indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
1009     Double_t width = 0, shift[2]={0};
1010     ix = 1;
1011     for (Int_t i=0; i<4; i++) xylim[i] = 999;
1012     Int_t nPix1 = nPix; nPix = 0;
1013     for (Int_t ipix=0; ipix<nPix1; ipix++) {
1014       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1015       if (nPix >= npadOK) { // too many pixels already
1016         fPixArray->RemoveAt(ipix); 
1017         delete pixPtr; 
1018         continue;
1019       }
1020       if (pixPtr->Charge() < pixMin) { // low charge
1021         fPixArray->RemoveAt(ipix); 
1022         delete pixPtr; 
1023         continue;
1024       }
1025       for (Int_t i=0; i<2; i++) {
1026         if (!i) {
1027           pixPtr->SetCharge(10);
1028           pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
1029           width = -pixPtr->Size(indx);
1030           pixPtr->Shift(indx, width);
1031           // Shift pixel position
1032           if (ix) {
1033             ix = 0;
1034             for (Int_t j=0; j<2; j++) {
1035               shift[j] = pixPtr->Coord(j) - xyCOG[j];
1036               shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
1037             }
1038             //cout << ipix << " " << i << " " << shift[0] << " " << shift[1] << endl;
1039           } // if (ix)
1040           pixPtr->Shift(0, -shift[0]);
1041           pixPtr->Shift(1, -shift[1]);
1042         } else {
1043           pixPtr = new AliMUONPixel(*pixPtr);
1044           pixPtr->Shift(indx, -2*width);
1045           fPixArray->Add((TObject*)pixPtr);
1046         } // else
1047         //pixPtr->Print();
1048         for (Int_t i=0; i<4; i++) 
1049           xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1050       } // for (Int_t i=0; i<2;
1051       nPix += 2;
1052     } // for (Int_t ipix=0;
1053
1054     fPixArray->Compress();
1055     nPix = fPixArray->GetEntriesFast();
1056
1057     // Remove excessive pixels
1058     if (nPix > npadOK) {
1059       for (Int_t ipix=npadOK; ipix<nPix; ipix++) { 
1060         pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1061         fPixArray->RemoveAt(ipix); 
1062         delete pixPtr;
1063       }
1064     } else {
1065       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(0);
1066       // add pixels if the maximum is at the limit of pixel area
1067       // start from Y-direction
1068       Int_t j = 0;
1069       for (Int_t i=3; i>-1; i--) {
1070         if (nPix < npadOK && 
1071             TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2)) {
1072           pixPtr = new AliMUONPixel(*pixPtr);
1073           pixPtr->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1074           j = TMath::Even (i/2);
1075           pixPtr->SetCoord(j, xyCOG[j]);
1076           fPixArray->Add((TObject*)pixPtr);
1077           nPix++;
1078         }
1079       }
1080     } // else    
1081
1082     fPixArray->Compress();
1083     nPix = fPixArray->GetEntriesFast();
1084     delete [] coef; delete [] probi; coef = 0; probi = 0;
1085   } // while (1)
1086
1087   // remove pixels with low signal or low visibility
1088   // Cuts are empirical !!!
1089   Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1090   thresh = TMath::Min (thresh,50.);
1091   Double_t cmax = -1, charge = 0;
1092   for (Int_t i=0; i<nPix; i++) cmax = TMath::Max (cmax,probi[i]); 
1093   //cout << thresh << " " << cmax << " " << cmax*0.9 << endl;
1094   // Mark pixels which should be removed
1095   for (Int_t i=0; i<nPix; i++) {
1096     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1097     charge = pixPtr->Charge();
1098     if (charge < thresh) pixPtr->SetCharge(-charge);
1099     //else if (cmax > 1.91) {
1100     //  if (probi[i] < 1.9) pixPtr->SetCharge(-charge);
1101     //}
1102     //AZ else if (probi[i] < cmax*0.9) pixPtr->SetCharge(-charge);
1103     //18-01-06 else if (probi[i] < cmax*0.8) pixPtr->SetCharge(-charge);
1104     //cout << i << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << " " << charge << " " << probi[i] << endl;
1105   }
1106   // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1107   Int_t near = 0;
1108   for (Int_t i=0; i<nPix; i++) {
1109     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1110     charge = pixPtr->Charge();
1111     if (charge > 0) continue;
1112     near = FindNearest(pixPtr);
1113     pixPtr->SetCharge(0);
1114     probi[i] = 0; // make it "invisible"
1115     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(near);
1116     pixPtr->SetCharge(pixPtr->Charge() + (-charge));
1117   }
1118   Mlem(coef,probi,2);
1119   // Update histogram
1120   for (Int_t i=0; i<nPix; i++) {
1121     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1122     ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1123     iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1124     mlem->SetBinContent(ix, iy, pixPtr->Charge());
1125   }
1126   if (fDraw) fDraw->DrawHist("c2", mlem);
1127
1128   // Try to split into clusters
1129   Bool_t ok = kTRUE;
1130   if (mlem->GetSum() < 1) ok = kFALSE;
1131   else Split(mlem, coef);
1132   delete [] coef; delete [] probi; coef = 0; probi = 0;
1133   fPixArray->Delete(); 
1134   return ok;
1135 }
1136
1137 //_____________________________________________________________________________
1138 void AliMUONClusterFinderAZ::Mlem(Double_t *coef, Double_t *probi, Int_t nIter)
1139 {
1140 /// Use MLEM to find pixel charges
1141   
1142   Int_t nPix = fPixArray->GetEntriesFast();
1143   Int_t npad = fnPads[0] + fnPads[1];
1144   Double_t *probi1 = new Double_t [nPix];
1145   Double_t probMax = 0;
1146   Int_t indx, indx1;
1147   AliMUONPixel *pixPtr;
1148
1149   for (Int_t ipix=0; ipix<nPix; ipix++) if (probi[ipix] > probMax) probMax = probi[ipix];
1150   for (Int_t iter=0; iter<nIter; iter++) {
1151     // Do iterations
1152     for (Int_t ipix=0; ipix<nPix; ipix++) {
1153       // Correct each pixel
1154       if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1155       Double_t sum = 0;
1156       //probi1[ipix] = probi[ipix];
1157       probi1[ipix] = probMax;
1158       for (Int_t j=0; j<npad; j++) {
1159         if (fPadIJ[1][j] < 0) continue; 
1160         Double_t sum1 = 0;
1161         indx1 = j*nPix;
1162         indx = indx1 + ipix;
1163         for (Int_t i=0; i<nPix; i++) {
1164           // Caculate expectation
1165           pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1166           sum1 += pixPtr->Charge()*coef[indx1+i];
1167         } // for (Int_t i=0;
1168         if (fXyq[2][j] > fgkSaturation-1 && sum1 > fXyq[2][j]) { probi1[ipix] -= coef[indx]; continue; } // correct for pad charge overflows
1169         //cout << sum1 << " " << fXyq[2][j] << " " << coef[j*nPix+ipix] << endl;
1170         if (coef[indx] > 1.e-6) sum += fXyq[2][j]*coef[indx]/sum1;
1171       } // for (Int_t j=0;
1172       pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1173       if (probi1[ipix] > 1.e-6) pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1174     } // for (Int_t ipix=0;
1175   } // for (Int_t iter=0;
1176   delete [] probi1;
1177   return;
1178 }
1179
1180 //_____________________________________________________________________________
1181 void AliMUONClusterFinderAZ::FindCOG(TH2D *mlem, Double_t *xyc)
1182 {
1183 /// Calculate position of the center-of-gravity around the maximum pixel
1184
1185   Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1186   Int_t i1 = -9, j1 = -9;
1187   mlem->GetMaximumBin(ixmax,iymax,ix);
1188   Int_t nx = mlem->GetNbinsX();
1189   Int_t ny = mlem->GetNbinsY();
1190   Double_t thresh = mlem->GetMaximum()/10;
1191   Double_t x, y, cont, xq=0, yq=0, qq=0;
1192
1193   for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1194     y = mlem->GetYaxis()->GetBinCenter(i);
1195     for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1196       cont = mlem->GetCellContent(j,i);
1197       if (cont < thresh) continue;
1198       if (i != i1) {i1 = i; nsumy++;}
1199       if (j != j1) {j1 = j; nsumx++;}
1200       x = mlem->GetXaxis()->GetBinCenter(j);
1201       xq += x*cont;
1202       yq += y*cont;
1203       qq += cont;
1204       nsum++;
1205     }
1206   }
1207
1208   Double_t cmax = 0;
1209   Int_t i2 = 0, j2 = 0;
1210   x = y = 0;
1211   if (nsumy == 1) {
1212     // one bin in Y - add one more (with the largest signal)
1213     for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1214       if (i == iymax) continue;
1215       for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1216         cont = mlem->GetCellContent(j,i);
1217         if (cont > cmax) {
1218           cmax = cont;
1219           x = mlem->GetXaxis()->GetBinCenter(j);
1220           y = mlem->GetYaxis()->GetBinCenter(i);
1221           i2 = i;
1222           j2 = j;
1223         }
1224       }
1225     }
1226     xq += x*cmax;
1227     yq += y*cmax;
1228     qq += cmax;
1229     if (i2 != i1) nsumy++;
1230     if (j2 != j1) nsumx++;
1231     nsum++;
1232   } // if (nsumy == 1)
1233
1234   if (nsumx == 1) {
1235     // one bin in X - add one more (with the largest signal)
1236     cmax = x = y = 0;
1237     for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1238       if (j == ixmax) continue;
1239       for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1240         cont = mlem->GetCellContent(j,i);
1241         if (cont > cmax) {
1242           cmax = cont;
1243           x = mlem->GetXaxis()->GetBinCenter(j);
1244           y = mlem->GetYaxis()->GetBinCenter(i);
1245           i2 = i;
1246           j2 = j;
1247         }
1248       }
1249     }
1250     xq += x*cmax;
1251     yq += y*cmax;
1252     qq += cmax;
1253     if (i2 != i1) nsumy++;
1254     if (j2 != j1) nsumx++;
1255     nsum++;
1256   } // if (nsumx == 1)
1257
1258   xyc[0] = xq/qq; xyc[1] = yq/qq;
1259   if (fDebug) cout << xyc[0] << " " << xyc[1] << " " << qq << " " << nsum << " " << nsumx << " " << nsumy << endl;
1260   return;
1261 }
1262
1263 //_____________________________________________________________________________
1264 Int_t AliMUONClusterFinderAZ::FindNearest(AliMUONPixel *pixPtr0)
1265 {
1266 /// Find the pixel nearest to the given one
1267 /// (algorithm may be not very efficient)
1268
1269   Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1270   Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1271   Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1272   AliMUONPixel *pixPtr;
1273
1274   for (Int_t i=0; i<nPix; i++) {
1275     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1276     if (pixPtr->Charge() < 0.5) continue;
1277     dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1278     dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1279     r = dx *dx + dy * dy;
1280     if (r < rmin) { rmin = r; imin = i; }
1281   }
1282   return imin;
1283 }
1284
1285 //_____________________________________________________________________________
1286 void AliMUONClusterFinderAZ::Split(TH2D *mlem, Double_t *coef)
1287 {
1288 /// The main steering function to work with clusters of pixels in anode
1289 /// plane (find clusters, decouple them from each other, merge them (if
1290 /// necessary), pick up coupled pads, call the fitting function)
1291   
1292   Int_t nx = mlem->GetNbinsX();
1293   Int_t ny = mlem->GetNbinsY();
1294   Int_t nPix = fPixArray->GetEntriesFast();
1295
1296   Bool_t *used = new Bool_t[ny*nx];
1297   Double_t cont;
1298   Int_t nclust = 0, indx, indx1;
1299
1300   for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE; 
1301
1302   TObjArray *clusters[200]={0};
1303   TObjArray *pix;
1304
1305   // Find clusters of histogram bins (easier to work in 2-D space)
1306   for (Int_t i=1; i<=ny; i++) {
1307     for (Int_t j=1; j<=nx; j++) {
1308       indx = (i-1)*nx + j - 1;
1309       if (used[indx]) continue;
1310       cont = mlem->GetCellContent(j,i);
1311       if (cont < 0.5) continue;
1312       pix = new TObjArray(20);
1313       used[indx] = 1;
1314       pix->Add(BinToPix(mlem,j,i));
1315       AddBin(mlem, i, j, 0, used, pix); // recursive call
1316       if (nclust >= 200) AliFatal(" Too many clusters !!!");
1317       clusters[nclust++] = pix;
1318     } // for (Int_t j=1; j<=nx; j++) {
1319   } // for (Int_t i=1; i<=ny;
1320   if (fDebug) cout << nclust << endl;
1321   delete [] used; used = 0;
1322   
1323   // Compute couplings between clusters and clusters to pads
1324   Int_t npad = fnPads[0] + fnPads[1];
1325
1326   // Write out some information for algorithm development
1327   Int_t cath=0, npadx[2]={0}, npady[2]={0};
1328   Double_t xlow[2]={9999,9999}, xhig[2]={-9999,-9999};
1329   Double_t ylow[2]={9999,9999}, yhig[2]={-9999,-9999};
1330   for (Int_t j=0; j<npad; j++) {
1331     if (fXyq[3][j] < 0) continue; // exclude virtual pads
1332     cath = fPadIJ[0][j];
1333     if (fXyq[0][j] < xlow[cath]-0.001) { 
1334       if (fXyq[0][j]+fXyq[3][j] <= xlow[cath] && npadx[cath]) npadx[cath]++;
1335       xlow[cath] = fXyq[0][j];
1336     }
1337     if (fXyq[0][j] > xhig[cath]+0.001) { 
1338       if (fXyq[0][j]-fXyq[3][j] >= xhig[cath]) npadx[cath]++; 
1339       xhig[cath] = fXyq[0][j]; 
1340     }
1341     if (fXyq[1][j] < ylow[cath]-0.001) { 
1342       if (fXyq[1][j]+fXyq[4][j] <= ylow[cath] && npady[cath]) npady[cath]++;
1343       ylow[cath] = fXyq[1][j];
1344     }
1345     if (fXyq[1][j] > yhig[cath]+0.001) { 
1346       if (fXyq[1][j]-fXyq[4][j] >= yhig[cath]) npady[cath]++;
1347       yhig[cath] = fXyq[1][j]; 
1348     }
1349   }
1350   //if (lun1) fprintf(lun1," %4d %2d %3d %3d %3d %3d \n",gAlice->GetHeader()->GetEvent(),AliMUONClusterInput::Instance()->Chamber(), npadx[0], npadx[1], npady[0], npady[1]);
1351
1352   // Exclude pads with overflows
1353   for (Int_t j=0; j<npad; j++) {
1354     if (fXyq[2][j] > fgkSaturation-1) fPadIJ[1][j] = -5;
1355     else fPadIJ[1][j] = 0;
1356   }
1357
1358   // Compute couplings of clusters to pads
1359   TMatrixD *aijclupad = new TMatrixD(nclust,npad);
1360   *aijclupad = 0;
1361   Int_t npxclu;
1362   for (Int_t iclust=0; iclust<nclust; iclust++) {
1363     pix = clusters[iclust];
1364     npxclu = pix->GetEntriesFast();
1365     for (Int_t i=0; i<npxclu; i++) {
1366       indx = fPixArray->IndexOf(pix->UncheckedAt(i));
1367       for (Int_t j=0; j<npad; j++) {
1368         if (fPadIJ[1][j] < 0 && fPadIJ[1][j] != -5) continue;
1369         if (coef[j*nPix+indx] < fgkCouplMin) continue;
1370         (*aijclupad)(iclust,j) += coef[j*nPix+indx];
1371       }
1372     }
1373   }
1374   // Compute couplings between clusters
1375   TMatrixD *aijcluclu = new TMatrixD(nclust,nclust);
1376   *aijcluclu = 0;
1377   for (Int_t iclust=0; iclust<nclust; iclust++) {
1378     for (Int_t j=0; j<npad; j++) {
1379       // Exclude overflows
1380       if (fPadIJ[1][j] < 0) continue;
1381       if ((*aijclupad)(iclust,j) < fgkCouplMin) continue;
1382       for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) {
1383         if ((*aijclupad)(iclust1,j) < fgkCouplMin) continue;
1384         (*aijcluclu)(iclust,iclust1) += 
1385           TMath::Sqrt ((*aijclupad)(iclust,j)*(*aijclupad)(iclust1,j));
1386       }
1387     }
1388   }
1389   for (Int_t iclust=0; iclust<nclust; iclust++) {
1390     for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) {
1391       (*aijcluclu)(iclust1,iclust) = (*aijcluclu)(iclust,iclust1);
1392     }
1393   }
1394
1395   if (fDebug && nclust > 1) aijcluclu->Print();
1396
1397   // Find groups of coupled clusters
1398   used = new Bool_t[nclust];
1399   for (Int_t i=0; i<nclust; i++) used[i] = kFALSE;
1400   Int_t *clustNumb = new Int_t[nclust];
1401   Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
1402   Double_t parOk[8];
1403
1404   for (Int_t igroup=0; igroup<nclust; igroup++) {
1405     if (used[igroup]) continue;
1406     used[igroup] = kTRUE;
1407     clustNumb[0] = igroup;
1408     nCoupled = 1;
1409     // Find group of coupled clusters
1410     AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
1411     if (fDebug) {
1412       cout << " nCoupled: " << nCoupled << endl;
1413       for (Int_t i=0; i<nCoupled; i++) cout << clustNumb[i] << " "; cout << endl; 
1414     }
1415     fnCoupled = nCoupled;
1416
1417     while (nCoupled > 0) {
1418
1419       if (nCoupled < 4) {
1420         nForFit = nCoupled;
1421         for (Int_t i=0; i<nCoupled; i++) clustFit[i] = clustNumb[i];
1422       } else {
1423         // Too many coupled clusters to fit - try to decouple them
1424         // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with 
1425         // all the others in the group 
1426         for (Int_t j=0; j<3; j++) minGroup[j] = -1;
1427         Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
1428
1429         // Flag clusters for fit
1430         nForFit = 0;
1431         while (minGroup[nForFit] >= 0 && nForFit < 3) {
1432           if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
1433           clustFit[nForFit] = clustNumb[minGroup[nForFit]];
1434           clustNumb[minGroup[nForFit]] -= 999;
1435           nForFit++;
1436         }
1437         if (fDebug) cout << nForFit << " " << coupl << endl;
1438       } // else
1439
1440       // Select pads for fit. 
1441       if (SelectPad(nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1) {
1442         // Deselect pads
1443         for (Int_t j=0; j<npad; j++) {
1444           if (TMath::Abs(fPadIJ[1][j]) == 1) fPadIJ[1][j] = 0;
1445           if (TMath::Abs(fPadIJ[1][j]) == -9) fPadIJ[1][j] = -5;
1446         }
1447         // Merge the failed cluster candidates (with too few pads to fit) with 
1448         // the one with the strongest coupling
1449         Merge(nForFit, nCoupled, clustNumb, clustFit, clusters, aijcluclu, aijclupad);
1450       } else {
1451         // Do the fit
1452         nfit = Fit(0, nForFit, clustFit, clusters, parOk);
1453       }
1454
1455       // Subtract the fitted charges from pads with strong coupling and/or
1456       // return pads for further use
1457       UpdatePads(nfit, parOk);
1458
1459       // Mark used pads
1460       for (Int_t j=0; j<npad; j++) {
1461         if (fPadIJ[1][j] == 1) fPadIJ[1][j] = -1;
1462         if (fPadIJ[1][j] == -9) fPadIJ[1][j] = -5;
1463       }
1464
1465       // Sort the clusters (move to the right the used ones)
1466       Int_t beg = 0, end = nCoupled - 1;
1467       while (beg < end) {
1468         if (clustNumb[beg] >= 0) { beg++; continue; }
1469         for (Int_t j=end; j>beg; j--) {
1470           if (clustNumb[j] < 0) continue;
1471           end = j - 1;
1472           indx = clustNumb[beg];
1473           clustNumb[beg] = clustNumb[j];
1474           clustNumb[j] = indx;
1475           break;
1476         }
1477         beg++;
1478       }
1479
1480       nCoupled -= nForFit;
1481       if (nCoupled > 3) {
1482         // Remove couplings of used clusters
1483         for (Int_t iclust=nCoupled; iclust<nCoupled+nForFit; iclust++) {
1484           indx = clustNumb[iclust] + 999;
1485           for (Int_t iclust1=0; iclust1<nCoupled; iclust1++) {
1486             indx1 = clustNumb[iclust1];
1487             (*aijcluclu)(indx,indx1) = (*aijcluclu)(indx1,indx) = 0;
1488           }
1489         }
1490
1491         // Update the remaining clusters couplings (exclude couplings from 
1492         // the used pads)
1493         for (Int_t j=0; j<npad; j++) {
1494           if (fPadIJ[1][j] != -1) continue;
1495           for (Int_t iclust=0; iclust<nCoupled; iclust++) {
1496             indx = clustNumb[iclust];
1497             if ((*aijclupad)(indx,j) < fgkCouplMin) continue;
1498             for (Int_t iclust1=iclust+1; iclust1<nCoupled; iclust1++) {
1499               indx1 = clustNumb[iclust1];
1500               if ((*aijclupad)(indx1,j) < fgkCouplMin) continue;
1501               // Check this
1502               (*aijcluclu)(indx,indx1) -= 
1503                 TMath::Sqrt ((*aijclupad)(indx,j)*(*aijclupad)(indx1,j));
1504               (*aijcluclu)(indx1,indx) = (*aijcluclu)(indx,indx1);
1505             }
1506           }
1507           fPadIJ[1][j] = -8;
1508         } // for (Int_t j=0; j<npad;
1509       } // if (nCoupled > 3)
1510     } // while (nCoupled > 0)
1511   } // for (Int_t igroup=0; igroup<nclust;
1512
1513   aijcluclu->Delete(); aijclupad->Delete();
1514   for (Int_t iclust=0; iclust<nclust; iclust++) {
1515     pix = clusters[iclust]; 
1516     pix->Clear();
1517     delete pix; pix = 0;
1518   }
1519   delete [] clustNumb; clustNumb = 0; delete [] used; used = 0;
1520 }
1521
1522 //_____________________________________________________________________________
1523 void AliMUONClusterFinderAZ::AddBin(TH2D *mlem, Int_t ic, Int_t jc, Int_t mode, Bool_t *used, TObjArray *pix)
1524 {
1525 /// Add a bin to the cluster
1526
1527   Int_t nx = mlem->GetNbinsX();
1528   Int_t ny = mlem->GetNbinsY();
1529   Double_t cont1, cont = mlem->GetCellContent(jc,ic);
1530   AliMUONPixel *pixPtr = 0;
1531
1532   for (Int_t i=TMath::Max(ic-1,1); i<=TMath::Min(ic+1,ny); i++) {
1533     for (Int_t j=TMath::Max(jc-1,1); j<=TMath::Min(jc+1,nx); j++) {
1534       if (i != ic && j != jc) continue;
1535       if (used[(i-1)*nx+j-1]) continue;
1536       cont1 = mlem->GetCellContent(j,i);
1537       if (mode && cont1 > cont) continue;
1538       used[(i-1)*nx+j-1] = kTRUE;
1539       if (cont1 < 0.5) continue;
1540       if (pix) pix->Add(BinToPix(mlem,j,i)); 
1541       else {
1542         pixPtr = new AliMUONPixel (mlem->GetXaxis()->GetBinCenter(j), 
1543                                    mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1544         fPixArray->Add((TObject*)pixPtr);
1545       }
1546       AddBin(mlem, i, j, mode, used, pix); // recursive call
1547     }
1548   }
1549 }
1550
1551 //_____________________________________________________________________________
1552 TObject* AliMUONClusterFinderAZ::BinToPix(TH2D *mlem, Int_t jc, Int_t ic)
1553 {
1554 /// Translate histogram bin to pixel 
1555   
1556   Double_t yc = mlem->GetYaxis()->GetBinCenter(ic);
1557   Double_t xc = mlem->GetXaxis()->GetBinCenter(jc);
1558   
1559   Int_t nPix = fPixArray->GetEntriesFast();
1560   AliMUONPixel *pixPtr = NULL;
1561
1562   // Compare pixel and bin positions
1563   for (Int_t i=0; i<nPix; i++) {
1564     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1565     if (pixPtr->Charge() < 0.5) continue;
1566     if (TMath::Abs(pixPtr->Coord(0)-xc)<1.e-4 && TMath::Abs(pixPtr->Coord(1)-yc)<1.e-4) return (TObject*) pixPtr;
1567   }
1568   AliError(Form(" Something wrong ??? %f %f ", xc, yc));
1569   return NULL;
1570 }
1571
1572 //_____________________________________________________________________________
1573 void AliMUONClusterFinderAZ::AddCluster(Int_t ic, Int_t nclust, TMatrixD *aijcluclu, Bool_t *used, Int_t *clustNumb, Int_t &nCoupled)
1574 {
1575 /// Add a cluster to the group of coupled clusters
1576
1577   for (Int_t i=0; i<nclust; i++) {
1578     if (used[i]) continue;
1579     if ((*aijcluclu)(i,ic) < fgkCouplMin) continue;
1580     used[i] = kTRUE;
1581     clustNumb[nCoupled++] = i;
1582     AddCluster(i, nclust, aijcluclu, used, clustNumb, nCoupled);
1583   }
1584 }
1585
1586 //_____________________________________________________________________________
1587 Double_t AliMUONClusterFinderAZ::MinGroupCoupl(Int_t nCoupled, Int_t *clustNumb, TMatrixD *aijcluclu, Int_t *minGroup)
1588 {
1589 /// Find group of clusters with minimum coupling to all the others
1590
1591   Int_t i123max = TMath::Min(3,nCoupled/2); 
1592   Int_t indx, indx1, indx2, indx3, nTot = 0;
1593   Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1594
1595   for (Int_t i123=1; i123<=i123max; i123++) {
1596
1597     if (i123 == 1) {
1598       coupl1 = new Double_t [nCoupled];
1599       for (Int_t i=0; i<nCoupled; i++) coupl1[i] = 0;
1600     }
1601     else if (i123 == 2) {
1602       nTot = nCoupled*nCoupled;
1603       coupl2 = new Double_t [nTot];
1604       for (Int_t i=0; i<nTot; i++) coupl2[i] = 9999;
1605     } else {
1606       nTot = nTot*nCoupled;
1607       coupl3 = new Double_t [nTot];
1608       for (Int_t i=0; i<nTot; i++) coupl3[i] = 9999;
1609     } // else
1610
1611     for (Int_t i=0; i<nCoupled; i++) {
1612       indx1 = clustNumb[i];
1613       for (Int_t j=i+1; j<nCoupled; j++) {
1614         indx2 = clustNumb[j];
1615         if (i123 == 1) {
1616           coupl1[i] += (*aijcluclu)(indx1,indx2);
1617           coupl1[j] += (*aijcluclu)(indx1,indx2);
1618         } 
1619         else if (i123 == 2) {
1620           indx = i*nCoupled + j;
1621           coupl2[indx] = coupl1[i] + coupl1[j];
1622           coupl2[indx] -= 2 * ((*aijcluclu)(indx1,indx2));
1623         } else {
1624           for (Int_t k=j+1; k<nCoupled; k++) {
1625             indx3 = clustNumb[k];
1626             indx = i*nCoupled*nCoupled + j*nCoupled + k;
1627             coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1628             coupl3[indx] -= 2 * ((*aijcluclu)(indx1,indx3)+(*aijcluclu)(indx2,indx3));
1629           }
1630         } // else
1631       } // for (Int_t j=i+1;
1632     } // for (Int_t i=0;
1633   } // for (Int_t i123=1;
1634
1635   // Find minimum coupling
1636   Double_t couplMin = 9999;
1637   Int_t locMin = 0;
1638
1639   for (Int_t i123=1; i123<=i123max; i123++) {
1640     if (i123 == 1) {
1641       locMin = TMath::LocMin(nCoupled, coupl1);
1642       couplMin = coupl1[locMin];
1643       minGroup[0] = locMin;
1644       delete [] coupl1; coupl1 = 0;
1645     } 
1646     else if (i123 == 2) {
1647       locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
1648       if (coupl2[locMin] < couplMin) {
1649         couplMin = coupl2[locMin];
1650         minGroup[0] = locMin/nCoupled;
1651         minGroup[1] = locMin%nCoupled;
1652       }
1653       delete [] coupl2; coupl2 = 0;
1654     } else {
1655       locMin = TMath::LocMin(nTot, coupl3);
1656       if (coupl3[locMin] < couplMin) {
1657         couplMin = coupl3[locMin];
1658         minGroup[0] = locMin/nCoupled/nCoupled;
1659         minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
1660         minGroup[2] = locMin%nCoupled;
1661       }
1662       delete [] coupl3; coupl3 = 0;
1663     } // else
1664   } // for (Int_t i123=1;
1665   return couplMin;
1666 }
1667
1668 //_____________________________________________________________________________
1669 Int_t AliMUONClusterFinderAZ::SelectPad(Int_t nCoupled, Int_t nForFit, Int_t *clustNumb, Int_t *clustFit, TMatrixD *aijclupad)
1670 {
1671 /// Select pads for fit. If too many coupled clusters, find pads giving 
1672 /// the strongest coupling with the rest of clusters and exclude them from the fit.
1673
1674   Int_t npad = fnPads[0] + fnPads[1];
1675   Double_t *padpix = 0;
1676
1677   if (nCoupled > 3) {
1678     padpix = new Double_t[npad];
1679     for (Int_t i=0; i<npad; i++) padpix[i] = 0; 
1680   }
1681
1682   Int_t nOK = 0, indx, indx1;
1683   for (Int_t iclust=0; iclust<nForFit; iclust++) {
1684     indx = clustFit[iclust];
1685     for (Int_t j=0; j<npad; j++) {
1686       if ((*aijclupad)(indx,j) < fgkCouplMin) continue;
1687       if (fPadIJ[1][j] == -5) fPadIJ[1][j] = -9; // flag overflow
1688       if (fPadIJ[1][j] < 0) continue; // exclude overflows and used pads
1689       if (!fPadIJ[1][j]) { fPadIJ[1][j] = 1; nOK++; } // pad to be used in fit
1690       if (nCoupled > 3) {
1691         // Check other clusters
1692         for (Int_t iclust1=0; iclust1<nCoupled; iclust1++) {
1693           indx1 = clustNumb[iclust1];
1694           if (indx1 < 0) continue;
1695           if ((*aijclupad)(indx1,j) < fgkCouplMin) continue;
1696           padpix[j] += (*aijclupad)(indx1,j);
1697         }
1698       } // if (nCoupled > 3)
1699     } // for (Int_t j=0; j<npad;
1700   } // for (Int_t iclust=0; iclust<nForFit
1701   if (nCoupled < 4) return nOK;
1702
1703   Double_t aaa = 0;
1704   for (Int_t j=0; j<npad; j++) {
1705     if (padpix[j] < fgkCouplMin) continue;
1706     if (fDebug) cout << j << " " << padpix[j] << " " << fXyq[0][j] << " " << fXyq[1][j] << endl;
1707     aaa += padpix[j];
1708     fPadIJ[1][j] = -1; // exclude pads with strong coupling to the other clusters
1709     nOK--;
1710   }
1711   delete [] padpix; padpix = 0;
1712   return nOK;
1713 }
1714   
1715 //_____________________________________________________________________________
1716 void AliMUONClusterFinderAZ::Merge(Int_t nForFit, Int_t nCoupled, Int_t *clustNumb, Int_t *clustFit, TObjArray **clusters, TMatrixD *aijcluclu, TMatrixD *aijclupad)
1717 {
1718 /// Merge the group of clusters with the one having the strongest coupling with them
1719
1720   Int_t indx, indx1, npxclu, npxclu1, imax=0;
1721   TObjArray *pix, *pix1;
1722   Double_t couplMax;
1723
1724   for (Int_t icl=0; icl<nForFit; icl++) {
1725     indx = clustFit[icl];
1726     pix = clusters[indx];
1727     npxclu = pix->GetEntriesFast();
1728     couplMax = -1;
1729     for (Int_t icl1=0; icl1<nCoupled; icl1++) {
1730       indx1 = clustNumb[icl1];
1731       if (indx1 < 0) continue;
1732       if ((*aijcluclu)(indx,indx1) > couplMax) {
1733         couplMax = (*aijcluclu)(indx,indx1);
1734         imax = indx1;
1735       }
1736     } // for (Int_t icl1=0;
1737     /*if (couplMax < fgkCouplMin) {
1738       cout << " Oops " << couplMax << endl;
1739       aijcluclu->Print();
1740       cout << icl << " " << indx << " " << npxclu << " " << nLinks << endl;
1741       ::exit(0);
1742       }*/
1743     // Add to it
1744     pix1 = clusters[imax];
1745     npxclu1 = pix1->GetEntriesFast();
1746     // Add pixels 
1747     for (Int_t i=0; i<npxclu; i++) { pix1->Add(pix->UncheckedAt(i)); pix->RemoveAt(i); }
1748     if (fDebug) cout << " New number of pixels: " << npxclu1 << " " << pix1->GetEntriesFast() << endl;
1749     //Add cluster-to-cluster couplings
1750     //aijcluclu->Print();
1751     for (Int_t icl1=0; icl1<nCoupled; icl1++) {
1752       indx1 = clustNumb[icl1];
1753       if (indx1 < 0 || indx1 == imax) continue;
1754       (*aijcluclu)(indx1,imax) += (*aijcluclu)(indx,indx1);
1755       (*aijcluclu)(imax,indx1) = (*aijcluclu)(indx1,imax);
1756     }
1757     (*aijcluclu)(indx,imax) = (*aijcluclu)(imax,indx) = 0;
1758     //aijcluclu->Print();
1759     //Add cluster-to-pad couplings
1760     for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
1761       if (fPadIJ[1][j] < 0 && fPadIJ[1][j] != -5) continue; // exclude used pads
1762       (*aijclupad)(imax,j) += (*aijclupad)(indx,j);
1763       (*aijclupad)(indx,j) = 0;
1764     }
1765   } // for (Int_t icl=0; icl<nForFit;
1766 }
1767
1768 //_____________________________________________________________________________
1769 Int_t AliMUONClusterFinderAZ::Fit(Int_t iSimple, Int_t nfit, Int_t *clustFit, TObjArray **clusters, Double_t *parOk)
1770 {
1771 /// Find selected clusters to selected pad charges
1772   
1773   TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
1774   Double_t xmin = mlem->GetXaxis()->GetXmin() - mlem->GetXaxis()->GetBinWidth(1);
1775   Double_t xmax = mlem->GetXaxis()->GetXmax() + mlem->GetXaxis()->GetBinWidth(1);
1776   Double_t ymin = mlem->GetYaxis()->GetXmin() - mlem->GetYaxis()->GetBinWidth(1);
1777   Double_t ymax = mlem->GetYaxis()->GetXmax() + mlem->GetYaxis()->GetBinWidth(1);
1778   Double_t step[3]={0.01,0.002,0.02}, xPad = 0, yPad = 99999;
1779
1780   // Number of pads to use and number of virtual pads
1781   Int_t npads = 0, nVirtual = 0, nfit0 = nfit;
1782   for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
1783     if (fXyq[3][i] < 0) nVirtual++;
1784     if (fPadIJ[1][i] != 1) continue;
1785     if (fXyq[3][i] > 0) {
1786       npads++;
1787       if (yPad > 9999) { 
1788         xPad = fXyq[0][i]; 
1789         yPad = fXyq[1][i]; 
1790       } else {
1791         if (fXyq[4][i] < fXyq[3][i]) yPad = fXyq[1][i]; 
1792         else xPad = fXyq[0][i]; 
1793       }
1794     }
1795   }
1796   if (fDebug) {
1797     for (Int_t i=0; i<nfit; i++) {cout << i+1 << " " << clustFit[i] << " ";}
1798     cout << nfit << endl;
1799     cout << " Number of pads to fit: " << npads << endl;
1800   }
1801   fNpar = 0;
1802   fQtot = 0;
1803   if (npads < 2) return 0; 
1804   
1805   Int_t digit = 0;
1806   AliMUONDigit *mdig = 0;
1807   Int_t tracks[3] = {-1, -1, -1};
1808   for (Int_t cath=0; cath<2; cath++) {  
1809     for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
1810       if (fPadIJ[0][i] != cath) continue;
1811       if (fPadIJ[1][i] != 1) continue;
1812       if (fXyq[3][i] < 0) continue; // exclude virtual pads
1813       digit = TMath::Nint (fXyq[5][i]);
1814       if (digit >= 0) mdig = fInput->Digit(cath,digit);
1815       else mdig = fInput->Digit(TMath::Even(cath),-digit-1);
1816       //if (!mdig) mdig = fInput->Digit(TMath::Even(cath),digit);
1817       if (!mdig) continue; // protection for cluster display
1818       if (mdig->Hit() >= 0) {
1819         if (tracks[0] < 0) {
1820           tracks[0] = mdig->Hit();
1821           tracks[1] = mdig->Track(0);
1822         } else if (mdig->Track(0) < tracks[1]) {
1823           tracks[0] = mdig->Hit();
1824           tracks[1] = mdig->Track(0);
1825         }
1826       }
1827       if (mdig->Track(1) >= 0 && mdig->Track(1) != tracks[1]) {
1828         if (tracks[2] < 0) tracks[2] = mdig->Track(1);
1829         else tracks[2] = TMath::Min (tracks[2], mdig->Track(1));
1830       }
1831       //if (!mdig) break;
1832       //cout << mdig->Hit() << " " << mdig->Track(0) << " " << mdig->Track(1) <<endl;
1833     } // for (Int_t i=0;
1834   } // for (Int_t cath=0;
1835   //cout << tracks[0] << " " << tracks[1] << " " << tracks[2] <<endl;
1836   
1837   // Get number of pads in X and Y 
1838   Int_t nInX = 0, nInY;
1839   PadsInXandY(nInX, nInY);
1840   //cout << " nInX and Y: " << nInX << " " << nInY << endl;
1841
1842   Int_t nfitMax = 3; 
1843   nfitMax = TMath::Min (nfitMax, (npads + 1) / 3);
1844   if (nfitMax > 1) {
1845     if (nInX < 3 && nInY < 3 || nInX == 3 && nInY < 3 || nInX < 3 && nInY == 3) nfitMax = 1; // not enough pads in each direction
1846   }
1847   if (nfit > nfitMax) nfit = nfitMax;
1848
1849   // Take cluster maxima as fitting seeds
1850   TObjArray *pix;
1851   AliMUONPixel *pixPtr;
1852   Int_t npxclu;
1853   Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8], qq = 0;
1854   Double_t xyseed[3][2], qseed[3], xyCand[3][2] = {{0},{0}}, sigCand[3][2] = {{0},{0}};
1855
1856   for (Int_t ifit=1; ifit<=nfit0; ifit++) {
1857     cmax = 0;
1858     pix = clusters[clustFit[ifit-1]];
1859     npxclu = pix->GetEntriesFast();
1860     //qq = 0;
1861     for (Int_t clu=0; clu<npxclu; clu++) {
1862       pixPtr = (AliMUONPixel*) pix->UncheckedAt(clu);
1863       cont = pixPtr->Charge();
1864       fQtot += cont;
1865       if (cont > cmax) { 
1866         cmax = cont; 
1867         xseed = pixPtr->Coord(0);
1868         yseed = pixPtr->Coord(1);
1869       }
1870       qq += cont;
1871       /*
1872       xyCand[ifit-1][0] += pixPtr->Coord(0) * cont;
1873       xyCand[ifit-1][1] += pixPtr->Coord(1) * cont;
1874       sigCand[ifit-1][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
1875       sigCand[ifit-1][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
1876       */
1877       xyCand[0][0] += pixPtr->Coord(0) * cont;
1878       xyCand[0][1] += pixPtr->Coord(1) * cont;
1879       sigCand[0][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
1880       sigCand[0][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
1881     }
1882     xyseed[ifit-1][0] = xseed;
1883     xyseed[ifit-1][1] = yseed;
1884     qseed[ifit-1] = cmax;
1885     /*
1886     xyCand[ifit-1][0] /= qq; // <x>
1887     xyCand[ifit-1][1] /= qq; // <y>
1888     sigCand[ifit-1][0] = sigCand[ifit-1][0]/qq - xyCand[ifit-1][0]*xyCand[ifit-1][0]; // <x^2> - <x>^2
1889     sigCand[ifit-1][0] = sigCand[ifit-1][0] > 0 ? TMath::Sqrt (sigCand[ifit-1][0]) : 0;
1890     sigCand[ifit-1][1] = sigCand[ifit-1][1]/qq - xyCand[ifit-1][1]*xyCand[ifit-1][1]; // <y^2> - <y>^2
1891     sigCand[ifit-1][1] = sigCand[ifit-1][1] > 0 ? TMath::Sqrt (sigCand[ifit-1][1]) : 0;
1892     cout << xyCand[ifit-1][0] << " " << xyCand[ifit-1][1] << " " << sigCand[ifit-1][0] << " " << sigCand[ifit-1][1] << endl;
1893     */
1894   } // for (Int_t ifit=1;
1895
1896   xyCand[0][0] /= qq; // <x>
1897   xyCand[0][1] /= qq; // <y>
1898   sigCand[0][0] = sigCand[0][0]/qq - xyCand[0][0]*xyCand[0][0]; // <x^2> - <x>^2
1899   sigCand[0][0] = sigCand[0][0] > 0 ? TMath::Sqrt (sigCand[0][0]) : 0;
1900   sigCand[0][1] = sigCand[0][1]/qq - xyCand[0][1]*xyCand[0][1]; // <y^2> - <y>^2
1901   sigCand[0][1] = sigCand[0][1] > 0 ? TMath::Sqrt (sigCand[0][1]) : 0;
1902   if (fDebug) cout << xyCand[0][0] << " " << xyCand[0][1] << " " << sigCand[0][0] << " " << sigCand[0][1] << endl;
1903
1904   Int_t nDof, maxSeed[3], nMax = 0;
1905   Double_t fmin, chi2o = 9999, chi2n;
1906
1907   TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
1908   /*
1909   Int_t itmp[100], localMax[100];
1910   Double_t maxVal[100];
1911   if (!iSimple && nfit < nfitMax) {
1912     // Try to split pixel cluster according to local maxima
1913     Int_t nfit1 = nfit;
1914     for (Int_t iclus = 0; iclus < nfit1; iclus++) {
1915       nMax = FindLocalMaxima (clusters[clustFit[maxSeed[iclus]]], localMax, maxVal);
1916       TH2D *hist = (TH2D*) gROOT->FindObject("anode1");
1917       if (nMax == 1) { hist->Delete(); continue; }
1918       // Add extra fitting seeds from local maxima
1919       Int_t ixseed = hist->GetXaxis()->FindBin(xyseed[maxSeed[iclus]][0]);
1920       Int_t iyseed = hist->GetYaxis()->FindBin(xyseed[maxSeed[iclus]][1]);
1921       Int_t nx = hist->GetNbinsX();
1922       TMath::Sort(nMax, maxVal, itmp, kTRUE); // in decreasing order
1923       for (Int_t j = 0; j < nMax; j++) {
1924         Int_t iyc = localMax[itmp[j]] / nx + 1;
1925         Int_t ixc = localMax[itmp[j]] % nx + 1;
1926         if (ixc == ixseed && iyc == iyseed) continue; // local max already taken for seeding
1927         xyseed[nfit][0] = hist->GetXaxis()->GetBinCenter(ixc);
1928         xyseed[nfit][1] = hist->GetYaxis()->GetBinCenter(iyc);
1929         qseed[nfit] = maxVal[itmp[j]];
1930         maxSeed[nfit] = nfit++;
1931         if (nfit >= nfitMax) break;
1932       }
1933       hist->Delete();
1934       if (nfit >= nfitMax) break;
1935     } // for (Int_t iclus = 0;
1936     //nfit0 = nfit;
1937     //TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
1938   } //if (!iSimple && nfit < nfitMax)
1939   */
1940
1941   Double_t *gin = 0, func0, func1, param[8], step0[8];
1942   Double_t param0[2][8]={{0},{0}}, deriv[2][8]={{0},{0}}; 
1943   Double_t shift[8], stepMax, derMax, parmin[8], parmax[8], func2[2], shift0;
1944   Double_t delta[8], scMax, dder[8], estim, shiftSave = 0;
1945   Int_t min, max, nCall = 0, memory[8] = {0}, nLoop, idMax = 0, iestMax = 0, nFail;
1946   Double_t rad, dist[3] = {0};
1947
1948   // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is 
1949   // lower, try 3-track (if number of pads is sufficient).
1950   for (Int_t iseed=0; iseed<nfit; iseed++) {
1951
1952     if (iseed) { for (Int_t j=0; j<fNpar; j++) param[j] = parOk[j]; } // for bounded params
1953     for (Int_t j=0; j<3; j++) step0[fNpar+j] = shift[fNpar+j] = step[j];
1954     if (nfit == 1) param[fNpar] = xyCand[0][0]; // take COG
1955     else param[fNpar] = xyseed[maxSeed[iseed]][0];
1956     parmin[fNpar] = xmin; 
1957     parmax[fNpar++] = xmax; 
1958     if (nfit == 1) param[fNpar] = xyCand[0][1]; // take COG
1959     else param[fNpar] = xyseed[maxSeed[iseed]][1];
1960     parmin[fNpar] = ymin; 
1961     parmax[fNpar++] = ymax; 
1962     if (fNpar > 2) {
1963       param[fNpar] = fNpar == 4 ? 0.5 : 0.3;
1964       parmin[fNpar] = 0; 
1965       parmax[fNpar++] = 1; 
1966     }
1967     if (iseed) { for (Int_t j=0; j<fNpar; j++) param0[1][j] = 0; }
1968
1969     // Try new algorithm
1970     min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
1971
1972     while (1) {
1973       max = !min;
1974       Fcn1(fNpar, gin, func0, param, 1); nCall++;
1975       //cout << " Func: " << func0 << endl;
1976
1977       func2[max] = func0;
1978       for (Int_t j=0; j<fNpar; j++) {
1979         param0[max][j] = param[j];
1980         delta[j] = step0[j];
1981         param[j] += delta[j] / 10;
1982         if (j > 0) param[j-1] -= delta[j-1] / 10;
1983         Fcn1(fNpar, gin, func1, param, 1); nCall++;
1984         deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
1985         //cout << j << " " << deriv[max][j] << endl;
1986         dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) / 
1987           (param0[0][j] - param0[1][j]) : 0; // second derivative
1988       }
1989       param[fNpar-1] -= delta[fNpar-1] / 10;
1990       if (nCall > 2000) break;
1991
1992       min = func2[0] < func2[1] ? 0 : 1;
1993       nFail = min == max ? 0 : nFail + 1;
1994
1995       stepMax = derMax = estim = 0;
1996       for (Int_t j=0; j<fNpar; j++) { 
1997         // Estimated distance to minimum
1998         shift0 = shift[j];
1999         if (nLoop == 1) shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
2000         else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3) shift[j] = 0;
2001         else if (deriv[min][j]*deriv[!min][j] > 0 && TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])
2002                  //|| TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3) {
2003           || TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3 || TMath::Abs(dder[j]) < 1.e-6) {
2004           shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
2005           if (min == max) { 
2006             if (memory[j] > 1) { shift[j] *= 2; } //cout << " Memory " << memory[j] << " " << shift[j] << endl; }
2007             memory[j]++;
2008           }
2009         } else {
2010           shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
2011           memory[j] = 0;
2012         }
2013         if (TMath::Abs(shift[j])/step0[j] > estim) { 
2014           estim = TMath::Abs(shift[j])/step0[j];
2015           iestMax = j;
2016         }
2017
2018         // Too big step
2019         if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; // 
2020
2021         // Failed to improve minimum
2022         if (min != max) {
2023           memory[j] = 0;
2024           param[j] = param0[min][j];
2025           if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j]) shift[j] = (shift[j] + shift0) / 2;
2026           else shift[j] /= -2;
2027         } 
2028
2029         // Too big step
2030         if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min]) 
2031           shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
2032
2033         // Introduce step relaxation factor
2034         if (memory[j] < 3) {
2035           scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
2036           if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax) 
2037             shift[j] = TMath::Sign (shift0*scMax, shift[j]);
2038         }
2039         param[j] += shift[j]; 
2040         //AZ Check parameter limits 27-12-2004
2041         if (param[j] < parmin[j]) { 
2042           shift[j] = parmin[j] - param[j]; 
2043           param[j] = parmin[j]; 
2044         } else if (param[j] > parmax[j]) {
2045           shift[j] = parmax[j] - param[j];
2046           param[j] = parmax[j];
2047         }
2048         //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
2049         stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
2050         if (TMath::Abs(deriv[min][j]) > derMax) {
2051           idMax = j;
2052           derMax = TMath::Abs (deriv[min][j]);
2053         }
2054       } // for (Int_t j=0; j<fNpar;
2055       //cout << max << " " << func2[min] << " " << derMax << " " << stepMax << " " << estim << " " << iestMax << " " << nCall << endl;
2056       if (estim < 1 && derMax < 2 || nLoop > 150) break; // minimum was found
2057
2058       nLoop++;
2059       // Check for small step
2060       if (shift[idMax] == 0) { shift[idMax] = step0[idMax]/10; param[idMax] += shift[idMax]; continue; }
2061       if (!memory[idMax] && derMax > 0.5 && nLoop > 10) {
2062         //cout << " ok " << deriv[min][idMax] << " " << deriv[!min][idMax] << " " << dder[idMax]*shift[idMax] << " " << shift[idMax] << endl;
2063         if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10) {
2064           if (min == max) dder[idMax] = -dder[idMax];
2065           shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10; 
2066           param[idMax] += shift[idMax];
2067           stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
2068           //cout << shift[idMax] << " " << param[idMax] << endl;
2069           if (min == max) shiftSave = shift[idMax];
2070         }
2071         if (nFail > 10) {
2072           param[idMax] -= shift[idMax];
2073           shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
2074           param[idMax] += shift[idMax];
2075           //cout << shift[idMax] << endl;
2076         }
2077       }      
2078     } // while (1)
2079     fmin = func2[min];
2080
2081     nDof = npads - fNpar + nVirtual;
2082     if (!nDof) nDof++;
2083     chi2n = fmin / nDof;
2084     if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
2085
2086     if (chi2n*1.2+1.e-6 > chi2o ) { fNpar -= 3; break; }
2087
2088     // Save parameters and errors
2089
2090     if (nInX == 1) {
2091       // One pad per direction 
2092       for (Int_t i=0; i<fNpar; i++) if (i == 0 || i == 2 || i == 5) param0[min][i] = xPad;
2093     }
2094     if (nInY == 1) {
2095       // One pad per direction 
2096       for (Int_t i=0; i<fNpar; i++) if (i == 1 || i == 3 || i == 6) param0[min][i] = yPad;
2097     }
2098
2099     /*
2100     if (iseed > 0) {
2101       // Find distance to the nearest neighbour
2102       dist[0] = dist[1] = TMath::Sqrt ((param0[min][0]-param0[min][2])*
2103                                        (param0[min][0]-param0[min][2])
2104                                       +(param0[min][1]-param0[min][3])*
2105                                        (param0[min][1]-param0[min][3]));
2106       if (iseed > 1) {
2107         dist[2] = TMath::Sqrt ((param0[min][0]-param0[min][5])*
2108                                (param0[min][0]-param0[min][5])
2109                               +(param0[min][1]-param0[min][6])*
2110                                (param0[min][1]-param0[min][6]));
2111         rad = TMath::Sqrt ((param0[min][2]-param0[min][5])*
2112                            (param0[min][2]-param0[min][5])
2113                           +(param0[min][3]-param0[min][6])*
2114                            (param0[min][3]-param0[min][6]));
2115         if (dist[2] < dist[0]) dist[0] = dist[2];
2116         if (rad < dist[1]) dist[1] = rad;
2117         if (rad < dist[2]) dist[2] = rad;
2118       }
2119       cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
2120       if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
2121     }
2122     */
2123
2124     for (Int_t i=0; i<fNpar; i++) {
2125       parOk[i] = param0[min][i];
2126       //errOk[i] = fmin;
2127       errOk[i] = chi2n;
2128       // Bounded params
2129       parOk[i] = TMath::Max (parOk[i], parmin[i]);
2130       parOk[i] = TMath::Min (parOk[i], parmax[i]);
2131     }
2132
2133     chi2o = chi2n;
2134     if (fmin < 0.1) break; // !!!???
2135   } // for (Int_t iseed=0; 
2136
2137   if (fDebug) {
2138     for (Int_t i=0; i<fNpar; i++) {
2139       if (i == 4 || i == 7) {
2140         if (i == 7 || i == 4 && fNpar < 7) cout << parOk[i] << endl;
2141         else cout << parOk[i] * (1-parOk[7]) << endl;
2142         continue;
2143       }
2144       cout << parOk[i] << " " << errOk[i] << endl;
2145     }
2146   }
2147   nfit = (fNpar + 1) / 3;
2148   dist[0] = dist[1] = dist[2] = 0;
2149
2150   if (nfit > 1) {
2151     // Find distance to the nearest neighbour
2152     dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
2153                                      (parOk[0]-parOk[2])
2154                                     +(parOk[1]-parOk[3])*
2155                                      (parOk[1]-parOk[3]));
2156     if (nfit > 2) {
2157       dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
2158                              (parOk[0]-parOk[5])
2159                             +(parOk[1]-parOk[6])*
2160                              (parOk[1]-parOk[6]));
2161       rad = TMath::Sqrt ((parOk[2]-parOk[5])*
2162                          (parOk[2]-parOk[5])
2163                         +(parOk[3]-parOk[6])*
2164                          (parOk[3]-parOk[6]));
2165       if (dist[2] < dist[0]) dist[0] = dist[2];
2166       if (rad < dist[1]) dist[1] = rad;
2167       if (rad < dist[2]) dist[2] = rad;
2168     }
2169   }
2170
2171   Int_t indx;
2172   fnPads[1] -= nVirtual;
2173   if (!fDraw) {
2174     Double_t coef = 0;
2175     if (iSimple) fnCoupled = 0;
2176     //for (Int_t j=0; j<nfit; j++) {
2177     for (Int_t j=nfit-1; j>=0; j--) {
2178       indx = j<2 ? j*2 : j*2+1;  
2179       if (nfit == 1) coef = 1;
2180       else coef = j==nfit-1 ? parOk[indx+2] : 1-coef;
2181       coef = TMath::Max (coef, 0.);
2182       if (nfit == 3 && j < 2) coef = j==1 ? coef*parOk[indx+2] : coef - parOk[7];
2183       coef = TMath::Max (coef, 0.);
2184       AddRawCluster (parOk[indx], parOk[indx+1], coef*fQtot, errOk[indx], nfit0+10*nfit+100*nMax+10000*fnCoupled, tracks,
2185                      //sigCand[maxSeed[j]][0], sigCand[maxSeed[j]][1]);
2186                      //sigCand[0][0], sigCand[0][1], dist[j]);
2187                      sigCand[0][0], sigCand[0][1], dist[TMath::LocMin(nfit,dist)]);
2188     }
2189   } else fDraw->FillMuon(nfit, parOk, errOk);
2190   return nfit;
2191 }  
2192
2193 //_____________________________________________________________________________
2194 void AliMUONClusterFinderAZ::Fcn1(Int_t & /*npar*/, Double_t * /*gin*/, Double_t &f, Double_t *par, Int_t /*iflag*/)
2195 {
2196 /// Fit for one track
2197 /// AZ for Muinuit AliMUONClusterFinderAZ& c = *(AliMUONClusterFinderAZ::fgClusterFinder);    
2198
2199   AliMUONClusterFinderAZ& c = *this; //AZ
2200   
2201   Int_t cath, ix, iy, indx, npads=0;
2202   Double_t charge, delta, coef=0, chi2=0, qTot = 0;
2203   for (Int_t j=0; j<c.fnPads[0]+c.fnPads[1]; j++) {
2204     if (c.fPadIJ[1][j] != 1) continue;
2205     cath = c.fPadIJ[0][j];
2206     if (c.fXyq[3][j] > 0) npads++; // exclude virtual pads
2207     qTot += c.fXyq[2][j];
2208     ix = c.fPadIJ[2][j];
2209     iy = c.fPadIJ[3][j];
2210     c.fSegmentation[cath]->SetPad(ix, iy);
2211     charge = 0;
2212     for (Int_t i=c.fNpar/3; i>=0; i--) { // sum over tracks
2213       indx = i<2 ? 2*i : 2*i+1;
2214       c.fSegmentation[cath]->SetHit(par[indx], par[indx+1], c.fZpad);
2215       if (c.fNpar == 2) coef = 1;
2216       else coef = i==c.fNpar/3 ? par[indx+2] : 1-coef;
2217       coef = TMath::Max (coef, 0.);
2218       if (c.fNpar == 8 && i < 2) coef = i==1 ? coef*par[indx+2] : coef - par[7];
2219       coef = TMath::Max (coef, 0.);
2220       charge += c.fInput->Mathieson()->IntXY(fInput->DetElemId(), c.fInput->Segmentation2(cath))*coef;
2221     }
2222     charge *= c.fQtot;
2223     delta = charge - c.fXyq[2][j];
2224     delta *= delta;
2225     delta /= c.fXyq[2][j];
2226     //if (cath) delta /= 5; // just for test
2227     chi2 += delta;
2228   } // for (Int_t j=0;
2229   f = chi2; 
2230   Double_t qAver = qTot/npads; //(c.fnPads[0]+c.fnPads[1]);
2231   f = chi2/qAver;
2232 }
2233
2234 //_____________________________________________________________________________
2235 void AliMUONClusterFinderAZ::UpdatePads(Int_t /*nfit*/, Double_t *par)
2236 {
2237 /// Subtract the fitted charges from pads with strong coupling
2238
2239   Int_t cath, ix, iy, indx;
2240   Double_t charge, coef=0;
2241   for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2242     if (fPadIJ[1][j] != -1) continue;
2243     if (fNpar != 0) {
2244       cath = fPadIJ[0][j];
2245       ix = fPadIJ[2][j];
2246       iy = fPadIJ[3][j];
2247       fSegmentation[cath]->SetPad(ix, iy);
2248       charge = 0;
2249       for (Int_t i=fNpar/3; i>=0; i--) { // sum over tracks
2250         indx = i<2 ? 2*i : 2*i+1;
2251         fSegmentation[cath]->SetHit(par[indx], par[indx+1], fZpad);
2252         if (fNpar == 2) coef = 1;
2253         else coef = i==fNpar/3 ? par[indx+2] : 1-coef;
2254         coef = TMath::Max (coef, 0.);
2255         if (fNpar == 8 && i < 2) coef = i==1 ? coef*par[indx+2] : coef - par[7];
2256         coef = TMath::Max (coef, 0.);
2257         charge += fInput->Mathieson()->IntXY(fInput->DetElemId(),fInput->Segmentation2(cath))*coef;
2258       }
2259       charge *= fQtot;
2260       fXyq[2][j] -= charge;
2261     } // if (fNpar != 0)
2262     if (fXyq[2][j] > fgkZeroSuppression) fPadIJ[1][j] = 0; // return pad for further using
2263   } // for (Int_t j=0;
2264 }  
2265
2266 //_____________________________________________________________________________
2267 Bool_t AliMUONClusterFinderAZ::TestTrack(Int_t /*t*/) const 
2268 {
2269 /// Test if track was user selected
2270
2271   return kTRUE;
2272   /*
2273     if (fTrack[0]==-1 || fTrack[1]==-1) {
2274         return kTRUE;
2275     } else if (t==fTrack[0] || t==fTrack[1]) {
2276         return kTRUE;
2277     } else {
2278         return kFALSE;
2279     }
2280   */
2281 }
2282
2283 //_____________________________________________________________________________
2284 void AliMUONClusterFinderAZ::AddRawCluster(Double_t x, Double_t y, Double_t qTot, Double_t fmin, Int_t nfit, Int_t *tracks, Double_t /*sigx*/, Double_t /*sigy*/, Double_t /*dist*/)
2285 {
2286 /// Add a raw cluster copy to the list
2287
2288   if (qTot <= 0.501) return; 
2289   AliMUONRawCluster cnew;
2290
2291   Int_t cath, npads[2] = {0}, nover[2] = {0};
2292   for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2293     cath = fPadIJ[0][j];
2294     // There was an overflow
2295     if (fPadIJ[1][j] == -9) nover[cath]++;
2296     if (fPadIJ[1][j] != 1 && fPadIJ[1][j] != -9) continue;
2297     cnew.SetMultiplicity(cath,cnew.GetMultiplicity(cath)+1);
2298     if (fXyq[2][j] > cnew.GetPeakSignal(cath)) cnew.SetPeakSignal(cath,fXyq[2][j]);
2299     //cnew.SetCharge(cath,cnew.GetCharge(cath) + TMath::Nint (fXyq[2][j]));
2300     cnew.SetContrib(npads[cath],cath,fXyq[2][j]);
2301     cnew.SetIndex(npads[cath],cath,TMath::Nint (fXyq[5][j]));
2302     cnew.SetDetElemId(fInput->DetElemId());
2303     npads[cath]++;
2304   }
2305
2306   cnew.SetClusterType(nover[0] + nover[1] * 100);
2307   for (Int_t j=0; j<3; j++) cnew.SetTrack(j,tracks[j]);
2308
2309   Double_t xg, yg, zg;
2310   for (cath=0; cath<2; cath++) {
2311     // Perform local-to-global transformation
2312     fInput->Segmentation2(cath)->GetTransformer()->Local2Global(fInput->DetElemId(), x, y, fZpad, xg, yg, zg);
2313     cnew.SetX(cath, xg);
2314     cnew.SetY(cath, yg);
2315     cnew.SetZ(cath, zg);
2316     cnew.SetCharge(cath, TMath::Nint(qTot));
2317     //cnew.SetPeakSignal(cath,20);
2318     //cnew.SetMultiplicity(cath, 5);
2319     cnew.SetNcluster(cath, nfit);
2320     cnew.SetChi2(cath, fmin); //0.;1
2321   } 
2322   // Evaluate measurement errors
2323   //AZ Errors(&cnew);
2324
2325   cnew.SetGhost(nfit); //cnew.SetX(1,sigx); cnew.SetY(1,sigy); cnew.SetZ(1,dist);
2326   //cnew.fClusterType=cnew.PhysicsContribution();
2327   new((*fRawClusters)[fNRawClusters++]) AliMUONRawCluster(cnew); 
2328   if (fDebug) cout << fNRawClusters << " " << fInput->Chamber() << endl;
2329   //fNPeaks++;
2330 }
2331
2332 //_____________________________________________________________________________
2333 Int_t AliMUONClusterFinderAZ::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
2334 {
2335 /// Find local maxima in pixel space for large preclusters in order to
2336 /// try to split them into smaller pieces (to speed up the MLEM procedure)
2337 /// or to find additional fitting seeds if clusters were not completely resolved  
2338
2339   TH2D *hist = NULL;
2340   //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
2341   //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
2342   //if (hist) hist->Delete();
2343
2344   Double_t xylim[4] = {999, 999, 999, 999};
2345   Int_t nPix = pixArray->GetEntriesFast();
2346   AliMUONPixel *pixPtr = 0;
2347   for (Int_t ipix=0; ipix<nPix; ipix++) {
2348     pixPtr = (AliMUONPixel*) pixArray->UncheckedAt(ipix);
2349     for (Int_t i=0; i<4; i++) 
2350          xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
2351   }
2352   for (Int_t i=0; i<4; i++) xylim[i] -= pixPtr->Size(i/2); 
2353
2354   Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
2355   Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
2356   if (pixArray == fPixArray) hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
2357   else hist = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
2358   for (Int_t ipix=0; ipix<nPix; ipix++) {
2359     pixPtr = (AliMUONPixel*) pixArray->UncheckedAt(ipix);
2360     hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
2361   }
2362   if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
2363
2364   Int_t nMax = 0, indx;
2365   Int_t *isLocalMax = new Int_t[ny*nx];
2366   for (Int_t i=0; i<ny*nx; i++) isLocalMax[i] = 0;
2367
2368   for (Int_t i=1; i<=ny; i++) {
2369     indx = (i-1) * nx;
2370     for (Int_t j=1; j<=nx; j++) {
2371       if (hist->GetCellContent(j,i) < 0.5) continue;
2372       //if (isLocalMax[indx+j-1] < 0) continue;
2373       if (isLocalMax[indx+j-1] != 0) continue;
2374       FlagLocalMax(hist, i, j, isLocalMax);
2375     }
2376   }
2377
2378   for (Int_t i=1; i<=ny; i++) {
2379     indx = (i-1) * nx;
2380     for (Int_t j=1; j<=nx; j++) {
2381       if (isLocalMax[indx+j-1] > 0) { 
2382         localMax[nMax] = indx + j - 1; 
2383         maxVal[nMax++] = hist->GetCellContent(j,i);
2384         if (nMax > 99) AliFatal(" Too many local maxima !!!");
2385       }
2386     }
2387   }
2388   if (fDebug) cout << " Local max: " << nMax << endl;
2389   delete [] isLocalMax; isLocalMax = 0;
2390   return nMax;
2391 }
2392
2393 //_____________________________________________________________________________
2394 void AliMUONClusterFinderAZ::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
2395 {
2396 /// Flag pixels (whether or not local maxima)
2397
2398   Int_t nx = hist->GetNbinsX();
2399   Int_t ny = hist->GetNbinsY();
2400   Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
2401   Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
2402
2403   for (Int_t i1=i-1; i1<i+2; i1++) {
2404     if (i1 < 1 || i1 > ny) continue;
2405     indx1 = (i1 - 1) * nx;
2406     for (Int_t j1=j-1; j1<j+2; j1++) {
2407       if (j1 < 1 || j1 > nx) continue;
2408       if (i == i1 && j == j1) continue;
2409       indx2 = indx1 + j1 - 1;
2410       cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
2411       if (cont < cont1) { isLocalMax[indx] = -1; return; }
2412       else if (cont > cont1) isLocalMax[indx2] = -1;
2413       else { // the same charge
2414         isLocalMax[indx] = 1; 
2415         if (isLocalMax[indx2] == 0) {
2416           FlagLocalMax(hist, i1, j1, isLocalMax);
2417           if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
2418           else isLocalMax[indx2] = -1;
2419         }
2420       } 
2421     }
2422   }
2423   isLocalMax[indx] = 1; // local maximum
2424 }
2425
2426 //_____________________________________________________________________________
2427 void AliMUONClusterFinderAZ::FindCluster(Int_t *localMax, Int_t iMax)
2428 {
2429 /// Find pixel cluster around local maximum \a iMax and pick up pads
2430 /// overlapping with it
2431
2432   TH2D *hist = (TH2D*) gROOT->FindObject("anode");
2433   Int_t nx = hist->GetNbinsX();
2434   Int_t ny = hist->GetNbinsY();
2435   Int_t ic = localMax[iMax] / nx + 1;
2436   Int_t jc = localMax[iMax] % nx + 1;
2437   Bool_t *used = new Bool_t[ny*nx];
2438   for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE;
2439
2440   // Drop all pixels from the array - pick up only the ones from the cluster
2441   fPixArray->Delete();
2442
2443   Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2; 
2444   Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;  
2445   Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
2446   Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
2447   Double_t cont = hist->GetCellContent(jc,ic);
2448   AliMUONPixel *pixPtr = new AliMUONPixel (xc, yc, wx, wy, cont);
2449   fPixArray->Add((TObject*)pixPtr);
2450   used[(ic-1)*nx+jc-1] = kTRUE;
2451   AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
2452
2453   Int_t nPix = fPixArray->GetEntriesFast(), npad = fnPads[0] + fnPads[1];
2454   for (Int_t i=0; i<nPix; i++) {
2455     ((AliMUONPixel*)fPixArray->UncheckedAt(i))->SetSize(0,wx); 
2456     ((AliMUONPixel*)fPixArray->UncheckedAt(i))->SetSize(1,wy); 
2457   }
2458   if (fDebug) cout << iMax << " " << nPix << endl;
2459
2460   Float_t xy[4], xy12[4];
2461   // Pick up pads which overlap with found pixels
2462   for (Int_t i=0; i<npad; i++) fPadIJ[1][i] = -1;
2463   for (Int_t i=0; i<nPix; i++) {
2464     pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
2465     for (Int_t j=0; j<4; j++) 
2466       xy[j] = pixPtr->Coord(j/2) + (j%2 ? 1 : -1)*pixPtr->Size(j/2);
2467     for (Int_t j=0; j<npad; j++) 
2468       if (Overlap(xy, j, xy12, 0)) fPadIJ[1][j] = 0; // flag for use
2469   }
2470
2471   delete [] used; used = 0;
2472 }
2473
2474 //_____________________________________________________________________________
2475 void AliMUONClusterFinderAZ::AddVirtualPad()
2476 {
2477 /// Add virtual pad (with small charge) to improve fit for some
2478 /// clusters (when pad with max charge is at the extreme of the cluster)
2479
2480   // Get number of pads in X and Y-directions
2481   Int_t nInX = -1, nInY;
2482   PadsInXandY(nInX, nInY);
2483   //return;
2484
2485   // Add virtual pad only if number of pads per direction == 2
2486   if (nInX != 2 && nInY != 2) return;
2487
2488   // Find pads with max charge
2489   Int_t maxpad[2][2] = {{-1, -1}, {-1, -1}}, cath;
2490   Double_t sigmax[2] = {0}, aamax[2] = {0};
2491   for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2492     if (fPadIJ[1][j] != 0) continue;
2493     cath = fPadIJ[0][j];
2494     if (fXyq[2][j] > sigmax[cath]) {
2495       maxpad[cath][1] = maxpad[cath][0];
2496       aamax[cath] = sigmax[cath];
2497       sigmax[cath] = fXyq[2][j];
2498       maxpad[cath][0] = j;
2499     }
2500   }
2501   if (maxpad[0][0] >= 0 && maxpad[0][1] < 0 || maxpad[1][0] >= 0 && maxpad[1][1] < 0) {
2502     for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2503       if (fPadIJ[1][j] != 0) continue;
2504       cath = fPadIJ[0][j];
2505       if (j == maxpad[cath][0] || j == maxpad[cath][1]) continue;
2506       if (fXyq[2][j] > aamax[cath]) {
2507         aamax[cath] = fXyq[2][j];
2508         maxpad[cath][1] = j;
2509       }
2510     }
2511   }
2512   // Check for mirrors (side X on cathode 0) 
2513   Bool_t mirror = kFALSE;
2514   if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0) {
2515     mirror = fXyq[3][maxpad[0][0]] < fXyq[4][maxpad[0][0]]; 
2516     if (!mirror && TMath::Abs(fXyq[3][maxpad[0][0]]-fXyq[3][maxpad[1][0]]) < 0.001) {
2517       // Special case when pads on both cathodes have the same size
2518       Int_t yud[2] = {0};
2519       for (Int_t j = 0; j < fnPads[0]+fnPads[1]; j++) {
2520         cath = fPadIJ[0][j];
2521         if (j == maxpad[cath][0]) continue;
2522         if (fPadIJ[2][j] != fPadIJ[2][maxpad[cath][0]]) continue;
2523         if (fPadIJ[3][j] + 1 == fPadIJ[3][maxpad[cath][0]] || 
2524             fPadIJ[3][j] - 1 == fPadIJ[3][maxpad[cath][0]]) yud[cath]++;
2525       }
2526       if (!yud[0]) mirror = kTRUE; // take the other cathode
2527     } // if (!mirror &&...
2528   } // if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
2529
2530   // Find neughbours of pads with max charges
2531   Int_t nn, xList[10], yList[10], ix0, iy0, ix, iy, neighb;
2532   for (cath=0; cath<2; cath++) {
2533     if (!cath && maxpad[0][0] < 0) continue; // one-sided cluster - cathode 1
2534     if (cath && maxpad[1][0] < 0) break; // one-sided cluster - cathode 0
2535     if (maxpad[1][0] >= 0) {
2536       if (!mirror) {
2537         if (!cath && nInY != 2) continue;
2538         if (cath && nInX != 2 && (maxpad[0][0] >= 0 || nInY != 2)) continue;
2539       } else {
2540         if (!cath && nInX != 2) continue;
2541         if (cath && nInY != 2 && (maxpad[0][0] >= 0 || nInX != 2)) continue;
2542       }
2543     }
2544
2545     Int_t iAddX = 0, iAddY = 0, ix1 = 0, iy1 = 0, iPad = 0;
2546     if (maxpad[0][0] < 0) iPad = 1;
2547
2548     for (iPad=0; iPad<2; iPad++) {
2549       if (maxpad[cath][iPad] < 0) continue;
2550       if (iPad && !iAddX && !iAddY) break;
2551       if (iPad && fXyq[2][maxpad[cath][1]] / sigmax[cath] < 0.5) break;
2552
2553       Int_t neighbx = 0, neighby = 0;
2554       ix0 = fPadIJ[2][maxpad[cath][iPad]];
2555       iy0 = fPadIJ[3][maxpad[cath][iPad]];
2556       fSegmentation[cath]->Neighbours(ix0, iy0, &nn, xList, yList); 
2557       Float_t zpad; 
2558       for (Int_t j=0; j<nn; j++) {
2559         if (TMath::Abs(xList[j]-ix0) == 1 || xList[j]*ix0 == -1) neighbx++;
2560         if (TMath::Abs(yList[j]-iy0) == 1 || yList[j]*iy0 == -1) neighby++;
2561       }
2562       if (!mirror) {
2563         if (cath) neighb = neighbx; 
2564         else neighb = neighby;
2565         if (maxpad[0][0] < 0) neighb += neighby;
2566         else if (maxpad[1][0] < 0) neighb += neighbx;
2567       } else {
2568         if (!cath) neighb = neighbx; 
2569         else neighb = neighby;
2570         if (maxpad[0][0] < 0) neighb += neighbx;
2571         else if (maxpad[1][0] < 0) neighb += neighby;
2572       }
2573
2574       for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2575         if (fPadIJ[0][j] != cath) continue;
2576         ix = fPadIJ[2][j];
2577         iy = fPadIJ[3][j];
2578         if (iy == iy0 && ix == ix0) continue; 
2579         for (Int_t k=0; k<nn; k++) {
2580           if (xList[k] != ix || yList[k] != iy) continue;
2581           if (!mirror) {
2582             if ((!cath || maxpad[0][0] < 0) && 
2583                 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
2584               if (!iPad && TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) ix1 = xList[k]; //19-12-05 
2585               xList[k] = yList[k] = 0; 
2586               neighb--;
2587               break;
2588             }
2589             if ((cath || maxpad[1][0] < 0) && 
2590                 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
2591               if (!iPad) ix1 = xList[k]; //19-12-05
2592               xList[k] = yList[k] = 0; 
2593               neighb--;
2594             }
2595           } else {
2596             if ((!cath || maxpad[0][0] < 0) && 
2597                 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
2598               if (!iPad) ix1 = xList[k]; //19-12-05
2599               xList[k] = yList[k] = 0; 
2600               neighb--;
2601               break;
2602             }
2603             if ((cath || maxpad[1][0] < 0) && 
2604                 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
2605               xList[k] = yList[k] = 0; 
2606               neighb--;
2607             }
2608           }
2609           break;
2610         } // for (Int_t k=0; k<nn;
2611         if (!neighb) break;
2612       } // for (Int_t j=0; j<fnPads[0]+fnPads[1];
2613       if (!neighb) continue;
2614       
2615       // Add virtual pad
2616       Int_t npads, isec;
2617       isec = 0;
2618       for (Int_t j=0; j<nn; j++) {
2619         if (xList[j] == 0 && yList[j] == 0) continue;
2620         npads = fnPads[0] + fnPads[1];
2621         fPadIJ[0][npads] = cath;
2622         fPadIJ[1][npads] = 0;
2623         ix = xList[j]; 
2624         iy = yList[j];
2625         if (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) {
2626           if (iy != iy0) continue; // new segmentation - check
2627           if (nInX != 2) continue; // new
2628           if (!mirror) {
2629             if (!cath && maxpad[1][0] >= 0) continue;
2630           } else {
2631             if (cath && maxpad[0][0] >= 0) continue;
2632           }
2633           if (iPad && !iAddX) continue;
2634           fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
2635           if (fXyq[0][npads] > 1.e+5) continue; // temporary fix
2636           if (ix == ix1) continue; //19-12-05
2637           if (ix1 == ix0) continue;
2638           if (maxpad[1][0] < 0 || mirror && maxpad[0][0] >= 0) {
2639             if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/100, 5.);
2640             else fXyq[2][npads] = TMath::Min (aamax[0]/100, 5.);
2641           }
2642           else {
2643             if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/100, 5.);
2644             else fXyq[2][npads] = TMath::Min (aamax[1]/100, 5.);
2645           }
2646           fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
2647           fXyq[3][npads] = -2; // flag
2648           fPadIJ[2][npads] = ix;
2649           fPadIJ[3][npads] = iy;
2650           fnPads[1]++;
2651           iAddX = npads;
2652           if (fDebug) printf(" ***** Add virtual pad in X ***** %f %f %f %3d %3d \n", fXyq[2][npads], 
2653                              fXyq[0][npads], fXyq[1][npads], ix, iy);
2654           ix1 = ix0;
2655           continue;
2656         } 
2657         if (nInY != 2) continue;
2658         if (!mirror && cath && maxpad[0][0] >= 0) continue;
2659         if (mirror && !cath && maxpad[1][0] >= 0) continue;
2660         if (TMath::Abs(iy-iy0) == 1 || TMath::Abs(iy*iy0) == 1) {
2661           if (ix != ix0) continue; // new segmentation - check
2662           if (iPad && !iAddY) continue;
2663           fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
2664           if (iy1 == iy0) continue;
2665           //if (iPad && iy1 == iy0) continue;
2666           if (maxpad[0][0] < 0 || mirror && maxpad[1][0] >= 0) {
2667             if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/15, fgkZeroSuppression);
2668             else fXyq[2][npads] = TMath::Min (aamax[1]/15, fgkZeroSuppression);
2669           }
2670           else {
2671             if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/15, fgkZeroSuppression);
2672             else fXyq[2][npads] = TMath::Min (aamax[0]/15, fgkZeroSuppression);
2673           }
2674           fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
2675           fXyq[3][npads] = -2; // flag
2676           fPadIJ[2][npads] = ix;
2677           fPadIJ[3][npads] = iy;
2678           fnPads[1]++;
2679           iAddY = npads;
2680           if (fDebug) printf(" ***** Add virtual pad in Y ***** %f %f %f %3d %3d \n", fXyq[2][npads], 
2681                              fXyq[0][npads], fXyq[1][npads], ix, iy);
2682           iy1 = iy0;
2683         }
2684       } // for (Int_t j=0; j<nn;
2685     } // for (Int_t iPad=0;
2686   } // for (cath=0; cath<2;
2687   return;
2688 }
2689
2690 //_____________________________________________________________________________
2691 void AliMUONClusterFinderAZ::PadsInXandY(Int_t &nInX, Int_t &nInY)
2692 {
2693 /// Find number of pads in X and Y-directions (excluding virtual ones and
2694 /// overflows)
2695
2696   static Int_t nXsaved = 0, nYsaved = 0;
2697   nXsaved = nYsaved = 0;
2698   //if (nInX >= 0) {nInX = nXsaved; nInY = nYsaved; return; }
2699   Float_t *xPad0 = NULL, *yPad0 = NULL, *xPad1 = NULL, *yPad1 = NULL;
2700   Float_t wMinX[2] = {99, 99}, wMinY[2] = {99, 99};
2701   Int_t *nPad0 = NULL, *nPad1 = NULL;
2702   Int_t nPads = fnPads[0] + fnPads[1];
2703   if (fnPads[0]) {
2704     xPad0 = new Float_t[nPads];
2705     yPad0 = new Float_t[nPads];
2706     nPad0 = new Int_t[nPads];
2707   }
2708   if (fnPads[1]) {
2709     xPad1 = new Float_t[nPads];
2710     yPad1 = new Float_t[nPads];
2711     nPad1 = new Int_t[nPads];
2712   }
2713   Int_t n0 = 0, n1 = 0, cath, npadx[2] = {1, 1}, npady[2] = {1, 1};
2714   for (Int_t j = 0; j < nPads; j++) {
2715     if (nInX < 0 && fPadIJ[1][j] != 0) continue; // before fit
2716     else if (nInX == 0 && fPadIJ[1][j] != 1) continue; // fit - exclude overflows
2717     else if (nInX > 0 && fPadIJ[1][j] != 1 && fPadIJ[1][j] != -9) continue; // exclude non-marked
2718     if (nInX <= 0 && fXyq[2][j] > fgkSaturation-1) continue; // skip overflows
2719     cath = fPadIJ[0][j];
2720     if (fXyq[3][j] > 0) { // exclude virtual pads
2721       wMinX[cath] = TMath::Min (wMinX[cath], fXyq[3][j]);
2722       wMinY[cath] = TMath::Min (wMinY[cath], fXyq[4][j]);
2723       //20-12-05 }
2724       if (cath) { xPad1[n1] = fXyq[0][j]; yPad1[n1++] = fXyq[1][j]; }
2725       else { xPad0[n0] = fXyq[0][j]; yPad0[n0++] = fXyq[1][j]; }
2726     }
2727   }
2728
2729   // Sort
2730   if (n0) { 
2731     TMath::Sort (n0, xPad0, nPad0); // in X
2732     for (Int_t i = 1; i < n0; i++) 
2733       if (xPad0[nPad0[i]] - xPad0[nPad0[i-1]] < -0.01) npadx[0]++;
2734     TMath::Sort (n0, yPad0, nPad0); // in Y
2735     for (Int_t i = 1; i < n0; i++) 
2736       if (yPad0[nPad0[i]] - yPad0[nPad0[i-1]] < -0.01) npady[0]++;
2737   }
2738   
2739   if (n1) { 
2740     TMath::Sort (n1, xPad1, nPad1); // in X
2741     for (Int_t i = 1; i < n1; i++) 
2742       if (xPad1[nPad1[i]] - xPad1[nPad1[i-1]] < -0.01) npadx[1]++;
2743     TMath::Sort (n1, yPad1, nPad1); // in Y
2744     for (Int_t i = 1; i < n1; i++) 
2745       if (yPad1[nPad1[i]] - yPad1[nPad1[i-1]] < -0.01) npady[1]++;
2746   }
2747   if (fnPads[0]) { delete [] xPad0; delete [] yPad0; delete [] nPad0; }
2748   if (fnPads[1]) { delete [] xPad1; delete [] yPad1; delete [] nPad1; }
2749   if (TMath::Abs (wMinY[0] - wMinY[1]) < 1.e-3) nInY = TMath::Max (npady[0], npady[1]);
2750   else nInY = wMinY[0] < wMinY[1] ? npady[0] : npady[1];
2751   if (TMath::Abs (wMinX[0] - wMinX[1]) < 1.e-3) nInX = TMath::Max (npadx[0], npadx[1]);
2752   else nInX = wMinX[0] < wMinX[1] ? npadx[0] : npadx[1];
2753 }
2754
2755 //_____________________________________________________________________________
2756 void AliMUONClusterFinderAZ::Simple()
2757 {
2758 /// Process simple cluster (small number of pads) without EM-procedure
2759
2760   Int_t nForFit = 1, clustFit[1] = {0}, nfit;
2761   Double_t parOk[3] = {0.}; 
2762   TObjArray *clusters[1]; 
2763   clusters[0] = fPixArray;
2764   for (Int_t i = 0; i < fnPads[0]+fnPads[1]; i++) {
2765     if (fXyq[2][i] > fgkSaturation-1) fPadIJ[1][i] = -9;
2766     else fPadIJ[1][i] = 1;
2767   }
2768   nfit = Fit(1, nForFit, clustFit, clusters, parOk);
2769 }
2770
2771 //_____________________________________________________________________________
2772 void AliMUONClusterFinderAZ::Errors(AliMUONRawCluster *clus)
2773 {
2774 /// Correct reconstructed coordinates for some clusters and evaluate errors
2775
2776   Double_t qTot = clus->GetCharge(0), fmin = clus->GetChi2(0);
2777   Double_t xreco = clus->GetX(0), yreco = clus->GetY(0), zreco = clus->GetZ(0);
2778   Double_t sigmax[2] = {0};
2779
2780   Int_t nInX = 1, nInY, maxdig[2] ={-1, -1}, digit, cath1, isec;
2781   PadsInXandY(nInX, nInY);
2782
2783   // Find pad with maximum signal
2784   for (Int_t cath = 0; cath < 2; cath++) {
2785     for (Int_t j = 0; j < clus->GetMultiplicity(cath); j++) {
2786       cath1 = cath;
2787       digit = clus->GetIndex(j, cath);
2788       if (digit < 0) { cath1 = TMath::Even(cath); digit = -digit - 1; } // from the other cathode
2789
2790       if (clus->GetContrib(j,cath) > sigmax[cath1]) {
2791         sigmax[cath1] = clus->GetContrib(j,cath);
2792         maxdig[cath1] = digit;
2793       }
2794     }
2795   }
2796
2797   // Size of pad with maximum signal and reco coordinate distance from the pad center
2798   AliMUONDigit *mdig = 0;
2799   Double_t wx[2], wy[2], dxc[2], dyc[2];
2800   Float_t xpad, ypad, zpad;
2801   Int_t ix, iy;
2802   for (Int_t cath = 0; cath < 2; cath++) {
2803     if (maxdig[cath] < 0) continue;
2804     mdig = fInput->Digit(cath,maxdig[cath]);
2805     isec = fSegmentation[cath]->Sector(mdig->PadX(), mdig->PadY());
2806     wx[cath] = fSegmentation[cath]->Dpx(isec);
2807     wy[cath] = fSegmentation[cath]->Dpy(isec);
2808     fSegmentation[cath]->GetPadI(xreco, yreco, zreco, ix, iy);
2809     isec = fSegmentation[cath]->Sector(ix, iy);
2810     if (isec > 0) {
2811       fSegmentation[cath]->GetPadC(ix, iy, xpad, ypad, zpad);
2812       dxc[cath] = xreco - xpad;
2813       dyc[cath] = yreco - ypad;
2814     }
2815   }
2816
2817   // Check if pad with max charge at the edge (number of neughbours)
2818   Int_t nn, xList[10], yList[10], neighbx[2][2] = {{0,0}, {0,0}}, neighby[2][2]= {{0,0}, {0,0}};
2819   for (Int_t cath = 0; cath < 2; cath++) {
2820     if (maxdig[cath] < 0) continue;
2821     mdig = fInput->Digit(cath,maxdig[cath]);
2822     fSegmentation[cath]->Neighbours(mdig->PadX(), mdig->PadY(), &nn, xList, yList); 
2823     isec = fSegmentation[cath]->Sector(mdig->PadX(), mdig->PadY());
2824     /*??
2825     Float_t sprX = fResponse->SigmaIntegration() * fResponse->ChargeSpreadX();
2826     Float_t sprY = fResponse->SigmaIntegration() * fResponse->ChargeSpreadY();
2827     //fSegmentation[cath]->FirstPad(fInput->DetElemId(),muons[ihit][1], muons[ihit][2], muons[ihit][3], sprX, sprY);  
2828     //fSegmentation[cath]->FirstPad(fInput->DetElemId(),xreco, yreco, zreco, sprX, sprY);  
2829     fSegmentation[cath]->FirstPad(xreco, yreco, zreco, sprX, sprY);  
2830     Int_t border = 0;
2831     //if (fSegmentation[cath]->Sector(fInput->DetElemId(),fSegmentation[cath]->Ix(),fSegmentation[cath]->Iy()) <= 0) {
2832     if (fSegmentation[cath]->Sector(fSegmentation[cath]->Ix(), fSegmentation[cath]->Iy()) <= 0) {
2833       //fSegmentation[cath]->NextPad(fInput->DetElemId());
2834       fSegmentation[cath]->NextPad();
2835       border = 1;
2836     } 
2837     */
2838     for (Int_t j=0; j<nn; j++) {
2839       //if (border && yList[j] < fSegmentation[cath]->Iy()) continue;
2840       fSegmentation[cath]->GetPadC(xList[j], yList[j], xpad, ypad, zpad);
2841       //cout << ch << " " << xList[j] << " " << yList[j] << " " << border << " " << x << " " << y << " " << xpad << " " << ypad << endl;
2842       if (TMath::Abs(xpad) < 1 && TMath::Abs(ypad) < 1) continue;
2843       if (xList[j] == mdig->PadX()-1 || mdig->PadX() == 1 && 
2844           xList[j] == -1) neighbx[cath][0] = 1;
2845       else if (xList[j] == mdig->PadX()+1 || mdig->PadX() == -1 && 
2846                xList[j] == 1) neighbx[cath][1] = 1;
2847       if (yList[j] == mdig->PadY()-1 || mdig->PadY() == 1 &&
2848           yList[j] == -1) neighby[cath][0] = 1;
2849       else if (yList[j] == mdig->PadY()+1 || mdig->PadY() == -1 &&
2850                yList[j] == 1) neighby[cath][1] = 1;
2851     } // for (Int_t j=0; j<nn;
2852     if (neighbx[cath][0] && neighbx[cath][1]) neighbx[cath][0] = 0;
2853     else if (neighbx[cath][1]) neighbx[cath][0] = -1;
2854     else neighbx[cath][0] = 1;
2855     if (neighby[cath][0] && neighby[cath][1]) neighby[cath][0] = 0;
2856     else if (neighby[cath][1]) neighby[cath][0] = -1;
2857     else neighby[cath][0] = 1;
2858   }
2859
2860   Int_t iOver = clus->GetClusterType();
2861   // One-sided cluster
2862   if (!clus->GetMultiplicity(0)) {
2863     neighby[0][0] = neighby[1][0];
2864     wy[0] = wy[1];
2865     if (iOver < 99) iOver += 100 * iOver;
2866     dyc[0] = dyc[1];
2867   } else if (!clus->GetMultiplicity(1)) {
2868     neighbx[1][0] = neighbx[0][0];
2869     wx[1] = wx[0];
2870     if (iOver < 99) iOver += 100 * iOver;
2871     dxc[1] = dxc[0];
2872   }
2873
2874   // Apply corrections and evaluate errors
2875   Double_t errY, errX;
2876   Errors(nInY, nInX, neighby[0][0],neighbx[1][0], fmin, wy[0]*10, wx[1]*10, iOver, 
2877          dyc[0], dxc[1], qTot, yreco, xreco, errY, errX);
2878   errY = TMath::Max (errY, 0.01);
2879   //errY = 0.01;
2880   //errX = TMath::Max (errX, 0.144);
2881   clus->SetX(0, xreco); clus->SetY(0, yreco);
2882   clus->SetErrX(errX); clus->SetErrY(errY);
2883 }
2884
2885 //_____________________________________________________________________________
2886 void AliMUONClusterFinderAZ::Errors(Int_t ny, Int_t nx, Int_t iby, Int_t ibx, Double_t fmin,
2887                                     Double_t wy, Double_t wx, Int_t iover, 
2888                                     Double_t dyc, Double_t /*dxc*/, Double_t qtot, 
2889                                     Double_t &yrec, Double_t &xrec, Double_t &erry, Double_t &errx)
2890 {
2891 /// Correct reconstructed coordinates for some clusters and evaluate errors
2892
2893     erry = 0.01;
2894     errx = 0.144;
2895     Int_t iovery = iover % 100;
2896     Double_t corr = 0;
2897
2898 /* ---> Ny = 1 */
2899     if (ny == 1) {
2900       if (iby != 0) {
2901         // edge effect 
2902         yrec += iby * (0.1823+0.2008)/2;
2903         erry = 0.04587;
2904       } else {
2905         // Find "effective pad width" 
2906         Double_t width = 0.218 / (1.31e-4 * TMath::Exp (2.688 * TMath::Log(qtot)) + 1) * 2;
2907         width = TMath::Min (width, 0.4);
2908         erry = width / TMath::Sqrt(12.);
2909         erry = TMath::Max (erry, 0.01293);
2910       }
2911       goto x; //return;
2912     }
2913
2914 /* ---> "Bad" fit */
2915     if (fmin > 0.4) {
2916       erry = 0.1556;
2917       if (ny == 5) erry = 0.06481;
2918       goto x; //return;
2919     }
2920
2921 /* ---> By != 0 */
2922     if (iby != 0) {
2923       if (ny > 2) {
2924         erry = 0.00417; //0.01010
2925       } else {
2926         // ny = 2 
2927         if (dyc * iby > -0.05) {
2928           Double_t dyc2 = dyc * dyc;
2929           if (iby < 0) {
2930             corr = 0.019 - 0.602 * dyc + 8.739 * dyc2 - 44.209 * dyc2 * dyc;
2931             corr = TMath::Min (corr, TMath::Abs(-0.25-dyc));
2932             yrec -= corr;
2933             //dyc -= corr;
2934             erry = 0.00814;
2935           } else {
2936             corr = 0.006 + 0.300 * dyc + 6.147 * dyc2 + 42.039 * dyc2 * dyc;
2937             corr = TMath::Min (corr, 0.25-dyc);
2938             yrec += corr;
2939             //dyc += corr;
2940             erry = 0.01582;
2941           }
2942         } else {
2943           erry = (0.00303 + 0.00296) / 2;
2944         }
2945       }
2946       goto x; //return;
2947     }
2948
2949 /* ---> Overflows */
2950     if (iovery != 0) {
2951       if (qtot < 3000) {
2952         erry = 0.0671;
2953       } else {
2954         if (iovery > 1) {
2955           erry = 0.09214;
2956         } else if (TMath::Abs(wy - 5) < 0.1) {
2957           erry = 0.061; //0.06622
2958         } else {
2959           erry = 0.00812; // 0.01073 
2960         }
2961       }
2962       goto x; //return;
2963     }
2964
2965 /* ---> "Good" but very high signal */
2966     if (qtot > 4000) {
2967       if (TMath::Abs(wy - 4) < 0.1) {
2968         erry = 0.00117;
2969       } else if (fmin < 0.03 && qtot < 6000) {
2970         erry = 0.01003;
2971       } else {
2972         erry = 0.1931;
2973       }
2974       goto x; //return;
2975     }
2976
2977 /* ---> "Good" clusters */
2978     if (ny > 3) {
2979       if (TMath::Abs(wy - 5) < 0.1) {
2980         erry = 0.0011; //0.00304 
2981       } else if (qtot < 400.) {
2982         erry = 0.0165;
2983       } else {
2984         erry = 0.00135; // 0.00358 
2985       }
2986     } else if (ny == 3) {
2987       if (TMath::Abs(wy - 4) < 0.1) {
2988         erry = 35.407 / (1 + TMath::Exp(5.511*TMath::Log(qtot/265.51))) + 11.564;
2989         //erry = 83.512 / (1 + TMath::Exp(3.344*TMath::Log(qtot/211.58))) + 12.260;
2990       } else {
2991         erry = 147.03 / (1 + TMath::Exp(1.713*TMath::Log(qtot/73.151))) + 9.575;
2992         //erry = 91.743 / (1 + TMath::Exp(2.332*TMath::Log(qtot/151.67))) + 11.453;
2993       }
2994       erry *= 1.e-4;
2995     } else {
2996       // ny = 2 
2997       if (TMath::Abs(wy - 4) < 0.1) {
2998         erry = 60.800 / (1 + TMath::Exp(3.305*TMath::Log(qtot/104.53))) + 11.702;
2999         //erry = 73.128 / (1 + TMath::Exp(5.676*TMath::Log(qtot/120.93))) + 17.839;
3000       } else {
3001         erry = 117.98 / (1 + TMath::Exp(2.005*TMath::Log(qtot/37.649))) + 21.431;
3002         //erry = 99.066 / (1 + TMath::Exp(4.900*TMath::Log(qtot/107.57))) + 25.315;
3003       }
3004       erry *= 1.e-4;
3005     }
3006     //return;
3007
3008  x:
3009 /* ---> X-coordinate */
3010 /* ---> Y-side */    
3011     if (wx > 11) { 
3012       errx = 0.0036;
3013       xrec -= 0.1385;
3014       return;
3015     }
3016 /* ---> Nx = 1 */
3017     if (nx == 1) {
3018       if (TMath::Abs(wx - 6) < 0.1) {
3019         if (qtot < 40) errx = 0.1693;
3020         else errx = 0.06241;
3021       } else if (TMath::Abs(wx - 7.5) < 0.1) {
3022         if (qtot < 40) errx = 0.2173;
3023         else errx = 0.07703;
3024       } else if (TMath::Abs(wx - 10) < 0.1) {
3025         if (ibx == 0) {
3026           if (qtot < 40) errx = 0.2316;
3027           else errx = 0.1426;
3028         } else {
3029           xrec += (0.2115 + 0.1942) / 2 * ibx;
3030           errx = 0.1921;
3031         }
3032       }
3033       return;
3034     }
3035 /* ---> "Bad" fit */
3036     if (fmin > 0.5) {
3037       errx = 0.1591;
3038       return;
3039     }
3040 /* ---> Bx != 0 */
3041     if (ibx != 0) {
3042       if (ibx > 0) { errx = 0.06761; xrec -= 0.03832; }
3043       else { errx = 0.06653; xrec += 0.02581; }
3044       return;
3045     }
3046 /* ---> Overflows */
3047     if (iover != 0) {
3048       if (TMath::Abs(wx - 6) < 0.1) errx = 0.06979;
3049       else if (TMath::Abs(wx - 7.5) < 0.1) errx = 0.1089;
3050       else if (TMath::Abs(wx - 10) < 0.1) errx = 0.09847;
3051       return;
3052     }
3053 /* ---> Good */
3054     if (TMath::Abs(wx - 6) < 0.1) errx = 0.06022;
3055     else if (TMath::Abs(wx - 7.5) < 0.1) errx = 0.07247;
3056     else if (TMath::Abs(wx - 10) < 0.1) errx = 0.07359;
3057 }
3058