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