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>
48 #include "AliCodeTimer.h"
51 ClassImp(AliMUONClusterFinderMLEM)
54 const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
55 //const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
56 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
57 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
58 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
60 // Status flags for pads
61 const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
62 const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
63 const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
64 const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
65 const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
66 const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
68 //_____________________________________________________________________________
69 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
70 : AliMUONVClusterFinder(),
71 fPreClusterFinder(clusterFinder),
79 fPixArray(new TObjArray(20)),
88 fSegmentation[1] = fSegmentation[0] = 0x0;
90 if (fPlot) fDebug = 1;
93 //_____________________________________________________________________________
94 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
97 delete fPixArray; fPixArray = 0;
99 delete fPreClusterFinder;
101 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
102 fNClusters,fNAddVirtualPads));
105 //_____________________________________________________________________________
107 AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
108 TClonesArray* pads[2],
109 const AliMpArea& area,
110 const AliMpVSegmentation* seg[2])
112 /// Prepare for clustering
113 // AliCodeTimerAuto("")
115 for ( Int_t i = 0; i < 2; ++i )
117 fSegmentation[i] = seg[i];
120 // Find out the DetElemId
121 fDetElemId = detElemId;
124 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
125 fSplitter->SetDebug(fDebug);
127 // find out current event number, and reset the cluster number
128 fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
130 fClusterList.Delete();
132 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
134 if ( fPreClusterFinder->NeedSegmentation() )
136 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
140 return fPreClusterFinder->Prepare(detElemId,pads,area);
144 //_____________________________________________________________________________
146 AliMUONClusterFinderMLEM::NextCluster()
148 /// Return next cluster
149 // AliCodeTimerAuto("")
151 // if the list of clusters is not void, pick one from there
152 TObject* o = fClusterList.At(++fClusterNumber);
153 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
155 //FIXME : at this point, must check whether we've used all the digits
156 //from precluster : if not, let the preclustering know about those unused
157 //digits, so it can reuse them
159 // if the cluster list is exhausted, we need to go to the next
160 // pre-cluster and treat it
162 fPreCluster = fPreClusterFinder->NextCluster();
170 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
171 fClusterNumber = -1; //AZ
175 // WorkOnPreCluster may have used only part of the pads, so we check that
176 // now, and let the unused pads be reused by the preclustering...
178 Int_t mult = fPreCluster->Multiplicity();
179 for ( Int_t i = 0; i < mult; ++i )
181 AliMUONPad* pad = fPreCluster->Pad(i);
182 if ( !pad->IsUsed() )
184 fPreClusterFinder->UsePad(*pad);
188 return NextCluster();
191 //_____________________________________________________________________________
193 AliMUONClusterFinderMLEM::WorkOnPreCluster()
195 /// Starting from a precluster, builds a pixel array, and then
196 /// extract clusters from this array
198 // AliCodeTimerAuto("")
201 cout << " *** Event # " << AliRunLoader::GetRunLoader()->GetEventNumber()
202 << " det. elem.: " << fDetElemId << endl;
203 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
204 AliMUONPad* pad = fPreCluster->Pad(j);
205 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
206 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
207 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
211 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
212 if (!cluster) return kFALSE;
214 BuildPixArray(*cluster);
216 if ( fPixArray->GetLast() < 0 )
218 AliDebug(1,"No pixel for the above cluster");
223 // Use MLEM for cluster finder
224 Int_t nMax = 1, localMax[100], maxPos[100];
225 Double_t maxVal[100];
227 Int_t iSimple = 0, nInX = -1, nInY;
229 PadsInXandY(*cluster,nInX, nInY);
231 if (nInX < 4 && nInY < 4)
233 iSimple = 1; // simple cluster
237 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
239 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
240 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
244 for (Int_t i = 0; i < nMax; ++i)
248 FindCluster(*cluster,localMax, maxPos[i]);
251 MainLoop(*cluster,iSimple);
255 Int_t mult = cluster->Multiplicity();
256 for (Int_t j = 0; j < mult; ++j)
258 AliMUONPad* pad = cluster->Pad(j);
259 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
260 if ( pad->Status() != fgkOver ) continue; // pad was not used
262 pad->SetStatus(fgkZero);
263 pad->RevertCharge(); // use backup charge value
266 } // for (Int_t i=0; i<nMax;
269 fHistMlem = fHistAnode = 0x0;
274 //_____________________________________________________________________________
276 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
278 /// Check if the pad and the pixel overlaps
280 // make a fake pad from the pixel
281 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
282 pixel.Coord(0),pixel.Coord(1),
283 pixel.Size(0),pixel.Size(1),0);
285 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
288 //_____________________________________________________________________________
290 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
292 /// Check precluster in order to attempt to simplify it (mostly for
293 /// two-cathode preclusters)
297 // Disregard small clusters (leftovers from splitting or noise)
298 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
299 origCluster.Charge(0)+origCluster.Charge(1) < 10)
304 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
306 AliDebug(2,"Start of CheckPreCluster=");
307 //StdoutToAliDebug(2,cluster->Print("full"));
309 AliMUONCluster* rv(0x0);
311 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
313 rv = CheckPreclusterTwoCathodes(cluster);
322 //_____________________________________________________________________________
324 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
326 /// Check two-cathode cluster
328 Int_t npad = cluster->Multiplicity();
329 Int_t* flags = new Int_t[npad];
330 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
332 // Check pad overlaps
333 for ( Int_t i = 0; i < npad; ++i)
335 AliMUONPad* padi = cluster->Pad(i);
336 if ( padi->Cathode() != 0 ) continue;
337 for (Int_t j = i+1; j < npad; ++j)
339 AliMUONPad* padj = cluster->Pad(j);
340 if ( padj->Cathode() != 1 ) continue;
341 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
342 flags[i] = flags[j] = 1; // mark overlapped pads
346 // Check if all pads overlap
348 for (Int_t i = 0; i < npad; ++i)
350 if (!flags[i]) ++nFlags;
355 // not all pads overlap.
356 if (fDebug) cout << " nFlags: " << nFlags << endl;
357 TObjArray toBeRemoved;
358 for (Int_t i = 0; i < npad; ++i)
360 AliMUONPad* pad = cluster->Pad(i);
361 if (flags[i]) continue;
362 Int_t cath = pad->Cathode();
363 Int_t cath1 = TMath::Even(cath);
364 // Check for edge effect (missing pads on the _other_ cathode)
365 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
366 if (!mpPad.IsValid()) continue;
367 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
368 if (nFlags == 1 && pad->Charge() < 20) continue;
369 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
370 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
371 toBeRemoved.AddLast(pad);
372 fPreCluster->Pad(i)->Release();
374 Int_t nRemove = toBeRemoved.GetEntriesFast();
375 for ( Int_t i = 0; i < nRemove; ++i )
377 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
381 // Check correlations of cathode charges
382 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
385 Int_t cathode = cluster->MaxRawChargeCathode();
391 // get min and max pad charges on the cathode opposite to the
392 // max pad (given by MaxRawChargeCathode())
394 Int_t mult = cluster->Multiplicity();
395 for ( Int_t i = 0; i < mult; ++i )
397 AliMUONPad* pad = cluster->Pad(i);
398 if ( pad->Cathode() != cathode || !pad->IsReal() )
400 // only consider pads in the opposite cathode, and
401 // only consider real pads (i.e. exclude the virtual ones)
404 if ( pad->Charge() < cmin )
406 cmin = pad->Charge();
413 else if ( pad->Charge() > cmax )
415 cmax = pad->Charge();
419 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
420 imin,imax,cmin,cmax));
422 // arrange pads according to their distance to the max, normalized
424 Double_t* dist = new Double_t[mult];
429 AliMUONPad* padmax = cluster->Pad(imax);
431 for ( Int_t i = 0; i < mult; ++i )
434 if ( i == imax) continue;
435 AliMUONPad* pad = cluster->Pad(i);
436 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
437 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
438 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
439 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
442 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
448 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
449 Double_t xmax(-1), distPrev(999);
450 TObjArray toBeRemoved;
452 for ( Int_t i = 0; i < mult; ++i )
454 Int_t indx = flags[i];
455 AliMUONPad* pad = cluster->Pad(indx);
456 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
457 if ( dist[indx] > dmin )
459 // farther than the minimum pad
460 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
461 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
464 if (dx >= 0 && dy >= 0) continue;
465 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
466 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
468 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
469 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
472 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
474 cmax = TMath::Max(pad->Charge(),cmax);
478 cmax = pad->Charge();
481 distPrev = dist[indx];
482 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
483 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
486 toBeRemoved.AddLast(pad);
487 fPreCluster->Pad(indx)->Release();
490 Int_t nRemove = toBeRemoved.GetEntriesFast();
491 for ( Int_t i = 0; i < nRemove; ++i )
493 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
500 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
501 //StdoutToAliDebug(2,cluster->Print("full"));
506 //_____________________________________________________________________________
508 AliMUONClusterFinderMLEM::CheckOverlaps()
510 /// For debug only : check if some pixels overlap...
512 Int_t nPix = fPixArray->GetLast()+1;
515 for ( Int_t i = 0; i < nPix; ++i )
517 AliMUONPad* pixelI = Pixel(i);
518 AliMUONPad pi(dummy,dummy,dummy,dummy,
519 pixelI->Coord(0),pixelI->Coord(1),
520 pixelI->Size(0),pixelI->Size(1),0.0);
522 for ( Int_t j = i+1; j < nPix; ++j )
524 AliMUONPad* pixelJ = Pixel(j);
525 AliMUONPad pj(dummy,dummy,dummy,dummy,
526 pixelJ->Coord(0),pixelJ->Coord(1),
527 pixelJ->Size(0),pixelJ->Size(1),0.0);
530 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
532 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
534 StdoutToAliInfo(pixelI->Print();
535 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
537 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
538 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
539 cout << "-------" << endl;
547 //_____________________________________________________________________________
548 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
550 /// Build pixel array for MLEM method
552 Int_t npad = cluster.Multiplicity();
555 AliWarning("Got no pad at all ?!");
559 BuildPixArrayOneCathode(cluster);
561 Int_t nPix = fPixArray->GetLast()+1;
563 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
567 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
569 // Too many pixels - sort and remove pixels with the lowest signal
571 for ( Int_t i = npad; i < nPix; ++i )
575 fPixArray->Compress();
576 } // if (nPix > npad)
578 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
579 // fPixArray->Print(););
580 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
583 //_____________________________________________________________________________
584 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
586 /// Build the pixel array
588 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
590 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
591 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
592 Int_t found[2] = {0}, mult = cluster.Multiplicity();
594 for ( Int_t i = 0; i < mult; ++i) {
595 AliMUONPad* pad = cluster.Pad(i);
596 for (Int_t j = 0; j < 2; ++j) {
597 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
598 xy0[j] = pad->Coord(j);
602 if (found[0] && found[1]) break;
605 Double_t min[2], max[2];
606 Int_t cath0 = 0, cath1 = 1;
607 if (cluster.Multiplicity(0) == 0) cath0 = 1;
608 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
610 TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
611 TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
612 min[0] = leftDown.X();
613 min[1] = leftDown.Y();
614 max[0] = rightUp.X();
615 max[1] = rightUp.Y();
616 if (cath1 != cath0) {
617 leftDown = cluster.Area(cath1).LeftDownCorner();
618 rightUp = cluster.Area(cath1).RightUpCorner();
619 min[0] = TMath::Max (min[0], leftDown.X());
620 min[1] = TMath::Max (min[1], leftDown.Y());
621 max[0] = TMath::Min (max[0], rightUp.X());
622 max[1] = TMath::Min (max[1], rightUp.Y());
626 //width[0] /= 2; width[1] /= 2; // just for check
628 for (Int_t i = 0; i < 2; ++i) {
629 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
630 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
631 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
632 + TMath::Sign(0.5,dist)) * width[i] * 2;
633 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
634 if (nbins[i] == 0) ++nbins[i];
635 max[i] = min[i] + nbins[i] * width[i] * 2;
636 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
640 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
641 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
642 TAxis *xaxis = hist1->GetXaxis();
643 TAxis *yaxis = hist1->GetYaxis();
646 for ( Int_t i = 0; i < mult; ++i) {
647 AliMUONPad* pad = cluster.Pad(i);
648 Int_t ix0 = xaxis->FindBin(pad->X());
649 Int_t iy0 = yaxis->FindBin(pad->Y());
650 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
654 for (Int_t i = 1; i <= nbins[0]; ++i) {
655 Double_t x = xaxis->GetBinCenter(i);
656 for (Int_t j = 1; j <= nbins[1]; ++j) {
657 if (hist2->GetCellContent(i,j) < 0.1) continue;
658 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
659 // cluster.Multiplicity(1)) continue;
660 if (cath0 != cath1) {
662 Double_t cont = hist2->GetCellContent(i,j);
663 if (cont < 999.) continue;
664 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
666 Double_t y = yaxis->GetBinCenter(j);
667 Double_t charge = hist1->GetCellContent(i,j);
668 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
669 fPixArray->Add(pixPtr);
673 if (fPixArray->GetEntriesFast() == 1) {
674 // Split pixel into 2
675 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
676 pixPtr->SetSize(0,width[0]/2.);
677 pixPtr->Shift(0,-width[0]/4.);
678 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
679 fPixArray->Add(pixPtr);
682 //fPixArray->Print();
687 //_____________________________________________________________________________
688 void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
689 TH2D *hist1, TH2D *hist2)
691 /// "Span" pad over histogram in the direction idir
693 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
694 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
695 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
697 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
699 for (Int_t i = 0; i < nbinPad; ++i) {
700 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
701 if (ixy > nbins) break;
702 Double_t lowEdge = axis->GetBinLowEdge(ixy);
703 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
704 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
707 Double_t cont = pad->Charge();
708 if (hist2->GetCellContent(ix0, ixy) > 0.1)
709 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
710 hist1->SetCellContent(ix0, ixy, cont);
711 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
712 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
716 for (Int_t i = -1; i > -nbinPad; --i) {
717 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
719 Double_t upEdge = axis->GetBinUpEdge(ixy);
720 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
721 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
724 Double_t cont = pad->Charge();
725 if (hist2->GetCellContent(ix0, ixy) > 0.1)
726 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
727 hist1->SetCellContent(ix0, ixy, cont);
728 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
729 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
734 //_____________________________________________________________________________
736 AliMUONClusterFinderMLEM::Plot(const char* basename)
738 /// Make a plot and save it as png
743 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
748 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
749 fDetElemId,fClusterNumber));
752 //_____________________________________________________________________________
754 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
758 /// Compute coefficients needed for MLEM algorithm
760 Int_t nPix = fPixArray->GetLast()+1;
762 //memset(probi,0,nPix*sizeof(Double_t));
763 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
765 Int_t mult = cluster.Multiplicity();
766 for ( Int_t j = 0; j < mult; ++j )
768 AliMUONPad* pad = cluster.Pad(j);
771 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
773 Int_t indx1 = indx + ipix;
774 //if (pad->Status() < 0)
775 if (pad->Status() != fgkZero)
780 AliMUONPad* pixPtr = Pixel(ipix);
781 // coef is the charge (given by Mathieson integral) on pad, assuming
782 // the Mathieson is center at pixel.
783 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
784 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
785 // pad->Ix(),pad->Iy(),
786 // pad->X(),pad->Y(),
787 // pad->DX(),pad->DY(),
788 // pixPtr->Coord(0),pixPtr->Coord(1),
789 // pixPtr->Size(0),pixPtr->Size(1),
792 probi[ipix] += coef[indx1];
797 //_____________________________________________________________________________
798 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
800 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
802 // AliCodeTimerAuto("")
804 Int_t nPix = fPixArray->GetLast()+1;
806 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
807 //StdoutToAliDebug(2,cluster.Print("full"););
811 AliDebug(1,"No pixels, why am I here then ?");
814 AddVirtualPad(cluster); // add virtual pads if necessary
816 Int_t npadTot = cluster.Multiplicity();
818 for (Int_t i = 0; i < npadTot; ++i)
820 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
821 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
825 Double_t* probi(0x0);
826 Int_t lc(0); // loop counter
828 //Plot("mlem.start");
829 AliMUONPad* pixPtr = Pixel(0);
830 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
837 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
838 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
839 //StdoutToAliDebug(2,fPixArray->Print("","full"));
841 coef = new Double_t [npadTot*nPix];
842 probi = new Double_t [nPix];
844 // Calculate coefficients and pixel visibilities
845 ComputeCoefficients(cluster,coef,probi);
847 for (Int_t ipix = 0; ipix < nPix; ++ipix)
849 if (probi[ipix] < 0.01)
851 AliMUONPad* pixel = Pixel(ipix);
852 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
853 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
854 pixel->SetCharge(0); // "invisible" pixel
859 Mlem(cluster,coef, probi, 15);
861 // Find histogram limits for the 1'st pass only - for others computed below
863 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
865 pixPtr = Pixel(ipix);
866 for ( Int_t i = 0; i < 2; ++i )
869 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
870 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
873 } else pixPtr = Pixel(0);
875 for (Int_t i = 0; i < 4; i++)
877 xylim[i] -= pixPtr->Size(i/2);
880 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
881 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
883 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
884 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
885 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
886 xylim[0],-xylim[1],xylim[2],-xylim[3]
889 fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
891 for (Int_t ipix = 0; ipix < nPix; ++ipix)
893 AliMUONPad* pixPtr = Pixel(ipix);
894 fHistMlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
897 // Check if the total charge of pixels is too low
899 for ( Int_t i = 0; i < nPix; ++i)
901 qTot += Pixel(i)->Charge();
904 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
906 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
910 for ( Int_t i = 0; i < npadTot; ++i)
912 AliMUONPad* pad = cluster.Pad(i);
913 //if ( pad->Status() == 0) pad->SetStatus(-1);
914 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
921 // Simple cluster - skip further passes thru EM-procedure
929 // Calculate position of the center-of-gravity around the maximum pixel
933 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
934 pixPtr->Size(0) > pixPtr->Size(1)) break;
936 // Sort pixels according to the charge
937 MaskPeaks(1); // mask local maxima
939 MaskPeaks(0); // unmask local maxima
940 Double_t pixMin = 0.01*Pixel(0)->Charge();
941 pixMin = TMath::Min(pixMin,50.);
943 // Decrease pixel size and shift pixels to make them centered at
945 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
948 Double_t shift[2] = { 0.0, 0.0 };
949 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
952 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
954 AliMUONPad* pixPtr = Pixel(ipix);
955 if ( nPix >= npadOK // too many pixels already
957 pixPtr->Charge() < pixMin && pixPtr->Status() != fgkMustKeep // too low charge
963 for (Int_t i = 0; i < 2; ++i)
967 pixPtr->SetCharge(10);
968 pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
969 width = -pixPtr->Size(indx);
970 pixPtr->Shift(indx, width);
971 // Shift pixel position
975 for (Int_t j = 0; j < 2; ++j)
977 shift[j] = pixPtr->Coord(j) - xyCOG[j];
978 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
981 pixPtr->Shift(0, -shift[0]);
982 pixPtr->Shift(1, -shift[1]);
985 else if (nPix < npadOK)
987 pixPtr = new AliMUONPad(*pixPtr);
988 pixPtr->Shift(indx, -2*width);
989 pixPtr->SetStatus(fgkZero);
990 fPixArray->Add(pixPtr);
993 else continue; // skip adjustment of histo limits
994 for (Int_t j = 0; j < 4; ++j)
996 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr->Coord(j/2));
998 } // for (Int_t i=0; i<2;
999 } // for (Int_t ipix=0;
1001 fPixArray->Compress();
1003 AliDebug(2,Form("After shift:"));
1004 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1005 //Plot(Form("mlem.lc%d",lc+1));
1007 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1010 xylim[2],xylim[3]));
1014 AliMUONPad* pixPtr = Pixel(0);
1015 // add pixels if the maximum is at the limit of pixel area:
1016 // start from Y-direction
1018 for (Int_t i = 3; i > -1; --i)
1020 if (nPix < npadOK &&
1021 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2))
1023 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1024 AliMUONPad* p = new AliMUONPad(*pixPtr);
1025 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1026 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1027 j = TMath::Even (i/2);
1028 p->SetCoord(j, xyCOG[j]);
1029 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1030 //StdoutToAliDebug(2,cout << " ---- ";
1031 // p->Print("corners"););
1037 nPix = fPixArray->GetEntriesFast();
1042 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1043 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1045 // remove pixels with low signal or low visibility
1046 // Cuts are empirical !!!
1047 Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,1.);
1048 thresh = TMath::Min (thresh,50.);
1049 Double_t charge = 0;
1051 // Mark pixels which should be removed
1052 for (Int_t i = 0; i < nPix; ++i)
1054 AliMUONPad* pixPtr = Pixel(i);
1055 charge = pixPtr->Charge();
1056 if (charge < thresh)
1058 pixPtr->SetCharge(-charge);
1062 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1064 for (Int_t i = 0; i < nPix; ++i)
1066 AliMUONPad* pixPtr = Pixel(i);
1067 charge = pixPtr->Charge();
1068 if (charge > 0) continue;
1069 near = FindNearest(pixPtr);
1070 pixPtr->SetCharge(0);
1071 probi[i] = 0; // make it "invisible"
1072 AliMUONPad* pnear = Pixel(near);
1073 pnear->SetCharge(pnear->Charge() + (-charge));
1075 Mlem(cluster,coef,probi,2);
1077 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1078 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1079 //Plot("mlem.beforesplit");
1082 for (Int_t i = 0; i < nPix; ++i)
1084 AliMUONPad* pixPtr = Pixel(i);
1085 Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1086 Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1087 fHistMlem->SetBinContent(ix, iy, pixPtr->Charge());
1090 // Try to split into clusters
1092 if (fHistMlem->GetSum() < 1)
1098 fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
1103 fPixArray->Delete();
1108 //_____________________________________________________________________________
1109 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1111 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1112 /// - to avoid loosing low charge pixels after sorting)
1114 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1115 AliMUONPad* pix = Pixel(i);
1116 if (pix->Status() == fgkMustKeep) {
1117 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1118 else pix->SetCharge(pix->Charge()-10000.);
1123 //_____________________________________________________________________________
1124 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1125 Double_t* coef, Double_t* probi,
1128 /// Use MLEM to find pixel charges
1130 Int_t nPix = fPixArray->GetEntriesFast();
1132 Int_t npad = cluster.Multiplicity();
1134 Double_t* probi1 = new Double_t[nPix];
1135 Double_t probMax = TMath::MaxElement(nPix,probi);
1137 for (Int_t iter = 0; iter < nIter; ++iter)
1140 for (Int_t ipix = 0; ipix < nPix; ++ipix)
1142 Pixel(ipix)->SetChargeBackup(0);
1143 // Correct each pixel
1144 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1146 probi1[ipix] = probMax;
1147 for (Int_t j = 0; j < npad; ++j)
1149 AliMUONPad* pad = cluster.Pad(j);
1150 //if (pad->Status() < 0) continue;
1151 if (pad->Status() != fgkZero) continue;
1153 Int_t indx1 = j*nPix;
1154 Int_t indx = indx1 + ipix;
1155 // Calculate expectation
1156 for (Int_t i = 0; i < nPix; ++i)
1158 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1159 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1161 if ( pad->IsSaturated() && sum1 > pad->Charge() )
1163 // correct for pad charge overflows
1164 probi1[ipix] -= coef[indx];
1166 //sum1 = pad->Charge();
1169 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1170 } // for (Int_t j=0;
1171 AliMUONPad* pixPtr = Pixel(ipix);
1172 if (probi1[ipix] > 1.e-6)
1174 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1175 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1177 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1178 } // for (Int_t ipix=0;
1180 for (Int_t i = 0; i < nPix; ++i) {
1181 AliMUONPad* pixPtr = Pixel(i);
1182 pixPtr->RevertCharge();
1183 qTot += pixPtr->Charge();
1186 // Can happen in clusters with large number of overflows - speeding up
1190 } // for (Int_t iter=0;
1194 //_____________________________________________________________________________
1195 void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
1197 /// Calculate position of the center-of-gravity around the maximum pixel
1199 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1200 Int_t i1 = -9, j1 = -9;
1201 fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1202 Int_t nx = fHistMlem->GetNbinsX();
1203 Int_t ny = fHistMlem->GetNbinsY();
1204 Double_t thresh = fHistMlem->GetMaximum()/10;
1205 Double_t x, y, cont, xq=0, yq=0, qq=0;
1207 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1208 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1209 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1210 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1211 cont = fHistMlem->GetCellContent(j,i);
1212 if (cont < thresh) continue;
1213 if (i != i1) {i1 = i; nsumy++;}
1214 if (j != j1) {j1 = j; nsumx++;}
1215 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1224 Int_t i2 = 0, j2 = 0;
1227 // one bin in Y - add one more (with the largest signal)
1228 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1229 if (i == iymax) continue;
1230 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1231 cont = fHistMlem->GetCellContent(j,i);
1234 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1235 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1244 if (i2 != i1) nsumy++;
1245 if (j2 != j1) nsumx++;
1247 } // if (nsumy == 1)
1250 // one bin in X - add one more (with the largest signal)
1252 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1253 if (j == ixmax) continue;
1254 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1255 cont = fHistMlem->GetCellContent(j,i);
1258 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1259 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1268 if (i2 != i1) nsumy++;
1269 if (j2 != j1) nsumx++;
1271 } // if (nsumx == 1)
1273 xyc[0] = xq/qq; xyc[1] = yq/qq;
1274 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1277 //_____________________________________________________________________________
1278 Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1280 /// Find the pixel nearest to the given one
1281 /// (algorithm may be not very efficient)
1283 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1284 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1285 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1288 for (Int_t i = 0; i < nPix; ++i) {
1289 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1290 if (pixPtr == pixPtr0 || pixPtr->Charge() < 0.5) continue;
1291 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1292 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1293 r = dx *dx + dy * dy;
1294 if (r < rmin) { rmin = r; imin = i; }
1299 //_____________________________________________________________________________
1301 AliMUONClusterFinderMLEM::Paint(Option_t*)
1303 /// Paint cluster and pixels
1305 AliMpArea area(fPreCluster->Area());
1307 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1309 gVirtualX->SetFillStyle(1001);
1310 gVirtualX->SetFillColor(3);
1311 gVirtualX->SetLineColor(3);
1315 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1317 AliMUONPad* pixel = Pixel(i);
1319 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1320 pixel->Coord(1)-pixel->Size(1)*s,
1321 pixel->Coord(0)+pixel->Size(0)*s,
1322 pixel->Coord(1)+pixel->Size(1)*s);
1324 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1326 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1327 // pixel->Coord(1) + sign*pixel->Size(1),
1328 // pixel->Coord(0) + pixel->Size(0),
1329 // pixel->Coord(1) - sign*pixel->Size(1)
1335 gVirtualX->SetFillStyle(0);
1337 fPreCluster->Paint();
1339 gVirtualX->SetLineColor(1);
1340 gVirtualX->SetLineWidth(2);
1341 gVirtualX->SetFillStyle(0);
1342 gVirtualX->SetTextColor(1);
1343 gVirtualX->SetTextAlign(22);
1345 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1347 AliMUONPad* pixel = Pixel(i);
1348 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1349 pixel->Coord(1)-pixel->Size(1),
1350 pixel->Coord(0)+pixel->Size(0),
1351 pixel->Coord(1)+pixel->Size(1));
1352 gVirtualX->SetTextSize(pixel->Size(0)*60);
1354 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1358 //_____________________________________________________________________________
1359 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1361 /// Find local maxima in pixel space for large preclusters in order to
1362 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1363 /// or to find additional fitting seeds if clusters were not completely resolved
1365 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1367 Double_t xylim[4] = {999, 999, 999, 999};
1369 Int_t nPix = pixArray->GetEntriesFast();
1370 AliMUONPad *pixPtr = 0;
1371 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1372 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1373 for (Int_t i = 0; i < 4; ++i)
1374 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1376 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
1378 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1379 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1380 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1381 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1382 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1383 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1384 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1386 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1388 Int_t nMax = 0, indx, nxy = ny * nx;
1389 Int_t *isLocalMax = new Int_t[nxy];
1390 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
1392 for (Int_t i = 1; i <= ny; ++i) {
1394 for (Int_t j = 1; j <= nx; ++j) {
1395 if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
1396 //if (isLocalMax[indx+j-1] < 0) continue;
1397 if (isLocalMax[indx+j-1] != 0) continue;
1398 FlagLocalMax(fHistAnode, i, j, isLocalMax);
1402 for (Int_t i = 1; i <= ny; ++i) {
1404 for (Int_t j = 1; j <= nx; ++j) {
1405 if (isLocalMax[indx+j-1] > 0) {
1406 localMax[nMax] = indx + j - 1;
1407 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
1408 ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
1409 if (nMax > 99) AliFatal(" Too many local maxima !!!");
1413 if (fDebug) cout << " Local max: " << nMax << endl;
1414 delete [] isLocalMax;
1418 //_____________________________________________________________________________
1419 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1421 /// Flag pixels (whether or not local maxima)
1423 Int_t nx = hist->GetNbinsX();
1424 Int_t ny = hist->GetNbinsY();
1425 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1426 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1428 Int_t ie = i + 2, je = j + 2;
1429 for (Int_t i1 = i-1; i1 < ie; ++i1) {
1430 if (i1 < 1 || i1 > ny) continue;
1431 indx1 = (i1 - 1) * nx;
1432 for (Int_t j1 = j-1; j1 < je; ++j1) {
1433 if (j1 < 1 || j1 > nx) continue;
1434 if (i == i1 && j == j1) continue;
1435 indx2 = indx1 + j1 - 1;
1436 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1437 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1438 else if (cont > cont1) isLocalMax[indx2] = -1;
1439 else { // the same charge
1440 isLocalMax[indx] = 1;
1441 if (isLocalMax[indx2] == 0) {
1442 FlagLocalMax(hist, i1, j1, isLocalMax);
1443 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1444 else isLocalMax[indx2] = -1;
1449 isLocalMax[indx] = 1; // local maximum
1452 //_____________________________________________________________________________
1453 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1454 Int_t *localMax, Int_t iMax)
1456 /// Find pixel cluster around local maximum \a iMax and pick up pads
1457 /// overlapping with it
1460 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1462 hist->Draw("lego1Fb"); // debug
1467 Int_t nx = fHistAnode->GetNbinsX();
1468 Int_t ny = fHistAnode->GetNbinsY();
1469 Int_t ic = localMax[iMax] / nx + 1;
1470 Int_t jc = localMax[iMax] % nx + 1;
1471 Int_t nxy = ny * nx;
1472 Bool_t *used = new Bool_t[nxy];
1473 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
1475 // Drop all pixels from the array - pick up only the ones from the cluster
1476 fPixArray->Delete();
1478 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1479 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1480 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1481 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1482 Double_t cont = fHistAnode->GetCellContent(jc,ic);
1483 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1484 used[(ic-1)*nx+jc-1] = kTRUE;
1485 AddBinSimple(fHistAnode, ic, jc);
1486 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1488 Int_t nPix = fPixArray->GetEntriesFast();
1489 Int_t npad = cluster.Multiplicity();
1491 for (Int_t i = 0; i < nPix; ++i)
1493 AliMUONPad* pixPtr = Pixel(i);
1494 pixPtr->SetSize(0,wx);
1495 pixPtr->SetSize(1,wy);
1498 // Pick up pads which overlap with found pixels
1499 for (Int_t i = 0; i < npad; ++i)
1501 //cluster.Pad(i)->SetStatus(-1);
1502 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1505 for (Int_t i = 0; i < nPix; ++i)
1507 AliMUONPad* pixPtr = Pixel(i);
1508 for (Int_t j = 0; j < npad; ++j)
1510 AliMUONPad* pad = cluster.Pad(j);
1511 //if (pad->Status() == 0) continue;
1512 if (pad->Status() == fgkZero) continue;
1513 if ( Overlap(*pad,*pixPtr) )
1515 //pad->SetStatus(0);
1516 pad->SetStatus(fgkZero);
1517 if (fDebug) { cout << j << " "; pad->Print("full"); }
1525 //_____________________________________________________________________________
1527 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
1529 /// Add adjacent bins (+-1 in X and Y) to the cluster
1531 Int_t nx = hist->GetNbinsX();
1532 Int_t ny = hist->GetNbinsY();
1533 Double_t cont1, cont = hist->GetCellContent(jc,ic);
1534 AliMUONPad *pixPtr = 0;
1536 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1537 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1538 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1539 cont1 = hist->GetCellContent(j,i);
1540 if (cont1 > cont) continue;
1541 if (cont1 < 0.5) continue;
1542 pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j),
1543 hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1544 fPixArray->Add(pixPtr);
1549 //_____________________________________________________________________________
1550 AliMUONClusterFinderMLEM&
1551 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1553 /// Protected assignement operator
1555 if (this == &rhs) return *this;
1557 AliFatal("Not implemented.");
1562 //_____________________________________________________________________________
1563 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1565 /// Add virtual pad (with small charge) to improve fit for clusters
1566 /// with number of pads == 2 per direction
1568 // Find out non-bending and bending planes
1569 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1571 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1572 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1573 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1578 Bool_t same = kFALSE;
1579 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1582 Bool_t check[2] = {kFALSE, kFALSE};
1584 nxy[0] = nxy[1] = 0;
1585 for (Int_t inb = 0; inb < 2; ++inb) {
1586 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1587 if (inb == 0 && cn.GetFirst() == 2) check[inb] = kTRUE; // check non-bending plane
1588 else if (inb == 1 && cn.GetSecond() == 2) check[inb] = kTRUE; // check bending plane
1590 nxy[0] = TMath::Max (nxy[0], cn.GetFirst());
1591 nxy[1] = TMath::Max (nxy[1], cn.GetSecond());
1592 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1593 else if (inb == 1 && cn.GetSecond() < 2) nonb[inb] = !nonb[inb];
1597 if (nxy[0] > 2) check[0] = kFALSE;
1598 if (nxy[1] > 2) check[1] = kFALSE;
1600 if (!check[0] && !check[1]) return;
1602 for (Int_t inb = 0; inb < 2; ++inb) {
1603 if (!check[inb]) continue;
1605 // Find pads with maximum and next to maximum charges
1606 Int_t maxPads[2] = {-1, -1};
1607 Double_t amax[2] = {0};
1608 Int_t mult = cluster.Multiplicity();
1609 for (Int_t j = 0; j < mult; ++j) {
1610 AliMUONPad *pad = cluster.Pad(j);
1611 if (pad->Cathode() != nonb[inb]) continue;
1612 for (Int_t j2 = 0; j2 < 2; ++j2) {
1613 if (pad->Charge() > amax[j2]) {
1614 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1615 amax[j2] = pad->Charge();
1622 // Find min and max dimensions of the cluster
1623 Double_t limits[2] = {9999, -9999};
1624 for (Int_t j = 0; j < mult; ++j) {
1625 AliMUONPad *pad = cluster.Pad(j);
1626 if (pad->Cathode() != nonb[inb]) continue;
1627 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1628 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1631 // Loop over max and next to max pads
1632 Bool_t add = kFALSE;
1634 for (Int_t j = 0; j < 2; ++j) {
1636 if (maxPads[j] < 0) continue;
1638 if (amax[1] / amax[0] < 0.5) break;
1640 // Check if pad at the cluster limit
1641 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1643 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1644 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1646 //cout << " *** Pad not at the cluster limit: " << j << endl;
1649 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1651 // Add pad (if it exists)
1653 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1654 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1655 //AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1656 AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kFALSE);
1657 if (!mppad.IsValid()) continue; // non-existing pad
1658 cn = mppad.GetIndices();
1659 AliMUONPad muonPad(fDetElemId, nonb[inb], cn.GetFirst(), cn.GetSecond(),
1660 mppad.Position().X(), mppad.Position().Y(),
1661 mppad.Dimensions().X(), mppad.Dimensions().Y(), 0);
1662 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, 5.));
1663 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1664 else muonPad.SetCharge(TMath::Min (amax[j]/15, 6.));
1665 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1666 muonPad.SetReal(kFALSE);
1667 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1668 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1669 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1670 cluster.AddPad(muonPad); // add pad to the cluster
1677 //_____________________________________________________________________________
1678 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1679 Int_t &nInX, Int_t &nInY) const
1681 /// Find number of pads in X and Y-directions (excluding virtual ones and
1684 //Int_t statusToTest = 1;
1685 Int_t statusToTest = fgkUseForFit;
1687 //if ( nInX < 0 ) statusToTest = 0;
1688 if ( nInX < 0 ) statusToTest = fgkZero;
1690 Bool_t mustMatch(kTRUE);
1692 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1694 nInX = cn.GetFirst();
1695 nInY = cn.GetSecond();
1698 //_____________________________________________________________________________
1699 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1701 /// Remove pixel at index i
1702 AliMUONPad* pixPtr = Pixel(i);
1703 fPixArray->RemoveAt(i);
1707 //_____________________________________________________________________________
1708 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1710 /// Process simple cluster (small number of pads) without EM-procedure
1712 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1713 Double_t parOk[3] = {0.};
1714 TObjArray *clusters[1];
1715 clusters[0] = fPixArray;
1717 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1719 Int_t mult = cluster.Multiplicity();
1720 for (Int_t i = 0; i < mult; ++i)
1722 AliMUONPad* pad = cluster.Pad(i);
1724 if ( pad->IsSaturated())
1733 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1735 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
1738 //_____________________________________________________________________________
1740 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1742 /// Returns pixel at index i
1743 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1746 //_____________________________________________________________________________
1748 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1751 TString swhat(what);
1753 if ( swhat.Contains("precluster") )
1755 if ( fPreCluster) fPreCluster->Print();