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