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 AliMUONClusterFinderMLEM
21 /// Clusterizer class based on the Expectation-Maximization algorithm
23 /// Pre-clustering is handled by AliMUONPreClusterFinder
24 /// From a precluster a pixel array is built, and from this pixel array
25 /// a list of clusters is output (using AliMUONClusterSplitterMLEM).
27 /// \author Laurent Aphecetche (for the "new" C++ structure) and
28 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
29 //-----------------------------------------------------------------------------
31 #include "AliMUONClusterFinderMLEM.h"
33 #include "AliMUONCluster.h"
34 #include "AliMUONClusterSplitterMLEM.h"
35 #include "AliMUONVDigit.h"
36 #include "AliMUONPad.h"
37 #include "AliMUONPreClusterFinder.h"
39 #include "AliMpVPadIterator.h"
40 #include "AliMpVSegmentation.h"
41 #include "AliRunLoader.h"
42 #include "AliMUONVDigitStore.h"
43 #include <Riostream.h>
49 #include "AliCodeTimer.h"
52 ClassImp(AliMUONClusterFinderMLEM)
55 const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
56 //const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
57 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
58 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
59 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
61 TMinuit* AliMUONClusterFinderMLEM::fgMinuit = 0x0;
62 // Status flags for pads
63 const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
64 const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
65 const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
66 const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
67 const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
68 const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
70 //_____________________________________________________________________________
71 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
72 : AliMUONVClusterFinder(),
73 fPreClusterFinder(clusterFinder),
80 fPixArray(new TObjArray(20)),
89 fSegmentation[1] = fSegmentation[0] = 0x0;
91 fPadBeg[0] = fPadBeg[1] = fCathBeg;
93 if (!fgMinuit) fgMinuit = new TMinuit(8);
95 if (fPlot) fDebug = 1;
98 //_____________________________________________________________________________
99 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
102 delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
104 delete fPreClusterFinder;
106 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
107 fNClusters,fNAddVirtualPads));
110 //_____________________________________________________________________________
112 AliMUONClusterFinderMLEM::Prepare(const AliMpVSegmentation* segmentations[2],
113 const AliMUONVDigitStore& digitStore)
115 /// Prepare for clustering
116 // AliCodeTimerAuto("")
118 for ( Int_t i = 0; i < 2; ++i )
120 fSegmentation[i] = segmentations[i];
123 // Find out the DetElemId
126 TIter next(digitStore.CreateIterator());
127 AliMUONVDigit* d = static_cast<AliMUONVDigit*>(next());
131 fDetElemId = d->DetElemId();
134 if ( fDetElemId < 0 )
136 AliWarning("Could not find DE. Probably no digits at all ?");
141 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
142 fSplitter->SetDebug(fDebug);
144 // find out current event number, and reset the cluster number
145 fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
147 fClusterList.Delete();
149 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
151 return fPreClusterFinder->Prepare(segmentations,digitStore);
154 //_____________________________________________________________________________
156 AliMUONClusterFinderMLEM::NextCluster()
158 /// Return next cluster
159 // AliCodeTimerAuto("")
161 // if the list of clusters is not void, pick one from there
162 TObject* o = fClusterList.At(++fClusterNumber);
163 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
165 //FIXME : at this point, must check whether we've used all the digits
166 //from precluster : if not, let the preclustering know about those unused
167 //digits, so it can reuse them
169 // if the cluster list is exhausted, we need to go to the next
170 // pre-cluster and treat it
172 fPreCluster = fPreClusterFinder->NextCluster();
180 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
181 fClusterNumber = -1; //AZ
185 // WorkOnPreCluster may have used only part of the pads, so we check that
186 // now, and let the unused pads be reused by the preclustering...
188 Int_t mult = fPreCluster->Multiplicity();
189 for ( Int_t i = 0; i < mult; ++i )
191 AliMUONPad* pad = fPreCluster->Pad(i);
192 if ( !pad->IsUsed() )
194 fPreClusterFinder->UsePad(*pad);
198 return NextCluster();
201 //_____________________________________________________________________________
203 AliMUONClusterFinderMLEM::WorkOnPreCluster()
205 /// Starting from a precluster, builds a pixel array, and then
206 /// extract clusters from this array
208 // AliCodeTimerAuto("")
211 cout << " *** Event # " << AliRunLoader::GetRunLoader()->GetEventNumber()
212 << " det. elem.: " << fDetElemId << endl;
213 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
214 AliMUONPad* pad = fPreCluster->Pad(j);
215 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
216 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
217 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
221 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
222 if (!cluster) return kFALSE;
224 BuildPixArray(*cluster);
226 if ( fPixArray->GetLast() < 0 )
228 AliDebug(1,"No pixel for the above cluster");
233 // Use MLEM for cluster finder
234 Int_t nMax = 1, localMax[100], maxPos[100];
235 Double_t maxVal[100];
237 Int_t iSimple = 0, nInX = -1, nInY;
239 PadsInXandY(*cluster,nInX, nInY);
241 if (nInX < 4 && nInY < 4)
243 iSimple = 1; // simple cluster
247 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
249 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
250 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
254 for (Int_t i = 0; i < nMax; ++i)
258 FindCluster(*cluster,localMax, maxPos[i]);
261 MainLoop(*cluster,iSimple);
265 Int_t mult = cluster->Multiplicity();
266 for (Int_t j = 0; j < mult; ++j)
268 AliMUONPad* pad = cluster->Pad(j);
269 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
270 if ( pad->Status() != fgkOver ) continue; // pad was not used
272 pad->SetStatus(fgkZero);
273 pad->RevertCharge(); // use backup charge value
276 } // for (Int_t i=0; i<nMax;
277 //if (nMax > 1) ((TH2D*) gROOT->FindObject("anode"))->Delete();
278 delete ((TH2D*) gROOT->FindObject("anode"));
279 //TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
280 //if (mlem) mlem->Delete();
285 //_____________________________________________________________________________
287 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
289 /// Check if the pad and the pixel overlaps
291 // make a fake pad from the pixel
292 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
293 pixel.Coord(0),pixel.Coord(1),
294 pixel.Size(0),pixel.Size(1),0);
296 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
299 //_____________________________________________________________________________
301 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
303 /// Check precluster in order to attempt to simplify it (mostly for
304 /// two-cathode preclusters)
306 // AliCodeTimerAuto("")
308 // Disregard small clusters (leftovers from splitting or noise)
309 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
310 origCluster.Charge(0)+origCluster.Charge(1) < 10)
315 //AliMUONCluster* cluster = static_cast<AliMUONCluster*>(origCluster.Clone());
316 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
318 AliDebug(2,"Start of CheckPreCluster=");
319 //StdoutToAliDebug(2,cluster->Print("full"));
321 AliMUONCluster* rv(0x0);
323 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
325 rv = CheckPreclusterTwoCathodes(cluster);
334 //_____________________________________________________________________________
336 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
338 /// Check two-cathode cluster
340 Int_t npad = cluster->Multiplicity();
341 Int_t* flags = new Int_t[npad];
342 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
344 // Check pad overlaps
345 for ( Int_t i = 0; i < npad; ++i)
347 AliMUONPad* padi = cluster->Pad(i);
348 if ( padi->Cathode() != 0 ) continue;
349 for (Int_t j = i+1; j < npad; ++j)
351 AliMUONPad* padj = cluster->Pad(j);
352 if ( padj->Cathode() != 1 ) continue;
353 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
354 flags[i] = flags[j] = 1; // mark overlapped pads
358 // Check if all pads overlap
360 for (Int_t i = 0; i < npad; ++i)
362 if (!flags[i]) ++nFlags;
367 // not all pads overlap.
368 if (fDebug) cout << " nFlags: " << nFlags << endl;
369 TObjArray toBeRemoved;
370 for (Int_t i = 0; i < npad; ++i)
372 AliMUONPad* pad = cluster->Pad(i);
373 if (flags[i]) continue;
374 Int_t cath = pad->Cathode();
375 Int_t cath1 = TMath::Even(cath);
376 // Check for edge effect (missing pads on the _other_ cathode)
377 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
378 if (!mpPad.IsValid()) continue;
379 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
380 if (nFlags == 1 && pad->Charge() < 20) continue;
381 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
382 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
383 toBeRemoved.AddLast(pad);
384 fPreCluster->Pad(i)->Release();
386 Int_t nRemove = toBeRemoved.GetEntriesFast();
387 for ( Int_t i = 0; i < nRemove; ++i )
389 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
393 // Check correlations of cathode charges
394 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
397 Int_t cathode = cluster->MaxRawChargeCathode();
403 // get min and max pad charges on the cathode opposite to the
404 // max pad (given by MaxRawChargeCathode())
406 Int_t mult = cluster->Multiplicity();
407 for ( Int_t i = 0; i < mult; ++i )
409 AliMUONPad* pad = cluster->Pad(i);
410 if ( pad->Cathode() != cathode || !pad->IsReal() )
412 // only consider pads in the opposite cathode, and
413 // only consider real pads (i.e. exclude the virtual ones)
416 if ( pad->Charge() < cmin )
418 cmin = pad->Charge();
425 else if ( pad->Charge() > cmax )
427 cmax = pad->Charge();
431 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
432 imin,imax,cmin,cmax));
434 // arrange pads according to their distance to the max, normalized
436 Double_t* dist = new Double_t[mult];
441 AliMUONPad* padmax = cluster->Pad(imax);
443 for ( Int_t i = 0; i < mult; ++i )
446 if ( i == imax) continue;
447 AliMUONPad* pad = cluster->Pad(i);
448 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
449 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
450 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
451 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
454 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
460 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
461 Double_t xmax(-1), distPrev(999);
462 TObjArray toBeRemoved;
464 for ( Int_t i = 0; i < mult; ++i )
466 Int_t indx = flags[i];
467 AliMUONPad* pad = cluster->Pad(indx);
468 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
469 if ( dist[indx] > dmin )
471 // farther than the minimum pad
472 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
473 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
476 if (dx >= 0 && dy >= 0) continue;
477 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
478 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
480 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
481 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
484 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
486 cmax = TMath::Max(pad->Charge(),cmax);
490 cmax = pad->Charge();
493 distPrev = dist[indx];
494 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
495 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
498 toBeRemoved.AddLast(pad);
499 fPreCluster->Pad(indx)->Release();
502 Int_t nRemove = toBeRemoved.GetEntriesFast();
503 for ( Int_t i = 0; i < nRemove; ++i )
505 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
512 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
513 //StdoutToAliDebug(2,cluster->Print("full"));
518 //_____________________________________________________________________________
520 AliMUONClusterFinderMLEM::CheckOverlaps()
522 /// For debug only : check if some pixels overlap...
524 Int_t nPix = fPixArray->GetLast()+1;
527 for ( Int_t i = 0; i < nPix; ++i )
529 AliMUONPad* pixelI = Pixel(i);
530 AliMUONPad pi(dummy,dummy,dummy,dummy,
531 pixelI->Coord(0),pixelI->Coord(1),
532 pixelI->Size(0),pixelI->Size(1),0.0);
534 for ( Int_t j = i+1; j < nPix; ++j )
536 AliMUONPad* pixelJ = Pixel(j);
537 AliMUONPad pj(dummy,dummy,dummy,dummy,
538 pixelJ->Coord(0),pixelJ->Coord(1),
539 pixelJ->Size(0),pixelJ->Size(1),0.0);
542 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
544 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
546 StdoutToAliInfo(pixelI->Print();
547 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
549 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
550 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
551 cout << "-------" << endl;
559 //_____________________________________________________________________________
560 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
562 /// Build pixel array for MLEM method
564 Int_t npad = cluster.Multiplicity();
567 AliWarning("Got no pad at all ?!");
571 BuildPixArrayOneCathode(cluster);
573 Int_t nPix = fPixArray->GetLast()+1;
575 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
579 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
581 // Too many pixels - sort and remove pixels with the lowest signal
583 for ( Int_t i = npad; i < nPix; ++i )
587 fPixArray->Compress();
588 } // if (nPix > npad)
590 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
591 // fPixArray->Print(););
592 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
595 //_____________________________________________________________________________
596 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
598 /// Build the pixel array
600 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
602 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
603 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
604 Int_t found[2] = {0}, mult = cluster.Multiplicity();
606 for ( Int_t i = 0; i < mult; ++i) {
607 AliMUONPad* pad = cluster.Pad(i);
608 for (Int_t j = 0; j < 2; ++j) {
609 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
610 xy0[j] = pad->Coord(j);
614 if (found[0] && found[1]) break;
617 Double_t min[2], max[2];
618 Int_t cath0 = 0, cath1 = 1;
619 if (cluster.Multiplicity(0) == 0) cath0 = 1;
620 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
622 TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
623 TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
624 min[0] = leftDown.X();
625 min[1] = leftDown.Y();
626 max[0] = rightUp.X();
627 max[1] = rightUp.Y();
628 if (cath1 != cath0) {
629 leftDown = cluster.Area(cath1).LeftDownCorner();
630 rightUp = cluster.Area(cath1).RightUpCorner();
631 min[0] = TMath::Max (min[0], leftDown.X());
632 min[1] = TMath::Max (min[1], leftDown.Y());
633 max[0] = TMath::Min (max[0], rightUp.X());
634 max[1] = TMath::Min (max[1], rightUp.Y());
638 //width[0] /= 2; width[1] /= 2; // just for check
640 for (Int_t i = 0; i < 2; ++i) {
641 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
642 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
643 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
644 + TMath::Sign(0.5,dist)) * width[i] * 2;
645 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
646 if (nbins[i] == 0) ++nbins[i];
647 max[i] = min[i] + nbins[i] * width[i] * 2;
648 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
652 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
653 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
654 TAxis *xaxis = hist1->GetXaxis();
655 TAxis *yaxis = hist1->GetYaxis();
658 for ( Int_t i = 0; i < mult; ++i) {
659 AliMUONPad* pad = cluster.Pad(i);
660 Int_t ix0 = xaxis->FindBin(pad->X());
661 Int_t iy0 = yaxis->FindBin(pad->Y());
662 PadOverHist(0, ix0, iy0, pad);
666 for (Int_t i = 1; i <= nbins[0]; ++i) {
667 Double_t x = xaxis->GetBinCenter(i);
668 for (Int_t j = 1; j <= nbins[1]; ++j) {
669 if (hist2->GetCellContent(i,j) < 0.1) continue;
670 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
671 // cluster.Multiplicity(1)) continue;
672 if (cath0 != cath1) {
674 Double_t cont = hist2->GetCellContent(i,j);
675 if (cont < 999.) continue;
676 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
678 Double_t y = yaxis->GetBinCenter(j);
679 Double_t charge = hist1->GetCellContent(i,j);
680 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
681 fPixArray->Add(pixPtr);
685 if (fPixArray->GetEntriesFast() == 1) {
686 // Split pixel into 2
687 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
688 pixPtr->SetSize(0,width[0]/2.);
689 pixPtr->Shift(0,-width[0]/4.);
690 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
691 fPixArray->Add(pixPtr);
694 //fPixArray->Print();
699 //_____________________________________________________________________________
700 void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad)
702 /// "Span" pad over histogram in the direction idir
704 TH2D *hist1 = static_cast<TH2D*> (gROOT->FindObject("Grid"));
705 TH2D *hist2 = static_cast<TH2D*> (gROOT->FindObject("Entries"));
706 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
707 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
708 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
710 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
712 for (Int_t i = 0; i < nbinPad; ++i) {
713 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
714 if (ixy > nbins) break;
715 Double_t lowEdge = axis->GetBinLowEdge(ixy);
716 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
717 if (idir == 0) PadOverHist(1, ixy, iy0, pad); // span in the other direction
720 Double_t cont = pad->Charge();
721 if (hist2->GetCellContent(ix0, ixy) > 0.1)
722 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
723 hist1->SetCellContent(ix0, ixy, cont);
724 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
725 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
729 for (Int_t i = -1; i > -nbinPad; --i) {
730 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
732 Double_t upEdge = axis->GetBinUpEdge(ixy);
733 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
734 if (idir == 0) PadOverHist(1, ixy, iy0, pad); // span in the other direction
737 Double_t cont = pad->Charge();
738 if (hist2->GetCellContent(ix0, ixy) > 0.1)
739 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
740 hist1->SetCellContent(ix0, ixy, cont);
741 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
742 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
747 //_____________________________________________________________________________
749 AliMUONClusterFinderMLEM::Plot(const char* basename)
751 /// Make a plot and save it as png
756 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
761 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
762 fDetElemId,fClusterNumber));
765 //_____________________________________________________________________________
767 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
771 /// Compute coefficients needed for MLEM algorithm
773 Int_t nPix = fPixArray->GetLast()+1;
775 //memset(probi,0,nPix*sizeof(Double_t));
776 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
778 Int_t mult = cluster.Multiplicity();
779 for ( Int_t j = 0; j < mult; ++j )
781 AliMUONPad* pad = cluster.Pad(j);
784 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
786 Int_t indx1 = indx + ipix;
787 //if (pad->Status() < 0)
788 if (pad->Status() != fgkZero)
793 AliMUONPad* pixPtr = Pixel(ipix);
794 // coef is the charge (given by Mathieson integral) on pad, assuming
795 // the Mathieson is center at pixel.
796 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
797 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
798 // pad->Ix(),pad->Iy(),
799 // pad->X(),pad->Y(),
800 // pad->DX(),pad->DY(),
801 // pixPtr->Coord(0),pixPtr->Coord(1),
802 // pixPtr->Size(0),pixPtr->Size(1),
805 probi[ipix] += coef[indx1];
810 //_____________________________________________________________________________
811 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
813 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
815 // AliCodeTimerAuto("")
817 Int_t nPix = fPixArray->GetLast()+1;
819 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
820 //StdoutToAliDebug(2,cluster.Print("full"););
824 AliDebug(1,"No pixels, why am I here then ?");
827 AddVirtualPad(cluster); // add virtual pads if necessary
829 Int_t npadTot = cluster.Multiplicity();
831 for (Int_t i = 0; i < npadTot; ++i)
833 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
834 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
839 Double_t* probi(0x0);
840 Int_t lc(0); // loop counter
842 //Plot("mlem.start");
843 AliMUONPad* pixPtr = Pixel(0);
844 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
849 mlem = (TH2D*) gROOT->FindObject("mlem");
852 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
853 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
854 //StdoutToAliDebug(2,fPixArray->Print("","full"));
856 coef = new Double_t [npadTot*nPix];
857 probi = new Double_t [nPix];
859 // Calculate coefficients and pixel visibilities
860 ComputeCoefficients(cluster,coef,probi);
862 for (Int_t ipix = 0; ipix < nPix; ++ipix)
864 if (probi[ipix] < 0.01)
866 AliMUONPad* pixel = Pixel(ipix);
867 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
868 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
869 pixel->SetCharge(0); // "invisible" pixel
874 Mlem(cluster,coef, probi, 15);
876 // Find histogram limits for the 1'st pass only - for others computed below
878 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
880 pixPtr = Pixel(ipix);
881 for ( Int_t i = 0; i < 2; ++i )
884 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
885 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
888 } else pixPtr = Pixel(0);
890 for (Int_t i = 0; i < 4; i++)
892 xylim[i] -= pixPtr->Size(i/2);
895 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
896 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
898 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
899 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
900 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
901 xylim[0],-xylim[1],xylim[2],-xylim[3]
904 mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
906 for (Int_t ipix = 0; ipix < nPix; ++ipix)
908 AliMUONPad* pixPtr = Pixel(ipix);
909 mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
912 // Check if the total charge of pixels is too low
914 for ( Int_t i = 0; i < nPix; ++i)
916 qTot += Pixel(i)->Charge();
919 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
921 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
925 for ( Int_t i = 0; i < npadTot; ++i)
927 AliMUONPad* pad = cluster.Pad(i);
928 //if ( pad->Status() == 0) pad->SetStatus(-1);
929 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
936 // Simple cluster - skip further passes thru EM-procedure
944 // Calculate position of the center-of-gravity around the maximum pixel
946 FindCOG(mlem, xyCOG);
948 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
949 pixPtr->Size(0) > pixPtr->Size(1)) break;
951 // Sort pixels according to the charge
952 MaskPeaks(1); // mask local maxima
954 MaskPeaks(0); // unmask local maxima
955 Double_t pixMin = 0.01*Pixel(0)->Charge();
956 pixMin = TMath::Min(pixMin,50.);
958 // Decrease pixel size and shift pixels to make them centered at
960 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
963 Double_t shift[2] = { 0.0, 0.0 };
964 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
967 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
969 AliMUONPad* pixPtr = Pixel(ipix);
970 if ( nPix >= npadOK // too many pixels already
972 pixPtr->Charge() < pixMin && pixPtr->Status() != fgkMustKeep // too low charge
978 for (Int_t i = 0; i < 2; ++i)
982 pixPtr->SetCharge(10);
983 pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
984 width = -pixPtr->Size(indx);
985 pixPtr->Shift(indx, width);
986 // Shift pixel position
990 for (Int_t j = 0; j < 2; ++j)
992 shift[j] = pixPtr->Coord(j) - xyCOG[j];
993 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
996 pixPtr->Shift(0, -shift[0]);
997 pixPtr->Shift(1, -shift[1]);
1000 else if (nPix < npadOK)
1002 pixPtr = new AliMUONPad(*pixPtr);
1003 pixPtr->Shift(indx, -2*width);
1004 pixPtr->SetStatus(fgkZero);
1005 fPixArray->Add(pixPtr);
1008 else continue; // skip adjustment of histo limits
1009 for (Int_t j = 0; j < 4; ++j)
1011 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr->Coord(j/2));
1013 } // for (Int_t i=0; i<2;
1014 } // for (Int_t ipix=0;
1016 fPixArray->Compress();
1018 AliDebug(2,Form("After shift:"));
1019 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1020 //Plot(Form("mlem.lc%d",lc+1));
1022 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1025 xylim[2],xylim[3]));
1029 AliMUONPad* pixPtr = Pixel(0);
1030 // add pixels if the maximum is at the limit of pixel area:
1031 // start from Y-direction
1033 for (Int_t i = 3; i > -1; --i)
1035 if (nPix < npadOK &&
1036 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2))
1038 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1039 AliMUONPad* p = new AliMUONPad(*pixPtr);
1040 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1041 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1042 j = TMath::Even (i/2);
1043 p->SetCoord(j, xyCOG[j]);
1044 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1045 //StdoutToAliDebug(2,cout << " ---- ";
1046 // p->Print("corners"););
1052 nPix = fPixArray->GetEntriesFast();
1057 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1058 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1060 // remove pixels with low signal or low visibility
1061 // Cuts are empirical !!!
1062 Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1063 thresh = TMath::Min (thresh,50.);
1064 Double_t charge = 0;
1066 // Mark pixels which should be removed
1067 for (Int_t i = 0; i < nPix; ++i)
1069 AliMUONPad* pixPtr = Pixel(i);
1070 charge = pixPtr->Charge();
1071 if (charge < thresh)
1073 pixPtr->SetCharge(-charge);
1077 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1079 for (Int_t i = 0; i < nPix; ++i)
1081 AliMUONPad* pixPtr = Pixel(i);
1082 charge = pixPtr->Charge();
1083 if (charge > 0) continue;
1084 near = FindNearest(pixPtr);
1085 pixPtr->SetCharge(0);
1086 probi[i] = 0; // make it "invisible"
1087 AliMUONPad* pnear = Pixel(near);
1088 pnear->SetCharge(pnear->Charge() + (-charge));
1090 Mlem(cluster,coef,probi,2);
1092 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1093 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1094 //Plot("mlem.beforesplit");
1097 for (Int_t i = 0; i < nPix; ++i)
1099 AliMUONPad* pixPtr = Pixel(i);
1100 Int_t ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1101 Int_t iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1102 mlem->SetBinContent(ix, iy, pixPtr->Charge());
1105 // Try to split into clusters
1107 if (mlem->GetSum() < 1)
1113 fSplitter->Split(cluster,mlem,coef,fClusterList);
1118 fPixArray->Delete();
1123 //_____________________________________________________________________________
1124 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1126 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1127 /// - to avoid loosing low charge pixels after sorting)
1129 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1130 AliMUONPad* pix = Pixel(i);
1131 if (pix->Status() == fgkMustKeep) {
1132 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1133 else pix->SetCharge(pix->Charge()-10000.);
1138 //_____________________________________________________________________________
1139 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1140 Double_t* coef, Double_t* probi,
1143 /// Use MLEM to find pixel charges
1145 Int_t nPix = fPixArray->GetEntriesFast();
1147 Int_t npad = cluster.Multiplicity();
1149 Double_t* probi1 = new Double_t[nPix];
1150 Double_t probMax = TMath::MaxElement(nPix,probi);
1152 for (Int_t iter = 0; iter < nIter; ++iter)
1155 for (Int_t ipix = 0; ipix < nPix; ++ipix)
1157 Pixel(ipix)->SetChargeBackup(0);
1158 // Correct each pixel
1159 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1161 probi1[ipix] = probMax;
1162 for (Int_t j = 0; j < npad; ++j)
1164 AliMUONPad* pad = cluster.Pad(j);
1165 //if (pad->Status() < 0) continue;
1166 if (pad->Status() != fgkZero) continue;
1168 Int_t indx1 = j*nPix;
1169 Int_t indx = indx1 + ipix;
1170 // Calculate expectation
1171 for (Int_t i = 0; i < nPix; ++i)
1173 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1174 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1176 if ( pad->IsSaturated() && sum1 > pad->Charge() )
1178 // correct for pad charge overflows
1179 probi1[ipix] -= coef[indx];
1181 //sum1 = pad->Charge();
1184 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1185 //if (coef[indx] > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1186 } // for (Int_t j=0;
1187 AliMUONPad* pixPtr = Pixel(ipix);
1188 if (probi1[ipix] > 1.e-6)
1190 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1191 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1193 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1194 } // for (Int_t ipix=0;
1196 for (Int_t i = 0; i < nPix; ++i) {
1197 AliMUONPad* pixPtr = Pixel(i);
1198 pixPtr->RevertCharge();
1199 qTot += pixPtr->Charge();
1202 // Can happen in clusters with large number of overflows - speeding up
1206 } // for (Int_t iter=0;
1210 //_____________________________________________________________________________
1211 void AliMUONClusterFinderMLEM::FindCOG(TH2D *mlem, Double_t *xyc)
1213 /// Calculate position of the center-of-gravity around the maximum pixel
1215 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1216 Int_t i1 = -9, j1 = -9;
1217 mlem->GetMaximumBin(ixmax,iymax,ix);
1218 Int_t nx = mlem->GetNbinsX();
1219 Int_t ny = mlem->GetNbinsY();
1220 Double_t thresh = mlem->GetMaximum()/10;
1221 Double_t x, y, cont, xq=0, yq=0, qq=0;
1223 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1224 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1225 y = mlem->GetYaxis()->GetBinCenter(i);
1226 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1227 cont = mlem->GetCellContent(j,i);
1228 if (cont < thresh) continue;
1229 if (i != i1) {i1 = i; nsumy++;}
1230 if (j != j1) {j1 = j; nsumx++;}
1231 x = mlem->GetXaxis()->GetBinCenter(j);
1240 Int_t i2 = 0, j2 = 0;
1243 // one bin in Y - add one more (with the largest signal)
1244 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1245 if (i == iymax) continue;
1246 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1247 cont = mlem->GetCellContent(j,i);
1250 x = mlem->GetXaxis()->GetBinCenter(j);
1251 y = mlem->GetYaxis()->GetBinCenter(i);
1260 if (i2 != i1) nsumy++;
1261 if (j2 != j1) nsumx++;
1263 } // if (nsumy == 1)
1266 // one bin in X - add one more (with the largest signal)
1268 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1269 if (j == ixmax) continue;
1270 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1271 cont = mlem->GetCellContent(j,i);
1274 x = mlem->GetXaxis()->GetBinCenter(j);
1275 y = mlem->GetYaxis()->GetBinCenter(i);
1284 if (i2 != i1) nsumy++;
1285 if (j2 != j1) nsumx++;
1287 } // if (nsumx == 1)
1289 xyc[0] = xq/qq; xyc[1] = yq/qq;
1290 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1293 //_____________________________________________________________________________
1294 Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1296 /// Find the pixel nearest to the given one
1297 /// (algorithm may be not very efficient)
1299 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1300 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1301 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1304 for (Int_t i = 0; i < nPix; ++i) {
1305 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1306 if (pixPtr == pixPtr0 || pixPtr->Charge() < 0.5) continue;
1307 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1308 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1309 r = dx *dx + dy * dy;
1310 if (r < rmin) { rmin = r; imin = i; }
1315 //_____________________________________________________________________________
1317 AliMUONClusterFinderMLEM::Paint(Option_t*)
1319 /// Paint cluster and pixels
1321 AliMpArea area(fPreCluster->Area());
1323 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1325 gVirtualX->SetFillStyle(1001);
1326 gVirtualX->SetFillColor(3);
1327 gVirtualX->SetLineColor(3);
1331 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1333 AliMUONPad* pixel = Pixel(i);
1335 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1336 pixel->Coord(1)-pixel->Size(1)*s,
1337 pixel->Coord(0)+pixel->Size(0)*s,
1338 pixel->Coord(1)+pixel->Size(1)*s);
1340 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1342 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1343 // pixel->Coord(1) + sign*pixel->Size(1),
1344 // pixel->Coord(0) + pixel->Size(0),
1345 // pixel->Coord(1) - sign*pixel->Size(1)
1351 gVirtualX->SetFillStyle(0);
1353 fPreCluster->Paint();
1355 gVirtualX->SetLineColor(1);
1356 gVirtualX->SetLineWidth(2);
1357 gVirtualX->SetFillStyle(0);
1358 gVirtualX->SetTextColor(1);
1359 gVirtualX->SetTextAlign(22);
1361 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1363 AliMUONPad* pixel = Pixel(i);
1364 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1365 pixel->Coord(1)-pixel->Size(1),
1366 pixel->Coord(0)+pixel->Size(0),
1367 pixel->Coord(1)+pixel->Size(1));
1368 gVirtualX->SetTextSize(pixel->Size(0)*60);
1370 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1374 //_____________________________________________________________________________
1375 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1377 /// Find local maxima in pixel space for large preclusters in order to
1378 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1379 /// or to find additional fitting seeds if clusters were not completely resolved
1381 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1384 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
1385 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
1386 //if (hist) hist->Delete();
1388 Double_t xylim[4] = {999, 999, 999, 999};
1390 Int_t nPix = pixArray->GetEntriesFast();
1391 AliMUONPad *pixPtr = 0;
1392 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1393 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1394 for (Int_t i = 0; i < 4; ++i)
1395 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1397 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
1399 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1400 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1401 if (pixArray == fPixArray) hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1402 else hist = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1403 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1404 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1405 hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1407 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1409 Int_t nMax = 0, indx, nxy = ny * nx;
1410 Int_t *isLocalMax = new Int_t[nxy];
1411 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
1413 for (Int_t i = 1; i <= ny; ++i) {
1415 for (Int_t j = 1; j <= nx; ++j) {
1416 if (hist->GetCellContent(j,i) < 0.5) continue;
1417 //if (isLocalMax[indx+j-1] < 0) continue;
1418 if (isLocalMax[indx+j-1] != 0) continue;
1419 FlagLocalMax(hist, i, j, isLocalMax);
1423 for (Int_t i = 1; i <= ny; ++i) {
1425 for (Int_t j = 1; j <= nx; ++j) {
1426 if (isLocalMax[indx+j-1] > 0) {
1427 localMax[nMax] = indx + j - 1;
1428 maxVal[nMax++] = hist->GetCellContent(j,i);
1429 ((AliMUONPad*)fSplitter->BinToPix(hist, j, i))->SetStatus(fgkMustKeep);
1430 if (nMax > 99) AliFatal(" Too many local maxima !!!");
1434 if (fDebug) cout << " Local max: " << nMax << endl;
1435 delete [] isLocalMax;
1439 //_____________________________________________________________________________
1440 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1442 /// Flag pixels (whether or not local maxima)
1444 Int_t nx = hist->GetNbinsX();
1445 Int_t ny = hist->GetNbinsY();
1446 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1447 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1449 Int_t ie = i + 2, je = j + 2;
1450 for (Int_t i1 = i-1; i1 < ie; ++i1) {
1451 if (i1 < 1 || i1 > ny) continue;
1452 indx1 = (i1 - 1) * nx;
1453 for (Int_t j1 = j-1; j1 < je; ++j1) {
1454 if (j1 < 1 || j1 > nx) continue;
1455 if (i == i1 && j == j1) continue;
1456 indx2 = indx1 + j1 - 1;
1457 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1458 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1459 else if (cont > cont1) isLocalMax[indx2] = -1;
1460 else { // the same charge
1461 isLocalMax[indx] = 1;
1462 if (isLocalMax[indx2] == 0) {
1463 FlagLocalMax(hist, i1, j1, isLocalMax);
1464 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1465 else isLocalMax[indx2] = -1;
1470 isLocalMax[indx] = 1; // local maximum
1473 //_____________________________________________________________________________
1474 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1475 Int_t *localMax, Int_t iMax)
1477 /// Find pixel cluster around local maximum \a iMax and pick up pads
1478 /// overlapping with it
1480 TH2D *hist = (TH2D*) gROOT->FindObject("anode");
1482 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1484 hist->Draw("lego1Fb"); // debug
1489 Int_t nx = hist->GetNbinsX();
1490 Int_t ny = hist->GetNbinsY();
1491 Int_t ic = localMax[iMax] / nx + 1;
1492 Int_t jc = localMax[iMax] % nx + 1;
1493 Int_t nxy = ny * nx;
1494 Bool_t *used = new Bool_t[nxy];
1495 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
1497 // Drop all pixels from the array - pick up only the ones from the cluster
1498 fPixArray->Delete();
1500 Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2;
1501 Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;
1502 Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
1503 Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
1504 Double_t cont = hist->GetCellContent(jc,ic);
1505 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1506 used[(ic-1)*nx+jc-1] = kTRUE;
1507 AddBinSimple(hist, ic, jc);
1508 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1510 Int_t nPix = fPixArray->GetEntriesFast();
1511 Int_t npad = cluster.Multiplicity();
1513 for (Int_t i = 0; i < nPix; ++i)
1515 AliMUONPad* pixPtr = Pixel(i);
1516 pixPtr->SetSize(0,wx);
1517 pixPtr->SetSize(1,wy);
1520 // Pick up pads which overlap with found pixels
1521 for (Int_t i = 0; i < npad; ++i)
1523 //cluster.Pad(i)->SetStatus(-1);
1524 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1527 for (Int_t i = 0; i < nPix; ++i)
1529 AliMUONPad* pixPtr = Pixel(i);
1530 for (Int_t j = 0; j < npad; ++j)
1532 AliMUONPad* pad = cluster.Pad(j);
1533 //if (pad->Status() == 0) continue;
1534 if (pad->Status() == fgkZero) continue;
1535 if ( Overlap(*pad,*pixPtr) )
1537 //pad->SetStatus(0);
1538 pad->SetStatus(fgkZero);
1539 if (fDebug) { cout << j << " "; pad->Print("full"); }
1547 //_____________________________________________________________________________
1549 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *mlem, Int_t ic, Int_t jc)
1551 /// Add adjacent bins (+-1 in X and Y) to the cluster
1553 Int_t nx = mlem->GetNbinsX();
1554 Int_t ny = mlem->GetNbinsY();
1555 Double_t cont1, cont = mlem->GetCellContent(jc,ic);
1556 AliMUONPad *pixPtr = 0;
1558 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1559 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1560 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1561 cont1 = mlem->GetCellContent(j,i);
1562 if (cont1 > cont) continue;
1563 if (cont1 < 0.5) continue;
1564 pixPtr = new AliMUONPad (mlem->GetXaxis()->GetBinCenter(j),
1565 mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1566 fPixArray->Add(pixPtr);
1571 //_____________________________________________________________________________
1572 AliMUONClusterFinderMLEM&
1573 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1575 /// Protected assignement operator
1577 if (this == &rhs) return *this;
1579 AliFatal("Not implemented.");
1584 //_____________________________________________________________________________
1585 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1587 /// Add virtual pad (with small charge) to improve fit for clusters
1588 /// with number of pads == 2 per direction
1590 // Find out non-bending and bending planes
1591 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1593 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1594 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1595 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1600 Bool_t same = kFALSE;
1601 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1604 Bool_t check[2] = {kFALSE, kFALSE};
1606 nxy[0] = nxy[1] = 0;
1607 for (Int_t inb = 0; inb < 2; ++inb) {
1608 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1609 if (inb == 0 && cn.GetFirst() == 2) check[inb] = kTRUE; // check non-bending plane
1610 else if (inb == 1 && cn.GetSecond() == 2) check[inb] = kTRUE; // check bending plane
1612 nxy[0] = TMath::Max (nxy[0], cn.GetFirst());
1613 nxy[1] = TMath::Max (nxy[1], cn.GetSecond());
1614 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1615 else if (inb == 1 && cn.GetSecond() < 2) nonb[inb] = !nonb[inb];
1619 if (nxy[0] > 2) check[0] = kFALSE;
1620 if (nxy[1] > 2) check[1] = kFALSE;
1622 if (!check[0] && !check[1]) return;
1624 for (Int_t inb = 0; inb < 2; ++inb) {
1625 if (!check[inb]) continue;
1627 // Find pads with maximum and next to maximum charges
1628 Int_t maxPads[2] = {-1, -1};
1629 Double_t amax[2] = {0};
1630 Int_t mult = cluster.Multiplicity();
1631 for (Int_t j = 0; j < mult; ++j) {
1632 AliMUONPad *pad = cluster.Pad(j);
1633 if (pad->Cathode() != nonb[inb]) continue;
1634 for (Int_t j2 = 0; j2 < 2; ++j2) {
1635 if (pad->Charge() > amax[j2]) {
1636 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1637 amax[j2] = pad->Charge();
1644 // Find min and max dimensions of the cluster
1645 Double_t limits[2] = {9999, -9999};
1646 for (Int_t j = 0; j < mult; ++j) {
1647 AliMUONPad *pad = cluster.Pad(j);
1648 if (pad->Cathode() != nonb[inb]) continue;
1649 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1650 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1653 // Loop over max and next to max pads
1654 Bool_t add = kFALSE;
1656 for (Int_t j = 0; j < 2; ++j) {
1658 if (maxPads[j] < 0) continue;
1660 if (amax[1] / amax[0] < 0.5) break;
1662 // Check if pad at the cluster limit
1663 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1665 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1666 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1668 //cout << " *** Pad not at the cluster limit: " << j << endl;
1671 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1673 // Add pad (if it exists)
1675 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1676 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1677 //AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1678 AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kFALSE);
1679 if (!mppad.IsValid()) continue; // non-existing pad
1680 cn = mppad.GetIndices();
1681 AliMUONPad muonPad(fDetElemId, nonb[inb], cn.GetFirst(), cn.GetSecond(),
1682 mppad.Position().X(), mppad.Position().Y(),
1683 mppad.Dimensions().X(), mppad.Dimensions().Y(), 0);
1684 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, 5.));
1685 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1686 else muonPad.SetCharge(TMath::Min (amax[j]/15, 6.));
1687 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1688 muonPad.SetReal(kFALSE);
1689 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1690 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1691 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1692 cluster.AddPad(muonPad); // add pad to the cluster
1699 //_____________________________________________________________________________
1700 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1701 Int_t &nInX, Int_t &nInY) const
1703 /// Find number of pads in X and Y-directions (excluding virtual ones and
1706 //Int_t statusToTest = 1;
1707 Int_t statusToTest = fgkUseForFit;
1709 //if ( nInX < 0 ) statusToTest = 0;
1710 if ( nInX < 0 ) statusToTest = fgkZero;
1712 Bool_t mustMatch(kTRUE);
1714 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1716 nInX = cn.GetFirst();
1717 nInY = cn.GetSecond();
1720 //_____________________________________________________________________________
1721 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1723 /// Remove pixel at index i
1724 AliMUONPad* pixPtr = Pixel(i);
1725 fPixArray->RemoveAt(i);
1729 //_____________________________________________________________________________
1730 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1732 /// Process simple cluster (small number of pads) without EM-procedure
1734 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1735 Double_t parOk[3] = {0.};
1736 TObjArray *clusters[1];
1737 clusters[0] = fPixArray;
1739 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1741 Int_t mult = cluster.Multiplicity();
1742 for (Int_t i = 0; i < mult; ++i)
1744 AliMUONPad* pad = cluster.Pad(i);
1746 if ( pad->IsSaturated())
1755 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1757 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList);
1760 //_____________________________________________________________________________
1762 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1764 /// Returns pixel at index i
1765 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1768 //_____________________________________________________________________________
1770 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1773 TString swhat(what);
1775 if ( swhat.Contains("precluster") )
1777 if ( fPreCluster) fPreCluster->Print();