]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterSplitterMLEM.cxx
Fixing NEGATIVE_RETURNS defect reported by Coverity
[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) qAver = qTot / npads;
237   f = chi2 / qAver;
238 }
239
240 //_____________________________________________________________________________
241 Double_t AliMUONClusterSplitterMLEM::Param2Coef(Int_t icand, Double_t coef, Double_t *par) const
242 {
243   /// Extract hit contribution scale factor from fit parameters
244   
245   if (fNpar == 2) return 1.;
246   if (fNpar == 5) return icand==0 ? par[2] : TMath::Max(1.-par[2],0.);
247   if (icand == 0) return par[2];
248   if (icand == 1) return TMath::Max((1.-par[2])*par[5], 0.);
249   return TMath::Max(1.-par[2]-coef,0.);
250 }
251
252 //_____________________________________________________________________________
253 Int_t 
254 AliMUONClusterSplitterMLEM::Fit(const AliMUONCluster& cluster,
255                                 Int_t iSimple, Int_t nfit, 
256                                 const Int_t *clustFit, TObjArray **clusters, 
257                                 Double_t *parOk,
258                                 TObjArray& clusterList, TH2 *mlem)
259 {
260   /// Steering function and fitting procedure for the fit of pad charge distribution
261   
262   //  AliDebug(2,Form("iSimple=%d nfit=%d",iSimple,nfit));
263   
264   Double_t xmin = mlem->GetXaxis()->GetXmin() - mlem->GetXaxis()->GetBinWidth(1);
265   Double_t xmax = mlem->GetXaxis()->GetXmax() + mlem->GetXaxis()->GetBinWidth(1);
266   Double_t ymin = mlem->GetYaxis()->GetXmin() - mlem->GetYaxis()->GetBinWidth(1);
267   Double_t ymax = mlem->GetYaxis()->GetXmax() + mlem->GetYaxis()->GetBinWidth(1);
268   Double_t xPad = 0, yPad = 99999;
269   
270   // Number of pads to use and number of virtual pads
271   Int_t npads = 0, nVirtual = 0, nfit0 = nfit;
272   //cluster.Print("full");
273   Int_t mult = cluster.Multiplicity();
274   for (Int_t i = 0; i < mult; ++i ) 
275   {
276     AliMUONPad* pad = cluster.Pad(i);
277     if ( !pad->IsReal() ) ++nVirtual;
278     //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
279     if ( pad->Status() != AliMUONClusterFinderMLEM::GetUseForFitFlag() ) continue;
280     if ( pad->IsReal() )
281     {
282       ++npads;
283       if (yPad > 9999) 
284       { 
285         xPad = pad->X();
286         yPad = pad->Y();
287       } 
288       else 
289       {
290         if (pad->DY() < pad->DX() ) 
291         {
292           yPad = pad->Y();
293         }
294         else 
295         {
296           xPad = pad->X();
297         }
298       }
299     }
300   }
301   
302   fNpar = 0;
303   fQtot = 0;
304   
305   if (npads < 2) return 0; 
306   
307   // FIXME : AliWarning("Reconnect the following code for hit/track passing ?");
308   
309   //  Int_t tracks[3] = {-1, -1, -1};
310   
311   /*
312    Int_t digit = 0;
313    AliMUONDigit *mdig = 0;
314    for (Int_t cath=0; cath<2; cath++) {  
315      for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
316        if (fPadIJ[0][i] != cath) continue;
317        if (fPadIJ[1][i] != 1) continue;
318        if (fXyq[3][i] < 0) continue; // exclude virtual pads
319        digit = TMath::Nint (fXyq[5][i]);
320        if (digit >= 0) mdig = fInput->Digit(cath,digit);
321        else mdig = fInput->Digit(TMath::Even(cath),-digit-1);
322        //if (!mdig) mdig = fInput->Digit(TMath::Even(cath),digit);
323        if (!mdig) continue; // protection for cluster display
324        if (mdig->Hit() >= 0) {
325          if (tracks[0] < 0) {
326            tracks[0] = mdig->Hit();
327            tracks[1] = mdig->Track(0);
328          } else if (mdig->Track(0) < tracks[1]) {
329            tracks[0] = mdig->Hit();
330            tracks[1] = mdig->Track(0);
331          }
332        }
333        if (mdig->Track(1) >= 0 && mdig->Track(1) != tracks[1]) {
334          if (tracks[2] < 0) tracks[2] = mdig->Track(1);
335          else tracks[2] = TMath::Min (tracks[2], mdig->Track(1));
336        }
337      } // for (Int_t i=0;
338   } // for (Int_t cath=0;
339    */
340   
341   // Get number of pads in X and Y 
342   //const Int_t kStatusToTest(1);
343   const Int_t kStatusToTest(AliMUONClusterFinderMLEM::GetUseForFitFlag());
344   
345   Long_t nofPads = cluster.NofPads(kStatusToTest);
346   Int_t nInX = AliMp::PairFirst(nofPads);
347   Int_t nInY = AliMp::PairSecond(nofPads);
348
349   if (fDebug) {
350     Int_t npadOK = 0;
351     for (Int_t j = 0; j < cluster.Multiplicity(); ++j) {
352       AliMUONPad *pad = cluster.Pad(j);
353       //if (pad->Status() == 1 && !pad->IsSaturated()) npadOK++;
354       if (pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() && !pad->IsSaturated()) npadOK++;
355     }
356     cout << " Number of pads to fit: " << npadOK << endl;
357     cout << " nInX and Y: " << nInX << " " << nInY << endl;
358   }
359   
360   Int_t nfitMax = 3; 
361   nfitMax = TMath::Min (nfitMax, (npads + 1) / 3);
362   if (nfitMax > 1) {
363     if (((nInX < 3) && (nInY < 3)) || ((nInX == 3) && (nInY < 3)) || ((nInX < 3) && (nInY == 3))) nfitMax = 1; // not enough pads in each direction
364   }
365   if (nfit > nfitMax) nfit = nfitMax;
366   
367   // Take cluster maxima as fitting seeds
368   TObjArray *pix;
369   AliMUONPad *pixPtr;
370   Int_t npxclu;
371   Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8], qq = 0;
372   Double_t xyseed[3][2], qseed[3], xyCand[3][2] = {{0},{0}}, sigCand[3][2] = {{0},{0}};
373   
374   for (Int_t ifit = 1; ifit <= nfit0; ++ifit) 
375   {
376     cmax = 0;
377     pix = clusters[clustFit[ifit-1]];
378     npxclu = pix->GetEntriesFast();
379     //qq = 0;
380     for (Int_t clu = 0; clu < npxclu; ++clu) 
381     {
382       pixPtr = (AliMUONPad*) pix->UncheckedAt(clu);
383       cont = pixPtr->Charge();
384       fQtot += cont;
385       if (cont > cmax) 
386       { 
387         cmax = cont; 
388         xseed = pixPtr->Coord(0);
389         yseed = pixPtr->Coord(1);
390       }
391       qq += cont;
392       xyCand[0][0] += pixPtr->Coord(0) * cont;
393       xyCand[0][1] += pixPtr->Coord(1) * cont;
394       sigCand[0][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
395       sigCand[0][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
396     }
397     xyseed[ifit-1][0] = xseed;
398     xyseed[ifit-1][1] = yseed;
399     qseed[ifit-1] = cmax;
400   } // for (Int_t ifit=1;
401   
402   xyCand[0][0] /= qq; // <x>
403   xyCand[0][1] /= qq; // <y>
404   sigCand[0][0] = sigCand[0][0]/qq - xyCand[0][0]*xyCand[0][0]; // <x^2> - <x>^2
405   sigCand[0][0] = sigCand[0][0] > 0 ? TMath::Sqrt (sigCand[0][0]) : 0;
406   sigCand[0][1] = sigCand[0][1]/qq - xyCand[0][1]*xyCand[0][1]; // <y^2> - <y>^2
407   sigCand[0][1] = sigCand[0][1] > 0 ? TMath::Sqrt (sigCand[0][1]) : 0;
408   if (fDebug) cout << xyCand[0][0] << " " << xyCand[0][1] << " " << sigCand[0][0] << " " << sigCand[0][1] << endl;
409   
410   Int_t nDof, maxSeed[3];//, nMax = 0;
411
412   if ( nfit0 < 0 || nfit0 > 3 ) {
413      AliErrorStream() << "Wrong nfit0 value: " << nfit0 << endl;
414      return nfit;
415   }   
416   TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
417     
418   Double_t step[3]={0.01,0.002,0.02}, fmin, chi2o = 9999, chi2n;
419   Double_t *gin = 0, func0, func1, param[8], step0[8];
420   Double_t param0[2][8]={{0},{0}}, deriv[2][8]={{0},{0}}; 
421   Double_t shift[8], stepMax, derMax, parmin[8], parmax[8], func2[2], shift0;
422   Double_t delta[8], scMax, dder[8], estim, shiftSave = 0;
423   Int_t min, max, nCall = 0, nLoop, idMax = 0, iestMax = 0, nFail;
424   Double_t rad, dist[3] = {0};
425     
426   // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is 
427   // lower, try 3-track (if number of pads is sufficient).
428   Int_t iflag = 0; // for the first call of fcn1
429   for (Int_t iseed = 0; iseed < nfit; ++iseed) 
430   {
431       
432     Int_t memory[8] = {0};
433     if (iseed) 
434     { 
435       for (Int_t j = 0; j < fNpar; ++j) 
436       {
437         param[j] = parOk[j]; 
438       }
439       param[fNpar] = 0.6;
440       parmin[fNpar] = 0; 
441       parmax[fNpar++] = 1; 
442     }
443       
444     if (nfit == 1) 
445     {
446       param[fNpar] = xyCand[0][0]; // take COG
447     }
448     else 
449     {
450       param[fNpar] = xyseed[maxSeed[iseed]][0];
451       //param[fNpar] = fNpar==0 ? -16.1651 : -15.2761; 
452     }
453     parmin[fNpar] = xmin; 
454     parmax[fNpar++] = xmax; 
455     if (nfit == 1) 
456     {
457       param[fNpar] = xyCand[0][1]; // take COG
458     }
459     else 
460     {
461       param[fNpar] = xyseed[maxSeed[iseed]][1];
462       //param[fNpar] = fNpar==1 ? -15.1737 : -15.8487;
463     }
464     parmin[fNpar] = ymin; 
465     parmax[fNpar++] = ymax; 
466
467     for (Int_t j = 0; j < fNpar; ++j) 
468     {
469       step0[j] = shift[j] = step[j%3];
470     }
471
472     if (iseed) 
473     { 
474       for (Int_t j = 0; j < fNpar; ++j) 
475       {
476         param0[1][j] = 0; 
477       }
478     }
479     if (fDebug) {
480       for (Int_t j = 0; j < fNpar; ++j) cout << param[j] << " "; 
481       cout << endl;
482     }
483       
484     // Try new algorithm
485     min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
486       
487     while (1) 
488     {
489       max = !min;
490       Fcn1(cluster,fNpar, gin, func0, param, iflag); nCall++;
491       iflag = 1;
492       //cout << " Func: " << func0 << endl;
493       
494       func2[max] = func0;
495       for (Int_t j = 0; j < fNpar; ++j) 
496       {
497         param0[max][j] = param[j];
498         delta[j] = step0[j];
499         param[j] += delta[j] / 10;
500         if (j > 0) param[j-1] -= delta[j-1] / 10;
501         Fcn1(cluster,fNpar, gin, func1, param, iflag); nCall++;
502         deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
503         //cout << j << " " << deriv[max][j] << endl;
504         dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) / 
505           (param0[0][j] - param0[1][j]) : 0; // second derivative
506       }
507       param[fNpar-1] -= delta[fNpar-1] / 10;
508       if (nCall > 2000) break;
509         
510       min = func2[0] < func2[1] ? 0 : 1;
511       nFail = min == max ? 0 : nFail + 1;
512         
513       stepMax = derMax = estim = 0;
514       for (Int_t j = 0; j < fNpar; ++j) 
515       { 
516         // Estimated distance to minimum
517         shift0 = shift[j];
518         if (nLoop == 1) 
519         {
520           shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
521         }
522         else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3) 
523         {
524           shift[j] = 0;
525         }
526         else if (((deriv[min][j]*deriv[!min][j] > 0) && (TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])))
527                  || (TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3) || (TMath::Abs(dder[j]) < 1.e-6)) 
528         {
529           shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
530           if (min == max) 
531           { 
532             if (memory[j] > 1) 
533             { 
534               shift[j] *= 2; 
535             } 
536             memory[j]++;
537           }
538         }
539         else 
540         {
541           shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
542           memory[j] = 0;
543         }
544           
545         Double_t es = TMath::Abs(shift[j]) / step0[j];
546         if (es > estim) 
547         { 
548           estim = es;
549           iestMax = j;
550         }
551           
552         // Too big step
553         if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; // 
554         
555         // Failed to improve minimum
556         if (min != max) 
557         {
558           memory[j] = 0;
559           param[j] = param0[min][j];
560           if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j]) 
561           {
562             shift[j] = (shift[j] + shift0) / 2;
563           }
564           else 
565           {
566             shift[j] /= -2;
567           }
568         } 
569           
570         // Too big step
571         if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min]) 
572         {
573           shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
574         }
575           
576         // Introduce step relaxation factor
577         if (memory[j] < 3) 
578         {
579           scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
580           if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax) 
581           {
582             shift[j] = TMath::Sign (shift0*scMax, shift[j]);
583           }
584         }
585         param[j] += shift[j]; 
586         // Check parameter limits
587         if (param[j] < parmin[j]) 
588         { 
589           shift[j] = parmin[j] - param[j]; 
590           param[j] = parmin[j]; 
591         } 
592         else if (param[j] > parmax[j]) 
593         {
594           shift[j] = parmax[j] - param[j];
595           param[j] = parmax[j];
596         }
597         //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
598         stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
599         if (TMath::Abs(deriv[min][j]) > derMax) 
600         {
601           idMax = j;
602           derMax = TMath::Abs (deriv[min][j]);
603         }
604       } // for (Int_t j=0; j<fNpar;
605         
606       if (((estim < 1) && (derMax < 2)) || nLoop > 150) break; // minimum was found
607         
608       nLoop++;
609         
610       // Check for small step
611       if (shift[idMax] == 0) 
612       { 
613         shift[idMax] = step0[idMax]/10; 
614         param[idMax] += shift[idMax]; 
615         continue; 
616       }
617         
618       if (!memory[idMax] && derMax > 0.5 && nLoop > 10) 
619       {
620         if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10) 
621         {
622           if (min == max) dder[idMax] = -dder[idMax];
623           shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10; 
624           param[idMax] += shift[idMax];
625           stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
626           if (min == max) shiftSave = shift[idMax];
627         }
628         if (nFail > 10) 
629         {
630           param[idMax] -= shift[idMax];
631           shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
632           param[idMax] += shift[idMax];
633         }
634       }      
635     } // while (1)
636       
637     fmin = func2[min];
638     
639     nDof = npads - fNpar + nVirtual;
640     if (!nDof) nDof++;
641     chi2n = fmin / nDof;
642     if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
643       
644     //if (fNpar > 2) cout << param0[min][fNpar-3] << " " << chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) << endl;
645     //if (chi2n*1.2+1.e-6 > chi2o ) 
646     if (fNpar > 2 && (chi2n > chi2o || ((iseed == nfit-1) 
647                                         && (chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) > chi2o)))) 
648       { fNpar -= 3; break; }
649       
650     // Save parameters and errors
651       
652     if (nInX == 1) {
653       // One pad per direction 
654       //for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5) param0[min][i] = xPad;
655       for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5) 
656         param0[min][i] = xyCand[0][0];
657     }
658     if (nInY == 1) {
659       // One pad per direction 
660       //for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6) param0[min][i] = yPad;
661       for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6) 
662         param0[min][i] = xyCand[0][1];
663     }
664       
665     /*
666       if (iseed > 0) {
667       // Find distance to the nearest neighbour
668       dist[0] = dist[1] = TMath::Sqrt ((param0[min][0]-param0[min][2])*
669       (param0[min][0]-param0[min][2])
670       +(param0[min][1]-param0[min][3])*
671       (param0[min][1]-param0[min][3]));
672       if (iseed > 1) {
673       dist[2] = TMath::Sqrt ((param0[min][0]-param0[min][5])*
674       (param0[min][0]-param0[min][5])
675       +(param0[min][1]-param0[min][6])*
676       (param0[min][1]-param0[min][6]));
677       rad = TMath::Sqrt ((param0[min][2]-param0[min][5])*
678       (param0[min][2]-param0[min][5])
679       +(param0[min][3]-param0[min][6])*
680       (param0[min][3]-param0[min][6]));
681       if (dist[2] < dist[0]) dist[0] = dist[2];
682       if (rad < dist[1]) dist[1] = rad;
683       if (rad < dist[2]) dist[2] = rad;
684       }
685       cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
686       if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
687       }
688     */
689       
690     for (Int_t i = 0; i < fNpar; ++i) {
691       parOk[i] = param0[min][i];
692       //errOk[i] = fmin;
693       errOk[i] = chi2n;
694       // Bounded params
695       parOk[i] = TMath::Max (parOk[i], parmin[i]);
696       parOk[i] = TMath::Min (parOk[i], parmax[i]);
697     }
698       
699     chi2o = chi2n;
700     if (fmin < 0.1) break; // !!!???
701   } // for (Int_t iseed=0; 
702    
703   if (fDebug) {
704     for (Int_t i=0; i<fNpar; ++i) {
705       if (i == 4 || i == 7) {
706         if ((i == 7) || ((i == 4) && (fNpar < 7))) cout << parOk[i] << endl;
707         else cout << parOk[i] * (1-parOk[7]) << endl;
708         continue;
709       }
710       cout << parOk[i] << " " << errOk[i] << endl;
711     }
712   }
713   nfit = (fNpar + 1) / 3;
714   dist[0] = dist[1] = dist[2] = 0;
715   
716   if (nfit > 1) {
717     // Find distance to the nearest neighbour
718     dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
719                                      (parOk[0]-parOk[2])
720                                      +(parOk[1]-parOk[3])*
721                                      (parOk[1]-parOk[3]));
722     if (nfit > 2) {
723       dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
724                              (parOk[0]-parOk[5])
725                              +(parOk[1]-parOk[6])*
726                              (parOk[1]-parOk[6]));
727       rad = TMath::Sqrt ((parOk[2]-parOk[5])*
728                          (parOk[2]-parOk[5])
729                          +(parOk[3]-parOk[6])*
730                          (parOk[3]-parOk[6]));
731       if (dist[2] < dist[0]) dist[0] = dist[2];
732       if (rad < dist[1]) dist[1] = rad;
733       if (rad < dist[2]) dist[2] = rad;
734     }
735   }
736     
737   Int_t indx;
738   
739   Double_t coef = 0;
740   if (iSimple) fnCoupled = 0;
741   for (Int_t j = 0; j < nfit; ++j) {
742     indx = 3 * j;
743     coef = Param2Coef(j, coef, parOk);
744       
745     //void AliMUONClusterFinderMLEM::AddRawCluster(Double_t x, Double_t y, 
746     //                                             Double_t qTot, Double_t fmin,
747     //                                             Int_t nfit, Int_t *tracks, 
748     //                                             Double_t /*sigx*/, 
749     //                                             Double_t /*sigy*/, 
750     //                                             Double_t /*dist*/)
751     
752     if ( coef*fQtot >= fLowestClusterCharge ) 
753     {
754       //AZ AliMUONCluster* cluster1 = new AliMUONCluster();
755       AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
756       
757       cluster1->SetCharge(coef*fQtot,coef*fQtot);
758       cluster1->SetPosition(TVector2(parOk[indx],parOk[indx+1]),TVector2(sigCand[0][0],sigCand[0][1]));
759       //cluster1->SetChi2(dist[TMath::LocMin(nfit,dist)]);
760       Int_t idx = TMath::LocMin(nfit,dist);
761       if ( idx < 0 || idx > 2 ) {
762         AliErrorStream() << "Wrong index value: " << idx << endl;
763         return nfit;
764       }  
765       cluster1->SetChi2(dist[idx]);
766       
767       // FIXME: we miss some information in this cluster, as compared to 
768       // the original AddRawCluster code.
769       
770       AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
771                       fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
772                       cluster1->Position().X(),cluster1->Position().Y()));
773         
774       clusterList.Add(cluster1);
775     }
776     //      AddRawCluster (parOk[indx], // double x
777     //                     parOk[indx+1], // double y
778     //                     coef*qTot, // double charge
779     //                     errOk[indx], // double fmin
780     //                     nfit0+10*nfit+100*nMax+10000*fnCoupled, // int nfit
781     //                     tracks, // int* tracks
782     //                     sigCand[0][0], // double sigx
783     //                     sigCand[0][1], // double sigy
784     //                     dist[TMath::LocMin(nfit,dist)] // double dist
785     //                     );
786   }
787   return nfit;
788 }  
789
790
791 //_____________________________________________________________________________
792 void
793 AliMUONClusterSplitterMLEM::Split(const AliMUONCluster& cluster,
794                                   TH2 *mlem, Double_t *coef,
795                                   TObjArray& clusterList)
796 {
797   /// The main steering function to work with clusters of pixels in anode
798   /// plane (find clusters, decouple them from each other, merge them (if
799   /// necessary), pick up coupled pads, call the fitting function)
800   
801   Int_t nx = mlem->GetNbinsX();
802   Int_t ny = mlem->GetNbinsY();
803   Int_t nPix = fPixArray->GetEntriesFast();
804   
805   Double_t cont;
806   Int_t nclust = 0, indx, indx1, nxy = ny * nx; 
807   Bool_t *used = new Bool_t[nxy];
808   
809   for (Int_t j = 0; j < nxy; ++j) used[j] = kFALSE; 
810   
811   TObjArray *clusters[200]={0};
812   TObjArray *pix;
813   
814   // Find clusters of histogram bins (easier to work in 2-D space)
815   for (Int_t i = 1; i <= ny; ++i) 
816   {
817     for (Int_t j = 1; j <= nx; ++j) 
818     {
819       indx = (i-1)*nx + j - 1;
820       if (used[indx]) continue;
821       cont = mlem->GetCellContent(j,i);
822       if (cont < fLowestPixelCharge) continue;
823       pix = new TObjArray(20);
824       used[indx] = 1;
825       pix->Add(BinToPix(mlem,j,i));
826       AddBin(mlem, i, j, 0, used, pix); // recursive call
827       if (nclust >= 200) AliFatal(" Too many clusters !!!");
828       clusters[nclust++] = pix;
829     } // for (Int_t j=1; j<=nx; j++) {
830   } // for (Int_t i=1; i<=ny;
831   if (fDebug) cout << nclust << endl;
832   delete [] used;
833   
834   // Compute couplings between clusters and clusters to pads
835   Int_t npad = cluster.Multiplicity();
836   
837   // Exclude pads with overflows
838   /*
839   for (Int_t j = 0; j < npad; ++j) 
840   {
841     AliMUONPad* pad = cluster.Pad(j);
842     if ( pad->IsSaturated() )
843     {
844       pad->SetStatus(-5); 
845     }
846     else 
847     {
848       pad->SetStatus(0);
849     }
850   }
851   */
852   
853   // Compute couplings of clusters to pads (including overflows)
854   TMatrixD aijclupad(nclust,npad);
855   aijclupad = 0;
856   Int_t npxclu;
857   for (Int_t iclust = 0; iclust < nclust; ++iclust) 
858   {
859     pix = clusters[iclust];
860     npxclu = pix->GetEntriesFast();
861     for (Int_t i = 0; i < npxclu; ++i) 
862     {
863       indx = fPixArray->IndexOf(pix->UncheckedAt(i));
864       for (Int_t j = 0; j < npad; ++j) 
865       {
866         //AliMUONPad* pad = cluster.Pad(j);
867         //if ( pad->Status() < 0 && pad->Status() != -5) continue;
868         if (coef[j*nPix+indx] < fgkCouplMin) continue;
869         aijclupad(iclust,j) += coef[j*nPix+indx];
870       }
871     }
872   }
873   
874   // Compute couplings between clusters (exclude overflows)
875   TMatrixD aijcluclu(nclust,nclust);
876   aijcluclu = 0;
877   for (Int_t iclust = 0; iclust < nclust; ++iclust) 
878   {
879     for (Int_t j = 0; j < npad; ++j) 
880     {
881       // Exclude overflows
882       //if ( cluster.Pad(j)->Status() < 0) continue;
883       if ( cluster.Pad(j)->IsSaturated()) continue;
884       if (aijclupad(iclust,j) < fgkCouplMin) continue;
885       for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) 
886       {
887         if (aijclupad(iclust1,j) < fgkCouplMin) continue;
888         aijcluclu(iclust,iclust1) += 
889           TMath::Sqrt (aijclupad(iclust,j)*aijclupad(iclust1,j));
890       }
891     }
892   }
893   for (Int_t iclust = 0; iclust < nclust; ++iclust) 
894   {
895     for (Int_t iclust1 = iclust+1; iclust1 < nclust; ++iclust1) 
896     {
897       aijcluclu(iclust1,iclust) = aijcluclu(iclust,iclust1);
898     }
899   }
900   
901   if (fDebug && nclust > 1) aijcluclu.Print();
902
903   // Find groups of coupled clusters
904   used = new Bool_t[nclust];
905   for (Int_t j = 0; j < nclust; ++j) used[j] = kFALSE;
906
907   Int_t *clustNumb = new Int_t[nclust];
908   Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
909   //Double_t parOk[8];
910   Double_t parOk[8] = {0}; //AZ
911   
912   for (Int_t igroup = 0; igroup < nclust; ++igroup) 
913   {
914     if (used[igroup]) continue;
915     used[igroup] = kTRUE;
916     clustNumb[0] = igroup;
917     nCoupled = 1;
918     // Find group of coupled clusters
919     AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
920     
921     if (fDebug) {                                                                      
922       cout << " nCoupled: " << nCoupled << endl;
923       for (Int_t i=0; i<nCoupled; ++i) cout << clustNumb[i] << " "; cout << endl;
924     }
925     
926     fnCoupled = nCoupled;
927     
928     while (nCoupled > 0) 
929     {
930       if (nCoupled < 4) 
931       {
932         nForFit = nCoupled;
933         for (Int_t i = 0; i < nCoupled; ++i) clustFit[i] = clustNumb[i];
934       } 
935       else 
936       {
937         // Too many coupled clusters to fit - try to decouple them
938         // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with 
939         // all the others in the group 
940         for (Int_t j = 0; j < 3; ++j) minGroup[j] = -1;
941         Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
942         
943         // Flag clusters for fit
944         nForFit = 0;
945         while (minGroup[nForFit] >= 0 && nForFit < 3)
946         {
947           if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
948           clustFit[nForFit] = clustNumb[minGroup[nForFit]];
949           clustNumb[minGroup[nForFit]] -= 999;
950           nForFit++;
951         }
952         if (fDebug) cout << " nForFit " << nForFit << " " << coupl << endl;
953       } // else
954       
955       // Select pads for fit. 
956       if (SelectPad(cluster,nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1) 
957       {
958         // Deselect pads
959         for (Int_t j = 0; j < npad; ++j)
960         {
961           AliMUONPad* pad = cluster.Pad(j);
962           //if ( pad->Status()==1 ) pad->SetStatus(0);
963           //if ( pad->Status()==-9) pad->SetStatus(-5);
964           if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() ||
965                pad->Status() == AliMUONClusterFinderMLEM::GetCoupledFlag()) 
966             pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
967         }
968         // Merge the failed cluster candidates (with too few pads to fit) with 
969         // the one with the strongest coupling
970         Merge(cluster,nForFit, nCoupled, clustNumb, clustFit, clusters, aijcluclu, aijclupad);
971       } 
972       else 
973       {
974         // Do the fit
975         nfit = Fit(cluster,0, nForFit, clustFit, clusters, parOk, clusterList, mlem);
976         if (nfit == 0) { 
977           //cout << " (nfit == 0) " << fNpar << " " << cluster.Multiplicity() << endl; 
978           fNpar = 0; // should be 0 by itself but just in case ...
979         }
980       }
981       
982       // Subtract the fitted charges from pads with strong coupling and/or
983       // return pads for further use
984       UpdatePads(cluster,nfit, parOk);
985       
986       // Mark used pads
987       for (Int_t j = 0; j < npad; ++j) 
988       {
989         AliMUONPad* pad = cluster.Pad(j);
990         //if ( pad->Status()==1 ) pad->SetStatus(-2);
991         //if ( pad->Status()==-9) pad->SetStatus(-5);
992         if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() ) 
993           pad->SetStatus(AliMUONClusterFinderMLEM::GetModifiedFlag());
994       }
995       
996       // Sort the clusters (move to the right the used ones)
997       Int_t beg = 0, end = nCoupled - 1;
998       while (beg < end) 
999       {
1000         if (clustNumb[beg] >= 0) { ++beg; continue; }
1001         for (Int_t j = end; j > beg; --j) 
1002         {
1003           if (clustNumb[j] < 0) continue;
1004           end = j - 1;
1005           indx = clustNumb[beg];
1006           clustNumb[beg] = clustNumb[j];
1007           clustNumb[j] = indx;
1008           break;
1009         }
1010         ++beg;
1011       }
1012       
1013       nCoupled -= nForFit;
1014       if (nCoupled > 3) 
1015       {
1016         // Remove couplings of used clusters
1017         for (Int_t iclust = nCoupled; iclust < nCoupled+nForFit; ++iclust) 
1018         {
1019           indx = clustNumb[iclust] + 999;
1020           for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1) 
1021           {
1022             indx1 = clustNumb[iclust1];
1023             aijcluclu(indx,indx1) = aijcluclu(indx1,indx) = 0;
1024           }
1025         }
1026         
1027         // Update the remaining clusters couplings (subtract couplings from 
1028         // the used pads) - overflows excluded
1029         for (Int_t j = 0; j < npad; ++j) 
1030         {
1031           AliMUONPad* pad = cluster.Pad(j);
1032           //if ( pad->Status() != -2) continue;
1033           if ( pad->Status() != AliMUONClusterFinderMLEM::GetModifiedFlag()) continue;
1034           for (Int_t iclust=0; iclust<nCoupled; ++iclust) 
1035           {
1036             indx = clustNumb[iclust];
1037             if (aijclupad(indx,j) < fgkCouplMin) continue;
1038             for (Int_t iclust1 = iclust+1; iclust1 < nCoupled; ++iclust1) 
1039             {
1040               indx1 = clustNumb[iclust1];
1041               if (aijclupad(indx1,j) < fgkCouplMin) continue;
1042               // Check this
1043               aijcluclu(indx,indx1) -= 
1044                 TMath::Sqrt (aijclupad(indx,j)*aijclupad(indx1,j));
1045               aijcluclu(indx1,indx) = aijcluclu(indx,indx1);
1046             }
1047           }
1048           //pad->SetStatus(-8);
1049           pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag());
1050         } // for (Int_t j=0; j<npad;
1051       } // if (nCoupled > 3)
1052     } // while (nCoupled > 0)
1053   } // for (Int_t igroup=0; igroup<nclust;
1054   
1055   for (Int_t iclust = 0; iclust < nclust; ++iclust)
1056   {
1057     pix = clusters[iclust]; 
1058     pix->Clear();
1059     delete pix; 
1060   }
1061   delete [] clustNumb; 
1062   delete [] used; 
1063
1064 }
1065
1066 //_____________________________________________________________________________
1067 void 
1068 AliMUONClusterSplitterMLEM::Merge(const AliMUONCluster& cluster,
1069                                      Int_t nForFit, Int_t nCoupled, 
1070                                      const Int_t *clustNumb, const Int_t *clustFit, 
1071                                      TObjArray **clusters, 
1072                                      TMatrixD& aijcluclu, TMatrixD& aijclupad)
1073 {
1074   /// Merge the group of clusters with the one having the strongest coupling with them
1075   
1076   Int_t indx, indx1, npxclu, npxclu1, imax=0;
1077   TObjArray *pix, *pix1;
1078   Double_t couplMax;
1079   
1080   for (Int_t icl = 0; icl < nForFit; ++icl) 
1081   {
1082     indx = clustFit[icl];
1083     pix = clusters[indx];
1084     npxclu = pix->GetEntriesFast();
1085     couplMax = -1;
1086     for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1) 
1087     {
1088       indx1 = clustNumb[icl1];
1089       if (indx1 < 0) continue;
1090       if ( aijcluclu(indx,indx1) > couplMax) 
1091       {
1092         couplMax = aijcluclu(indx,indx1);
1093         imax = indx1;
1094       }
1095     } // for (Int_t icl1=0;
1096       // Add to it
1097     pix1 = clusters[imax];
1098     npxclu1 = pix1->GetEntriesFast();
1099     // Add pixels 
1100     for (Int_t i = 0; i < npxclu; ++i) 
1101     { 
1102       pix1->Add(pix->UncheckedAt(i)); 
1103       pix->RemoveAt(i); 
1104     }
1105     
1106     //Add cluster-to-cluster couplings
1107     for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1) 
1108     {
1109       indx1 = clustNumb[icl1];
1110       if (indx1 < 0 || indx1 == imax) continue;
1111       aijcluclu(indx1,imax) += aijcluclu(indx,indx1);
1112       aijcluclu(imax,indx1) = aijcluclu(indx1,imax);
1113     }
1114     aijcluclu(indx,imax) = aijcluclu(imax,indx) = 0;
1115     
1116     //Add cluster-to-pad couplings
1117     Int_t mult = cluster.Multiplicity();
1118     for (Int_t j = 0; j < mult; ++j) 
1119     {
1120       AliMUONPad* pad = cluster.Pad(j);
1121       //if ( pad->Status() < 0 && pad->Status() != -5 ) continue;// exclude used pads
1122       if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()) continue;// exclude used pads
1123         aijclupad(imax,j) += aijclupad(indx,j);
1124         aijclupad(indx,j) = 0;
1125     }
1126   } // for (Int_t icl=0; icl<nForFit;
1127 }
1128
1129
1130 //_____________________________________________________________________________
1131 Double_t 
1132 AliMUONClusterSplitterMLEM::MinGroupCoupl(Int_t nCoupled, const Int_t *clustNumb, 
1133                                           const TMatrixD& aijcluclu, Int_t *minGroup)
1134 {
1135   /// Find group of clusters with minimum coupling to all the others
1136   
1137   Int_t i123max = TMath::Min(3,nCoupled/2); 
1138   Int_t indx, indx1, indx2, indx3, nTot = 0;
1139   Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1140   
1141   for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1142     
1143     if (i123 == 1) {
1144       coupl1 = new Double_t [nCoupled];
1145       for (Int_t i = 0; i < nCoupled; ++i) coupl1[i] = 0;
1146     }
1147     else if (i123 == 2) {
1148       nTot = nCoupled*nCoupled;
1149       coupl2 = new Double_t [nTot];
1150       for (Int_t i = 0; i < nTot; ++i) coupl2[i] = 9999;
1151     } else {
1152       nTot = nTot*nCoupled;
1153       coupl3 = new Double_t [nTot];
1154       for (Int_t i = 0; i < nTot; ++i) coupl3[i] = 9999;
1155     } // else
1156     
1157     for (Int_t i = 0; i < nCoupled; ++i) {
1158       indx1 = clustNumb[i];
1159       for (Int_t j = i+1; j < nCoupled; ++j) {
1160         indx2 = clustNumb[j];
1161         if (i123 == 1) {
1162           coupl1[i] += aijcluclu(indx1,indx2);
1163           coupl1[j] += aijcluclu(indx1,indx2);
1164         } 
1165         else if (i123 == 2) {
1166           indx = i*nCoupled + j;
1167           coupl2[indx] = coupl1[i] + coupl1[j];
1168           coupl2[indx] -= 2 * (aijcluclu(indx1,indx2));
1169         } else {
1170           for (Int_t k = j+1; k < nCoupled; ++k) {
1171             indx3 = clustNumb[k];
1172             indx = i*nCoupled*nCoupled + j*nCoupled + k;
1173             coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1174             coupl3[indx] -= 2 * (aijcluclu(indx1,indx3)+aijcluclu(indx2,indx3));
1175           }
1176         } // else
1177       } // for (Int_t j=i+1;
1178     } // for (Int_t i=0;
1179   } // for (Int_t i123=1;
1180   
1181   // Find minimum coupling
1182   Double_t couplMin = 9999;
1183   Int_t locMin = 0;
1184   
1185   for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1186     if (i123 == 1) {
1187       locMin = TMath::LocMin(nCoupled, coupl1);
1188       couplMin = coupl1[locMin];
1189       minGroup[0] = locMin;
1190       delete [] coupl1;
1191     } 
1192     else if (i123 == 2) {
1193       locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
1194       if (coupl2[locMin] < couplMin) {
1195         couplMin = coupl2[locMin];
1196         minGroup[0] = locMin/nCoupled;
1197         minGroup[1] = locMin%nCoupled;
1198       }
1199       delete [] coupl2;
1200     } else {
1201       locMin = TMath::LocMin(nTot, coupl3);
1202       if (coupl3[locMin] < couplMin) {
1203         couplMin = coupl3[locMin];
1204         minGroup[0] = locMin/nCoupled/nCoupled;
1205         minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
1206         minGroup[2] = locMin%nCoupled;
1207       }
1208       delete [] coupl3; 
1209     } // else
1210   } // for (Int_t i123=1;
1211   return couplMin;
1212 }
1213
1214 //_____________________________________________________________________________
1215 Int_t 
1216 AliMUONClusterSplitterMLEM::SelectPad(const AliMUONCluster& cluster,
1217                                           Int_t nCoupled, Int_t nForFit, 
1218                                           const Int_t *clustNumb, const Int_t *clustFit, 
1219                                           const TMatrixD& aijclupad)
1220 {
1221   /// Select pads for fit. If too many coupled clusters, find pads giving 
1222   /// the strongest coupling with the rest of clusters and exclude them from the fit.
1223   
1224   Int_t npad = cluster.Multiplicity();
1225   Double_t *padpix = 0;
1226   
1227   if (nCoupled > 3) 
1228   {
1229     padpix = new Double_t[npad];
1230     for (Int_t i = 0; i < npad; ++i) padpix[i] = 0.;
1231   }
1232   
1233   Int_t nOK = 0, indx, indx1;
1234   for (Int_t iclust = 0; iclust < nForFit; ++iclust)
1235   {
1236     indx = clustFit[iclust];
1237     for (Int_t j = 0; j < npad; ++j) 
1238     {
1239       if ( aijclupad(indx,j) < fgkCouplMin) continue;
1240       AliMUONPad* pad = cluster.Pad(j);
1241       /*
1242       if ( pad->Status() == -5 ) pad->SetStatus(-9); // flag overflow
1243       if ( pad->Status() < 0 ) continue; // exclude overflows and used pads
1244       if ( !pad->Status() ) 
1245       {
1246         pad->SetStatus(1);
1247         ++nOK; // pad to be used in fit
1248       }      
1249       */
1250       if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag() 
1251            || pad->IsSaturated() ) continue; // used pads and overflows
1252       pad->SetStatus(AliMUONClusterFinderMLEM::GetUseForFitFlag());
1253       ++nOK; // pad to be used in fit
1254
1255       if (nCoupled > 3) 
1256       {
1257         // Check other clusters
1258         for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1) 
1259         {
1260           indx1 = clustNumb[iclust1];
1261           if (indx1 < 0) continue;
1262           if ( aijclupad(indx1,j) < fgkCouplMin ) continue;
1263           padpix[j] += aijclupad(indx1,j);
1264         }
1265       } // if (nCoupled > 3)
1266     } // for (Int_t j=0; j<npad;
1267   } // for (Int_t iclust=0; iclust<nForFit
1268   if (nCoupled < 4) return nOK;
1269   
1270   Double_t aaa = 0;
1271   for (Int_t j = 0; j < npad; ++j) 
1272   {
1273     if (padpix[j] < fgkCouplMin) continue;
1274     aaa += padpix[j];
1275     //cluster.Pad(j)->SetStatus(-1); // exclude pads with strong coupling to the other clusters
1276     cluster.Pad(j)->SetStatus(AliMUONClusterFinderMLEM::GetCoupledFlag()); // exclude pads with strong coupling to the other clusters
1277     nOK--;
1278   }
1279   delete [] padpix; 
1280   return nOK;
1281 }
1282
1283 //_____________________________________________________________________________
1284 void 
1285 AliMUONClusterSplitterMLEM::UpdatePads(const AliMUONCluster& cluster,
1286                                           Int_t /*nfit*/, Double_t *par)
1287 {
1288   /// Subtract the fitted charges from pads with strong coupling
1289   
1290   Int_t indx, mult = cluster.Multiplicity(), iend = fNpar/3;
1291   Double_t charge, coef=0;
1292   
1293   for (Int_t j = 0; j < mult; ++j) 
1294   {
1295     AliMUONPad* pad = cluster.Pad(j);
1296     //if ( pad->Status() != -1 ) continue;
1297     if ( pad->Status() != AliMUONClusterFinderMLEM::GetCoupledFlag() ) continue;
1298     if (fNpar != 0) 
1299     {
1300       charge = 0;
1301       for (Int_t i = 0; i <= iend; ++i) 
1302       { 
1303         // sum over hits
1304         indx = 3 * i;
1305         coef = Param2Coef(i, coef, par);
1306         charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
1307       }
1308       charge *= fQtot;
1309       pad->SetCharge(pad->Charge()-charge);
1310     } // if (fNpar != 0)
1311     
1312     //if (pad->Charge() > 6 /*fgkZeroSuppression*/) pad->SetStatus(0); 
1313     if (pad->Charge() > fLowestPadCharge) pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
1314     // return pad for further using // FIXME: remove usage of zerosuppression here
1315     else pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag()); // do not use anymore
1316     
1317   } // for (Int_t j=0;
1318 }  
1319
1320