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