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