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