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;
373 for ( int i = 0; i < 8; ++i ) errOk[i]=0.0;
375 Double_t xyseed[3][2], qseed[3], xyCand[3][2] = {{0},{0}}, sigCand[3][2] = {{0},{0}};
377 for (Int_t ifit = 1; ifit <= nfit0; ++ifit)
380 pix = clusters[clustFit[ifit-1]];
381 npxclu = pix->GetEntriesFast();
383 for (Int_t clu = 0; clu < npxclu; ++clu)
385 pixPtr = (AliMUONPad*) pix->UncheckedAt(clu);
386 cont = pixPtr->Charge();
391 xseed = pixPtr->Coord(0);
392 yseed = pixPtr->Coord(1);
395 xyCand[0][0] += pixPtr->Coord(0) * cont;
396 xyCand[0][1] += pixPtr->Coord(1) * cont;
397 sigCand[0][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
398 sigCand[0][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
400 xyseed[ifit-1][0] = xseed;
401 xyseed[ifit-1][1] = yseed;
402 qseed[ifit-1] = cmax;
403 } // for (Int_t ifit=1;
405 xyCand[0][0] /= qq; // <x>
406 xyCand[0][1] /= qq; // <y>
407 sigCand[0][0] = sigCand[0][0]/qq - xyCand[0][0]*xyCand[0][0]; // <x^2> - <x>^2
408 sigCand[0][0] = sigCand[0][0] > 0 ? TMath::Sqrt (sigCand[0][0]) : 0;
409 sigCand[0][1] = sigCand[0][1]/qq - xyCand[0][1]*xyCand[0][1]; // <y^2> - <y>^2
410 sigCand[0][1] = sigCand[0][1] > 0 ? TMath::Sqrt (sigCand[0][1]) : 0;
411 if (fDebug) cout << xyCand[0][0] << " " << xyCand[0][1] << " " << sigCand[0][0] << " " << sigCand[0][1] << endl;
413 Int_t nDof, maxSeed[3];//, nMax = 0;
415 if ( nfit0 < 0 || nfit0 > 3 ) {
416 AliErrorStream() << "Wrong nfit0 value: " << nfit0 << endl;
419 TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
421 Double_t step[3]={0.01,0.002,0.02}, fmin, chi2o = 9999, chi2n;
422 Double_t *gin = 0, func0, func1, param[8], step0[8];
423 Double_t param0[2][8]={{0},{0}}, deriv[2][8]={{0},{0}};
424 Double_t shift[8], stepMax, derMax, parmin[8], parmax[8], func2[2], shift0;
425 Double_t delta[8], scMax, dder[8], estim, shiftSave = 0;
426 Int_t min, max, nCall = 0, nLoop, idMax = 0, iestMax = 0, nFail;
427 Double_t rad, dist[3] = {0};
429 // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is
430 // lower, try 3-track (if number of pads is sufficient).
431 Int_t iflag = 0; // for the first call of fcn1
432 for (Int_t iseed = 0; iseed < nfit; ++iseed)
435 Int_t memory[8] = {0};
438 for (Int_t j = 0; j < fNpar; ++j)
449 param[fNpar] = xyCand[0][0]; // take COG
453 param[fNpar] = xyseed[maxSeed[iseed]][0];
454 //param[fNpar] = fNpar==0 ? -16.1651 : -15.2761;
456 parmin[fNpar] = xmin;
457 parmax[fNpar++] = xmax;
460 param[fNpar] = xyCand[0][1]; // take COG
464 param[fNpar] = xyseed[maxSeed[iseed]][1];
465 //param[fNpar] = fNpar==1 ? -15.1737 : -15.8487;
467 parmin[fNpar] = ymin;
468 parmax[fNpar++] = ymax;
470 for (Int_t j = 0; j < fNpar; ++j)
472 step0[j] = shift[j] = step[j%3];
477 for (Int_t j = 0; j < fNpar; ++j)
483 for (Int_t j = 0; j < fNpar; ++j) cout << param[j] << " ";
488 min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
493 Fcn1(cluster,fNpar, gin, func0, param, iflag); nCall++;
495 //cout << " Func: " << func0 << endl;
498 for (Int_t j = 0; j < fNpar; ++j)
500 param0[max][j] = param[j];
502 param[j] += delta[j] / 10;
503 if (j > 0) param[j-1] -= delta[j-1] / 10;
504 Fcn1(cluster,fNpar, gin, func1, param, iflag); nCall++;
505 deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
506 //cout << j << " " << deriv[max][j] << endl;
507 dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) /
508 (param0[0][j] - param0[1][j]) : 0; // second derivative
510 param[fNpar-1] -= delta[fNpar-1] / 10;
511 if (nCall > 2000) break;
513 min = func2[0] < func2[1] ? 0 : 1;
514 nFail = min == max ? 0 : nFail + 1;
516 stepMax = derMax = estim = 0;
517 for (Int_t j = 0; j < fNpar; ++j)
519 // Estimated distance to minimum
523 shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
525 else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3)
529 else if (((deriv[min][j]*deriv[!min][j] > 0) && (TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])))
530 || (TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3) || (TMath::Abs(dder[j]) < 1.e-6))
532 shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
544 shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
548 Double_t es = TMath::Abs(shift[j]) / step0[j];
556 if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; //
558 // Failed to improve minimum
562 param[j] = param0[min][j];
563 if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j])
565 shift[j] = (shift[j] + shift0) / 2;
574 if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min])
576 shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
579 // Introduce step relaxation factor
582 scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
583 if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax)
585 shift[j] = TMath::Sign (shift0*scMax, shift[j]);
588 param[j] += shift[j];
589 // Check parameter limits
590 if (param[j] < parmin[j])
592 shift[j] = parmin[j] - param[j];
593 param[j] = parmin[j];
595 else if (param[j] > parmax[j])
597 shift[j] = parmax[j] - param[j];
598 param[j] = parmax[j];
600 //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
601 stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
602 if (TMath::Abs(deriv[min][j]) > derMax)
605 derMax = TMath::Abs (deriv[min][j]);
607 } // for (Int_t j=0; j<fNpar;
609 if (((estim < 1) && (derMax < 2)) || nLoop > 150) break; // minimum was found
613 // Check for small step
614 if (shift[idMax] == 0)
616 shift[idMax] = step0[idMax]/10;
617 param[idMax] += shift[idMax];
621 if (!memory[idMax] && derMax > 0.5 && nLoop > 10)
623 if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10)
625 if (min == max) dder[idMax] = -dder[idMax];
626 shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10;
627 param[idMax] += shift[idMax];
628 stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
629 if (min == max) shiftSave = shift[idMax];
633 param[idMax] -= shift[idMax];
634 shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
635 param[idMax] += shift[idMax];
642 nDof = npads - fNpar + nVirtual;
645 if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
647 //if (fNpar > 2) cout << param0[min][fNpar-3] << " " << chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) << endl;
648 //if (chi2n*1.2+1.e-6 > chi2o )
649 if (fNpar > 2 && (chi2n > chi2o || ((iseed == nfit-1)
650 && (chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) > chi2o))))
651 { fNpar -= 3; break; }
653 // Save parameters and errors
656 // One pad per direction
657 //for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5) param0[min][i] = xPad;
658 for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5)
659 param0[min][i] = xyCand[0][0];
662 // One pad per direction
663 //for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6) param0[min][i] = yPad;
664 for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6)
665 param0[min][i] = xyCand[0][1];
670 // Find distance to the nearest neighbour
671 dist[0] = dist[1] = TMath::Sqrt ((param0[min][0]-param0[min][2])*
672 (param0[min][0]-param0[min][2])
673 +(param0[min][1]-param0[min][3])*
674 (param0[min][1]-param0[min][3]));
676 dist[2] = TMath::Sqrt ((param0[min][0]-param0[min][5])*
677 (param0[min][0]-param0[min][5])
678 +(param0[min][1]-param0[min][6])*
679 (param0[min][1]-param0[min][6]));
680 rad = TMath::Sqrt ((param0[min][2]-param0[min][5])*
681 (param0[min][2]-param0[min][5])
682 +(param0[min][3]-param0[min][6])*
683 (param0[min][3]-param0[min][6]));
684 if (dist[2] < dist[0]) dist[0] = dist[2];
685 if (rad < dist[1]) dist[1] = rad;
686 if (rad < dist[2]) dist[2] = rad;
688 cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
689 if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
693 for (Int_t i = 0; i < fNpar; ++i) {
694 parOk[i] = param0[min][i];
698 parOk[i] = TMath::Max (parOk[i], parmin[i]);
699 parOk[i] = TMath::Min (parOk[i], parmax[i]);
703 if (fmin < 0.1) break; // !!!???
704 } // for (Int_t iseed=0;
707 for (Int_t i=0; i<fNpar; ++i) {
708 if (i == 4 || i == 7) {
709 if ((i == 7) || ((i == 4) && (fNpar < 7))) cout << parOk[i] << endl;
710 else cout << parOk[i] * (1-parOk[7]) << endl;
713 cout << parOk[i] << " " << errOk[i] << endl;
716 nfit = (fNpar + 1) / 3;
717 dist[0] = dist[1] = dist[2] = 0;
720 // Find distance to the nearest neighbour
721 dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
723 +(parOk[1]-parOk[3])*
724 (parOk[1]-parOk[3]));
726 dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
728 +(parOk[1]-parOk[6])*
729 (parOk[1]-parOk[6]));
730 rad = TMath::Sqrt ((parOk[2]-parOk[5])*
732 +(parOk[3]-parOk[6])*
733 (parOk[3]-parOk[6]));
734 if (dist[2] < dist[0]) dist[0] = dist[2];
735 if (rad < dist[1]) dist[1] = rad;
736 if (rad < dist[2]) dist[2] = rad;
743 if (iSimple) fnCoupled = 0;
744 for (Int_t j = 0; j < nfit; ++j) {
746 coef = Param2Coef(j, coef, parOk);
748 //void AliMUONClusterFinderMLEM::AddRawCluster(Double_t x, Double_t y,
749 // Double_t qTot, Double_t fmin,
750 // Int_t nfit, Int_t *tracks,
751 // Double_t /*sigx*/,
752 // Double_t /*sigy*/,
753 // Double_t /*dist*/)
755 if ( coef*fQtot >= fLowestClusterCharge )
757 //AZ AliMUONCluster* cluster1 = new AliMUONCluster();
758 AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
760 cluster1->SetCharge(coef*fQtot,coef*fQtot);
761 cluster1->SetPosition(TVector2(parOk[indx],parOk[indx+1]),TVector2(sigCand[0][0],sigCand[0][1]));
762 //cluster1->SetChi2(dist[TMath::LocMin(nfit,dist)]);
763 Int_t idx = TMath::LocMin(nfit,dist);
764 if ( idx < 0 || idx > 2 ) {
765 AliErrorStream() << "Wrong index value: " << idx << endl;
768 cluster1->SetChi2(dist[idx]);
770 // FIXME: we miss some information in this cluster, as compared to
771 // the original AddRawCluster code.
773 AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
774 fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
775 cluster1->Position().X(),cluster1->Position().Y()));
777 clusterList.Add(cluster1);
779 // AddRawCluster (parOk[indx], // double x
780 // parOk[indx+1], // double y
781 // coef*qTot, // double charge
782 // errOk[indx], // double fmin
783 // nfit0+10*nfit+100*nMax+10000*fnCoupled, // int nfit
784 // tracks, // int* tracks
785 // sigCand[0][0], // double sigx
786 // sigCand[0][1], // double sigy
787 // dist[TMath::LocMin(nfit,dist)] // double dist
794 //_____________________________________________________________________________
796 AliMUONClusterSplitterMLEM::Split(const AliMUONCluster& cluster,
797 TH2 *mlem, Double_t *coef,
798 TObjArray& clusterList)
800 /// The main steering function to work with clusters of pixels in anode
801 /// plane (find clusters, decouple them from each other, merge them (if
802 /// necessary), pick up coupled pads, call the fitting function)
804 Int_t nx = mlem->GetNbinsX();
805 Int_t ny = mlem->GetNbinsY();
806 Int_t nPix = fPixArray->GetEntriesFast();
809 Int_t nclust = 0, indx, indx1, nxy = ny * nx;
810 Bool_t *used = new Bool_t[nxy];
812 for (Int_t j = 0; j < nxy; ++j) used[j] = kFALSE;
814 TObjArray *clusters[200]={0};
817 // Find clusters of histogram bins (easier to work in 2-D space)
818 for (Int_t i = 1; i <= ny; ++i)
820 for (Int_t j = 1; j <= nx; ++j)
822 indx = (i-1)*nx + j - 1;
823 if (used[indx]) continue;
824 cont = mlem->GetCellContent(j,i);
825 if (cont < fLowestPixelCharge) continue;
826 pix = new TObjArray(20);
828 pix->Add(BinToPix(mlem,j,i));
829 AddBin(mlem, i, j, 0, used, pix); // recursive call
830 if (nclust >= 200) AliFatal(" Too many clusters !!!");
831 clusters[nclust++] = pix;
832 } // for (Int_t j=1; j<=nx; j++) {
833 } // for (Int_t i=1; i<=ny;
834 if (fDebug) cout << nclust << endl;
837 // Compute couplings between clusters and clusters to pads
838 Int_t npad = cluster.Multiplicity();
840 // Exclude pads with overflows
842 for (Int_t j = 0; j < npad; ++j)
844 AliMUONPad* pad = cluster.Pad(j);
845 if ( pad->IsSaturated() )
856 // Compute couplings of clusters to pads (including overflows)
857 TMatrixD aijclupad(nclust,npad);
860 for (Int_t iclust = 0; iclust < nclust; ++iclust)
862 pix = clusters[iclust];
863 npxclu = pix->GetEntriesFast();
864 for (Int_t i = 0; i < npxclu; ++i)
866 indx = fPixArray->IndexOf(pix->UncheckedAt(i));
867 for (Int_t j = 0; j < npad; ++j)
869 //AliMUONPad* pad = cluster.Pad(j);
870 //if ( pad->Status() < 0 && pad->Status() != -5) continue;
871 if (coef[j*nPix+indx] < fgkCouplMin) continue;
872 aijclupad(iclust,j) += coef[j*nPix+indx];
877 // Compute couplings between clusters (exclude overflows)
878 TMatrixD aijcluclu(nclust,nclust);
880 for (Int_t iclust = 0; iclust < nclust; ++iclust)
882 for (Int_t j = 0; j < npad; ++j)
885 //if ( cluster.Pad(j)->Status() < 0) continue;
886 if ( cluster.Pad(j)->IsSaturated()) continue;
887 if (aijclupad(iclust,j) < fgkCouplMin) continue;
888 for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++)
890 if (aijclupad(iclust1,j) < fgkCouplMin) continue;
891 aijcluclu(iclust,iclust1) +=
892 TMath::Sqrt (aijclupad(iclust,j)*aijclupad(iclust1,j));
896 for (Int_t iclust = 0; iclust < nclust; ++iclust)
898 for (Int_t iclust1 = iclust+1; iclust1 < nclust; ++iclust1)
900 aijcluclu(iclust1,iclust) = aijcluclu(iclust,iclust1);
904 if (fDebug && nclust > 1) aijcluclu.Print();
906 // Find groups of coupled clusters
907 used = new Bool_t[nclust];
908 for (Int_t j = 0; j < nclust; ++j) used[j] = kFALSE;
910 Int_t *clustNumb = new Int_t[nclust];
911 Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
913 Double_t parOk[8] = {0}; //AZ
915 for (Int_t igroup = 0; igroup < nclust; ++igroup)
917 if (used[igroup]) continue;
918 used[igroup] = kTRUE;
919 clustNumb[0] = igroup;
921 // Find group of coupled clusters
922 AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
925 cout << " nCoupled: " << nCoupled << endl;
926 for (Int_t i=0; i<nCoupled; ++i) cout << clustNumb[i] << " "; cout << endl;
929 fnCoupled = nCoupled;
936 for (Int_t i = 0; i < nCoupled; ++i) clustFit[i] = clustNumb[i];
940 // Too many coupled clusters to fit - try to decouple them
941 // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with
942 // all the others in the group
943 for (Int_t j = 0; j < 3; ++j) minGroup[j] = -1;
944 Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
946 // Flag clusters for fit
948 while (minGroup[nForFit] >= 0 && nForFit < 3)
950 if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
951 clustFit[nForFit] = clustNumb[minGroup[nForFit]];
952 clustNumb[minGroup[nForFit]] -= 999;
955 if (fDebug) cout << " nForFit " << nForFit << " " << coupl << endl;
958 // Select pads for fit.
959 if (SelectPad(cluster,nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1)
962 for (Int_t j = 0; j < npad; ++j)
964 AliMUONPad* pad = cluster.Pad(j);
965 //if ( pad->Status()==1 ) pad->SetStatus(0);
966 //if ( pad->Status()==-9) pad->SetStatus(-5);
967 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() ||
968 pad->Status() == AliMUONClusterFinderMLEM::GetCoupledFlag())
969 pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
971 // Merge the failed cluster candidates (with too few pads to fit) with
972 // the one with the strongest coupling
973 Merge(cluster,nForFit, nCoupled, clustNumb, clustFit, clusters, aijcluclu, aijclupad);
978 nfit = Fit(cluster,0, nForFit, clustFit, clusters, parOk, clusterList, mlem);
980 //cout << " (nfit == 0) " << fNpar << " " << cluster.Multiplicity() << endl;
981 fNpar = 0; // should be 0 by itself but just in case ...
985 // Subtract the fitted charges from pads with strong coupling and/or
986 // return pads for further use
987 UpdatePads(cluster,nfit, parOk);
990 for (Int_t j = 0; j < npad; ++j)
992 AliMUONPad* pad = cluster.Pad(j);
993 //if ( pad->Status()==1 ) pad->SetStatus(-2);
994 //if ( pad->Status()==-9) pad->SetStatus(-5);
995 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() )
996 pad->SetStatus(AliMUONClusterFinderMLEM::GetModifiedFlag());
999 // Sort the clusters (move to the right the used ones)
1000 Int_t beg = 0, end = nCoupled - 1;
1003 if (clustNumb[beg] >= 0) { ++beg; continue; }
1004 for (Int_t j = end; j > beg; --j)
1006 if (clustNumb[j] < 0) continue;
1008 indx = clustNumb[beg];
1009 clustNumb[beg] = clustNumb[j];
1010 clustNumb[j] = indx;
1016 nCoupled -= nForFit;
1019 // Remove couplings of used clusters
1020 for (Int_t iclust = nCoupled; iclust < nCoupled+nForFit; ++iclust)
1022 indx = clustNumb[iclust] + 999;
1023 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1025 indx1 = clustNumb[iclust1];
1026 aijcluclu(indx,indx1) = aijcluclu(indx1,indx) = 0;
1030 // Update the remaining clusters couplings (subtract couplings from
1031 // the used pads) - overflows excluded
1032 for (Int_t j = 0; j < npad; ++j)
1034 AliMUONPad* pad = cluster.Pad(j);
1035 //if ( pad->Status() != -2) continue;
1036 if ( pad->Status() != AliMUONClusterFinderMLEM::GetModifiedFlag()) continue;
1037 for (Int_t iclust=0; iclust<nCoupled; ++iclust)
1039 indx = clustNumb[iclust];
1040 if (aijclupad(indx,j) < fgkCouplMin) continue;
1041 for (Int_t iclust1 = iclust+1; iclust1 < nCoupled; ++iclust1)
1043 indx1 = clustNumb[iclust1];
1044 if (aijclupad(indx1,j) < fgkCouplMin) continue;
1046 aijcluclu(indx,indx1) -=
1047 TMath::Sqrt (aijclupad(indx,j)*aijclupad(indx1,j));
1048 aijcluclu(indx1,indx) = aijcluclu(indx,indx1);
1051 //pad->SetStatus(-8);
1052 pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag());
1053 } // for (Int_t j=0; j<npad;
1054 } // if (nCoupled > 3)
1055 } // while (nCoupled > 0)
1056 } // for (Int_t igroup=0; igroup<nclust;
1058 for (Int_t iclust = 0; iclust < nclust; ++iclust)
1060 pix = clusters[iclust];
1064 delete [] clustNumb;
1069 //_____________________________________________________________________________
1071 AliMUONClusterSplitterMLEM::Merge(const AliMUONCluster& cluster,
1072 Int_t nForFit, Int_t nCoupled,
1073 const Int_t *clustNumb, const Int_t *clustFit,
1074 TObjArray **clusters,
1075 TMatrixD& aijcluclu, TMatrixD& aijclupad)
1077 /// Merge the group of clusters with the one having the strongest coupling with them
1079 Int_t indx, indx1, npxclu, npxclu1, imax=0;
1080 TObjArray *pix, *pix1;
1083 for (Int_t icl = 0; icl < nForFit; ++icl)
1085 indx = clustFit[icl];
1086 pix = clusters[indx];
1087 npxclu = pix->GetEntriesFast();
1089 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
1091 indx1 = clustNumb[icl1];
1092 if (indx1 < 0) continue;
1093 if ( aijcluclu(indx,indx1) > couplMax)
1095 couplMax = aijcluclu(indx,indx1);
1098 } // for (Int_t icl1=0;
1100 pix1 = clusters[imax];
1101 npxclu1 = pix1->GetEntriesFast();
1103 for (Int_t i = 0; i < npxclu; ++i)
1105 pix1->Add(pix->UncheckedAt(i));
1109 //Add cluster-to-cluster couplings
1110 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
1112 indx1 = clustNumb[icl1];
1113 if (indx1 < 0 || indx1 == imax) continue;
1114 aijcluclu(indx1,imax) += aijcluclu(indx,indx1);
1115 aijcluclu(imax,indx1) = aijcluclu(indx1,imax);
1117 aijcluclu(indx,imax) = aijcluclu(imax,indx) = 0;
1119 //Add cluster-to-pad couplings
1120 Int_t mult = cluster.Multiplicity();
1121 for (Int_t j = 0; j < mult; ++j)
1123 AliMUONPad* pad = cluster.Pad(j);
1124 //if ( pad->Status() < 0 && pad->Status() != -5 ) continue;// exclude used pads
1125 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()) continue;// exclude used pads
1126 aijclupad(imax,j) += aijclupad(indx,j);
1127 aijclupad(indx,j) = 0;
1129 } // for (Int_t icl=0; icl<nForFit;
1133 //_____________________________________________________________________________
1135 AliMUONClusterSplitterMLEM::MinGroupCoupl(Int_t nCoupled, const Int_t *clustNumb,
1136 const TMatrixD& aijcluclu, Int_t *minGroup)
1138 /// Find group of clusters with minimum coupling to all the others
1140 Int_t i123max = TMath::Min(3,nCoupled/2);
1141 Int_t indx, indx1, indx2, indx3, nTot = 0;
1142 Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1144 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1147 coupl1 = new Double_t [nCoupled];
1148 for (Int_t i = 0; i < nCoupled; ++i) coupl1[i] = 0;
1150 else if (i123 == 2) {
1151 nTot = nCoupled*nCoupled;
1152 coupl2 = new Double_t [nTot];
1153 for (Int_t i = 0; i < nTot; ++i) coupl2[i] = 9999;
1155 nTot = nTot*nCoupled;
1156 coupl3 = new Double_t [nTot];
1157 for (Int_t i = 0; i < nTot; ++i) coupl3[i] = 9999;
1160 for (Int_t i = 0; i < nCoupled; ++i) {
1161 indx1 = clustNumb[i];
1162 for (Int_t j = i+1; j < nCoupled; ++j) {
1163 indx2 = clustNumb[j];
1165 coupl1[i] += aijcluclu(indx1,indx2);
1166 coupl1[j] += aijcluclu(indx1,indx2);
1168 else if (i123 == 2) {
1169 indx = i*nCoupled + j;
1170 coupl2[indx] = coupl1[i] + coupl1[j];
1171 coupl2[indx] -= 2 * (aijcluclu(indx1,indx2));
1173 for (Int_t k = j+1; k < nCoupled; ++k) {
1174 indx3 = clustNumb[k];
1175 indx = i*nCoupled*nCoupled + j*nCoupled + k;
1176 coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1177 coupl3[indx] -= 2 * (aijcluclu(indx1,indx3)+aijcluclu(indx2,indx3));
1180 } // for (Int_t j=i+1;
1181 } // for (Int_t i=0;
1182 } // for (Int_t i123=1;
1184 // Find minimum coupling
1185 Double_t couplMin = 9999;
1188 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1190 locMin = TMath::LocMin(nCoupled, coupl1);
1191 couplMin = coupl1[locMin];
1192 minGroup[0] = locMin;
1195 else if (i123 == 2) {
1196 locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
1197 if (coupl2[locMin] < couplMin) {
1198 couplMin = coupl2[locMin];
1199 minGroup[0] = locMin/nCoupled;
1200 minGroup[1] = locMin%nCoupled;
1204 locMin = TMath::LocMin(nTot, coupl3);
1205 if (coupl3[locMin] < couplMin) {
1206 couplMin = coupl3[locMin];
1207 minGroup[0] = locMin/nCoupled/nCoupled;
1208 minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
1209 minGroup[2] = locMin%nCoupled;
1213 } // for (Int_t i123=1;
1217 //_____________________________________________________________________________
1219 AliMUONClusterSplitterMLEM::SelectPad(const AliMUONCluster& cluster,
1220 Int_t nCoupled, Int_t nForFit,
1221 const Int_t *clustNumb, const Int_t *clustFit,
1222 const TMatrixD& aijclupad)
1224 /// Select pads for fit. If too many coupled clusters, find pads giving
1225 /// the strongest coupling with the rest of clusters and exclude them from the fit.
1227 Int_t npad = cluster.Multiplicity();
1228 Double_t *padpix = 0;
1232 padpix = new Double_t[npad];
1233 for (Int_t i = 0; i < npad; ++i) padpix[i] = 0.;
1236 Int_t nOK = 0, indx, indx1;
1237 for (Int_t iclust = 0; iclust < nForFit; ++iclust)
1239 indx = clustFit[iclust];
1240 for (Int_t j = 0; j < npad; ++j)
1242 if ( aijclupad(indx,j) < fgkCouplMin) continue;
1243 AliMUONPad* pad = cluster.Pad(j);
1245 if ( pad->Status() == -5 ) pad->SetStatus(-9); // flag overflow
1246 if ( pad->Status() < 0 ) continue; // exclude overflows and used pads
1247 if ( !pad->Status() )
1250 ++nOK; // pad to be used in fit
1253 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()
1254 || pad->IsSaturated() ) continue; // used pads and overflows
1255 pad->SetStatus(AliMUONClusterFinderMLEM::GetUseForFitFlag());
1256 ++nOK; // pad to be used in fit
1260 // Check other clusters
1261 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1263 indx1 = clustNumb[iclust1];
1264 if (indx1 < 0) continue;
1265 if ( aijclupad(indx1,j) < fgkCouplMin ) continue;
1266 padpix[j] += aijclupad(indx1,j);
1268 } // if (nCoupled > 3)
1269 } // for (Int_t j=0; j<npad;
1270 } // for (Int_t iclust=0; iclust<nForFit
1271 if (nCoupled < 4) return nOK;
1274 for (Int_t j = 0; j < npad; ++j)
1276 if (padpix[j] < fgkCouplMin) continue;
1278 //cluster.Pad(j)->SetStatus(-1); // exclude pads with strong coupling to the other clusters
1279 cluster.Pad(j)->SetStatus(AliMUONClusterFinderMLEM::GetCoupledFlag()); // exclude pads with strong coupling to the other clusters
1286 //_____________________________________________________________________________
1288 AliMUONClusterSplitterMLEM::UpdatePads(const AliMUONCluster& cluster,
1289 Int_t /*nfit*/, Double_t *par)
1291 /// Subtract the fitted charges from pads with strong coupling
1293 Int_t indx, mult = cluster.Multiplicity(), iend = fNpar/3;
1294 Double_t charge, coef=0;
1296 for (Int_t j = 0; j < mult; ++j)
1298 AliMUONPad* pad = cluster.Pad(j);
1299 //if ( pad->Status() != -1 ) continue;
1300 if ( pad->Status() != AliMUONClusterFinderMLEM::GetCoupledFlag() ) continue;
1304 for (Int_t i = 0; i <= iend; ++i)
1308 coef = Param2Coef(i, coef, par);
1309 charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
1312 pad->SetCharge(pad->Charge()-charge);
1313 } // if (fNpar != 0)
1315 //if (pad->Charge() > 6 /*fgkZeroSuppression*/) pad->SetStatus(0);
1316 if (pad->Charge() > fLowestPadCharge) pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
1317 // return pad for further using // FIXME: remove usage of zerosuppression here
1318 else pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag()); // do not use anymore
1320 } // for (Int_t j=0;