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->GetCellContent(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->GetCellContent(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;
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);
274 Double_t xPad = 0, yPad = 99999;
276 // Number of pads to use and number of virtual pads
277 Int_t npads = 0, nVirtual = 0, nfit0 = nfit;
278 //cluster.Print("full");
279 Int_t mult = cluster.Multiplicity();
280 for (Int_t i = 0; i < mult; ++i )
282 AliMUONPad* pad = cluster.Pad(i);
283 if ( !pad->IsReal() ) ++nVirtual;
284 //if ( pad->Status() !=1 || pad->IsSaturated() ) continue;
285 if ( pad->Status() != AliMUONClusterFinderMLEM::GetUseForFitFlag() ) continue;
296 if (pad->DY() < pad->DX() )
311 if (npads < 2) return 0;
313 // FIXME : AliWarning("Reconnect the following code for hit/track passing ?");
315 // Int_t tracks[3] = {-1, -1, -1};
319 AliMUONDigit *mdig = 0;
320 for (Int_t cath=0; cath<2; cath++) {
321 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
322 if (fPadIJ[0][i] != cath) continue;
323 if (fPadIJ[1][i] != 1) continue;
324 if (fXyq[3][i] < 0) continue; // exclude virtual pads
325 digit = TMath::Nint (fXyq[5][i]);
326 if (digit >= 0) mdig = fInput->Digit(cath,digit);
327 else mdig = fInput->Digit(TMath::Even(cath),-digit-1);
328 //if (!mdig) mdig = fInput->Digit(TMath::Even(cath),digit);
329 if (!mdig) continue; // protection for cluster display
330 if (mdig->Hit() >= 0) {
332 tracks[0] = mdig->Hit();
333 tracks[1] = mdig->Track(0);
334 } else if (mdig->Track(0) < tracks[1]) {
335 tracks[0] = mdig->Hit();
336 tracks[1] = mdig->Track(0);
339 if (mdig->Track(1) >= 0 && mdig->Track(1) != tracks[1]) {
340 if (tracks[2] < 0) tracks[2] = mdig->Track(1);
341 else tracks[2] = TMath::Min (tracks[2], mdig->Track(1));
344 } // for (Int_t cath=0;
347 // Get number of pads in X and Y
348 //const Int_t kStatusToTest(1);
349 const Int_t kStatusToTest(AliMUONClusterFinderMLEM::GetUseForFitFlag());
351 Long_t nofPads = cluster.NofPads(kStatusToTest);
352 Int_t nInX = AliMp::PairFirst(nofPads);
353 Int_t nInY = AliMp::PairSecond(nofPads);
357 for (Int_t j = 0; j < cluster.Multiplicity(); ++j) {
358 AliMUONPad *pad = cluster.Pad(j);
359 //if (pad->Status() == 1 && !pad->IsSaturated()) npadOK++;
360 if (pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() && !pad->IsSaturated()) npadOK++;
362 cout << " Number of pads to fit: " << npadOK << endl;
363 cout << " nInX and Y: " << nInX << " " << nInY << endl;
367 nfitMax = TMath::Min (nfitMax, (npads + 1) / 3);
369 if (((nInX < 3) && (nInY < 3)) || ((nInX == 3) && (nInY < 3)) || ((nInX < 3) && (nInY == 3))) nfitMax = 1; // not enough pads in each direction
371 if (nfit > nfitMax) nfit = nfitMax;
373 // Take cluster maxima as fitting seeds
377 Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8], qq = 0;
379 for ( int i = 0; i < 8; ++i ) errOk[i]=0.0;
381 Double_t xyseed[3][2], qseed[3], xyCand[3][2] = {{0},{0}}, sigCand[3][2] = {{0},{0}};
383 for (Int_t ifit = 1; ifit <= nfit0; ++ifit)
386 pix = clusters[clustFit[ifit-1]];
387 npxclu = pix->GetEntriesFast();
389 for (Int_t clu = 0; clu < npxclu; ++clu)
391 pixPtr = (AliMUONPad*) pix->UncheckedAt(clu);
392 cont = pixPtr->Charge();
397 xseed = pixPtr->Coord(0);
398 yseed = pixPtr->Coord(1);
401 xyCand[0][0] += pixPtr->Coord(0) * cont;
402 xyCand[0][1] += pixPtr->Coord(1) * cont;
403 sigCand[0][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
404 sigCand[0][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
406 xyseed[ifit-1][0] = xseed;
407 xyseed[ifit-1][1] = yseed;
408 qseed[ifit-1] = cmax;
409 } // for (Int_t ifit=1;
411 xyCand[0][0] /= qq; // <x>
412 xyCand[0][1] /= qq; // <y>
413 sigCand[0][0] = sigCand[0][0]/qq - xyCand[0][0]*xyCand[0][0]; // <x^2> - <x>^2
414 sigCand[0][0] = sigCand[0][0] > 0 ? TMath::Sqrt (sigCand[0][0]) : 0;
415 sigCand[0][1] = sigCand[0][1]/qq - xyCand[0][1]*xyCand[0][1]; // <y^2> - <y>^2
416 sigCand[0][1] = sigCand[0][1] > 0 ? TMath::Sqrt (sigCand[0][1]) : 0;
417 if (fDebug) cout << xyCand[0][0] << " " << xyCand[0][1] << " " << sigCand[0][0] << " " << sigCand[0][1] << endl;
419 Int_t nDof, maxSeed[3];//, nMax = 0;
421 if ( nfit0 < 0 || nfit0 > 3 ) {
422 AliErrorStream() << "Wrong nfit0 value: " << nfit0 << endl;
425 TMath::Sort(nfit0, qseed, maxSeed, kTRUE); // in decreasing order
427 Double_t step[3]={0.01,0.002,0.02}, fmin, chi2o = 9999, chi2n;
428 Double_t *gin = 0, func0, func1, param[8]={0}, step0[8]={0};
429 Double_t param0[2][8]={{0},{0}}, deriv[2][8]={{0},{0}};
430 Double_t shift[8]={0}, stepMax, derMax, parmin[8]={0}, parmax[8]={0}, func2[2]={0}, shift0;
431 Double_t delta[8]={0}, scMax, dder[8], estim, shiftSave = 0;
432 Int_t min, max, nCall = 0, nLoop, idMax = 0, iestMax = 0, nFail;
433 Double_t rad, dist[3] = {0};
435 // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is
436 // lower, try 3-track (if number of pads is sufficient).
437 Int_t iflag = 0; // for the first call of fcn1
438 for (Int_t iseed = 0; iseed < nfit; ++iseed)
441 Int_t memory[8] = {0};
444 for (Int_t j = 0; j < fNpar; ++j)
449 parmin[fNpar] = 1E-9;
455 param[fNpar] = xyCand[0][0]; // take COG
459 param[fNpar] = xyseed[maxSeed[iseed]][0];
460 //param[fNpar] = fNpar==0 ? -16.1651 : -15.2761;
462 parmin[fNpar] = xmin;
463 parmax[fNpar++] = xmax;
466 param[fNpar] = xyCand[0][1]; // take COG
470 param[fNpar] = xyseed[maxSeed[iseed]][1];
471 //param[fNpar] = fNpar==1 ? -15.1737 : -15.8487;
473 parmin[fNpar] = ymin;
474 parmax[fNpar++] = ymax;
476 for (Int_t j = 0; j < fNpar; ++j)
478 step0[j] = shift[j] = step[j%3];
483 for (Int_t j = 0; j < fNpar; ++j)
489 for (Int_t j = 0; j < fNpar; ++j) cout << param[j] << " ";
494 min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
499 Fcn1(cluster,fNpar, gin, func0, param, iflag); nCall++;
501 //cout << " Func: " << func0 << endl;
504 for (Int_t j = 0; j < fNpar; ++j)
506 param0[max][j] = param[j];
508 param[j] += delta[j] / 10;
509 if (j > 0) param[j-1] -= delta[j-1] / 10;
510 Fcn1(cluster,fNpar, gin, func1, param, iflag); nCall++;
511 deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
512 //cout << j << " " << deriv[max][j] << endl;
513 dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) /
514 (param0[0][j] - param0[1][j]) : 0; // second derivative
516 param[fNpar-1] -= delta[fNpar-1] / 10;
517 if (nCall > 2000) break;
519 min = func2[0] < func2[1] ? 0 : 1;
520 nFail = min == max ? 0 : nFail + 1;
522 stepMax = derMax = estim = 0;
523 for (Int_t j = 0; j < fNpar; ++j)
525 // Estimated distance to minimum
529 shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
531 else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3)
535 else if (((deriv[min][j]*deriv[!min][j] > 0) && (TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])))
536 || (TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3) || (TMath::Abs(dder[j]) < 1.e-6))
538 shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
550 shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
554 Double_t es = TMath::Abs(shift[j]) / step0[j];
562 if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; //
564 // Failed to improve minimum
568 param[j] = param0[min][j];
569 if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j])
571 shift[j] = (shift[j] + shift0) / 2;
580 if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min])
582 shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
585 // Introduce step relaxation factor
588 scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
589 if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax)
591 shift[j] = TMath::Sign (shift0*scMax, shift[j]);
594 param[j] += shift[j];
595 // Check parameter limits
596 if (param[j] < parmin[j])
598 shift[j] = parmin[j] - param[j];
599 param[j] = parmin[j];
601 else if (param[j] > parmax[j])
603 shift[j] = parmax[j] - param[j];
604 param[j] = parmax[j];
606 //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
607 stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
608 if (TMath::Abs(deriv[min][j]) > derMax)
611 derMax = TMath::Abs (deriv[min][j]);
613 } // for (Int_t j=0; j<fNpar;
615 if (((estim < 1) && (derMax < 2)) || nLoop > 150) break; // minimum was found
619 // Check for small step
620 if (shift[idMax] == 0)
622 shift[idMax] = step0[idMax]/10;
623 param[idMax] += shift[idMax];
627 if (!memory[idMax] && derMax > 0.5 && nLoop > 10)
629 if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10)
631 if (min == max) dder[idMax] = -dder[idMax];
632 shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10;
633 param[idMax] += shift[idMax];
634 stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
635 if (min == max) shiftSave = shift[idMax];
639 param[idMax] -= shift[idMax];
640 shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
641 param[idMax] += shift[idMax];
648 nDof = npads - fNpar + nVirtual;
651 if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
653 //if (fNpar > 2) cout << param0[min][fNpar-3] << " " << chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) << endl;
654 //if (chi2n*1.2+1.e-6 > chi2o )
655 if (fNpar > 2 && (chi2n > chi2o || ((iseed == nfit-1)
656 && (chi2n * (1+TMath::Min(1-param0[min][fNpar-3],0.25)) > chi2o))))
657 { fNpar -= 3; break; }
659 // Save parameters and errors
662 // One pad per direction
663 //for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5) param0[min][i] = xPad;
664 for (Int_t i=0; i<fNpar; ++i) if (i == 0 || i == 2 || i == 5)
665 param0[min][i] = xyCand[0][0];
668 // One pad per direction
669 //for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6) param0[min][i] = yPad;
670 for (Int_t i=0; i<fNpar; ++i) if (i == 1 || i == 3 || i == 6)
671 param0[min][i] = xyCand[0][1];
676 // Find distance to the nearest neighbour
677 dist[0] = dist[1] = TMath::Sqrt ((param0[min][0]-param0[min][2])*
678 (param0[min][0]-param0[min][2])
679 +(param0[min][1]-param0[min][3])*
680 (param0[min][1]-param0[min][3]));
682 dist[2] = TMath::Sqrt ((param0[min][0]-param0[min][5])*
683 (param0[min][0]-param0[min][5])
684 +(param0[min][1]-param0[min][6])*
685 (param0[min][1]-param0[min][6]));
686 rad = TMath::Sqrt ((param0[min][2]-param0[min][5])*
687 (param0[min][2]-param0[min][5])
688 +(param0[min][3]-param0[min][6])*
689 (param0[min][3]-param0[min][6]));
690 if (dist[2] < dist[0]) dist[0] = dist[2];
691 if (rad < dist[1]) dist[1] = rad;
692 if (rad < dist[2]) dist[2] = rad;
694 cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
695 if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
699 for (Int_t i = 0; i < fNpar; ++i) {
700 parOk[i] = param0[min][i];
704 parOk[i] = TMath::Max (parOk[i], parmin[i]);
705 parOk[i] = TMath::Min (parOk[i], parmax[i]);
709 if (fmin < 0.1) break; // !!!???
710 } // for (Int_t iseed=0;
713 for (Int_t i=0; i<fNpar; ++i) {
714 if (i == 4 || i == 7) {
715 if ((i == 7) || ((i == 4) && (fNpar < 7))) cout << parOk[i] << endl;
716 else cout << parOk[i] * (1-parOk[7]) << endl;
719 cout << parOk[i] << " " << errOk[i] << endl;
722 nfit = (fNpar + 1) / 3;
723 dist[0] = dist[1] = dist[2] = 0;
726 // Find distance to the nearest neighbour
727 dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
729 +(parOk[1]-parOk[3])*
730 (parOk[1]-parOk[3]));
732 dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
734 +(parOk[1]-parOk[6])*
735 (parOk[1]-parOk[6]));
736 rad = TMath::Sqrt ((parOk[2]-parOk[5])*
738 +(parOk[3]-parOk[6])*
739 (parOk[3]-parOk[6]));
740 if (dist[2] < dist[0]) dist[0] = dist[2];
741 if (rad < dist[1]) dist[1] = rad;
742 if (rad < dist[2]) dist[2] = rad;
749 if (iSimple) fnCoupled = 0;
750 for (Int_t j = 0; j < nfit; ++j) {
752 coef = Param2Coef(j, coef, parOk);
754 //void AliMUONClusterFinderMLEM::AddRawCluster(Double_t x, Double_t y,
755 // Double_t qTot, Double_t fmin,
756 // Int_t nfit, Int_t *tracks,
757 // Double_t /*sigx*/,
758 // Double_t /*sigy*/,
759 // Double_t /*dist*/)
761 if ( coef*fQtot >= fLowestClusterCharge )
763 //AZ AliMUONCluster* cluster1 = new AliMUONCluster();
764 AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
766 cluster1->SetCharge(coef*fQtot,coef*fQtot);
767 cluster1->SetPosition(TVector2(parOk[indx],parOk[indx+1]),TVector2(sigCand[0][0],sigCand[0][1]));
768 //cluster1->SetChi2(dist[TMath::LocMin(nfit,dist)]);
769 Int_t idx = TMath::LocMin(nfit,dist);
770 if ( idx < 0 || idx > 2 ) {
771 AliErrorStream() << "Wrong index value: " << idx << endl;
774 cluster1->SetChi2(dist[idx]);
776 // FIXME: we miss some information in this cluster, as compared to
777 // the original AddRawCluster code.
779 AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
780 fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
781 cluster1->Position().X(),cluster1->Position().Y()));
783 clusterList.Add(cluster1);
785 // AddRawCluster (parOk[indx], // double x
786 // parOk[indx+1], // double y
787 // coef*qTot, // double charge
788 // errOk[indx], // double fmin
789 // nfit0+10*nfit+100*nMax+10000*fnCoupled, // int nfit
790 // tracks, // int* tracks
791 // sigCand[0][0], // double sigx
792 // sigCand[0][1], // double sigy
793 // dist[TMath::LocMin(nfit,dist)] // double dist
800 //_____________________________________________________________________________
802 AliMUONClusterSplitterMLEM::Split(const AliMUONCluster& cluster,
803 TH2 *mlem, Double_t *coef,
804 TObjArray& clusterList)
806 /// The main steering function to work with clusters of pixels in anode
807 /// plane (find clusters, decouple them from each other, merge them (if
808 /// necessary), pick up coupled pads, call the fitting function)
810 Int_t nx = mlem->GetNbinsX();
811 Int_t ny = mlem->GetNbinsY();
812 Int_t nPix = fPixArray->GetEntriesFast();
815 Int_t nclust = 0, indx, indx1, nxy = ny * nx;
816 Bool_t *used = new Bool_t[nxy];
818 for (Int_t j = 0; j < nxy; ++j) used[j] = kFALSE;
820 TObjArray *clusters[200]={0};
823 // Find clusters of histogram bins (easier to work in 2-D space)
824 for (Int_t i = 1; i <= ny; ++i)
826 for (Int_t j = 1; j <= nx; ++j)
828 indx = (i-1)*nx + j - 1;
829 if (used[indx]) continue;
830 cont = mlem->GetCellContent(j,i);
831 if (cont < fLowestPixelCharge) continue;
832 pix = new TObjArray(20);
834 pix->Add(BinToPix(mlem,j,i));
835 AddBin(mlem, i, j, 0, used, pix); // recursive call
836 if (nclust >= 200) AliFatal(" Too many clusters !!!");
837 clusters[nclust++] = pix;
838 } // for (Int_t j=1; j<=nx; j++) {
839 } // for (Int_t i=1; i<=ny;
840 if (fDebug) cout << nclust << endl;
843 // Compute couplings between clusters and clusters to pads
844 Int_t npad = cluster.Multiplicity();
846 // Exclude pads with overflows
848 for (Int_t j = 0; j < npad; ++j)
850 AliMUONPad* pad = cluster.Pad(j);
851 if ( pad->IsSaturated() )
862 // Compute couplings of clusters to pads (including overflows)
863 TMatrixD aijclupad(nclust,npad);
866 for (Int_t iclust = 0; iclust < nclust; ++iclust)
868 pix = clusters[iclust];
869 npxclu = pix->GetEntriesFast();
870 for (Int_t i = 0; i < npxclu; ++i)
872 indx = fPixArray->IndexOf(pix->UncheckedAt(i));
873 for (Int_t j = 0; j < npad; ++j)
875 //AliMUONPad* pad = cluster.Pad(j);
876 //if ( pad->Status() < 0 && pad->Status() != -5) continue;
877 if (coef[j*nPix+indx] < fgkCouplMin) continue;
878 aijclupad(iclust,j) += coef[j*nPix+indx];
883 // Compute couplings between clusters (exclude overflows)
884 TMatrixD aijcluclu(nclust,nclust);
886 for (Int_t iclust = 0; iclust < nclust; ++iclust)
888 for (Int_t j = 0; j < npad; ++j)
891 //if ( cluster.Pad(j)->Status() < 0) continue;
892 if ( cluster.Pad(j)->IsSaturated()) continue;
893 if (aijclupad(iclust,j) < fgkCouplMin) continue;
894 for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++)
896 if (aijclupad(iclust1,j) < fgkCouplMin) continue;
897 aijcluclu(iclust,iclust1) +=
898 TMath::Sqrt (aijclupad(iclust,j)*aijclupad(iclust1,j));
902 for (Int_t iclust = 0; iclust < nclust; ++iclust)
904 for (Int_t iclust1 = iclust+1; iclust1 < nclust; ++iclust1)
906 aijcluclu(iclust1,iclust) = aijcluclu(iclust,iclust1);
910 if (fDebug && nclust > 1) aijcluclu.Print();
912 // Find groups of coupled clusters
913 used = new Bool_t[nclust];
914 for (Int_t j = 0; j < nclust; ++j) used[j] = kFALSE;
916 Int_t *clustNumb = new Int_t[nclust];
917 Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
919 Double_t parOk[8] = {0}; //AZ
921 for (Int_t igroup = 0; igroup < nclust; ++igroup)
923 if (used[igroup]) continue;
924 used[igroup] = kTRUE;
925 clustNumb[0] = igroup;
927 // Find group of coupled clusters
928 AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
931 cout << " nCoupled: " << nCoupled << endl;
932 for (Int_t i=0; i<nCoupled; ++i) cout << clustNumb[i] << " "; cout << endl;
935 fnCoupled = nCoupled;
942 for (Int_t i = 0; i < nCoupled; ++i) clustFit[i] = clustNumb[i];
946 // Too many coupled clusters to fit - try to decouple them
947 // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with
948 // all the others in the group
949 for (Int_t j = 0; j < 3; ++j) minGroup[j] = -1;
950 Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
952 // Flag clusters for fit
954 while (nForFit < 3 && minGroup[nForFit] >= 0)
956 if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
957 clustFit[nForFit] = clustNumb[minGroup[nForFit]];
958 clustNumb[minGroup[nForFit]] -= 999;
961 if (fDebug) cout << " nForFit " << nForFit << " " << coupl << endl;
964 // Select pads for fit.
965 if (SelectPad(cluster,nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1)
968 for (Int_t j = 0; j < npad; ++j)
970 AliMUONPad* pad = cluster.Pad(j);
971 //if ( pad->Status()==1 ) pad->SetStatus(0);
972 //if ( pad->Status()==-9) pad->SetStatus(-5);
973 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() ||
974 pad->Status() == AliMUONClusterFinderMLEM::GetCoupledFlag())
975 pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
977 // Merge the failed cluster candidates (with too few pads to fit) with
978 // the one with the strongest coupling
979 Merge(cluster,nForFit, nCoupled, clustNumb, clustFit, clusters, aijcluclu, aijclupad);
984 nfit = Fit(cluster,0, nForFit, clustFit, clusters, parOk, clusterList, mlem);
986 //cout << " (nfit == 0) " << fNpar << " " << cluster.Multiplicity() << endl;
987 fNpar = 0; // should be 0 by itself but just in case ...
991 // Subtract the fitted charges from pads with strong coupling and/or
992 // return pads for further use
993 UpdatePads(cluster,nfit, parOk);
996 for (Int_t j = 0; j < npad; ++j)
998 AliMUONPad* pad = cluster.Pad(j);
999 //if ( pad->Status()==1 ) pad->SetStatus(-2);
1000 //if ( pad->Status()==-9) pad->SetStatus(-5);
1001 if ( pad->Status() == AliMUONClusterFinderMLEM::GetUseForFitFlag() )
1002 pad->SetStatus(AliMUONClusterFinderMLEM::GetModifiedFlag());
1005 // Sort the clusters (move to the right the used ones)
1006 Int_t beg = 0, end = nCoupled - 1;
1009 if (clustNumb[beg] >= 0) { ++beg; continue; }
1010 for (Int_t j = end; j > beg; --j)
1012 if (clustNumb[j] < 0) continue;
1014 indx = clustNumb[beg];
1015 clustNumb[beg] = clustNumb[j];
1016 clustNumb[j] = indx;
1022 nCoupled -= nForFit;
1025 // Remove couplings of used clusters
1026 for (Int_t iclust = nCoupled; iclust < nCoupled+nForFit; ++iclust)
1028 indx = clustNumb[iclust] + 999;
1029 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1031 indx1 = clustNumb[iclust1];
1032 aijcluclu(indx,indx1) = aijcluclu(indx1,indx) = 0;
1036 // Update the remaining clusters couplings (subtract couplings from
1037 // the used pads) - overflows excluded
1038 for (Int_t j = 0; j < npad; ++j)
1040 AliMUONPad* pad = cluster.Pad(j);
1041 //if ( pad->Status() != -2) continue;
1042 if ( pad->Status() != AliMUONClusterFinderMLEM::GetModifiedFlag()) continue;
1043 for (Int_t iclust=0; iclust<nCoupled; ++iclust)
1045 indx = clustNumb[iclust];
1046 if (aijclupad(indx,j) < fgkCouplMin) continue;
1047 for (Int_t iclust1 = iclust+1; iclust1 < nCoupled; ++iclust1)
1049 indx1 = clustNumb[iclust1];
1050 if (aijclupad(indx1,j) < fgkCouplMin) continue;
1052 aijcluclu(indx,indx1) -=
1053 TMath::Sqrt (aijclupad(indx,j)*aijclupad(indx1,j));
1054 aijcluclu(indx1,indx) = aijcluclu(indx,indx1);
1057 //pad->SetStatus(-8);
1058 pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag());
1059 } // for (Int_t j=0; j<npad;
1060 } // if (nCoupled > 3)
1061 } // while (nCoupled > 0)
1062 } // for (Int_t igroup=0; igroup<nclust;
1064 for (Int_t iclust = 0; iclust < nclust; ++iclust)
1066 pix = clusters[iclust];
1070 delete [] clustNumb;
1075 //_____________________________________________________________________________
1077 AliMUONClusterSplitterMLEM::Merge(const AliMUONCluster& cluster,
1078 Int_t nForFit, Int_t nCoupled,
1079 const Int_t *clustNumb, const Int_t *clustFit,
1080 TObjArray **clusters,
1081 TMatrixD& aijcluclu, TMatrixD& aijclupad)
1083 /// Merge the group of clusters with the one having the strongest coupling with them
1085 Int_t indx, indx1, npxclu, npxclu1, imax=0;
1086 TObjArray *pix, *pix1;
1089 for (Int_t icl = 0; icl < nForFit; ++icl)
1091 indx = clustFit[icl];
1092 pix = clusters[indx];
1093 npxclu = pix->GetEntriesFast();
1095 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
1097 indx1 = clustNumb[icl1];
1098 if (indx1 < 0) continue;
1099 if ( aijcluclu(indx,indx1) > couplMax)
1101 couplMax = aijcluclu(indx,indx1);
1104 } // for (Int_t icl1=0;
1106 pix1 = clusters[imax];
1107 npxclu1 = pix1->GetEntriesFast();
1109 for (Int_t i = 0; i < npxclu; ++i)
1111 pix1->Add(pix->UncheckedAt(i));
1115 //Add cluster-to-cluster couplings
1116 for (Int_t icl1 = 0; icl1 < nCoupled; ++icl1)
1118 indx1 = clustNumb[icl1];
1119 if (indx1 < 0 || indx1 == imax) continue;
1120 aijcluclu(indx1,imax) += aijcluclu(indx,indx1);
1121 aijcluclu(imax,indx1) = aijcluclu(indx1,imax);
1123 aijcluclu(indx,imax) = aijcluclu(imax,indx) = 0;
1125 //Add cluster-to-pad couplings
1126 Int_t mult = cluster.Multiplicity();
1127 for (Int_t j = 0; j < mult; ++j)
1129 AliMUONPad* pad = cluster.Pad(j);
1130 //if ( pad->Status() < 0 && pad->Status() != -5 ) continue;// exclude used pads
1131 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()) continue;// exclude used pads
1132 aijclupad(imax,j) += aijclupad(indx,j);
1133 aijclupad(indx,j) = 0;
1135 } // for (Int_t icl=0; icl<nForFit;
1139 //_____________________________________________________________________________
1141 AliMUONClusterSplitterMLEM::MinGroupCoupl(Int_t nCoupled, const Int_t *clustNumb,
1142 const TMatrixD& aijcluclu, Int_t *minGroup)
1144 /// Find group of clusters with minimum coupling to all the others
1146 Int_t i123max = TMath::Min(3,nCoupled/2);
1147 Int_t indx, indx1, indx2, indx3, nTot = 0;
1148 Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1150 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1153 coupl1 = new Double_t [nCoupled];
1154 for (Int_t i = 0; i < nCoupled; ++i) coupl1[i] = 0;
1156 else if (i123 == 2) {
1157 nTot = nCoupled*nCoupled;
1158 coupl2 = new Double_t [nTot];
1159 for (Int_t i = 0; i < nTot; ++i) coupl2[i] = 9999;
1161 nTot = nTot*nCoupled;
1162 coupl3 = new Double_t [nTot];
1163 for (Int_t i = 0; i < nTot; ++i) coupl3[i] = 9999;
1166 for (Int_t i = 0; i < nCoupled; ++i) {
1167 indx1 = clustNumb[i];
1168 for (Int_t j = i+1; j < nCoupled; ++j) {
1169 indx2 = clustNumb[j];
1171 coupl1[i] += aijcluclu(indx1,indx2);
1172 coupl1[j] += aijcluclu(indx1,indx2);
1174 else if (i123 == 2) {
1175 indx = i*nCoupled + j;
1176 coupl2[indx] = coupl1[i] + coupl1[j];
1177 coupl2[indx] -= 2 * (aijcluclu(indx1,indx2));
1179 for (Int_t k = j+1; k < nCoupled; ++k) {
1180 indx3 = clustNumb[k];
1181 indx = i*nCoupled*nCoupled + j*nCoupled + k;
1182 coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1183 coupl3[indx] -= 2 * (aijcluclu(indx1,indx3)+aijcluclu(indx2,indx3));
1186 } // for (Int_t j=i+1;
1187 } // for (Int_t i=0;
1188 } // for (Int_t i123=1;
1190 // Find minimum coupling
1191 Double_t couplMin = 9999;
1194 for (Int_t i123 = 1; i123 <= i123max; ++i123) {
1196 locMin = TMath::LocMin(nCoupled, coupl1);
1197 couplMin = coupl1[locMin];
1198 minGroup[0] = locMin;
1201 else if (i123 == 2) {
1202 locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
1203 if (coupl2[locMin] < couplMin) {
1204 couplMin = coupl2[locMin];
1205 minGroup[0] = locMin/nCoupled;
1206 minGroup[1] = locMin%nCoupled;
1210 locMin = TMath::LocMin(nTot, coupl3);
1211 if (coupl3[locMin] < couplMin) {
1212 couplMin = coupl3[locMin];
1213 minGroup[0] = locMin/nCoupled/nCoupled;
1214 minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
1215 minGroup[2] = locMin%nCoupled;
1219 } // for (Int_t i123=1;
1223 //_____________________________________________________________________________
1225 AliMUONClusterSplitterMLEM::SelectPad(const AliMUONCluster& cluster,
1226 Int_t nCoupled, Int_t nForFit,
1227 const Int_t *clustNumb, const Int_t *clustFit,
1228 const TMatrixD& aijclupad)
1230 /// Select pads for fit. If too many coupled clusters, find pads giving
1231 /// the strongest coupling with the rest of clusters and exclude them from the fit.
1233 Int_t npad = cluster.Multiplicity();
1234 Double_t *padpix = 0;
1238 padpix = new Double_t[npad];
1239 for (Int_t i = 0; i < npad; ++i) padpix[i] = 0.;
1242 Int_t nOK = 0, indx, indx1;
1243 for (Int_t iclust = 0; iclust < nForFit; ++iclust)
1245 indx = clustFit[iclust];
1246 for (Int_t j = 0; j < npad; ++j)
1248 if ( aijclupad(indx,j) < fgkCouplMin) continue;
1249 AliMUONPad* pad = cluster.Pad(j);
1251 if ( pad->Status() == -5 ) pad->SetStatus(-9); // flag overflow
1252 if ( pad->Status() < 0 ) continue; // exclude overflows and used pads
1253 if ( !pad->Status() )
1256 ++nOK; // pad to be used in fit
1259 if ( pad->Status() != AliMUONClusterFinderMLEM::GetZeroFlag()
1260 || pad->IsSaturated() ) continue; // used pads and overflows
1261 pad->SetStatus(AliMUONClusterFinderMLEM::GetUseForFitFlag());
1262 ++nOK; // pad to be used in fit
1266 // Check other clusters
1267 for (Int_t iclust1 = 0; iclust1 < nCoupled; ++iclust1)
1269 indx1 = clustNumb[iclust1];
1270 if (indx1 < 0) continue;
1271 if ( aijclupad(indx1,j) < fgkCouplMin ) continue;
1272 padpix[j] += aijclupad(indx1,j);
1274 } // if (nCoupled > 3)
1275 } // for (Int_t j=0; j<npad;
1276 } // for (Int_t iclust=0; iclust<nForFit
1277 if (nCoupled < 4) return nOK;
1280 for (Int_t j = 0; j < npad; ++j)
1282 if (padpix[j] < fgkCouplMin) continue;
1284 //cluster.Pad(j)->SetStatus(-1); // exclude pads with strong coupling to the other clusters
1285 cluster.Pad(j)->SetStatus(AliMUONClusterFinderMLEM::GetCoupledFlag()); // exclude pads with strong coupling to the other clusters
1292 //_____________________________________________________________________________
1294 AliMUONClusterSplitterMLEM::UpdatePads(const AliMUONCluster& cluster,
1295 Int_t /*nfit*/, Double_t *par)
1297 /// Subtract the fitted charges from pads with strong coupling
1299 Int_t indx, mult = cluster.Multiplicity(), iend = fNpar/3;
1300 Double_t charge, coef=0;
1302 for (Int_t j = 0; j < mult; ++j)
1304 AliMUONPad* pad = cluster.Pad(j);
1305 //if ( pad->Status() != -1 ) continue;
1306 if ( pad->Status() != AliMUONClusterFinderMLEM::GetCoupledFlag() ) continue;
1310 for (Int_t i = 0; i <= iend; ++i)
1314 coef = Param2Coef(i, coef, par);
1315 charge += ChargeIntegration(par[indx],par[indx+1],*pad) * coef;
1318 pad->SetCharge(pad->Charge()-charge);
1319 } // if (fNpar != 0)
1321 //if (pad->Charge() > 6 /*fgkZeroSuppression*/) pad->SetStatus(0);
1322 if (pad->Charge() > fLowestPadCharge) pad->SetStatus(AliMUONClusterFinderMLEM::GetZeroFlag());
1323 // return pad for further using // FIXME: remove usage of zerosuppression here
1324 else pad->SetStatus(AliMUONClusterFinderMLEM::GetOverFlag()); // do not use anymore
1326 } // for (Int_t j=0;