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