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