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