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