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