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