1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONClusterSplitterMLEM
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.
25 /// \author Laurent Aphecetche (for the "new" C++ structure) and
26 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
27 //-----------------------------------------------------------------------------
29 #include "AliMUONClusterSplitterMLEM.h"
30 #include "AliMUONClusterFinderMLEM.h" // for status flag constants
32 #include "AliMUONCluster.h"
33 #include "AliMUONPad.h"
34 #include "AliMUONPad.h"
35 #include "AliMUONConstants.h"
36 #include "AliMpDEManager.h"
37 #include "AliMUONMathieson.h"
39 #include "AliMpEncodePair.h"
43 #include <TClonesArray.h>
47 #include <TObjArray.h>
49 #include <Riostream.h>
52 ClassImp(AliMUONClusterSplitterMLEM)
55 //const Double_t AliMUONClusterSplitterMLEM::fgkCouplMin = 1.e-3; // threshold on coupling
56 const Double_t AliMUONClusterSplitterMLEM::fgkCouplMin = 1.e-2; // threshold on coupling
58 //_____________________________________________________________________________
59 AliMUONClusterSplitterMLEM::AliMUONClusterSplitterMLEM(Int_t detElemId,
61 Double_t lowestPixelCharge,
62 Double_t lowestPadCharge,
63 Double_t lowestClusterCharge)
67 fDetElemId(detElemId),
72 fLowestPixelCharge(lowestPixelCharge),
73 fLowestPadCharge(lowestPadCharge),
74 fLowestClusterCharge(lowestClusterCharge)
78 AliMq::Station12Type stationType = AliMpDEManager::GetStation12Type(fDetElemId);
80 Float_t kx3 = AliMUONConstants::SqrtKx3();
81 Float_t ky3 = AliMUONConstants::SqrtKy3();
82 Float_t pitch = AliMUONConstants::Pitch();
84 if ( stationType == AliMq::kStation1 )
86 kx3 = AliMUONConstants::SqrtKx3St1();
87 ky3 = AliMUONConstants::SqrtKy3St1();
88 pitch = AliMUONConstants::PitchSt1();
91 fMathieson = new AliMUONMathieson;
93 fMathieson->SetPitch(pitch);
94 fMathieson->SetSqrtKx3AndDeriveKx2Kx4(kx3);
95 fMathieson->SetSqrtKy3AndDeriveKy2Ky4(ky3);
99 //_____________________________________________________________________________
100 AliMUONClusterSplitterMLEM::~AliMUONClusterSplitterMLEM()
107 //_____________________________________________________________________________
109 AliMUONClusterSplitterMLEM::AddBin(TH2 *mlem,
110 Int_t ic, Int_t jc, Int_t mode,
111 Bool_t *used, TObjArray *pix)
113 /// Add a bin to the cluster
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;
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) {
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;
128 if (cont1 < fLowestPixelCharge) continue;
129 if (pix) pix->Add(BinToPix(mlem,j,i));
131 pixPtr = new AliMUONPad (mlem->GetXaxis()->GetBinCenter(j),
132 mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
133 fPixArray->Add(pixPtr);
135 AddBin(mlem, i, j, mode, used, pix); // recursive call
140 //_____________________________________________________________________________
142 AliMUONClusterSplitterMLEM::AddCluster(Int_t ic, Int_t nclust,
144 Bool_t *used, Int_t *clustNumb, Int_t &nCoupled)
146 /// Add a cluster to the group of coupled clusters
148 for (Int_t i = 0; i < nclust; ++i) {
149 if (used[i]) continue;
150 if (aijcluclu(i,ic) < fgkCouplMin) continue;
152 clustNumb[nCoupled++] = i;
153 AddCluster(i, nclust, aijcluclu, used, clustNumb, nCoupled);
157 //_____________________________________________________________________________
159 AliMUONClusterSplitterMLEM::BinToPix(TH2 *mlem,
162 /// Translate histogram bin to pixel
164 Double_t yc = mlem->GetYaxis()->GetBinCenter(ic);
165 Double_t xc = mlem->GetXaxis()->GetBinCenter(jc);
167 Int_t nPix = fPixArray->GetEntriesFast();
168 AliMUONPad *pixPtr = NULL;
170 // Compare pixel and bin positions
171 for (Int_t i = 0; i < nPix; ++i) {
172 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
173 if (pixPtr->Charge() < fLowestPixelCharge) continue;
174 if (TMath::Abs(pixPtr->Coord(0)-xc)<1.e-4 && TMath::Abs(pixPtr->Coord(1)-yc)<1.e-4)
176 //return (TObject*) pixPtr;
180 AliError(Form(" Something wrong ??? %f %f ", xc, yc));
184 //_____________________________________________________________________________
186 AliMUONClusterSplitterMLEM::ChargeIntegration(Double_t x, Double_t y,
187 const AliMUONPad& pad)
189 /// Compute the Mathieson integral on pad area, assuming the center
190 /// of the Mathieson is at (x,y)
192 TVector2 lowerLeft(TVector2(x,y)-pad.Position()-pad.Dimensions());
193 TVector2 upperRight(lowerLeft + pad.Dimensions()*2.0);
195 return fMathieson->IntXY(lowerLeft.X(),lowerLeft.Y(),
196 upperRight.X(),upperRight.Y());
199 //_____________________________________________________________________________
201 AliMUONClusterSplitterMLEM::Fcn1(const AliMUONCluster& cluster,
202 Int_t & /*fNpar*/, Double_t * /*gin*/,
203 Double_t &f, Double_t *par, Int_t iflag)
205 /// Computes the functional to be minimized
208 Double_t charge, delta, coef=0, chi2=0, qTot = 0;
209 static Double_t qAver = 0;
211 Int_t mult = cluster.Multiplicity(), iend = fNpar / 3;
212 for (Int_t j = 0; j < mult; ++j)
214 AliMUONPad* pad = cluster.Pad(j);
215 //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
216 if ( pad->Status() != AliMUONClusterFinderMLEM::GetUseForFitFlag() ||
217 pad->Charge() == 0 ) continue;
219 if ( pad->IsReal() ) npads++; // exclude virtual pads
220 qTot += pad->Charge();
223 for (Int_t i = 0; i <= iend; ++i)
227 coef = Param2Coef(i, coef, par);
228 charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
231 delta = charge - pad->Charge();
233 delta /= pad->Charge();
236 if (iflag == 0) qAver = qTot / npads;
240 //_____________________________________________________________________________
241 Double_t AliMUONClusterSplitterMLEM::Param2Coef(Int_t icand, Double_t coef, Double_t *par) const
243 /// Extract hit contribution scale factor from fit parameters
245 if (fNpar == 2) return 1.;
246 if (fNpar == 5) return icand==0 ? par[2] : TMath::Max(1.-par[2],0.);
247 if (icand == 0) return par[2];
248 if (icand == 1) return TMath::Max((1.-par[2])*par[5], 0.);
249 return TMath::Max(1.-par[2]-coef,0.);
252 //_____________________________________________________________________________
254 AliMUONClusterSplitterMLEM::Fit(const AliMUONCluster& cluster,
255 Int_t iSimple, Int_t nfit,
256 const Int_t *clustFit, TObjArray **clusters,
258 TObjArray& clusterList, TH2 *mlem)
260 /// Steering function and fitting procedure for the fit of pad charge distribution
262 // AliDebug(2,Form("iSimple=%d nfit=%d",iSimple,nfit));
264 Double_t xmin = mlem->GetXaxis()->GetXmin() - mlem->GetXaxis()->GetBinWidth(1);
265 Double_t xmax = mlem->GetXaxis()->GetXmax() + mlem->GetXaxis()->GetBinWidth(1);
266 Double_t ymin = mlem->GetYaxis()->GetXmin() - mlem->GetYaxis()->GetBinWidth(1);
267 Double_t ymax = mlem->GetYaxis()->GetXmax() + mlem->GetYaxis()->GetBinWidth(1);
268 Double_t xPad = 0, yPad = 99999;
270 // Number of pads to use and number of virtual pads
271 Int_t npads = 0, nVirtual = 0, nfit0 = nfit;
272 //cluster.Print("full");
273 Int_t mult = cluster.Multiplicity();
274 for (Int_t i = 0; i < mult; ++i )
276 AliMUONPad* pad = cluster.Pad(i);
277 if ( !pad->IsReal() ) ++nVirtual;
278 //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
279 if ( pad->Status() != AliMUONClusterFinderMLEM::GetUseForFitFlag() ) continue;
290 if (pad->DY() < pad->DX() )
305 if (npads < 2) return 0;
307 // FIXME : AliWarning("Reconnect the following code for hit/track passing ?");
309 // Int_t tracks[3] = {-1, -1, -1};
313 AliMUONDigit *mdig = 0;
314 for (Int_t cath=0; cath<2; cath++) {
315 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
316 if (fPadIJ[0][i] != cath) continue;
317 if (fPadIJ[1][i] != 1) continue;
318 if (fXyq[3][i] < 0) continue; // exclude virtual pads
319 digit = TMath::Nint (fXyq[5][i]);
320 if (digit >= 0) mdig = fInput->Digit(cath,digit);
321 else mdig = fInput->Digit(TMath::Even(cath),-digit-1);
322 //if (!mdig) mdig = fInput->Digit(TMath::Even(cath),digit);
323 if (!mdig) continue; // protection for cluster display
324 if (mdig->Hit() >= 0) {
326 tracks[0] = mdig->Hit();
327 tracks[1] = mdig->Track(0);
328 } else if (mdig->Track(0) < tracks[1]) {
329 tracks[0] = mdig->Hit();
330 tracks[1] = mdig->Track(0);
333 if (mdig->Track(1) >= 0 && mdig->Track(1) != tracks[1]) {
334 if (tracks[2] < 0) tracks[2] = mdig->Track(1);
335 else tracks[2] = TMath::Min (tracks[2], mdig->Track(1));
338 } // for (Int_t cath=0;
341 // Get number of pads in X and Y
342 //const Int_t kStatusToTest(1);
343 const Int_t kStatusToTest(AliMUONClusterFinderMLEM::GetUseForFitFlag());
345 Long_t nofPads = cluster.NofPads(kStatusToTest);
346 Int_t nInX = AliMp::PairFirst(nofPads);
347 Int_t nInY = AliMp::PairSecond(nofPads);
351 for (Int_t j = 0; j < cluster.Multiplicity(); ++j) {
352 AliMUONPad *pad = cluster.Pad(j);
353 //if (pad->Status() == 1 && !pad->IsSaturated()) npadOK++;
354 if (pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() && !pad->IsSaturated()) npadOK++;
356 cout << " Number of pads to fit: " << npadOK << endl;
357 cout << " nInX and Y: " << nInX << " " << nInY << endl;
361 nfitMax = TMath::Min (nfitMax, (npads + 1) / 3);
363 if (((nInX < 3) && (nInY < 3)) || ((nInX == 3) && (nInY < 3)) || ((nInX < 3) && (nInY == 3))) nfitMax = 1; // not enough pads in each direction
365 if (nfit > nfitMax) nfit = nfitMax;
367 // Take cluster maxima as fitting seeds
371 Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8], qq = 0;
372 Double_t xyseed[3][2], qseed[3], xyCand[3][2] = {{0},{0}}, sigCand[3][2] = {{0},{0}};
374 for (Int_t ifit = 1; ifit <= nfit0; ++ifit)
377 pix = clusters[clustFit[ifit-1]];
378 npxclu = pix->GetEntriesFast();
380 for (Int_t clu = 0; clu < npxclu; ++clu)
382 pixPtr = (AliMUONPad*) pix->UncheckedAt(clu);
383 cont = pixPtr->Charge();
388 xseed = pixPtr->Coord(0);
389 yseed = pixPtr->Coord(1);
392 xyCand[0][0] += pixPtr->Coord(0) * cont;
393 xyCand[0][1] += pixPtr->Coord(1) * cont;
394 sigCand[0][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
395 sigCand[0][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
397 xyseed[ifit-1][0] = xseed;
398 xyseed[ifit-1][1] = yseed;
399 qseed[ifit-1] = cmax;
400 } // for (Int_t ifit=1;
402 xyCand[0][0] /= qq; // <x>
403 xyCand[0][1] /= qq; // <y>
404 sigCand[0][0] = sigCand[0][0]/qq - xyCand[0][0]*xyCand[0][0]; // <x^2> - <x>^2
405 sigCand[0][0] = sigCand[0][0] > 0 ? TMath::Sqrt (sigCand[0][0]) : 0;
406 sigCand[0][1] = sigCand[0][1]/qq - xyCand[0][1]*xyCand[0][1]; // <y^2> - <y>^2
407 sigCand[0][1] = sigCand[0][1] > 0 ? TMath::Sqrt (sigCand[0][1]) : 0;
408 if (fDebug) cout << xyCand[0][0] << " " << xyCand[0][1] << " " << sigCand[0][0] << " " << sigCand[0][1] << endl;
410 Int_t nDof, maxSeed[3];//, nMax = 0;
412 if ( nfit0 < 0 || nfit0 > 3 ) {
413 AliErrorStream() << "Wrong nfit0 value: " << nfit0 << endl;
416 TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
418 Double_t step[3]={0.01,0.002,0.02}, fmin, chi2o = 9999, chi2n;
419 Double_t *gin = 0, func0, func1, param[8], step0[8];
420 Double_t param0[2][8]={{0},{0}}, deriv[2][8]={{0},{0}};
421 Double_t shift[8], stepMax, derMax, parmin[8], parmax[8], func2[2], shift0;
422 Double_t delta[8], scMax, dder[8], estim, shiftSave = 0;
423 Int_t min, max, nCall = 0, nLoop, idMax = 0, iestMax = 0, nFail;
424 Double_t rad, dist[3] = {0};
426 // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is
427 // lower, try 3-track (if number of pads is sufficient).
428 Int_t iflag = 0; // for the first call of fcn1
429 for (Int_t iseed = 0; iseed < nfit; ++iseed)
432 Int_t memory[8] = {0};
435 for (Int_t j = 0; j < fNpar; ++j)
446 param[fNpar] = xyCand[0][0]; // take COG
450 param[fNpar] = xyseed[maxSeed[iseed]][0];
451 //param[fNpar] = fNpar==0 ? -16.1651 : -15.2761;
453 parmin[fNpar] = xmin;
454 parmax[fNpar++] = xmax;
457 param[fNpar] = xyCand[0][1]; // take COG
461 param[fNpar] = xyseed[maxSeed[iseed]][1];
462 //param[fNpar] = fNpar==1 ? -15.1737 : -15.8487;
464 parmin[fNpar] = ymin;
465 parmax[fNpar++] = ymax;
467 for (Int_t j = 0; j < fNpar; ++j)
469 step0[j] = shift[j] = step[j%3];
474 for (Int_t j = 0; j < fNpar; ++j)
480 for (Int_t j = 0; j < fNpar; ++j) cout << param[j] << " ";
485 min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
490 Fcn1(cluster,fNpar, gin, func0, param, iflag); nCall++;
492 //cout << " Func: " << func0 << endl;
495 for (Int_t j = 0; j < fNpar; ++j)
497 param0[max][j] = param[j];
499 param[j] += delta[j] / 10;
500 if (j > 0) param[j-1] -= delta[j-1] / 10;
501 Fcn1(cluster,fNpar, gin, func1, param, iflag); nCall++;
502 deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
503 //cout << j << " " << deriv[max][j] << endl;
504 dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) /
505 (param0[0][j] - param0[1][j]) : 0; // second derivative
507 param[fNpar-1] -= delta[fNpar-1] / 10;
508 if (nCall > 2000) break;
510 min = func2[0] < func2[1] ? 0 : 1;
511 nFail = min == max ? 0 : nFail + 1;
513 stepMax = derMax = estim = 0;
514 for (Int_t j = 0; j < fNpar; ++j)
516 // Estimated distance to minimum
520 shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
522 else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3)
526 else if (((deriv[min][j]*deriv[!min][j] > 0) && (TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])))
527 || (TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3) || (TMath::Abs(dder[j]) < 1.e-6))
529 shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
541 shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
545 Double_t es = TMath::Abs(shift[j]) / step0[j];
553 if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; //
555 // Failed to improve minimum
559 param[j] = param0[min][j];
560 if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j])
562 shift[j] = (shift[j] + shift0) / 2;
571 if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min])
573 shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
576 // Introduce step relaxation factor
579 scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
580 if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax)
582 shift[j] = TMath::Sign (shift0*scMax, shift[j]);
585 param[j] += shift[j];
586 // Check parameter limits
587 if (param[j] < parmin[j])
589 shift[j] = parmin[j] - param[j];
590 param[j] = parmin[j];
592 else if (param[j] > parmax[j])
594 shift[j] = parmax[j] - param[j];
595 param[j] = parmax[j];
597 //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
598 stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
599 if (TMath::Abs(deriv[min][j]) > derMax)
602 derMax = TMath::Abs (deriv[min][j]);
604 } // for (Int_t j=0; j<fNpar;
606 if (((estim < 1) && (derMax < 2)) || nLoop > 150) break; // minimum was found
610 // Check for small step
611 if (shift[idMax] == 0)
613 shift[idMax] = step0[idMax]/10;
614 param[idMax] += shift[idMax];
618 if (!memory[idMax] && derMax > 0.5 && nLoop > 10)
620 if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10)
622 if (min == max) dder[idMax] = -dder[idMax];
623 shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10;
624 param[idMax] += shift[idMax];
625 stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
626 if (min == max) shiftSave = shift[idMax];
630 param[idMax] -= shift[idMax];
631 shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
632 param[idMax] += shift[idMax];
639 nDof = npads - fNpar + nVirtual;
642 if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
644 //if (fNpar > 2) cout << param0[min][fNpar-3] << " " << chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) << endl;
645 //if (chi2n*1.2+1.e-6 > chi2o )
646 if (fNpar > 2 && (chi2n > chi2o || ((iseed == nfit-1)
647 && (chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) > chi2o))))
648 { fNpar -= 3; break; }
650 // Save parameters and errors
653 // One pad per direction
654 //for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5) param0[min][i] = xPad;
655 for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5)
656 param0[min][i] = xyCand[0][0];
659 // One pad per direction
660 //for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6) param0[min][i] = yPad;
661 for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6)
662 param0[min][i] = xyCand[0][1];
667 // Find distance to the nearest neighbour
668 dist[0] = dist[1] = TMath::Sqrt ((param0[min][0]-param0[min][2])*
669 (param0[min][0]-param0[min][2])
670 +(param0[min][1]-param0[min][3])*
671 (param0[min][1]-param0[min][3]));
673 dist[2] = TMath::Sqrt ((param0[min][0]-param0[min][5])*
674 (param0[min][0]-param0[min][5])
675 +(param0[min][1]-param0[min][6])*
676 (param0[min][1]-param0[min][6]));
677 rad = TMath::Sqrt ((param0[min][2]-param0[min][5])*
678 (param0[min][2]-param0[min][5])
679 +(param0[min][3]-param0[min][6])*
680 (param0[min][3]-param0[min][6]));
681 if (dist[2] < dist[0]) dist[0] = dist[2];
682 if (rad < dist[1]) dist[1] = rad;
683 if (rad < dist[2]) dist[2] = rad;
685 cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
686 if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
690 for (Int_t i = 0; i < fNpar; ++i) {
691 parOk[i] = param0[min][i];
695 parOk[i] = TMath::Max (parOk[i], parmin[i]);
696 parOk[i] = TMath::Min (parOk[i], parmax[i]);
700 if (fmin < 0.1) break; // !!!???
701 } // for (Int_t iseed=0;
704 for (Int_t i=0; i<fNpar; ++i) {
705 if (i == 4 || i == 7) {
706 if ((i == 7) || ((i == 4) && (fNpar < 7))) cout << parOk[i] << endl;
707 else cout << parOk[i] * (1-parOk[7]) << endl;
710 cout << parOk[i] << " " << errOk[i] << endl;
713 nfit = (fNpar + 1) / 3;
714 dist[0] = dist[1] = dist[2] = 0;
717 // Find distance to the nearest neighbour
718 dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
720 +(parOk[1]-parOk[3])*
721 (parOk[1]-parOk[3]));
723 dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
725 +(parOk[1]-parOk[6])*
726 (parOk[1]-parOk[6]));
727 rad = TMath::Sqrt ((parOk[2]-parOk[5])*
729 +(parOk[3]-parOk[6])*
730 (parOk[3]-parOk[6]));
731 if (dist[2] < dist[0]) dist[0] = dist[2];
732 if (rad < dist[1]) dist[1] = rad;
733 if (rad < dist[2]) dist[2] = rad;
740 if (iSimple) fnCoupled = 0;
741 for (Int_t j = 0; j < nfit; ++j) {
743 coef = Param2Coef(j, coef, parOk);
745 //void AliMUONClusterFinderMLEM::AddRawCluster(Double_t x, Double_t y,
746 // Double_t qTot, Double_t fmin,
747 // Int_t nfit, Int_t *tracks,
748 // Double_t /*sigx*/,
749 // Double_t /*sigy*/,
750 // Double_t /*dist*/)
752 if ( coef*fQtot >= fLowestClusterCharge )
754 //AZ AliMUONCluster* cluster1 = new AliMUONCluster();
755 AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
757 cluster1->SetCharge(coef*fQtot,coef*fQtot);
758 cluster1->SetPosition(TVector2(parOk[indx],parOk[indx+1]),TVector2(sigCand[0][0],sigCand[0][1]));
759 cluster1->SetChi2(dist[TMath::LocMin(nfit,dist)]);
761 // FIXME: we miss some information in this cluster, as compared to
762 // the original AddRawCluster code.
764 AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
765 fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
766 cluster1->Position().X(),cluster1->Position().Y()));
768 clusterList.Add(cluster1);
770 // AddRawCluster (parOk[indx], // double x
771 // parOk[indx+1], // double y
772 // coef*qTot, // double charge
773 // errOk[indx], // double fmin
774 // nfit0+10*nfit+100*nMax+10000*fnCoupled, // int nfit
775 // tracks, // int* tracks
776 // sigCand[0][0], // double sigx
777 // sigCand[0][1], // double sigy
778 // dist[TMath::LocMin(nfit,dist)] // double dist
785 //_____________________________________________________________________________
787 AliMUONClusterSplitterMLEM::Split(const AliMUONCluster& cluster,
788 TH2 *mlem, Double_t *coef,
789 TObjArray& clusterList)
791 /// The main steering function to work with clusters of pixels in anode
792 /// plane (find clusters, decouple them from each other, merge them (if
793 /// necessary), pick up coupled pads, call the fitting function)
795 Int_t nx = mlem->GetNbinsX();
796 Int_t ny = mlem->GetNbinsY();
797 Int_t nPix = fPixArray->GetEntriesFast();
800 Int_t nclust = 0, indx, indx1, nxy = ny * nx;
801 Bool_t *used = new Bool_t[nxy];
803 for (Int_t j = 0; j < nxy; ++j) used[j] = kFALSE;
805 TObjArray *clusters[200]={0};
808 // Find clusters of histogram bins (easier to work in 2-D space)
809 for (Int_t i = 1; i <= ny; ++i)
811 for (Int_t j = 1; j <= nx; ++j)
813 indx = (i-1)*nx + j - 1;
814 if (used[indx]) continue;
815 cont = mlem->GetCellContent(j,i);
816 if (cont < fLowestPixelCharge) continue;
817 pix = new TObjArray(20);
819 pix->Add(BinToPix(mlem,j,i));
820 AddBin(mlem, i, j, 0, used, pix); // recursive call
821 if (nclust >= 200) AliFatal(" Too many clusters !!!");
822 clusters[nclust++] = pix;
823 } // for (Int_t j=1; j<=nx; j++) {
824 } // for (Int_t i=1; i<=ny;
825 if (fDebug) cout << nclust << endl;
828 // Compute couplings between clusters and clusters to pads
829 Int_t npad = cluster.Multiplicity();
831 // Exclude pads with overflows
833 for (Int_t j = 0; j < npad; ++j)
835 AliMUONPad* pad = cluster.Pad(j);
836 if ( pad->IsSaturated() )
847 // Compute couplings of clusters to pads (including overflows)
848 TMatrixD aijclupad(nclust,npad);
851 for (Int_t iclust = 0; iclust < nclust; ++iclust)
853 pix = clusters[iclust];
854 npxclu = pix->GetEntriesFast();
855 for (Int_t i = 0; i < npxclu; ++i)
857 indx = fPixArray->IndexOf(pix->UncheckedAt(i));
858 for (Int_t j = 0; j < npad; ++j)
860 //AliMUONPad* pad = cluster.Pad(j);
861 //if ( pad->Status() < 0 && pad->Status() != -5) continue;
862 if (coef[j*nPix+indx] < fgkCouplMin) continue;
863 aijclupad(iclust,j) += coef[j*nPix+indx];
868 // Compute couplings between clusters (exclude overflows)
869 TMatrixD aijcluclu(nclust,nclust);
871 for (Int_t iclust = 0; iclust < nclust; ++iclust)
873 for (Int_t j = 0; j < npad; ++j)
876 //if ( cluster.Pad(j)->Status() < 0) continue;
877 if ( cluster.Pad(j)->IsSaturated()) continue;
878 if (aijclupad(iclust,j) < fgkCouplMin) continue;
879 for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++)
881 if (aijclupad(iclust1,j) < fgkCouplMin) continue;
882 aijcluclu(iclust,iclust1) +=
883 TMath::Sqrt (aijclupad(iclust,j)*aijclupad(iclust1,j));
887 for (Int_t iclust = 0; iclust < nclust; ++iclust)
889 for (Int_t iclust1 = iclust+1; iclust1 < nclust; ++iclust1)
891 aijcluclu(iclust1,iclust) = aijcluclu(iclust,iclust1);
895 if (fDebug && nclust > 1) aijcluclu.Print();
897 // Find groups of coupled clusters
898 used = new Bool_t[nclust];
899 for (Int_t j = 0; j < nclust; ++j) used[j] = kFALSE;
901 Int_t *clustNumb = new Int_t[nclust];
902 Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
904 Double_t parOk[8] = {0}; //AZ
906 for (Int_t igroup = 0; igroup < nclust; ++igroup)
908 if (used[igroup]) continue;
909 used[igroup] = kTRUE;
910 clustNumb[0] = igroup;
912 // Find group of coupled clusters
913 AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
916 cout << " nCoupled: " << nCoupled << endl;
917 for (Int_t i=0; i<nCoupled; ++i) cout << clustNumb[i] << " "; cout << endl;
920 fnCoupled = nCoupled;
927 for (Int_t i = 0; i < nCoupled; ++i) clustFit[i] = clustNumb[i];
931 // Too many coupled clusters to fit - try to decouple them
932 // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with
933 // all the others in the group
934 for (Int_t j = 0; j < 3; ++j) minGroup[j] = -1;
935 Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
937 // Flag clusters for fit
939 while (minGroup[nForFit] >= 0 && nForFit < 3)
941 if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
942 clustFit[nForFit] = clustNumb[minGroup[nForFit]];
943 clustNumb[minGroup[nForFit]] -= 999;
946 if (fDebug) cout << " nForFit " << nForFit << " " << coupl << endl;
949 // Select pads for fit.
950 if (SelectPad(cluster,nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1)
953 for (Int_t j = 0; j < npad; ++j)
955 AliMUONPad* pad = cluster.Pad(j);
956 //if ( pad->Status()==1 ) pad->SetStatus(0);
957 //if ( pad->Status()==-9) pad->SetStatus(-5);
958 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() ||
959 pad->Status() == AliMUONClusterFinderMLEM::GetCoupledFlag())
960 pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
962 // Merge the failed cluster candidates (with too few pads to fit) with
963 // the one with the strongest coupling
964 Merge(cluster,nForFit, nCoupled, clustNumb, clustFit, clusters, aijcluclu, aijclupad);
969 nfit = Fit(cluster,0, nForFit, clustFit, clusters, parOk, clusterList, mlem);
971 //cout << " (nfit == 0) " << fNpar << " " << cluster.Multiplicity() << endl;
972 fNpar = 0; // should be 0 by itself but just in case ...
976 // Subtract the fitted charges from pads with strong coupling and/or
977 // return pads for further use
978 UpdatePads(cluster,nfit, parOk);
981 for (Int_t j = 0; j < npad; ++j)
983 AliMUONPad* pad = cluster.Pad(j);
984 //if ( pad->Status()==1 ) pad->SetStatus(-2);
985 //if ( pad->Status()==-9) pad->SetStatus(-5);
986 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() )
987 pad->SetStatus(AliMUONClusterFinderMLEM::GetModifiedFlag());
990 // Sort the clusters (move to the right the used ones)
991 Int_t beg = 0, end = nCoupled - 1;
994 if (clustNumb[beg] >= 0) { ++beg; continue; }
995 for (Int_t j = end; j > beg; --j)
997 if (clustNumb[j] < 0) continue;
999 indx = clustNumb[beg];
1000 clustNumb[beg] = clustNumb[j];
1001 clustNumb[j] = indx;
1007 nCoupled -= nForFit;
1010 // Remove couplings of used clusters
1011 for (Int_t iclust = nCoupled; iclust < nCoupled+nForFit; ++iclust)
1013 indx = clustNumb[iclust] + 999;
1014 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1016 indx1 = clustNumb[iclust1];
1017 aijcluclu(indx,indx1) = aijcluclu(indx1,indx) = 0;
1021 // Update the remaining clusters couplings (subtract couplings from
1022 // the used pads) - overflows excluded
1023 for (Int_t j = 0; j < npad; ++j)
1025 AliMUONPad* pad = cluster.Pad(j);
1026 //if ( pad->Status() != -2) continue;
1027 if ( pad->Status() != AliMUONClusterFinderMLEM::GetModifiedFlag()) continue;
1028 for (Int_t iclust=0; iclust<nCoupled; ++iclust)
1030 indx = clustNumb[iclust];
1031 if (aijclupad(indx,j) < fgkCouplMin) continue;
1032 for (Int_t iclust1 = iclust+1; iclust1 < nCoupled; ++iclust1)
1034 indx1 = clustNumb[iclust1];
1035 if (aijclupad(indx1,j) < fgkCouplMin) continue;
1037 aijcluclu(indx,indx1) -=
1038 TMath::Sqrt (aijclupad(indx,j)*aijclupad(indx1,j));
1039 aijcluclu(indx1,indx) = aijcluclu(indx,indx1);
1042 //pad->SetStatus(-8);
1043 pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag());
1044 } // for (Int_t j=0; j<npad;
1045 } // if (nCoupled > 3)
1046 } // while (nCoupled > 0)
1047 } // for (Int_t igroup=0; igroup<nclust;
1049 for (Int_t iclust = 0; iclust < nclust; ++iclust)
1051 pix = clusters[iclust];
1055 delete [] clustNumb;
1060 //_____________________________________________________________________________
1062 AliMUONClusterSplitterMLEM::Merge(const AliMUONCluster& cluster,
1063 Int_t nForFit, Int_t nCoupled,
1064 const Int_t *clustNumb, const Int_t *clustFit,
1065 TObjArray **clusters,
1066 TMatrixD& aijcluclu, TMatrixD& aijclupad)
1068 /// Merge the group of clusters with the one having the strongest coupling with them
1070 Int_t indx, indx1, npxclu, npxclu1, imax=0;
1071 TObjArray *pix, *pix1;
1074 for (Int_t icl = 0; icl < nForFit; ++icl)
1076 indx = clustFit[icl];
1077 pix = clusters[indx];
1078 npxclu = pix->GetEntriesFast();
1080 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
1082 indx1 = clustNumb[icl1];
1083 if (indx1 < 0) continue;
1084 if ( aijcluclu(indx,indx1) > couplMax)
1086 couplMax = aijcluclu(indx,indx1);
1089 } // for (Int_t icl1=0;
1091 pix1 = clusters[imax];
1092 npxclu1 = pix1->GetEntriesFast();
1094 for (Int_t i = 0; i < npxclu; ++i)
1096 pix1->Add(pix->UncheckedAt(i));
1100 //Add cluster-to-cluster couplings
1101 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
1103 indx1 = clustNumb[icl1];
1104 if (indx1 < 0 || indx1 == imax) continue;
1105 aijcluclu(indx1,imax) += aijcluclu(indx,indx1);
1106 aijcluclu(imax,indx1) = aijcluclu(indx1,imax);
1108 aijcluclu(indx,imax) = aijcluclu(imax,indx) = 0;
1110 //Add cluster-to-pad couplings
1111 Int_t mult = cluster.Multiplicity();
1112 for (Int_t j = 0; j < mult; ++j)
1114 AliMUONPad* pad = cluster.Pad(j);
1115 //if ( pad->Status() < 0 && pad->Status() != -5 ) continue;// exclude used pads
1116 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()) continue;// exclude used pads
1117 aijclupad(imax,j) += aijclupad(indx,j);
1118 aijclupad(indx,j) = 0;
1120 } // for (Int_t icl=0; icl<nForFit;
1124 //_____________________________________________________________________________
1126 AliMUONClusterSplitterMLEM::MinGroupCoupl(Int_t nCoupled, const Int_t *clustNumb,
1127 const TMatrixD& aijcluclu, Int_t *minGroup)
1129 /// Find group of clusters with minimum coupling to all the others
1131 Int_t i123max = TMath::Min(3,nCoupled/2);
1132 Int_t indx, indx1, indx2, indx3, nTot = 0;
1133 Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1135 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1138 coupl1 = new Double_t [nCoupled];
1139 for (Int_t i = 0; i < nCoupled; ++i) coupl1[i] = 0;
1141 else if (i123 == 2) {
1142 nTot = nCoupled*nCoupled;
1143 coupl2 = new Double_t [nTot];
1144 for (Int_t i = 0; i < nTot; ++i) coupl2[i] = 9999;
1146 nTot = nTot*nCoupled;
1147 coupl3 = new Double_t [nTot];
1148 for (Int_t i = 0; i < nTot; ++i) coupl3[i] = 9999;
1151 for (Int_t i = 0; i < nCoupled; ++i) {
1152 indx1 = clustNumb[i];
1153 for (Int_t j = i+1; j < nCoupled; ++j) {
1154 indx2 = clustNumb[j];
1156 coupl1[i] += aijcluclu(indx1,indx2);
1157 coupl1[j] += aijcluclu(indx1,indx2);
1159 else if (i123 == 2) {
1160 indx = i*nCoupled + j;
1161 coupl2[indx] = coupl1[i] + coupl1[j];
1162 coupl2[indx] -= 2 * (aijcluclu(indx1,indx2));
1164 for (Int_t k = j+1; k < nCoupled; ++k) {
1165 indx3 = clustNumb[k];
1166 indx = i*nCoupled*nCoupled + j*nCoupled + k;
1167 coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1168 coupl3[indx] -= 2 * (aijcluclu(indx1,indx3)+aijcluclu(indx2,indx3));
1171 } // for (Int_t j=i+1;
1172 } // for (Int_t i=0;
1173 } // for (Int_t i123=1;
1175 // Find minimum coupling
1176 Double_t couplMin = 9999;
1179 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1181 locMin = TMath::LocMin(nCoupled, coupl1);
1182 couplMin = coupl1[locMin];
1183 minGroup[0] = locMin;
1186 else if (i123 == 2) {
1187 locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
1188 if (coupl2[locMin] < couplMin) {
1189 couplMin = coupl2[locMin];
1190 minGroup[0] = locMin/nCoupled;
1191 minGroup[1] = locMin%nCoupled;
1195 locMin = TMath::LocMin(nTot, coupl3);
1196 if (coupl3[locMin] < couplMin) {
1197 couplMin = coupl3[locMin];
1198 minGroup[0] = locMin/nCoupled/nCoupled;
1199 minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
1200 minGroup[2] = locMin%nCoupled;
1204 } // for (Int_t i123=1;
1208 //_____________________________________________________________________________
1210 AliMUONClusterSplitterMLEM::SelectPad(const AliMUONCluster& cluster,
1211 Int_t nCoupled, Int_t nForFit,
1212 const Int_t *clustNumb, const Int_t *clustFit,
1213 const TMatrixD& aijclupad)
1215 /// Select pads for fit. If too many coupled clusters, find pads giving
1216 /// the strongest coupling with the rest of clusters and exclude them from the fit.
1218 Int_t npad = cluster.Multiplicity();
1219 Double_t *padpix = 0;
1223 padpix = new Double_t[npad];
1224 for (Int_t i = 0; i < npad; ++i) padpix[i] = 0.;
1227 Int_t nOK = 0, indx, indx1;
1228 for (Int_t iclust = 0; iclust < nForFit; ++iclust)
1230 indx = clustFit[iclust];
1231 for (Int_t j = 0; j < npad; ++j)
1233 if ( aijclupad(indx,j) < fgkCouplMin) continue;
1234 AliMUONPad* pad = cluster.Pad(j);
1236 if ( pad->Status() == -5 ) pad->SetStatus(-9); // flag overflow
1237 if ( pad->Status() < 0 ) continue; // exclude overflows and used pads
1238 if ( !pad->Status() )
1241 ++nOK; // pad to be used in fit
1244 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()
1245 || pad->IsSaturated() ) continue; // used pads and overflows
1246 pad->SetStatus(AliMUONClusterFinderMLEM::GetUseForFitFlag());
1247 ++nOK; // pad to be used in fit
1251 // Check other clusters
1252 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1254 indx1 = clustNumb[iclust1];
1255 if (indx1 < 0) continue;
1256 if ( aijclupad(indx1,j) < fgkCouplMin ) continue;
1257 padpix[j] += aijclupad(indx1,j);
1259 } // if (nCoupled > 3)
1260 } // for (Int_t j=0; j<npad;
1261 } // for (Int_t iclust=0; iclust<nForFit
1262 if (nCoupled < 4) return nOK;
1265 for (Int_t j = 0; j < npad; ++j)
1267 if (padpix[j] < fgkCouplMin) continue;
1269 //cluster.Pad(j)->SetStatus(-1); // exclude pads with strong coupling to the other clusters
1270 cluster.Pad(j)->SetStatus(AliMUONClusterFinderMLEM::GetCoupledFlag()); // exclude pads with strong coupling to the other clusters
1277 //_____________________________________________________________________________
1279 AliMUONClusterSplitterMLEM::UpdatePads(const AliMUONCluster& cluster,
1280 Int_t /*nfit*/, Double_t *par)
1282 /// Subtract the fitted charges from pads with strong coupling
1284 Int_t indx, mult = cluster.Multiplicity(), iend = fNpar/3;
1285 Double_t charge, coef=0;
1287 for (Int_t j = 0; j < mult; ++j)
1289 AliMUONPad* pad = cluster.Pad(j);
1290 //if ( pad->Status() != -1 ) continue;
1291 if ( pad->Status() != AliMUONClusterFinderMLEM::GetCoupledFlag() ) continue;
1295 for (Int_t i = 0; i <= iend; ++i)
1299 coef = Param2Coef(i, coef, par);
1300 charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
1303 pad->SetCharge(pad->Charge()-charge);
1304 } // if (fNpar != 0)
1306 //if (pad->Charge() > 6 /*fgkZeroSuppression*/) pad->SetStatus(0);
1307 if (pad->Charge() > fLowestPadCharge) pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
1308 // return pad for further using // FIXME: remove usage of zerosuppression here
1309 else pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag()); // do not use anymore
1311 } // for (Int_t j=0;