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