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