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