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