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>
54 ClassImp(AliMUONClusterSplitterMLEM)
57 //const Double_t AliMUONClusterSplitterMLEM::fgkCouplMin = 1.e-3; // threshold on coupling
58 const Double_t AliMUONClusterSplitterMLEM::fgkCouplMin = 1.e-2; // threshold on coupling
60 //_____________________________________________________________________________
61 AliMUONClusterSplitterMLEM::AliMUONClusterSplitterMLEM(Int_t detElemId,
63 Double_t lowestPixelCharge,
64 Double_t lowestPadCharge,
65 Double_t lowestClusterCharge)
69 fDetElemId(detElemId),
74 fLowestPixelCharge(lowestPixelCharge),
75 fLowestPadCharge(lowestPadCharge),
76 fLowestClusterCharge(lowestClusterCharge)
80 AliMq::Station12Type stationType = AliMpDEManager::GetStation12Type(fDetElemId);
82 Float_t kx3 = AliMUONConstants::SqrtKx3();
83 Float_t ky3 = AliMUONConstants::SqrtKy3();
84 Float_t pitch = AliMUONConstants::Pitch();
86 if ( stationType == AliMq::kStation1 )
88 kx3 = AliMUONConstants::SqrtKx3St1();
89 ky3 = AliMUONConstants::SqrtKy3St1();
90 pitch = AliMUONConstants::PitchSt1();
93 fMathieson = new AliMUONMathieson;
95 fMathieson->SetPitch(pitch);
96 fMathieson->SetSqrtKx3AndDeriveKx2Kx4(kx3);
97 fMathieson->SetSqrtKy3AndDeriveKy2Ky4(ky3);
101 //_____________________________________________________________________________
102 AliMUONClusterSplitterMLEM::~AliMUONClusterSplitterMLEM()
109 //_____________________________________________________________________________
111 AliMUONClusterSplitterMLEM::AddBin(TH2 *mlem,
112 Int_t ic, Int_t jc, Int_t mode,
113 Bool_t *used, TObjArray *pix)
115 /// Add a bin to the cluster
117 Int_t nx = mlem->GetNbinsX();
118 Int_t ny = mlem->GetNbinsY();
119 Double_t cont1, cont = mlem->GetBinContent(mlem->GetBin(jc,ic));
120 AliMUONPad *pixPtr = 0;
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) {
125 if (i != ic && j != jc) continue;
126 if (used[(i-1)*nx+j-1]) continue;
127 cont1 = mlem->GetBinContent(mlem->GetBin(j,i));
128 if (mode && cont1 > cont) continue;
129 used[(i-1)*nx+j-1] = kTRUE;
130 if (cont1 < fLowestPixelCharge) continue;
131 if (pix) pix->Add(BinToPix(mlem,j,i));
133 pixPtr = new AliMUONPad (mlem->GetXaxis()->GetBinCenter(j),
134 mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
135 fPixArray->Add(pixPtr);
137 AddBin(mlem, i, j, mode, used, pix); // recursive call
142 //_____________________________________________________________________________
144 AliMUONClusterSplitterMLEM::AddCluster(Int_t ic, Int_t nclust,
146 Bool_t *used, Int_t *clustNumb, Int_t &nCoupled)
148 /// Add a cluster to the group of coupled clusters
150 for (Int_t i = 0; i < nclust; ++i) {
151 if (used[i]) continue;
152 if (aijcluclu(i,ic) < fgkCouplMin) continue;
154 clustNumb[nCoupled++] = i;
155 AddCluster(i, nclust, aijcluclu, used, clustNumb, nCoupled);
159 //_____________________________________________________________________________
161 AliMUONClusterSplitterMLEM::BinToPix(TH2 *mlem,
164 /// Translate histogram bin to pixel
166 Double_t yc = mlem->GetYaxis()->GetBinCenter(ic);
167 Double_t xc = mlem->GetXaxis()->GetBinCenter(jc);
169 Int_t nPix = fPixArray->GetEntriesFast();
170 AliMUONPad *pixPtr = NULL;
172 // Compare pixel and bin positions
173 for (Int_t i = 0; i < nPix; ++i) {
174 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
175 if (pixPtr->Charge() < fLowestPixelCharge) continue;
176 if (TMath::Abs(pixPtr->Coord(0)-xc)<1.e-4 && TMath::Abs(pixPtr->Coord(1)-yc)<1.e-4)
178 //return (TObject*) pixPtr;
182 AliError(Form(" Something wrong ??? %f %f ", xc, yc));
186 //_____________________________________________________________________________
188 AliMUONClusterSplitterMLEM::ChargeIntegration(Double_t x, Double_t y,
189 const AliMUONPad& pad)
191 /// Compute the Mathieson integral on pad area, assuming the center
192 /// of the Mathieson is at (x,y)
194 TVector2 lowerLeft(TVector2(x,y)-pad.Position()-pad.Dimensions());
195 TVector2 upperRight(lowerLeft + pad.Dimensions()*2.0);
197 return fMathieson->IntXY(lowerLeft.X(),lowerLeft.Y(),
198 upperRight.X(),upperRight.Y());
201 //_____________________________________________________________________________
203 AliMUONClusterSplitterMLEM::Fcn1(const AliMUONCluster& cluster,
204 Int_t & /*fNpar*/, Double_t * /*gin*/,
205 Double_t &f, Double_t *par, Int_t iflag)
207 /// Computes the functional to be minimized
210 Double_t charge, delta, coef=0, chi2=0, qTot = 0;
211 static Double_t qAver = 0;
213 Int_t mult = cluster.Multiplicity(), iend = fNpar / 3;
214 for (Int_t j = 0; j < mult; ++j)
216 AliMUONPad* pad = cluster.Pad(j);
217 //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
218 if ( pad->Status() != AliMUONClusterFinderMLEM::GetUseForFitFlag() ||
219 pad->Charge() == 0 ) continue;
221 if ( pad->IsReal() ) npads++; // exclude virtual pads
222 qTot += pad->Charge();
225 for (Int_t i = 0; i <= iend; ++i)
229 coef = Param2Coef(i, coef, par);
230 charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
233 delta = charge - pad->Charge();
235 delta /= pad->Charge();
238 if (iflag == 0 && npads) qAver = qTot / npads;
239 if (!npads && iflag==0)
241 AliError(Form("Got npads=0. Please check"));
246 //_____________________________________________________________________________
247 Double_t AliMUONClusterSplitterMLEM::Param2Coef(Int_t icand, Double_t coef, Double_t *par) const
249 /// Extract hit contribution scale factor from fit parameters
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.);
258 //_____________________________________________________________________________
260 AliMUONClusterSplitterMLEM::Fit(const AliMUONCluster& cluster,
261 Int_t iSimple, Int_t nfit,
262 const Int_t *clustFit, TObjArray **clusters,
264 TObjArray& clusterList, TH2 *mlem)
266 /// Steering function and fitting procedure for the fit of pad charge distribution
268 // AliDebug(2,Form("iSimple=%d nfit=%d",iSimple,nfit));
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);
275 // Number of pads to use and number of virtual pads
276 Int_t npads = 0, nVirtual = 0, nfit0 = nfit;
277 //cluster.Print("full");
278 Int_t mult = cluster.Multiplicity();
279 for (Int_t i = 0; i < mult; ++i )
281 AliMUONPad* pad = cluster.Pad(i);
282 if ( !pad->IsReal() ) ++nVirtual;
283 //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
284 if ( pad->Status() != AliMUONClusterFinderMLEM::GetUseForFitFlag() ) continue;
294 if (npads < 2) return 0;
296 // FIXME : AliWarning("Reconnect the following code for hit/track passing ?");
298 // Int_t tracks[3] = {-1, -1, -1};
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) {
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);
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));
327 } // for (Int_t cath=0;
330 // Get number of pads in X and Y
331 //const Int_t kStatusToTest(1);
332 const Int_t kStatusToTest(AliMUONClusterFinderMLEM::GetUseForFitFlag());
334 Long_t nofPads = cluster.NofPads(kStatusToTest);
335 Int_t nInX = AliMp::PairFirst(nofPads);
336 Int_t nInY = AliMp::PairSecond(nofPads);
340 for (Int_t j = 0; j < cluster.Multiplicity(); ++j) {
341 AliMUONPad *pad = cluster.Pad(j);
342 //if (pad->Status() == 1 && !pad->IsSaturated()) npadOK++;
343 if (pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() && !pad->IsSaturated()) npadOK++;
345 cout << " Number of pads to fit: " << npadOK << endl;
346 cout << " nInX and Y: " << nInX << " " << nInY << endl;
350 nfitMax = TMath::Min (nfitMax, (npads + 1) / 3);
352 if (((nInX < 3) && (nInY < 3)) || ((nInX == 3) && (nInY < 3)) || ((nInX < 3) && (nInY == 3))) nfitMax = 1; // not enough pads in each direction
354 if (nfit > nfitMax) nfit = nfitMax;
356 // Take cluster maxima as fitting seeds
360 Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8], qq = 0;
362 for ( int i = 0; i < 8; ++i ) errOk[i]=0.0;
364 Double_t xyseed[3][2], qseed[3], xyCand[3][2] = {{0},{0}}, sigCand[3][2] = {{0},{0}};
366 for (Int_t ifit = 1; ifit <= nfit0; ++ifit)
369 pix = clusters[clustFit[ifit-1]];
370 npxclu = pix->GetEntriesFast();
372 for (Int_t clu = 0; clu < npxclu; ++clu)
374 pixPtr = (AliMUONPad*) pix->UncheckedAt(clu);
375 cont = pixPtr->Charge();
380 xseed = pixPtr->Coord(0);
381 yseed = pixPtr->Coord(1);
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;
389 xyseed[ifit-1][0] = xseed;
390 xyseed[ifit-1][1] = yseed;
391 qseed[ifit-1] = cmax;
392 } // for (Int_t ifit=1;
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;
400 if (fDebug) cout << xyCand[0][0] << " " << xyCand[0][1] << " " << sigCand[0][0] << " " << sigCand[0][1] << endl;
402 Int_t nDof, maxSeed[3];//, nMax = 0;
404 if ( nfit0 < 0 || nfit0 > 3 ) {
405 AliErrorStream() << "Wrong nfit0 value: " << nfit0 << endl;
408 TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
410 Double_t step[3]={0.01,0.002,0.02}, fmin, chi2o = 9999, chi2n;
411 Double_t *gin = 0, func0, func1, param[8]={0}, step0[8]={0};
412 Double_t param0[2][8]={{0},{0}}, deriv[2][8]={{0},{0}};
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;
415 Int_t min, max, nCall = 0, nLoop, idMax = 0, nFail;
416 Double_t rad, dist[3] = {0};
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)
424 Int_t memory[8] = {0};
427 for (Int_t j = 0; j < fNpar; ++j)
432 parmin[fNpar] = 1E-9;
438 param[fNpar] = xyCand[0][0]; // take COG
442 param[fNpar] = xyseed[maxSeed[iseed]][0];
443 //param[fNpar] = fNpar==0 ? -16.1651 : -15.2761;
445 parmin[fNpar] = xmin;
446 parmax[fNpar++] = xmax;
449 param[fNpar] = xyCand[0][1]; // take COG
453 param[fNpar] = xyseed[maxSeed[iseed]][1];
454 //param[fNpar] = fNpar==1 ? -15.1737 : -15.8487;
456 parmin[fNpar] = ymin;
457 parmax[fNpar++] = ymax;
459 for (Int_t j = 0; j < fNpar; ++j)
461 step0[j] = shift[j] = step[j%3];
466 for (Int_t j = 0; j < fNpar; ++j)
472 for (Int_t j = 0; j < fNpar; ++j) cout << param[j] << " ";
477 min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
482 Fcn1(cluster,fNpar, gin, func0, param, iflag); nCall++;
484 //cout << " Func: " << func0 << endl;
487 for (Int_t j = 0; j < fNpar; ++j)
489 param0[max][j] = param[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
499 param[fNpar-1] -= delta[fNpar-1] / 10;
500 if (nCall > 2000) break;
502 min = func2[0] < func2[1] ? 0 : 1;
503 nFail = min == max ? 0 : nFail + 1;
505 stepMax = derMax = estim = 0;
506 for (Int_t j = 0; j < fNpar; ++j)
508 // Estimated distance to minimum
512 shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
514 else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3)
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))
521 shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
533 shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
537 Double_t es = TMath::Abs(shift[j]) / step0[j];
544 if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; //
546 // Failed to improve minimum
550 param[j] = param0[min][j];
551 if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j])
553 shift[j] = (shift[j] + shift0) / 2;
562 if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min])
564 shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
567 // Introduce step relaxation factor
570 scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
571 if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax)
573 shift[j] = TMath::Sign (shift0*scMax, shift[j]);
576 param[j] += shift[j];
577 // Check parameter limits
578 if (param[j] < parmin[j])
580 shift[j] = parmin[j] - param[j];
581 param[j] = parmin[j];
583 else if (param[j] > parmax[j])
585 shift[j] = parmax[j] - param[j];
586 param[j] = parmax[j];
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)
593 derMax = TMath::Abs (deriv[min][j]);
595 } // for (Int_t j=0; j<fNpar;
597 if (((estim < 1) && (derMax < 2)) || nLoop > 150) break; // minimum was found
601 // Check for small step
602 if (shift[idMax] == 0)
604 shift[idMax] = step0[idMax]/10;
605 param[idMax] += shift[idMax];
609 if (!memory[idMax] && derMax > 0.5 && nLoop > 10)
611 if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10)
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];
621 param[idMax] -= shift[idMax];
622 shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
623 param[idMax] += shift[idMax];
630 nDof = npads - fNpar + nVirtual;
633 if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
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 )
637 if (fNpar > 2 && (chi2n > chi2o || ((iseed == nfit-1)
638 && (chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) > chi2o))))
639 { fNpar -= 3; break; }
641 // Save parameters and errors
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];
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];
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]));
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;
676 cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
677 if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
681 for (Int_t i = 0; i < fNpar; ++i) {
682 parOk[i] = param0[min][i];
686 parOk[i] = TMath::Max (parOk[i], parmin[i]);
687 parOk[i] = TMath::Min (parOk[i], parmax[i]);
691 if (fmin < 0.1) break; // !!!???
692 } // for (Int_t iseed=0;
695 for (Int_t i=0; i<fNpar; ++i) {
696 if (i == 4 || i == 7) {
697 if ((i == 7) || ((i == 4) && (fNpar < 7))) cout << parOk[i] << endl;
698 else cout << parOk[i] * (1-parOk[7]) << endl;
701 cout << parOk[i] << " " << errOk[i] << endl;
704 nfit = (fNpar + 1) / 3;
705 dist[0] = dist[1] = dist[2] = 0;
708 // Find distance to the nearest neighbour
709 dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
711 +(parOk[1]-parOk[3])*
712 (parOk[1]-parOk[3]));
714 dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
716 +(parOk[1]-parOk[6])*
717 (parOk[1]-parOk[6]));
718 rad = TMath::Sqrt ((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;
731 if (iSimple) fnCoupled = 0;
732 for (Int_t j = 0; j < nfit; ++j) {
734 coef = Param2Coef(j, coef, parOk);
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*/)
743 if ( coef*fQtot >= fLowestClusterCharge )
745 //AZ AliMUONCluster* cluster1 = new AliMUONCluster();
746 AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
748 cluster1->SetCharge(coef*fQtot,coef*fQtot);
749 cluster1->SetPosition(TVector2(parOk[indx],parOk[indx+1]),TVector2(sigCand[0][0],sigCand[0][1]));
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;
756 cluster1->SetChi2(dist[idx]);
758 // FIXME: we miss some information in this cluster, as compared to
759 // the original AddRawCluster code.
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()));
765 clusterList.Add(cluster1);
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
782 //_____________________________________________________________________________
784 AliMUONClusterSplitterMLEM::Split(const AliMUONCluster& cluster,
785 TH2 *mlem, Double_t *coef,
786 TObjArray& clusterList)
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)
792 Int_t nx = mlem->GetNbinsX();
793 Int_t ny = mlem->GetNbinsY();
794 Int_t nPix = fPixArray->GetEntriesFast();
797 Int_t nclust = 0, indx, indx1, nxy = ny * nx;
798 Bool_t *used = new Bool_t[nxy];
800 for (Int_t j = 0; j < nxy; ++j) used[j] = kFALSE;
802 TObjArray *clusters[200]={0};
805 // Find clusters of histogram bins (easier to work in 2-D space)
806 for (Int_t i = 1; i <= ny; ++i)
808 for (Int_t j = 1; j <= nx; ++j)
810 indx = (i-1)*nx + j - 1;
811 if (used[indx]) continue;
812 cont = mlem->GetBinContent(mlem->GetBin(j,i));
813 if (cont < fLowestPixelCharge) continue;
814 pix = new TObjArray(20);
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++) {
821 } // for (Int_t i=1; i<=ny;
822 if (fDebug) cout << nclust << endl;
825 // Compute couplings between clusters and clusters to pads
826 Int_t npad = cluster.Multiplicity();
828 // Exclude pads with overflows
830 for (Int_t j = 0; j < npad; ++j)
832 AliMUONPad* pad = cluster.Pad(j);
833 if ( pad->IsSaturated() )
844 // Compute couplings of clusters to pads (including overflows)
845 TMatrixD aijclupad(nclust,npad);
848 for (Int_t iclust = 0; iclust < nclust; ++iclust)
850 pix = clusters[iclust];
851 npxclu = pix->GetEntriesFast();
852 for (Int_t i = 0; i < npxclu; ++i)
854 indx = fPixArray->IndexOf(pix->UncheckedAt(i));
855 for (Int_t j = 0; j < npad; ++j)
857 //AliMUONPad* pad = cluster.Pad(j);
858 //if ( pad->Status() < 0 && pad->Status() != -5) continue;
859 if (coef[j*nPix+indx] < fgkCouplMin) continue;
860 aijclupad(iclust,j) += coef[j*nPix+indx];
865 // Compute couplings between clusters (exclude overflows)
866 TMatrixD aijcluclu(nclust,nclust);
868 for (Int_t iclust = 0; iclust < nclust; ++iclust)
870 for (Int_t j = 0; j < npad; ++j)
873 //if ( cluster.Pad(j)->Status() < 0) continue;
874 if ( cluster.Pad(j)->IsSaturated()) continue;
875 if (aijclupad(iclust,j) < fgkCouplMin) continue;
876 for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++)
878 if (aijclupad(iclust1,j) < fgkCouplMin) continue;
879 aijcluclu(iclust,iclust1) +=
880 TMath::Sqrt (aijclupad(iclust,j)*aijclupad(iclust1,j));
884 for (Int_t iclust = 0; iclust < nclust; ++iclust)
886 for (Int_t iclust1 = iclust+1; iclust1 < nclust; ++iclust1)
888 aijcluclu(iclust1,iclust) = aijcluclu(iclust,iclust1);
892 if (fDebug && nclust > 1) aijcluclu.Print();
894 // Find groups of coupled clusters
895 used = new Bool_t[nclust];
896 for (Int_t j = 0; j < nclust; ++j) used[j] = kFALSE;
898 Int_t *clustNumb = new Int_t[nclust];
899 Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
901 Double_t parOk[8] = {0}; //AZ
903 for (Int_t igroup = 0; igroup < nclust; ++igroup)
905 if (used[igroup]) continue;
906 used[igroup] = kTRUE;
907 clustNumb[0] = igroup;
909 // Find group of coupled clusters
910 AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
913 cout << " nCoupled: " << nCoupled << endl;
914 for (Int_t i=0; i<nCoupled; ++i) cout << clustNumb[i] << " "; cout << endl;
917 fnCoupled = nCoupled;
924 for (Int_t i = 0; i < nCoupled; ++i) clustFit[i] = clustNumb[i];
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
931 for (Int_t j = 0; j < 3; ++j) minGroup[j] = -1;
932 Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
934 // Flag clusters for fit
936 while (nForFit < 3 && minGroup[nForFit] >= 0)
938 if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
939 clustFit[nForFit] = clustNumb[minGroup[nForFit]];
940 clustNumb[minGroup[nForFit]] -= 999;
943 if (fDebug) cout << " nForFit " << nForFit << " " << coupl << endl;
946 // Select pads for fit.
947 if (SelectPad(cluster,nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1)
950 for (Int_t j = 0; j < npad; ++j)
952 AliMUONPad* pad = cluster.Pad(j);
953 //if ( pad->Status()==1 ) pad->SetStatus(0);
954 //if ( pad->Status()==-9) pad->SetStatus(-5);
955 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() ||
956 pad->Status() == AliMUONClusterFinderMLEM::GetCoupledFlag())
957 pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
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);
966 nfit = Fit(cluster,0, nForFit, clustFit, clusters, parOk, clusterList, mlem);
968 //cout << " (nfit == 0) " << fNpar << " " << cluster.Multiplicity() << endl;
969 fNpar = 0; // should be 0 by itself but just in case ...
973 // Subtract the fitted charges from pads with strong coupling and/or
974 // return pads for further use
975 UpdatePads(cluster,nfit, parOk);
978 for (Int_t j = 0; j < npad; ++j)
980 AliMUONPad* pad = cluster.Pad(j);
981 //if ( pad->Status()==1 ) pad->SetStatus(-2);
982 //if ( pad->Status()==-9) pad->SetStatus(-5);
983 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() )
984 pad->SetStatus(AliMUONClusterFinderMLEM::GetModifiedFlag());
987 // Sort the clusters (move to the right the used ones)
988 Int_t beg = 0, end = nCoupled - 1;
991 if (clustNumb[beg] >= 0) { ++beg; continue; }
992 for (Int_t j = end; j > beg; --j)
994 if (clustNumb[j] < 0) continue;
996 indx = clustNumb[beg];
997 clustNumb[beg] = clustNumb[j];
1004 nCoupled -= nForFit;
1007 // Remove couplings of used clusters
1008 for (Int_t iclust = nCoupled; iclust < nCoupled+nForFit; ++iclust)
1010 indx = clustNumb[iclust] + 999;
1011 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1013 indx1 = clustNumb[iclust1];
1014 aijcluclu(indx,indx1) = aijcluclu(indx1,indx) = 0;
1018 // Update the remaining clusters couplings (subtract couplings from
1019 // the used pads) - overflows excluded
1020 for (Int_t j = 0; j < npad; ++j)
1022 AliMUONPad* pad = cluster.Pad(j);
1023 //if ( pad->Status() != -2) continue;
1024 if ( pad->Status() != AliMUONClusterFinderMLEM::GetModifiedFlag()) continue;
1025 for (Int_t iclust=0; iclust<nCoupled; ++iclust)
1027 indx = clustNumb[iclust];
1028 if (aijclupad(indx,j) < fgkCouplMin) continue;
1029 for (Int_t iclust1 = iclust+1; iclust1 < nCoupled; ++iclust1)
1031 indx1 = clustNumb[iclust1];
1032 if (aijclupad(indx1,j) < fgkCouplMin) continue;
1034 aijcluclu(indx,indx1) -=
1035 TMath::Sqrt (aijclupad(indx,j)*aijclupad(indx1,j));
1036 aijcluclu(indx1,indx) = aijcluclu(indx,indx1);
1039 //pad->SetStatus(-8);
1040 pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag());
1041 } // for (Int_t j=0; j<npad;
1042 } // if (nCoupled > 3)
1043 } // while (nCoupled > 0)
1044 } // for (Int_t igroup=0; igroup<nclust;
1046 for (Int_t iclust = 0; iclust < nclust; ++iclust)
1048 pix = clusters[iclust];
1052 delete [] clustNumb;
1057 //_____________________________________________________________________________
1059 AliMUONClusterSplitterMLEM::Merge(const AliMUONCluster& cluster,
1060 Int_t nForFit, Int_t nCoupled,
1061 const Int_t *clustNumb, const Int_t *clustFit,
1062 TObjArray **clusters,
1063 TMatrixD& aijcluclu, TMatrixD& aijclupad)
1065 /// Merge the group of clusters with the one having the strongest coupling with them
1067 Int_t indx, indx1, npxclu, imax=0;
1068 TObjArray *pix, *pix1;
1071 for (Int_t icl = 0; icl < nForFit; ++icl)
1073 indx = clustFit[icl];
1074 pix = clusters[indx];
1075 npxclu = pix->GetEntriesFast();
1077 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
1079 indx1 = clustNumb[icl1];
1080 if (indx1 < 0) continue;
1081 if ( aijcluclu(indx,indx1) > couplMax)
1083 couplMax = aijcluclu(indx,indx1);
1086 } // for (Int_t icl1=0;
1088 pix1 = clusters[imax];
1090 for (Int_t i = 0; i < npxclu; ++i)
1092 pix1->Add(pix->UncheckedAt(i));
1096 //Add cluster-to-cluster couplings
1097 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
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);
1104 aijcluclu(indx,imax) = aijcluclu(imax,indx) = 0;
1106 //Add cluster-to-pad couplings
1107 Int_t mult = cluster.Multiplicity();
1108 for (Int_t j = 0; j < mult; ++j)
1110 AliMUONPad* pad = cluster.Pad(j);
1111 //if ( pad->Status() < 0 && pad->Status() != -5 ) continue;// exclude used pads
1112 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()) continue;// exclude used pads
1113 aijclupad(imax,j) += aijclupad(indx,j);
1114 aijclupad(indx,j) = 0;
1116 } // for (Int_t icl=0; icl<nForFit;
1120 //_____________________________________________________________________________
1122 AliMUONClusterSplitterMLEM::MinGroupCoupl(Int_t nCoupled, const Int_t *clustNumb,
1123 const TMatrixD& aijcluclu, Int_t *minGroup)
1125 /// Find group of clusters with minimum coupling to all the others
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;
1131 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1134 coupl1 = new Double_t [nCoupled];
1135 for (Int_t i = 0; i < nCoupled; ++i) coupl1[i] = 0;
1137 else if (i123 == 2) {
1138 nTot = nCoupled*nCoupled;
1139 coupl2 = new Double_t [nTot];
1140 for (Int_t i = 0; i < nTot; ++i) coupl2[i] = 9999;
1142 nTot = nTot*nCoupled;
1143 coupl3 = new Double_t [nTot];
1144 for (Int_t i = 0; i < nTot; ++i) coupl3[i] = 9999;
1147 for (Int_t i = 0; i < nCoupled; ++i) {
1148 indx1 = clustNumb[i];
1149 for (Int_t j = i+1; j < nCoupled; ++j) {
1150 indx2 = clustNumb[j];
1152 coupl1[i] += aijcluclu(indx1,indx2);
1153 coupl1[j] += aijcluclu(indx1,indx2);
1155 else if (i123 == 2) {
1156 indx = i*nCoupled + j;
1157 coupl2[indx] = coupl1[i] + coupl1[j];
1158 coupl2[indx] -= 2 * (aijcluclu(indx1,indx2));
1160 for (Int_t k = j+1; k < nCoupled; ++k) {
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));
1167 } // for (Int_t j=i+1;
1168 } // for (Int_t i=0;
1169 } // for (Int_t i123=1;
1171 // Find minimum coupling
1172 Double_t couplMin = 9999;
1175 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1177 locMin = TMath::LocMin(nCoupled, coupl1);
1178 couplMin = coupl1[locMin];
1179 minGroup[0] = locMin;
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;
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;
1200 } // for (Int_t i123=1;
1204 //_____________________________________________________________________________
1206 AliMUONClusterSplitterMLEM::SelectPad(const AliMUONCluster& cluster,
1207 Int_t nCoupled, Int_t nForFit,
1208 const Int_t *clustNumb, const Int_t *clustFit,
1209 const TMatrixD& aijclupad)
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.
1214 Int_t npad = cluster.Multiplicity();
1215 Double_t *padpix = 0;
1219 padpix = new Double_t[npad];
1220 for (Int_t i = 0; i < npad; ++i) padpix[i] = 0.;
1223 Int_t nOK = 0, indx, indx1;
1224 for (Int_t iclust = 0; iclust < nForFit; ++iclust)
1226 indx = clustFit[iclust];
1227 for (Int_t j = 0; j < npad; ++j)
1229 if ( aijclupad(indx,j) < fgkCouplMin) continue;
1230 AliMUONPad* pad = cluster.Pad(j);
1232 if ( pad->Status() == -5 ) pad->SetStatus(-9); // flag overflow
1233 if ( pad->Status() < 0 ) continue; // exclude overflows and used pads
1234 if ( !pad->Status() )
1237 ++nOK; // pad to be used in fit
1240 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()
1241 || pad->IsSaturated() ) continue; // used pads and overflows
1242 pad->SetStatus(AliMUONClusterFinderMLEM::GetUseForFitFlag());
1243 ++nOK; // pad to be used in fit
1247 // Check other clusters
1248 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1250 indx1 = clustNumb[iclust1];
1251 if (indx1 < 0) continue;
1252 if ( aijclupad(indx1,j) < fgkCouplMin ) continue;
1253 padpix[j] += aijclupad(indx1,j);
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;
1261 for (Int_t j = 0; j < npad; ++j)
1263 if (padpix[j] < fgkCouplMin) continue;
1265 //cluster.Pad(j)->SetStatus(-1); // exclude pads with strong coupling to the other clusters
1266 cluster.Pad(j)->SetStatus(AliMUONClusterFinderMLEM::GetCoupledFlag()); // exclude pads with strong coupling to the other clusters
1273 //_____________________________________________________________________________
1275 AliMUONClusterSplitterMLEM::UpdatePads(const AliMUONCluster& cluster,
1276 Int_t /*nfit*/, Double_t *par)
1278 /// Subtract the fitted charges from pads with strong coupling
1280 Int_t indx, mult = cluster.Multiplicity(), iend = fNpar/3;
1281 Double_t charge, coef=0;
1283 for (Int_t j = 0; j < mult; ++j)
1285 AliMUONPad* pad = cluster.Pad(j);
1286 //if ( pad->Status() != -1 ) continue;
1287 if ( pad->Status() != AliMUONClusterFinderMLEM::GetCoupledFlag() ) continue;
1291 for (Int_t i = 0; i <= iend; ++i)
1295 coef = Param2Coef(i, coef, par);
1296 charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
1299 pad->SetCharge(pad->Charge()-charge);
1300 } // if (fNpar != 0)
1302 //if (pad->Charge() > 6 /*fgkZeroSuppression*/) pad->SetStatus(0);
1303 if (pad->Charge() > fLowestPadCharge) pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
1304 // return pad for further using // FIXME: remove usage of zerosuppression here
1305 else pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag()); // do not use anymore
1307 } // for (Int_t j=0;