]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterSplitterMLEM.cxx
dfc92a7e102575e2adab1edf2e320e60e9b7843c
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterSplitterMLEM.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 AliMUONClusterSplitterMLEM
20 /// 
21 /// Splitter class for the MLEM algorithm. Performs fitting procedure
22 /// with up to 3 hit candidates and tries to split clusters if the number
23 /// of candidates exceeds 3.
24 ///
25 /// \author Laurent Aphecetche (for the "new" C++ structure) and 
26 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
27 //-----------------------------------------------------------------------------
28
29 #include "AliMUONClusterSplitterMLEM.h"
30 #include "AliMUONClusterFinderMLEM.h" // for status flag constants
31
32 #include "AliMUONCluster.h"
33 #include "AliMUONPad.h"
34 #include "AliMUONPad.h"
35 #include "AliMpStationType.h"
36 #include "AliMUONConstants.h"
37 #include "AliMpDEManager.h"
38 #include "AliMUONMathieson.h"
39
40 #include "AliLog.h"
41
42 #include <TClonesArray.h>
43 #include <TH2.h>
44 #include <TMath.h>
45 #include <TMatrixD.h>
46 #include <TObjArray.h>
47 #include <TROOT.h>
48 #include <TRandom.h>
49 #include <Riostream.h>
50
51 /// \cond CLASSIMP
52 ClassImp(AliMUONClusterSplitterMLEM)
53 /// \endcond
54
55 //const Double_t AliMUONClusterSplitterMLEM::fgkCouplMin = 1.e-3; // threshold on coupling 
56 const Double_t AliMUONClusterSplitterMLEM::fgkCouplMin = 1.e-2; // threshold on coupling 
57
58 //_____________________________________________________________________________
59 AliMUONClusterSplitterMLEM::AliMUONClusterSplitterMLEM(Int_t detElemId, 
60                                                        TObjArray* fPixArray) 
61 : TObject(),
62 fPixArray(fPixArray),
63 fMathieson(0x0),
64 fDetElemId(detElemId),
65 fNpar(0),
66 fQtot(0),
67 fnCoupled(0),
68 fDebug(0)
69 {
70   /// Constructor
71   
72   AliMp::StationType stationType = AliMpDEManager::GetStationType(fDetElemId);
73   
74   Float_t kx3 = AliMUONConstants::SqrtKx3();
75   Float_t ky3 = AliMUONConstants::SqrtKy3();
76   Float_t pitch = AliMUONConstants::Pitch();
77   
78   if ( stationType == AliMp::kStation1 )
79   {
80     kx3 = AliMUONConstants::SqrtKx3St1();
81     ky3 = AliMUONConstants::SqrtKy3St1();
82     pitch = AliMUONConstants::PitchSt1();
83   }
84   
85   fMathieson = new AliMUONMathieson;
86   
87   fMathieson->SetPitch(pitch);
88   fMathieson->SetSqrtKx3AndDeriveKx2Kx4(kx3);
89   fMathieson->SetSqrtKy3AndDeriveKy2Ky4(ky3);
90   
91 }
92
93 //_____________________________________________________________________________
94 AliMUONClusterSplitterMLEM::~AliMUONClusterSplitterMLEM()
95 {
96   /// Destructor
97   
98   delete fMathieson;
99 }
100
101 //_____________________________________________________________________________
102 void 
103 AliMUONClusterSplitterMLEM::AddBin(TH2 *mlem, 
104                                    Int_t ic, Int_t jc, Int_t mode, 
105                                    Bool_t *used, TObjArray *pix)
106 {
107   /// Add a bin to the cluster
108   
109   Int_t nx = mlem->GetNbinsX();
110   Int_t ny = mlem->GetNbinsY();
111   Double_t cont1, cont = mlem->GetCellContent(jc,ic);
112   AliMUONPad *pixPtr = 0;
113   
114   Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
115   for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
116     for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
117       if (i != ic && j != jc) continue;
118       if (used[(i-1)*nx+j-1]) continue;
119       cont1 = mlem->GetCellContent(j,i);
120       if (mode && cont1 > cont) continue;
121       used[(i-1)*nx+j-1] = kTRUE;
122       if (cont1 < 0.5) continue;
123       if (pix) pix->Add(BinToPix(mlem,j,i)); 
124       else {
125         pixPtr = new AliMUONPad (mlem->GetXaxis()->GetBinCenter(j), 
126                                  mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
127         fPixArray->Add(pixPtr);
128       }
129       AddBin(mlem, i, j, mode, used, pix); // recursive call
130     }
131   }
132 }
133
134 //_____________________________________________________________________________
135 void 
136 AliMUONClusterSplitterMLEM::AddCluster(Int_t ic, Int_t nclust, 
137                                        TMatrixD& aijcluclu, 
138                                        Bool_t *used, Int_t *clustNumb, Int_t &nCoupled)
139 {
140   /// Add a cluster to the group of coupled clusters
141   
142   for (Int_t i = 0; i < nclust; ++i) {
143     if (used[i]) continue;
144     if (aijcluclu(i,ic) < fgkCouplMin) continue;
145     used[i] = kTRUE;
146     clustNumb[nCoupled++] = i;
147     AddCluster(i, nclust, aijcluclu, used, clustNumb, nCoupled);
148   }
149 }
150
151 //_____________________________________________________________________________
152 TObject* 
153 AliMUONClusterSplitterMLEM::BinToPix(TH2 *mlem,
154                                      Int_t jc, Int_t ic)
155 {
156   /// Translate histogram bin to pixel 
157   
158   Double_t yc = mlem->GetYaxis()->GetBinCenter(ic);
159   Double_t xc = mlem->GetXaxis()->GetBinCenter(jc);
160   
161   Int_t nPix = fPixArray->GetEntriesFast();
162   AliMUONPad *pixPtr = NULL;
163   
164   // Compare pixel and bin positions
165   for (Int_t i = 0; i < nPix; ++i) {
166     pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
167     if (pixPtr->Charge() < 0.5) continue;
168     if (TMath::Abs(pixPtr->Coord(0)-xc)<1.e-4 && TMath::Abs(pixPtr->Coord(1)-yc)<1.e-4) 
169     {
170       //return (TObject*) pixPtr;
171       return pixPtr;
172     }
173   }
174   AliError(Form(" Something wrong ??? %f %f ", xc, yc));
175   return NULL;
176 }
177
178 //_____________________________________________________________________________
179 Float_t
180 AliMUONClusterSplitterMLEM::ChargeIntegration(Double_t x, Double_t y,
181                                               const AliMUONPad& pad)
182 {
183   /// Compute the Mathieson integral on pad area, assuming the center
184   /// of the Mathieson is at (x,y)
185   
186   TVector2 lowerLeft(TVector2(x,y)-pad.Position()-pad.Dimensions());
187   TVector2 upperRight(lowerLeft + pad.Dimensions()*2.0);
188   
189         return fMathieson->IntXY(lowerLeft.X(),lowerLeft.Y(),
190                            upperRight.X(),upperRight.Y());
191 }
192
193 //_____________________________________________________________________________
194 void 
195 AliMUONClusterSplitterMLEM::Fcn1(const AliMUONCluster& cluster, 
196                                     Int_t & /*fNpar*/, Double_t * /*gin*/, 
197                                     Double_t &f, Double_t *par, Int_t iflag)
198 {
199   /// Computes the functional to be minimized
200   
201   Int_t indx, npads=0;
202   Double_t charge, delta, coef=0, chi2=0, qTot = 0;
203   static Double_t qAver = 0;
204   
205   Int_t mult = cluster.Multiplicity(), iend = fNpar / 3;
206   for (Int_t j = 0; j < mult; ++j) 
207   {
208     AliMUONPad* pad = cluster.Pad(j);
209     //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
210     if ( pad->Status() != AliMUONClusterFinderMLEM::fgkUseForFit) continue;
211     if (iflag == 0) {
212       if ( pad->IsReal() ) npads++; // exclude virtual pads
213       qTot += pad->Charge(); 
214     }
215     charge = 0;
216     for (Int_t i = 0; i <= iend; ++i)
217     { 
218       // sum over hits
219       indx = 3 * i;
220       coef = Param2Coef(i, coef, par);
221       charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
222     }
223     charge *= fQtot;
224     delta = charge - pad->Charge(); 
225     delta *= delta;
226     delta /= pad->Charge(); 
227     chi2 += delta;
228   } // for (Int_t j=0;
229   if (iflag == 0) qAver = qTot / npads;
230   f = chi2 / qAver;
231 }
232
233 //_____________________________________________________________________________
234 Double_t AliMUONClusterSplitterMLEM::Param2Coef(Int_t icand, Double_t coef, Double_t *par)
235 {
236   /// Extract hit contribution scale factor from fit parameters
237   
238   if (fNpar == 2) return 1.;
239   if (fNpar == 5) return icand==0 ? par[2] : TMath::Max(1.-par[2],0.);
240   if (icand == 0) return par[2];
241   if (icand == 1) return TMath::Max((1.-par[2])*par[5], 0.);
242   return TMath::Max(1.-par[2]-coef,0.);
243 }
244
245 //_____________________________________________________________________________
246 Int_t 
247 AliMUONClusterSplitterMLEM::Fit(const AliMUONCluster& cluster,
248                                 Int_t iSimple, Int_t nfit, 
249                                 Int_t *clustFit, TObjArray **clusters, 
250                                 Double_t *parOk,
251                                 TObjArray& clusterList)
252 {
253   /// Steering function and fitting procedure for the fit of pad charge distribution
254   
255   //  AliDebug(2,Form("iSimple=%d nfit=%d",iSimple,nfit));
256   
257   TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
258   Double_t xmin = mlem->GetXaxis()->GetXmin() - mlem->GetXaxis()->GetBinWidth(1);
259   Double_t xmax = mlem->GetXaxis()->GetXmax() + mlem->GetXaxis()->GetBinWidth(1);
260   Double_t ymin = mlem->GetYaxis()->GetXmin() - mlem->GetYaxis()->GetBinWidth(1);
261   Double_t ymax = mlem->GetYaxis()->GetXmax() + mlem->GetYaxis()->GetBinWidth(1);
262   Double_t xPad = 0, yPad = 99999;
263   
264   // Number of pads to use and number of virtual pads
265   Int_t npads = 0, nVirtual = 0, nfit0 = nfit;
266   //cluster.Print("full");
267   Int_t mult = cluster.Multiplicity();
268   for (Int_t i = 0; i < mult; ++i ) 
269   {
270     AliMUONPad* pad = cluster.Pad(i);
271     if ( !pad->IsReal() ) ++nVirtual;
272     //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
273     if ( pad->Status() != AliMUONClusterFinderMLEM::fgkUseForFit ) continue;
274     if ( pad->IsReal() )
275     {
276       ++npads;
277       if (yPad > 9999) 
278       { 
279         xPad = pad->X();
280         yPad = pad->Y();
281       } 
282       else 
283       {
284         if (pad->DY() < pad->DX() ) 
285         {
286           yPad = pad->Y();
287         }
288         else 
289         {
290           xPad = pad->X();
291         }
292       }
293     }
294   }
295   
296   fNpar = 0;
297   fQtot = 0;
298   
299   if (npads < 2) return 0; 
300   
301   // FIXME : AliWarning("Reconnect the following code for hit/track passing ?");
302   
303   //  Int_t tracks[3] = {-1, -1, -1};
304   
305   /*
306    Int_t digit = 0;
307    AliMUONDigit *mdig = 0;
308    for (Int_t cath=0; cath<2; cath++) {  
309      for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
310        if (fPadIJ[0][i] != cath) continue;
311        if (fPadIJ[1][i] != 1) continue;
312        if (fXyq[3][i] < 0) continue; // exclude virtual pads
313        digit = TMath::Nint (fXyq[5][i]);
314        if (digit >= 0) mdig = fInput->Digit(cath,digit);
315        else mdig = fInput->Digit(TMath::Even(cath),-digit-1);
316        //if (!mdig) mdig = fInput->Digit(TMath::Even(cath),digit);
317        if (!mdig) continue; // protection for cluster display
318        if (mdig->Hit() >= 0) {
319          if (tracks[0] < 0) {
320            tracks[0] = mdig->Hit();
321            tracks[1] = mdig->Track(0);
322          } else if (mdig->Track(0) < tracks[1]) {
323            tracks[0] = mdig->Hit();
324            tracks[1] = mdig->Track(0);
325          }
326        }
327        if (mdig->Track(1) >= 0 && mdig->Track(1) != tracks[1]) {
328          if (tracks[2] < 0) tracks[2] = mdig->Track(1);
329          else tracks[2] = TMath::Min (tracks[2], mdig->Track(1));
330        }
331      } // for (Int_t i=0;
332   } // for (Int_t cath=0;
333    */
334   
335   // Get number of pads in X and Y 
336   //const Int_t kStatusToTest(1);
337   const Int_t kStatusToTest(AliMUONClusterFinderMLEM::fgkUseForFit);
338   
339   AliMpIntPair nofPads = cluster.NofPads(kStatusToTest);
340   Int_t nInX = nofPads.GetFirst();
341   Int_t nInY = nofPads.GetSecond();
342
343   if (fDebug) {
344     Int_t npadOK = 0;
345     for (Int_t j = 0; j < cluster.Multiplicity(); ++j) {
346       AliMUONPad *pad = cluster.Pad(j);
347       //if (pad->Status() == 1 && !pad->IsSaturated()) npadOK++;
348       if (pad->Status() == AliMUONClusterFinderMLEM::fgkUseForFit && !pad->IsSaturated()) npadOK++;
349     }
350     cout << " Number of pads to fit: " << npadOK << endl;
351     cout << " nInX and Y: " << nInX << " " << nInY << endl;
352   }
353   
354   Int_t nfitMax = 3; 
355   nfitMax = TMath::Min (nfitMax, (npads + 1) / 3);
356   if (nfitMax > 1) {
357     if (nInX < 3 && nInY < 3 || nInX == 3 && nInY < 3 || nInX < 3 && nInY == 3) nfitMax = 1; // not enough pads in each direction
358   }
359   if (nfit > nfitMax) nfit = nfitMax;
360   
361   // Take cluster maxima as fitting seeds
362   TObjArray *pix;
363   AliMUONPad *pixPtr;
364   Int_t npxclu;
365   Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8], qq = 0;
366   Double_t xyseed[3][2], qseed[3], xyCand[3][2] = {{0},{0}}, sigCand[3][2] = {{0},{0}};
367   
368   for (Int_t ifit = 1; ifit <= nfit0; ++ifit) 
369   {
370     cmax = 0;
371     pix = clusters[clustFit[ifit-1]];
372     npxclu = pix->GetEntriesFast();
373     //qq = 0;
374     for (Int_t clu = 0; clu < npxclu; ++clu) 
375     {
376       pixPtr = (AliMUONPad*) pix->UncheckedAt(clu);
377       cont = pixPtr->Charge();
378       fQtot += cont;
379       if (cont > cmax) 
380       { 
381         cmax = cont; 
382         xseed = pixPtr->Coord(0);
383         yseed = pixPtr->Coord(1);
384       }
385       qq += cont;
386       xyCand[0][0] += pixPtr->Coord(0) * cont;
387       xyCand[0][1] += pixPtr->Coord(1) * cont;
388       sigCand[0][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
389       sigCand[0][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
390     }
391     xyseed[ifit-1][0] = xseed;
392     xyseed[ifit-1][1] = yseed;
393     qseed[ifit-1] = cmax;
394   } // for (Int_t ifit=1;
395   
396   xyCand[0][0] /= qq; // <x>
397   xyCand[0][1] /= qq; // <y>
398   sigCand[0][0] = sigCand[0][0]/qq - xyCand[0][0]*xyCand[0][0]; // <x^2> - <x>^2
399   sigCand[0][0] = sigCand[0][0] > 0 ? TMath::Sqrt (sigCand[0][0]) : 0;
400   sigCand[0][1] = sigCand[0][1]/qq - xyCand[0][1]*xyCand[0][1]; // <y^2> - <y>^2
401   sigCand[0][1] = sigCand[0][1] > 0 ? TMath::Sqrt (sigCand[0][1]) : 0;
402   if (fDebug) cout << xyCand[0][0] << " " << xyCand[0][1] << " " << sigCand[0][0] << " " << sigCand[0][1] << endl;
403   
404   Int_t nDof, maxSeed[3];//, nMax = 0;
405     
406   TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
407     
408   Double_t step[3]={0.01,0.002,0.02}, fmin, chi2o = 9999, chi2n;
409   Double_t *gin = 0, func0, func1, param[8], step0[8];
410   Double_t param0[2][8]={{0},{0}}, deriv[2][8]={{0},{0}}; 
411   Double_t shift[8], stepMax, derMax, parmin[8], parmax[8], func2[2], shift0;
412   Double_t delta[8], scMax, dder[8], estim, shiftSave = 0;
413   Int_t min, max, nCall = 0, nLoop, idMax = 0, iestMax = 0, nFail;
414   Double_t rad, dist[3] = {0};
415     
416   // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is 
417   // lower, try 3-track (if number of pads is sufficient).
418   Int_t iflag = 0; // for the first call of fcn1
419   for (Int_t iseed = 0; iseed < nfit; ++iseed) 
420   {
421       
422     Int_t memory[8] = {0};
423     if (iseed) 
424     { 
425       for (Int_t j = 0; j < fNpar; ++j) 
426       {
427         param[j] = parOk[j]; 
428       }
429       param[fNpar] = 0.6;
430       parmin[fNpar] = 0; 
431       parmax[fNpar++] = 1; 
432     }
433       
434     if (nfit == 1) 
435     {
436       param[fNpar] = xyCand[0][0]; // take COG
437     }
438     else 
439     {
440       param[fNpar] = xyseed[maxSeed[iseed]][0];
441       //param[fNpar] = fNpar==0 ? -16.1651 : -15.2761; 
442     }
443     parmin[fNpar] = xmin; 
444     parmax[fNpar++] = xmax; 
445     if (nfit == 1) 
446     {
447       param[fNpar] = xyCand[0][1]; // take COG
448     }
449     else 
450     {
451       param[fNpar] = xyseed[maxSeed[iseed]][1];
452       //param[fNpar] = fNpar==1 ? -15.1737 : -15.8487;
453     }
454     parmin[fNpar] = ymin; 
455     parmax[fNpar++] = ymax; 
456
457     for (Int_t j = 0; j < fNpar; ++j) 
458     {
459       step0[j] = shift[j] = step[j%3];
460     }
461
462     if (iseed) 
463     { 
464       for (Int_t j = 0; j < fNpar; ++j) 
465       {
466         param0[1][j] = 0; 
467       }
468     }
469     if (fDebug) {
470       for (Int_t j = 0; j < fNpar; ++j) cout << param[j] << " "; 
471       cout << endl;
472     }
473       
474     // Try new algorithm
475     min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
476       
477     while (1) 
478     {
479       max = !min;
480       Fcn1(cluster,fNpar, gin, func0, param, iflag); nCall++;
481       iflag = 1;
482       //cout << " Func: " << func0 << endl;
483       
484       func2[max] = func0;
485       for (Int_t j = 0; j < fNpar; ++j) 
486       {
487         param0[max][j] = param[j];
488         delta[j] = step0[j];
489         param[j] += delta[j] / 10;
490         if (j > 0) param[j-1] -= delta[j-1] / 10;
491         Fcn1(cluster,fNpar, gin, func1, param, iflag); nCall++;
492         deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
493         //cout << j << " " << deriv[max][j] << endl;
494         dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) / 
495           (param0[0][j] - param0[1][j]) : 0; // second derivative
496       }
497       param[fNpar-1] -= delta[fNpar-1] / 10;
498       if (nCall > 2000) break;
499         
500       min = func2[0] < func2[1] ? 0 : 1;
501       nFail = min == max ? 0 : nFail + 1;
502         
503       stepMax = derMax = estim = 0;
504       for (Int_t j = 0; j < fNpar; ++j) 
505       { 
506         // Estimated distance to minimum
507         shift0 = shift[j];
508         if (nLoop == 1) 
509         {
510           shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
511         }
512         else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3) 
513         {
514           shift[j] = 0;
515         }
516         else if (deriv[min][j]*deriv[!min][j] > 0 && TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])
517                  || TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3 || TMath::Abs(dder[j]) < 1.e-6) 
518         {
519           shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
520           if (min == max) 
521           { 
522             if (memory[j] > 1) 
523             { 
524               shift[j] *= 2; 
525             } 
526             memory[j]++;
527           }
528         }
529         else 
530         {
531           shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
532           memory[j] = 0;
533         }
534           
535         Double_t es = TMath::Abs(shift[j]) / step0[j];
536         if (es > estim) 
537         { 
538           estim = es;
539           iestMax = j;
540         }
541           
542         // Too big step
543         if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; // 
544         
545         // Failed to improve minimum
546         if (min != max) 
547         {
548           memory[j] = 0;
549           param[j] = param0[min][j];
550           if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j]) 
551           {
552             shift[j] = (shift[j] + shift0) / 2;
553           }
554           else 
555           {
556             shift[j] /= -2;
557           }
558         } 
559           
560         // Too big step
561         if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min]) 
562         {
563           shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
564         }
565           
566         // Introduce step relaxation factor
567         if (memory[j] < 3) 
568         {
569           scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
570           if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax) 
571           {
572             shift[j] = TMath::Sign (shift0*scMax, shift[j]);
573           }
574         }
575         param[j] += shift[j]; 
576         // Check parameter limits
577         if (param[j] < parmin[j]) 
578         { 
579           shift[j] = parmin[j] - param[j]; 
580           param[j] = parmin[j]; 
581         } 
582         else if (param[j] > parmax[j]) 
583         {
584           shift[j] = parmax[j] - param[j];
585           param[j] = parmax[j];
586         }
587         //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
588         stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
589         if (TMath::Abs(deriv[min][j]) > derMax) 
590         {
591           idMax = j;
592           derMax = TMath::Abs (deriv[min][j]);
593         }
594       } // for (Int_t j=0; j<fNpar;
595         
596       if (estim < 1 && derMax < 2 || nLoop > 150) break; // minimum was found
597         
598       nLoop++;
599         
600       // Check for small step
601       if (shift[idMax] == 0) 
602       { 
603         shift[idMax] = step0[idMax]/10; 
604         param[idMax] += shift[idMax]; 
605         continue; 
606       }
607         
608       if (!memory[idMax] && derMax > 0.5 && nLoop > 10) 
609       {
610         if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10) 
611         {
612           if (min == max) dder[idMax] = -dder[idMax];
613           shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10; 
614           param[idMax] += shift[idMax];
615           stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
616           if (min == max) shiftSave = shift[idMax];
617         }
618         if (nFail > 10) 
619         {
620           param[idMax] -= shift[idMax];
621           shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
622           param[idMax] += shift[idMax];
623         }
624       }      
625     } // while (1)
626       
627     fmin = func2[min];
628     
629     nDof = npads - fNpar + nVirtual;
630     if (!nDof) nDof++;
631     chi2n = fmin / nDof;
632     if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
633       
634     //if (fNpar > 2) cout << param0[min][fNpar-3] << " " << chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) << endl;
635     //if (chi2n*1.2+1.e-6 > chi2o ) 
636     if (fNpar > 2 && (chi2n > chi2o || iseed == nfit-1 
637                       && chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) > chi2o)) 
638       { fNpar -= 3; break; }
639       
640     // Save parameters and errors
641       
642     if (nInX == 1) {
643       // One pad per direction 
644       //for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5) param0[min][i] = xPad;
645       for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5) 
646         param0[min][i] = xyCand[0][0];
647     }
648     if (nInY == 1) {
649       // One pad per direction 
650       //for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6) param0[min][i] = yPad;
651       for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6) 
652         param0[min][i] = xyCand[0][1];
653     }
654       
655     /*
656       if (iseed > 0) {
657       // Find distance to the nearest neighbour
658       dist[0] = dist[1] = TMath::Sqrt ((param0[min][0]-param0[min][2])*
659       (param0[min][0]-param0[min][2])
660       +(param0[min][1]-param0[min][3])*
661       (param0[min][1]-param0[min][3]));
662       if (iseed > 1) {
663       dist[2] = TMath::Sqrt ((param0[min][0]-param0[min][5])*
664       (param0[min][0]-param0[min][5])
665       +(param0[min][1]-param0[min][6])*
666       (param0[min][1]-param0[min][6]));
667       rad = TMath::Sqrt ((param0[min][2]-param0[min][5])*
668       (param0[min][2]-param0[min][5])
669       +(param0[min][3]-param0[min][6])*
670       (param0[min][3]-param0[min][6]));
671       if (dist[2] < dist[0]) dist[0] = dist[2];
672       if (rad < dist[1]) dist[1] = rad;
673       if (rad < dist[2]) dist[2] = rad;
674       }
675       cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
676       if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
677       }
678     */
679       
680     for (Int_t i = 0; i < fNpar; ++i) {
681       parOk[i] = param0[min][i];
682       //errOk[i] = fmin;
683       errOk[i] = chi2n;
684       // Bounded params
685       parOk[i] = TMath::Max (parOk[i], parmin[i]);
686       parOk[i] = TMath::Min (parOk[i], parmax[i]);
687     }
688       
689     chi2o = chi2n;
690     if (fmin < 0.1) break; // !!!???
691   } // for (Int_t iseed=0; 
692    
693   if (fDebug) {
694     for (Int_t i=0; i<fNpar; ++i) {
695       if (i == 4 || i == 7) {
696         if (i == 7 || i == 4 && fNpar < 7) cout << parOk[i] << endl;
697         else cout << parOk[i] * (1-parOk[7]) << endl;
698         continue;
699       }
700       cout << parOk[i] << " " << errOk[i] << endl;
701     }
702   }
703   nfit = (fNpar + 1) / 3;
704   dist[0] = dist[1] = dist[2] = 0;
705   
706   if (nfit > 1) {
707     // Find distance to the nearest neighbour
708     dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
709                                      (parOk[0]-parOk[2])
710                                      +(parOk[1]-parOk[3])*
711                                      (parOk[1]-parOk[3]));
712     if (nfit > 2) {
713       dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
714                              (parOk[0]-parOk[5])
715                              +(parOk[1]-parOk[6])*
716                              (parOk[1]-parOk[6]));
717       rad = TMath::Sqrt ((parOk[2]-parOk[5])*
718                          (parOk[2]-parOk[5])
719                          +(parOk[3]-parOk[6])*
720                          (parOk[3]-parOk[6]));
721       if (dist[2] < dist[0]) dist[0] = dist[2];
722       if (rad < dist[1]) dist[1] = rad;
723       if (rad < dist[2]) dist[2] = rad;
724     }
725   }
726     
727   Int_t indx;
728   
729   Double_t coef = 0;
730   if (iSimple) fnCoupled = 0;
731   for (Int_t j = 0; j < nfit; ++j) {
732     indx = 3 * j;
733     coef = Param2Coef(j, coef, parOk);
734       
735     //void AliMUONClusterFinderMLEM::AddRawCluster(Double_t x, Double_t y, 
736     //                                             Double_t qTot, Double_t fmin,
737     //                                             Int_t nfit, Int_t *tracks, 
738     //                                             Double_t /*sigx*/, 
739     //                                             Double_t /*sigy*/, 
740     //                                             Double_t /*dist*/)
741     
742     if ( coef*fQtot >= 14 )
743     {
744       //AZ AliMUONCluster* cluster1 = new AliMUONCluster();
745       AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
746       
747       cluster1->SetCharge(coef*fQtot,coef*fQtot);
748       cluster1->SetPosition(TVector2(parOk[indx],parOk[indx+1]),TVector2(sigCand[0][0],sigCand[0][1]));
749       cluster1->SetChi2(dist[TMath::LocMin(nfit,dist)]);
750       
751       // FIXME: we miss some information in this cluster, as compared to 
752       // the original AddRawCluster code.
753       
754       AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
755                       fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
756                       cluster1->Position().X(),cluster1->Position().Y()));
757         
758       clusterList.Add(cluster1);
759     }
760     //      AddRawCluster (parOk[indx], // double x
761     //                     parOk[indx+1], // double y
762     //                     coef*qTot, // double charge
763     //                     errOk[indx], // double fmin
764     //                     nfit0+10*nfit+100*nMax+10000*fnCoupled, // int nfit
765     //                     tracks, // int* tracks
766     //                     sigCand[0][0], // double sigx
767     //                     sigCand[0][1], // double sigy
768     //                     dist[TMath::LocMin(nfit,dist)] // double dist
769     //                     );
770   }
771   return nfit;
772 }  
773
774
775 //_____________________________________________________________________________
776 void
777 AliMUONClusterSplitterMLEM::Split(const AliMUONCluster& cluster,
778                                   TH2 *mlem, Double_t *coef,
779                                   TObjArray& clusterList)
780 {
781   /// The main steering function to work with clusters of pixels in anode
782   /// plane (find clusters, decouple them from each other, merge them (if
783   /// necessary), pick up coupled pads, call the fitting function)
784   
785   Int_t nx = mlem->GetNbinsX();
786   Int_t ny = mlem->GetNbinsY();
787   Int_t nPix = fPixArray->GetEntriesFast();
788   
789   Double_t cont;
790   Int_t nclust = 0, indx, indx1, nxy = ny * nx; 
791   Bool_t *used = new Bool_t[nxy];
792   
793   for (Int_t j = 0; j < nxy; ++j) used[j] = kFALSE; 
794   
795   TObjArray *clusters[200]={0};
796   TObjArray *pix;
797   
798   // Find clusters of histogram bins (easier to work in 2-D space)
799   for (Int_t i = 1; i <= ny; ++i) 
800   {
801     for (Int_t j = 1; j <= nx; ++j) 
802     {
803       indx = (i-1)*nx + j - 1;
804       if (used[indx]) continue;
805       cont = mlem->GetCellContent(j,i);
806       if (cont < 0.5) continue;
807       pix = new TObjArray(20);
808       used[indx] = 1;
809       pix->Add(BinToPix(mlem,j,i));
810       AddBin(mlem, i, j, 0, used, pix); // recursive call
811       if (nclust >= 200) AliFatal(" Too many clusters !!!");
812       clusters[nclust++] = pix;
813     } // for (Int_t j=1; j<=nx; j++) {
814   } // for (Int_t i=1; i<=ny;
815   if (fDebug) cout << nclust << endl;
816   delete [] used;
817   
818   // Compute couplings between clusters and clusters to pads
819   Int_t npad = cluster.Multiplicity();
820   
821   // Exclude pads with overflows
822   /*
823   for (Int_t j = 0; j < npad; ++j) 
824   {
825     AliMUONPad* pad = cluster.Pad(j);
826     if ( pad->IsSaturated() )
827     {
828       pad->SetStatus(-5); 
829     }
830     else 
831     {
832       pad->SetStatus(0);
833     }
834   }
835   */
836   
837   // Compute couplings of clusters to pads (including overflows)
838   TMatrixD aijclupad(nclust,npad);
839   aijclupad = 0;
840   Int_t npxclu;
841   for (Int_t iclust = 0; iclust < nclust; ++iclust) 
842   {
843     pix = clusters[iclust];
844     npxclu = pix->GetEntriesFast();
845     for (Int_t i = 0; i < npxclu; ++i) 
846     {
847       indx = fPixArray->IndexOf(pix->UncheckedAt(i));
848       for (Int_t j = 0; j < npad; ++j) 
849       {
850         //AliMUONPad* pad = cluster.Pad(j);
851         //if ( pad->Status() < 0 && pad->Status() != -5) continue;
852         if (coef[j*nPix+indx] < fgkCouplMin) continue;
853         aijclupad(iclust,j) += coef[j*nPix+indx];
854       }
855     }
856   }
857   
858   // Compute couplings between clusters (exclude overflows)
859   TMatrixD aijcluclu(nclust,nclust);
860   aijcluclu = 0;
861   for (Int_t iclust = 0; iclust < nclust; ++iclust) 
862   {
863     for (Int_t j = 0; j < npad; ++j) 
864     {
865       // Exclude overflows
866       //if ( cluster.Pad(j)->Status() < 0) continue;
867       if ( cluster.Pad(j)->IsSaturated()) continue;
868       if (aijclupad(iclust,j) < fgkCouplMin) continue;
869       for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) 
870       {
871         if (aijclupad(iclust1,j) < fgkCouplMin) continue;
872         aijcluclu(iclust,iclust1) += 
873           TMath::Sqrt (aijclupad(iclust,j)*aijclupad(iclust1,j));
874       }
875     }
876   }
877   for (Int_t iclust = 0; iclust < nclust; ++iclust) 
878   {
879     for (Int_t iclust1 = iclust+1; iclust1 < nclust; ++iclust1) 
880     {
881       aijcluclu(iclust1,iclust) = aijcluclu(iclust,iclust1);
882     }
883   }
884   
885   if (fDebug && nclust > 1) aijcluclu.Print();
886
887   // Find groups of coupled clusters
888   used = new Bool_t[nclust];
889   for (Int_t j = 0; j < nclust; ++j) used[j] = kFALSE;
890
891   Int_t *clustNumb = new Int_t[nclust];
892   Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
893   //Double_t parOk[8];
894   Double_t parOk[8] = {0}; //AZ
895   
896   for (Int_t igroup = 0; igroup < nclust; ++igroup) 
897   {
898     if (used[igroup]) continue;
899     used[igroup] = kTRUE;
900     clustNumb[0] = igroup;
901     nCoupled = 1;
902     // Find group of coupled clusters
903     AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
904     
905     if (fDebug) {                                                                      
906       cout << " nCoupled: " << nCoupled << endl;
907       for (Int_t i=0; i<nCoupled; ++i) cout << clustNumb[i] << " "; cout << endl;
908     }
909     
910     fnCoupled = nCoupled;
911     
912     while (nCoupled > 0) 
913     {
914       if (nCoupled < 4) 
915       {
916         nForFit = nCoupled;
917         for (Int_t i = 0; i < nCoupled; ++i) clustFit[i] = clustNumb[i];
918       } 
919       else 
920       {
921         // Too many coupled clusters to fit - try to decouple them
922         // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with 
923         // all the others in the group 
924         for (Int_t j = 0; j < 3; ++j) minGroup[j] = -1;
925         Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
926         
927         // Flag clusters for fit
928         nForFit = 0;
929         while (minGroup[nForFit] >= 0 && nForFit < 3)
930         {
931           if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
932           clustFit[nForFit] = clustNumb[minGroup[nForFit]];
933           clustNumb[minGroup[nForFit]] -= 999;
934           nForFit++;
935         }
936         if (fDebug) cout << " nForFit " << nForFit << " " << coupl << endl;
937       } // else
938       
939       // Select pads for fit. 
940       if (SelectPad(cluster,nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1) 
941       {
942         // Deselect pads
943         for (Int_t j = 0; j < npad; ++j)
944         {
945           AliMUONPad* pad = cluster.Pad(j);
946           //if ( pad->Status()==1 ) pad->SetStatus(0);
947           //if ( pad->Status()==-9) pad->SetStatus(-5);
948           if ( pad->Status() == AliMUONClusterFinderMLEM::fgkUseForFit ||
949                pad->Status() == AliMUONClusterFinderMLEM::fgkCoupled) 
950             pad->SetStatus(AliMUONClusterFinderMLEM::fgkZero);
951         }
952         // Merge the failed cluster candidates (with too few pads to fit) with 
953         // the one with the strongest coupling
954         Merge(cluster,nForFit, nCoupled, clustNumb, clustFit, clusters, aijcluclu, aijclupad);
955       } 
956       else 
957       {
958         // Do the fit
959         nfit = Fit(cluster,0, nForFit, clustFit, clusters, parOk, clusterList);
960         if (nfit == 0) { 
961           //cout << " (nfit == 0) " << fNpar << " " << cluster.Multiplicity() << endl; 
962           fNpar = 0; // should be 0 by itself but just in case ...
963         }
964       }
965       
966       // Subtract the fitted charges from pads with strong coupling and/or
967       // return pads for further use
968       UpdatePads(cluster,nfit, parOk);
969       
970       // Mark used pads
971       for (Int_t j = 0; j < npad; ++j) 
972       {
973         AliMUONPad* pad = cluster.Pad(j);
974         //if ( pad->Status()==1 ) pad->SetStatus(-2);
975         //if ( pad->Status()==-9) pad->SetStatus(-5);
976         if ( pad->Status() == AliMUONClusterFinderMLEM::fgkUseForFit ) 
977           pad->SetStatus(AliMUONClusterFinderMLEM::fgkModified);
978       }
979       
980       // Sort the clusters (move to the right the used ones)
981       Int_t beg = 0, end = nCoupled - 1;
982       while (beg < end) 
983       {
984         if (clustNumb[beg] >= 0) { ++beg; continue; }
985         for (Int_t j = end; j > beg; --j) 
986         {
987           if (clustNumb[j] < 0) continue;
988           end = j - 1;
989           indx = clustNumb[beg];
990           clustNumb[beg] = clustNumb[j];
991           clustNumb[j] = indx;
992           break;
993         }
994         ++beg;
995       }
996       
997       nCoupled -= nForFit;
998       if (nCoupled > 3) 
999       {
1000         // Remove couplings of used clusters
1001         for (Int_t iclust = nCoupled; iclust < nCoupled+nForFit; ++iclust) 
1002         {
1003           indx = clustNumb[iclust] + 999;
1004           for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1) 
1005           {
1006             indx1 = clustNumb[iclust1];
1007             aijcluclu(indx,indx1) = aijcluclu(indx1,indx) = 0;
1008           }
1009         }
1010         
1011         // Update the remaining clusters couplings (subtract couplings from 
1012         // the used pads) - overflows excluded
1013         for (Int_t j = 0; j < npad; ++j) 
1014         {
1015           AliMUONPad* pad = cluster.Pad(j);
1016           //if ( pad->Status() != -2) continue;
1017           if ( pad->Status() != AliMUONClusterFinderMLEM::fgkModified) continue;
1018           for (Int_t iclust=0; iclust<nCoupled; ++iclust) 
1019           {
1020             indx = clustNumb[iclust];
1021             if (aijclupad(indx,j) < fgkCouplMin) continue;
1022             for (Int_t iclust1 = iclust+1; iclust1 < nCoupled; ++iclust1) 
1023             {
1024               indx1 = clustNumb[iclust1];
1025               if (aijclupad(indx1,j) < fgkCouplMin) continue;
1026               // Check this
1027               aijcluclu(indx,indx1) -= 
1028                 TMath::Sqrt (aijclupad(indx,j)*aijclupad(indx1,j));
1029               aijcluclu(indx1,indx) = aijcluclu(indx,indx1);
1030             }
1031           }
1032           //pad->SetStatus(-8);
1033           pad->SetStatus(AliMUONClusterFinderMLEM::fgkOver);
1034         } // for (Int_t j=0; j<npad;
1035       } // if (nCoupled > 3)
1036     } // while (nCoupled > 0)
1037   } // for (Int_t igroup=0; igroup<nclust;
1038   
1039   for (Int_t iclust = 0; iclust < nclust; ++iclust)
1040   {
1041     pix = clusters[iclust]; 
1042     pix->Clear();
1043     delete pix; 
1044   }
1045   delete [] clustNumb; 
1046   delete [] used; 
1047
1048 }
1049
1050 //_____________________________________________________________________________
1051 void 
1052 AliMUONClusterSplitterMLEM::Merge(const AliMUONCluster& cluster,
1053                                      Int_t nForFit, Int_t nCoupled, 
1054                                      Int_t *clustNumb, Int_t *clustFit, 
1055                                      TObjArray **clusters, 
1056                                      TMatrixD& aijcluclu, TMatrixD& aijclupad)
1057 {
1058   /// Merge the group of clusters with the one having the strongest coupling with them
1059   
1060   Int_t indx, indx1, npxclu, npxclu1, imax=0;
1061   TObjArray *pix, *pix1;
1062   Double_t couplMax;
1063   
1064   for (Int_t icl = 0; icl < nForFit; ++icl) 
1065   {
1066     indx = clustFit[icl];
1067     pix = clusters[indx];
1068     npxclu = pix->GetEntriesFast();
1069     couplMax = -1;
1070     for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1) 
1071     {
1072       indx1 = clustNumb[icl1];
1073       if (indx1 < 0) continue;
1074       if ( aijcluclu(indx,indx1) > couplMax) 
1075       {
1076         couplMax = aijcluclu(indx,indx1);
1077         imax = indx1;
1078       }
1079     } // for (Int_t icl1=0;
1080       // Add to it
1081     pix1 = clusters[imax];
1082     npxclu1 = pix1->GetEntriesFast();
1083     // Add pixels 
1084     for (Int_t i = 0; i < npxclu; ++i) 
1085     { 
1086       pix1->Add(pix->UncheckedAt(i)); 
1087       pix->RemoveAt(i); 
1088     }
1089     
1090     //Add cluster-to-cluster couplings
1091     for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1) 
1092     {
1093       indx1 = clustNumb[icl1];
1094       if (indx1 < 0 || indx1 == imax) continue;
1095       aijcluclu(indx1,imax) += aijcluclu(indx,indx1);
1096       aijcluclu(imax,indx1) = aijcluclu(indx1,imax);
1097     }
1098     aijcluclu(indx,imax) = aijcluclu(imax,indx) = 0;
1099     
1100     //Add cluster-to-pad couplings
1101     Int_t mult = cluster.Multiplicity();
1102     for (Int_t j = 0; j < mult; ++j) 
1103     {
1104       AliMUONPad* pad = cluster.Pad(j);
1105       //if ( pad->Status() < 0 && pad->Status() != -5 ) continue;// exclude used pads
1106       if ( pad->Status() != AliMUONClusterFinderMLEM::fgkZero) continue;// exclude used pads
1107         aijclupad(imax,j) += aijclupad(indx,j);
1108         aijclupad(indx,j) = 0;
1109     }
1110   } // for (Int_t icl=0; icl<nForFit;
1111 }
1112
1113
1114 //_____________________________________________________________________________
1115 Double_t 
1116 AliMUONClusterSplitterMLEM::MinGroupCoupl(Int_t nCoupled, Int_t *clustNumb, 
1117                                           TMatrixD& aijcluclu, Int_t *minGroup)
1118 {
1119   /// Find group of clusters with minimum coupling to all the others
1120   
1121   Int_t i123max = TMath::Min(3,nCoupled/2); 
1122   Int_t indx, indx1, indx2, indx3, nTot = 0;
1123   Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1124   
1125   for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1126     
1127     if (i123 == 1) {
1128       coupl1 = new Double_t [nCoupled];
1129       for (Int_t i = 0; i < nCoupled; ++i) coupl1[i] = 0;
1130     }
1131     else if (i123 == 2) {
1132       nTot = nCoupled*nCoupled;
1133       coupl2 = new Double_t [nTot];
1134       for (Int_t i = 0; i < nTot; ++i) coupl2[i] = 9999;
1135     } else {
1136       nTot = nTot*nCoupled;
1137       coupl3 = new Double_t [nTot];
1138       for (Int_t i = 0; i < nTot; ++i) coupl3[i] = 9999;
1139     } // else
1140     
1141     for (Int_t i = 0; i < nCoupled; ++i) {
1142       indx1 = clustNumb[i];
1143       for (Int_t j = i+1; j < nCoupled; ++j) {
1144         indx2 = clustNumb[j];
1145         if (i123 == 1) {
1146           coupl1[i] += aijcluclu(indx1,indx2);
1147           coupl1[j] += aijcluclu(indx1,indx2);
1148         } 
1149         else if (i123 == 2) {
1150           indx = i*nCoupled + j;
1151           coupl2[indx] = coupl1[i] + coupl1[j];
1152           coupl2[indx] -= 2 * (aijcluclu(indx1,indx2));
1153         } else {
1154           for (Int_t k = j+1; k < nCoupled; ++k) {
1155             indx3 = clustNumb[k];
1156             indx = i*nCoupled*nCoupled + j*nCoupled + k;
1157             coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1158             coupl3[indx] -= 2 * (aijcluclu(indx1,indx3)+aijcluclu(indx2,indx3));
1159           }
1160         } // else
1161       } // for (Int_t j=i+1;
1162     } // for (Int_t i=0;
1163   } // for (Int_t i123=1;
1164   
1165   // Find minimum coupling
1166   Double_t couplMin = 9999;
1167   Int_t locMin = 0;
1168   
1169   for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1170     if (i123 == 1) {
1171       locMin = TMath::LocMin(nCoupled, coupl1);
1172       couplMin = coupl1[locMin];
1173       minGroup[0] = locMin;
1174       delete [] coupl1;
1175     } 
1176     else if (i123 == 2) {
1177       locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
1178       if (coupl2[locMin] < couplMin) {
1179         couplMin = coupl2[locMin];
1180         minGroup[0] = locMin/nCoupled;
1181         minGroup[1] = locMin%nCoupled;
1182       }
1183       delete [] coupl2;
1184     } else {
1185       locMin = TMath::LocMin(nTot, coupl3);
1186       if (coupl3[locMin] < couplMin) {
1187         couplMin = coupl3[locMin];
1188         minGroup[0] = locMin/nCoupled/nCoupled;
1189         minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
1190         minGroup[2] = locMin%nCoupled;
1191       }
1192       delete [] coupl3; 
1193     } // else
1194   } // for (Int_t i123=1;
1195   return couplMin;
1196 }
1197
1198 //_____________________________________________________________________________
1199 Int_t 
1200 AliMUONClusterSplitterMLEM::SelectPad(const AliMUONCluster& cluster,
1201                                           Int_t nCoupled, Int_t nForFit, 
1202                                           Int_t *clustNumb, Int_t *clustFit, 
1203                                           TMatrixD& aijclupad)
1204 {
1205   /// Select pads for fit. If too many coupled clusters, find pads giving 
1206   /// the strongest coupling with the rest of clusters and exclude them from the fit.
1207   
1208   Int_t npad = cluster.Multiplicity();
1209   Double_t *padpix = 0;
1210   
1211   if (nCoupled > 3) 
1212   {
1213     padpix = new Double_t[npad];
1214     for (Int_t i = 0; i < npad; ++i) padpix[i] = 0.;
1215   }
1216   
1217   Int_t nOK = 0, indx, indx1;
1218   for (Int_t iclust = 0; iclust < nForFit; ++iclust)
1219   {
1220     indx = clustFit[iclust];
1221     for (Int_t j = 0; j < npad; ++j) 
1222     {
1223       if ( aijclupad(indx,j) < fgkCouplMin) continue;
1224       AliMUONPad* pad = cluster.Pad(j);
1225       /*
1226       if ( pad->Status() == -5 ) pad->SetStatus(-9); // flag overflow
1227       if ( pad->Status() < 0 ) continue; // exclude overflows and used pads
1228       if ( !pad->Status() ) 
1229       {
1230         pad->SetStatus(1);
1231         ++nOK; // pad to be used in fit
1232       }      
1233       */
1234       if ( pad->Status() != AliMUONClusterFinderMLEM::fgkZero 
1235            || pad->IsSaturated() ) continue; // used pads and overflows
1236       pad->SetStatus(AliMUONClusterFinderMLEM::fgkUseForFit);
1237       ++nOK; // pad to be used in fit
1238
1239       if (nCoupled > 3) 
1240       {
1241         // Check other clusters
1242         for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1) 
1243         {
1244           indx1 = clustNumb[iclust1];
1245           if (indx1 < 0) continue;
1246           if ( aijclupad(indx1,j) < fgkCouplMin ) continue;
1247           padpix[j] += aijclupad(indx1,j);
1248         }
1249       } // if (nCoupled > 3)
1250     } // for (Int_t j=0; j<npad;
1251   } // for (Int_t iclust=0; iclust<nForFit
1252   if (nCoupled < 4) return nOK;
1253   
1254   Double_t aaa = 0;
1255   for (Int_t j = 0; j < npad; ++j) 
1256   {
1257     if (padpix[j] < fgkCouplMin) continue;
1258     aaa += padpix[j];
1259     //cluster.Pad(j)->SetStatus(-1); // exclude pads with strong coupling to the other clusters
1260     cluster.Pad(j)->SetStatus(AliMUONClusterFinderMLEM::fgkCoupled); // exclude pads with strong coupling to the other clusters
1261     nOK--;
1262   }
1263   delete [] padpix; 
1264   return nOK;
1265 }
1266
1267 //_____________________________________________________________________________
1268 void 
1269 AliMUONClusterSplitterMLEM::UpdatePads(const AliMUONCluster& cluster,
1270                                           Int_t /*nfit*/, Double_t *par)
1271 {
1272   /// Subtract the fitted charges from pads with strong coupling
1273   
1274   Int_t indx, mult = cluster.Multiplicity(), iend = fNpar/3;
1275   Double_t charge, coef=0;
1276   
1277   for (Int_t j = 0; j < mult; ++j) 
1278   {
1279     AliMUONPad* pad = cluster.Pad(j);
1280     //if ( pad->Status() != -1 ) continue;
1281     if ( pad->Status() != AliMUONClusterFinderMLEM::fgkCoupled ) continue;
1282     if (fNpar != 0) 
1283     {
1284       charge = 0;
1285       for (Int_t i = 0; i <= iend; ++i) 
1286       { 
1287         // sum over hits
1288         indx = 3 * i;
1289         coef = Param2Coef(i, coef, par);
1290         charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
1291       }
1292       charge *= fQtot;
1293       pad->SetCharge(pad->Charge()-charge);
1294     } // if (fNpar != 0)
1295     
1296     //if (pad->Charge() > 6 /*fgkZeroSuppression*/) pad->SetStatus(0); 
1297     if (pad->Charge() > 6 /*fgkZeroSuppression*/) pad->SetStatus(AliMUONClusterFinderMLEM::fgkZero); 
1298     // return pad for further using // FIXME: remove usage of zerosuppression here
1299     else pad->SetStatus(AliMUONClusterFinderMLEM::fgkOver); // do not use anymore
1300     
1301   } // for (Int_t j=0;
1302 }  
1303
1304