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::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
55 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
56 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
58 // Status flags for pads
59 const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
60 const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
61 const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
62 const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
63 const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
64 const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
66 //_____________________________________________________________________________
67 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
68 : AliMUONVClusterFinder(),
69 fPreClusterFinder(clusterFinder),
77 fPixArray(new TObjArray(20)),
83 fLowestPixelCharge(0),
85 fLowestClusterCharge(0)
89 fkSegmentation[1] = fkSegmentation[0] = 0x0;
91 if (fPlot) fDebug = 1;
94 //_____________________________________________________________________________
95 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
98 delete fPixArray; fPixArray = 0;
100 delete fPreClusterFinder;
102 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
103 fNClusters,fNAddVirtualPads));
106 //_____________________________________________________________________________
108 AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
109 TClonesArray* pads[2],
110 const AliMpArea& area,
111 const AliMpVSegmentation* seg[2])
113 /// Prepare for clustering
114 // AliCodeTimerAuto("",0)
116 for ( Int_t i = 0; i < 2; ++i )
118 fkSegmentation[i] = seg[i];
121 // Find out the DetElemId
122 fDetElemId = detElemId;
125 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,
129 fLowestClusterCharge);
130 fSplitter->SetDebug(fDebug);
132 // find out current event number, and reset the cluster number
133 AliRunLoader *runLoader = AliRunLoader::Instance();
134 fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
136 fClusterList.Delete();
139 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
141 if ( fPreClusterFinder->NeedSegmentation() )
143 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
147 return fPreClusterFinder->Prepare(detElemId,pads,area);
151 //_____________________________________________________________________________
153 AliMUONClusterFinderMLEM::NextCluster()
155 /// Return next cluster
156 // AliCodeTimerAuto("",0)
158 // if the list of clusters is not void, pick one from there
159 TObject* o = fClusterList.At(++fClusterNumber);
160 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
162 //FIXME : at this point, must check whether we've used all the digits
163 //from precluster : if not, let the preclustering know about those unused
164 //digits, so it can reuse them
166 // if the cluster list is exhausted, we need to go to the next
167 // pre-cluster and treat it
169 fPreCluster = fPreClusterFinder->NextCluster();
172 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
173 fClusterNumber = -1; //AZ
183 // WorkOnPreCluster may have used only part of the pads, so we check that
184 // now, and let the unused pads be reused by the preclustering...
186 Int_t mult = fPreCluster->Multiplicity();
187 for ( Int_t i = 0; i < mult; ++i )
189 AliMUONPad* pad = fPreCluster->Pad(i);
190 if ( !pad->IsUsed() )
192 fPreClusterFinder->UsePad(*pad);
196 return NextCluster();
199 //_____________________________________________________________________________
201 AliMUONClusterFinderMLEM::WorkOnPreCluster()
203 /// Starting from a precluster, builds a pixel array, and then
204 /// extract clusters from this array
206 // AliCodeTimerAuto("",0)
209 cout << " *** Event # " << fEventNumber
210 << " det. elem.: " << fDetElemId << endl;
211 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
212 AliMUONPad* pad = fPreCluster->Pad(j);
213 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
214 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
215 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
219 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
220 if (!cluster) return kFALSE;
222 BuildPixArray(*cluster);
224 if ( fPixArray->GetLast() < 0 )
226 AliDebug(1,"No pixel for the above cluster");
231 // Use MLEM for cluster finder
232 Int_t nMax = 1, localMax[100], maxPos[100];
233 Double_t maxVal[100];
235 Int_t iSimple = 0, nInX = -1, nInY;
237 PadsInXandY(*cluster,nInX, nInY);
239 if (nInX < 4 && nInY < 4)
241 iSimple = 1; // simple cluster
245 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
247 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
248 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
252 for (Int_t i = 0; i < nMax; ++i)
256 FindCluster(*cluster,localMax, maxPos[i]);
259 MainLoop(*cluster,iSimple);
263 Int_t mult = cluster->Multiplicity();
264 for (Int_t j = 0; j < mult; ++j)
266 AliMUONPad* pad = cluster->Pad(j);
267 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
268 if ( pad->Status() != fgkOver ) continue; // pad was not used
270 pad->SetStatus(fgkZero);
271 pad->RevertCharge(); // use backup charge value
274 } // for (Int_t i=0; i<nMax;
277 fHistMlem = fHistAnode = 0x0;
282 //_____________________________________________________________________________
284 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
286 /// Check if the pad and the pixel overlaps
288 // make a fake pad from the pixel
289 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
290 pixel.Coord(0),pixel.Coord(1),
291 pixel.Size(0),pixel.Size(1),0);
293 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
296 //_____________________________________________________________________________
298 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
300 /// Check precluster in order to attempt to simplify it (mostly for
301 /// two-cathode preclusters)
303 AliCodeTimerAuto("",0)
305 // Disregard small clusters (leftovers from splitting or noise)
306 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
307 origCluster.Charge(0)+origCluster.Charge(1) < fLowestClusterCharge )
312 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
314 AliDebug(2,"Start of CheckPreCluster=");
315 //StdoutToAliDebug(2,cluster->Print("full"));
317 AliMUONCluster* rv(0x0);
319 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
321 rv = CheckPreclusterTwoCathodes(cluster);
330 //_____________________________________________________________________________
332 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
334 /// Check two-cathode cluster
336 Int_t npad = cluster->Multiplicity();
337 Int_t* flags = new Int_t[npad];
338 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
340 // Check pad overlaps
341 for ( Int_t i = 0; i < npad; ++i)
343 AliMUONPad* padi = cluster->Pad(i);
344 if ( padi->Cathode() != 0 ) continue;
345 for (Int_t j = i+1; j < npad; ++j)
347 AliMUONPad* padj = cluster->Pad(j);
348 if ( padj->Cathode() != 1 ) continue;
349 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
350 flags[i] = flags[j] = 1; // mark overlapped pads
354 // Check if all pads overlap
356 for (Int_t i = 0; i < npad; ++i)
358 if (!flags[i]) ++nFlags;
363 // not all pads overlap.
364 if (fDebug) cout << " nFlags: " << nFlags << endl;
365 TObjArray toBeRemoved;
366 for (Int_t i = 0; i < npad; ++i)
368 AliMUONPad* pad = cluster->Pad(i);
369 if (flags[i]) continue;
370 Int_t cath = pad->Cathode();
371 Int_t cath1 = TMath::Even(cath);
372 // Check for edge effect (missing pads on the _other_ cathode)
374 fkSegmentation[cath1]->PadByPosition(pad->Position().X(),
375 pad->Position().Y(),kFALSE);
376 if (!mpPad.IsValid()) continue;
377 if (nFlags == 1 && pad->Charge() < fLowestPadCharge) continue;
378 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
379 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
380 toBeRemoved.AddLast(pad);
381 fPreCluster->Pad(i)->Release();
383 Int_t nRemove = toBeRemoved.GetEntriesFast();
384 for ( Int_t i = 0; i < nRemove; ++i )
386 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
390 // Check correlations of cathode charges
391 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
394 Int_t cathode = cluster->MaxRawChargeCathode();
400 // get min and max pad charges on the cathode opposite to the
401 // max pad (given by MaxRawChargeCathode())
403 Int_t mult = cluster->Multiplicity();
404 for ( Int_t i = 0; i < mult; ++i )
406 AliMUONPad* pad = cluster->Pad(i);
407 if ( pad->Cathode() != cathode || !pad->IsReal() )
409 // only consider pads in the opposite cathode, and
410 // only consider real pads (i.e. exclude the virtual ones)
413 if ( pad->Charge() < cmin )
415 cmin = pad->Charge();
422 else if ( pad->Charge() > cmax )
424 cmax = pad->Charge();
428 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
429 imin,imax,cmin,cmax));
431 // arrange pads according to their distance to the max, normalized
433 Double_t* dist = new Double_t[mult];
438 AliMUONPad* padmax = cluster->Pad(imax);
440 for ( Int_t i = 0; i < mult; ++i )
443 if ( i == imax) continue;
444 AliMUONPad* pad = cluster->Pad(i);
445 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
446 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
447 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
448 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
451 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
457 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
458 Double_t xmax(-1), distPrev(999);
459 TObjArray toBeRemoved;
461 for ( Int_t i = 0; i < mult; ++i )
463 Int_t indx = flags[i];
464 AliMUONPad* pad = cluster->Pad(indx);
465 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
466 if ( dist[indx] > dmin )
468 // farther than the minimum pad
469 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
470 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
473 if (dx >= 0 && dy >= 0) continue;
474 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
475 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
477 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
478 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
481 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
483 cmax = TMath::Max(pad->Charge(),cmax);
487 cmax = pad->Charge();
490 distPrev = dist[indx];
491 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
492 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
495 toBeRemoved.AddLast(pad);
496 fPreCluster->Pad(indx)->Release();
499 Int_t nRemove = toBeRemoved.GetEntriesFast();
500 for ( Int_t i = 0; i < nRemove; ++i )
502 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
509 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
510 //StdoutToAliDebug(2,cluster->Print("full"));
515 //_____________________________________________________________________________
517 AliMUONClusterFinderMLEM::CheckOverlaps()
519 /// For debug only : check if some pixels overlap...
521 Int_t nPix = fPixArray->GetLast()+1;
524 for ( Int_t i = 0; i < nPix; ++i )
526 AliMUONPad* pixelI = Pixel(i);
527 AliMUONPad pi(dummy,dummy,dummy,dummy,
528 pixelI->Coord(0),pixelI->Coord(1),
529 pixelI->Size(0),pixelI->Size(1),0.0);
531 for ( Int_t j = i+1; j < nPix; ++j )
533 AliMUONPad* pixelJ = Pixel(j);
534 AliMUONPad pj(dummy,dummy,dummy,dummy,
535 pixelJ->Coord(0),pixelJ->Coord(1),
536 pixelJ->Size(0),pixelJ->Size(1),0.0);
539 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
541 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
543 StdoutToAliInfo(pixelI->Print();
544 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
546 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
547 cout << " Area surface = " << area.GetDimensionX()*area.GetDimensionY()*4 << endl;
548 cout << "-------" << endl;
556 //_____________________________________________________________________________
557 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
559 /// Build pixel array for MLEM method
561 Int_t npad = cluster.Multiplicity();
564 AliWarning("Got no pad at all ?!");
568 BuildPixArrayOneCathode(cluster);
570 Int_t nPix = fPixArray->GetLast()+1;
572 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
576 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
578 // Too many pixels - sort and remove pixels with the lowest signal
580 for ( Int_t i = npad; i < nPix; ++i )
584 fPixArray->Compress();
585 } // if (nPix > npad)
587 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
588 // fPixArray->Print(););
589 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
592 //_____________________________________________________________________________
593 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
595 /// Build the pixel array
597 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
599 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
600 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2]={99999,99999};
601 Int_t found[2] = {0,0}, mult = cluster.Multiplicity();
603 for ( Int_t i = 0; i < mult; ++i) {
604 AliMUONPad* pad = cluster.Pad(i);
605 for (Int_t j = 0; j < 2; ++j) {
606 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
607 xy0[j] = pad->Coord(j);
611 if (found[0] && found[1]) break;
614 Double_t min[2], max[2];
615 Int_t cath0 = 0, cath1 = 1;
616 if (cluster.Multiplicity(0) == 0) cath0 = 1;
617 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
620 Double_t leftDownX, leftDownY;
621 cluster.Area(cath0).LeftDownCorner(leftDownX, leftDownY);
622 Double_t rightUpX, rightUpY;
623 cluster.Area(cath0).RightUpCorner(rightUpX, rightUpY);
628 if (cath1 != cath0) {
629 cluster.Area(cath1).LeftDownCorner(leftDownX, leftDownY);
630 cluster.Area(cath1).RightUpCorner(rightUpX, rightUpY);
631 min[0] = TMath::Max (min[0], leftDownX);
632 min[1] = TMath::Max (min[1], leftDownY);
633 max[0] = TMath::Min (max[0], rightUpX);
634 max[1] = TMath::Min (max[1], rightUpY);
638 //width[0] /= 2; width[1] /= 2; // just for check
639 Int_t nbins[2]={0,0};
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, hist1, hist2);
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,
701 TH2D *hist1, TH2D *hist2)
703 /// "Span" pad over histogram in the direction idir
705 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
706 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
707 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
709 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
711 for (Int_t i = 0; i < nbinPad; ++i) {
712 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
713 if (ixy > nbins) break;
714 Double_t lowEdge = axis->GetBinLowEdge(ixy);
715 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
716 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
719 Double_t cont = pad->Charge();
720 if (hist2->GetCellContent(ix0, ixy) > 0.1)
721 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
722 hist1->SetCellContent(ix0, ixy, cont);
723 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
724 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
728 for (Int_t i = -1; i > -nbinPad; --i) {
729 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
731 Double_t upEdge = axis->GetBinUpEdge(ixy);
732 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
733 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
736 Double_t cont = pad->Charge();
737 if (hist2->GetCellContent(ix0, ixy) > 0.1)
738 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
739 hist1->SetCellContent(ix0, ixy, cont);
740 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
741 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
746 //_____________________________________________________________________________
748 AliMUONClusterFinderMLEM::Plot(const char* basename)
750 /// Make a plot and save it as png
755 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
760 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
761 fDetElemId,fClusterNumber));
764 //_____________________________________________________________________________
766 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
770 /// Compute coefficients needed for MLEM algorithm
772 Int_t npadTot = cluster.Multiplicity();
773 Int_t nPix = fPixArray->GetLast()+1;
775 //memset(probi,0,nPix*sizeof(Double_t));
776 for (Int_t j = 0; j < npadTot*nPix; ++j) coef[j] = 0.;
777 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
779 Int_t mult = cluster.Multiplicity();
780 for ( Int_t j = 0; j < mult; ++j )
782 AliMUONPad* pad = cluster.Pad(j);
785 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
787 Int_t indx1 = indx + ipix;
788 //if (pad->Status() < 0)
789 if (pad->Status() != fgkZero)
794 AliMUONPad* pixPtr = Pixel(ipix);
795 // coef is the charge (given by Mathieson integral) on pad, assuming
796 // the Mathieson is center at pixel.
797 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
798 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
799 // pad->Ix(),pad->Iy(),
800 // pad->X(),pad->Y(),
801 // pad->DX(),pad->DY(),
802 // pixPtr->Coord(0),pixPtr->Coord(1),
803 // pixPtr->Size(0),pixPtr->Size(1),
806 probi[ipix] += coef[indx1];
811 //_____________________________________________________________________________
812 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
814 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
816 // AliCodeTimerAuto("",0)
818 Int_t nPix = fPixArray->GetLast()+1;
820 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
821 //StdoutToAliDebug(2,cluster.Print("full"););
825 AliDebug(1,"No pixels, why am I here then ?");
828 AddVirtualPad(cluster); // add virtual pads if necessary
830 Int_t npadTot = cluster.Multiplicity();
832 for (Int_t i = 0; i < npadTot; ++i)
834 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
835 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()};
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 AliDebug(2,Form("LowestPadCharge=%e",fLowestPadCharge));
908 fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
910 for (Int_t ipix = 0; ipix < nPix; ++ipix)
912 AliMUONPad* pixPtr2 = Pixel(ipix);
913 fHistMlem->Fill(pixPtr2->Coord(0),pixPtr2->Coord(1),pixPtr2->Charge());
916 // Check if the total charge of pixels is too low
918 for ( Int_t i = 0; i < nPix; ++i)
920 qTot += Pixel(i)->Charge();
923 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < fLowestClusterCharge ) )
925 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
929 for ( Int_t i = 0; i < npadTot; ++i)
931 AliMUONPad* pad = cluster.Pad(i);
932 //if ( pad->Status() == 0) pad->SetStatus(-1);
933 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
940 // Simple cluster - skip further passes thru EM-procedure
948 // Calculate position of the center-of-gravity around the maximum pixel
952 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
953 pixPtr->Size(0) > pixPtr->Size(1)) break;
955 // Sort pixels according to the charge
956 MaskPeaks(1); // mask local maxima
958 MaskPeaks(0); // unmask local maxima
959 Double_t pixMin = 0.01*Pixel(0)->Charge();
960 pixMin = TMath::Min(pixMin,100*fLowestPixelCharge);
962 // Decrease pixel size and shift pixels to make them centered at
964 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
967 Double_t shift[2] = { 0.0, 0.0 };
968 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
971 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
973 AliMUONPad* pixPtr2 = Pixel(ipix);
974 if ( nPix >= npadOK // too many pixels already
976 ((pixPtr2->Charge() < pixMin) && (pixPtr2->Status() != fgkMustKeep)) // too low charge
982 for (Int_t i = 0; i < 2; ++i)
986 pixPtr2->SetCharge(fLowestPadCharge);
987 pixPtr2->SetSize(indx, pixPtr2->Size(indx)/2);
988 width = -pixPtr2->Size(indx);
989 pixPtr2->Shift(indx, width);
990 // Shift pixel position
994 for (Int_t j = 0; j < 2; ++j)
996 shift[j] = pixPtr2->Coord(j) - xyCOG[j];
997 shift[j] -= ((Int_t)(shift[j]/pixPtr2->Size(j)/2))*pixPtr2->Size(j)*2;
1000 pixPtr2->Shift(0, -shift[0]);
1001 pixPtr2->Shift(1, -shift[1]);
1004 else if (nPix < npadOK)
1006 pixPtr2 = new AliMUONPad(*pixPtr2);
1007 pixPtr2->Shift(indx, -2*width);
1008 pixPtr2->SetStatus(fgkZero);
1009 fPixArray->Add(pixPtr2);
1012 else continue; // skip adjustment of histo limits
1013 for (Int_t j = 0; j < 4; ++j)
1015 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr2->Coord(j/2));
1017 } // for (Int_t i=0; i<2;
1018 } // for (Int_t ipix=0;
1020 fPixArray->Compress();
1022 AliDebug(2,Form("After shift:"));
1023 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1024 //Plot(Form("mlem.lc%d",lc+1));
1026 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1029 xylim[2],xylim[3]));
1033 AliMUONPad* pixPtr2 = Pixel(0);
1034 // add pixels if the maximum is at the limit of pixel area:
1035 // start from Y-direction
1037 for (Int_t i = 3; i > -1; --i)
1039 if (nPix < npadOK &&
1040 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr2->Size(i/2))
1042 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1043 AliMUONPad* p = new AliMUONPad(*pixPtr2);
1044 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr2->Size(i/2));
1045 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1046 j = TMath::Even (i/2);
1047 p->SetCoord(j, xyCOG[j]);
1048 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1049 //StdoutToAliDebug(2,cout << " ---- ";
1050 // p->Print("corners"););
1056 nPix = fPixArray->GetEntriesFast();
1061 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1062 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1064 // remove pixels with low signal or low visibility
1065 // Cuts are empirical !!!
1066 Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,2.0*fLowestPixelCharge);
1067 thresh = TMath::Min (thresh,100.0*fLowestPixelCharge);
1068 Double_t charge = 0;
1070 // Mark pixels which should be removed
1071 for (Int_t i = 0; i < nPix; ++i)
1073 AliMUONPad* pixPtr2 = Pixel(i);
1074 charge = pixPtr2->Charge();
1075 if (charge < thresh)
1077 pixPtr2->SetCharge(-charge);
1081 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1083 for (Int_t i = 0; i < nPix; ++i)
1085 AliMUONPad* pixPtr2 = Pixel(i);
1086 charge = pixPtr2->Charge();
1087 if (charge > 0) continue;
1088 near = FindNearest(pixPtr2);
1089 pixPtr2->SetCharge(0);
1090 probi[i] = 0; // make it "invisible"
1091 AliMUONPad* pnear = Pixel(near);
1092 pnear->SetCharge(pnear->Charge() + (-charge));
1094 Mlem(cluster,coef,probi,2);
1096 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1097 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1098 //Plot("mlem.beforesplit");
1101 for (Int_t i = 0; i < nPix; ++i)
1103 AliMUONPad* pixPtr2 = Pixel(i);
1104 Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr2->Coord(0));
1105 Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr2->Coord(1));
1106 fHistMlem->SetBinContent(ix, iy, pixPtr2->Charge());
1109 // Try to split into clusters
1111 if (fHistMlem->GetSum() < 2.0*fLowestPixelCharge)
1117 fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
1122 fPixArray->Delete();
1127 //_____________________________________________________________________________
1128 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1130 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1131 /// - to avoid loosing low charge pixels after sorting)
1133 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1134 AliMUONPad* pix = Pixel(i);
1135 if (pix->Status() == fgkMustKeep) {
1136 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1137 else pix->SetCharge(pix->Charge()-10000.);
1142 //_____________________________________________________________________________
1143 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1144 const Double_t* coef, Double_t* probi,
1147 /// Use MLEM to find pixel charges
1149 Int_t nPix = fPixArray->GetEntriesFast();
1151 Int_t npad = cluster.Multiplicity();
1153 Double_t* probi1 = new Double_t[nPix];
1154 Double_t probMax = TMath::MaxElement(nPix,probi);
1156 for (Int_t iter = 0; iter < nIter; ++iter)
1159 for (Int_t ipix = 0; ipix < nPix; ++ipix)
1161 Pixel(ipix)->SetChargeBackup(0);
1162 // Correct each pixel
1164 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1166 probi1[ipix] = probMax;
1167 for (Int_t j = 0; j < npad; ++j)
1169 AliMUONPad* pad = cluster.Pad(j);
1170 //if (pad->Status() < 0) continue;
1171 if (pad->Status() != fgkZero) continue;
1173 Int_t indx1 = j*nPix;
1174 Int_t indx = indx1 + ipix;
1175 // Calculate expectation
1176 for (Int_t i = 0; i < nPix; ++i)
1178 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1179 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1181 if ( pad->IsSaturated() && sum1 > pad->Charge() )
1183 // correct for pad charge overflows
1184 probi1[ipix] -= coef[indx];
1186 //sum1 = pad->Charge();
1189 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1190 } // for (Int_t j=0;
1191 AliMUONPad* pixPtr = Pixel(ipix);
1192 if (probi1[ipix] > 1.e-6)
1194 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1195 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1197 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1198 } // for (Int_t ipix=0;
1200 for (Int_t i = 0; i < nPix; ++i) {
1201 AliMUONPad* pixPtr = Pixel(i);
1202 pixPtr->RevertCharge();
1203 qTot += pixPtr->Charge();
1206 // Can happen in clusters with large number of overflows - speeding up
1210 } // for (Int_t iter=0;
1214 //_____________________________________________________________________________
1215 void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
1217 /// Calculate position of the center-of-gravity around the maximum pixel
1219 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1220 Int_t i1 = -9, j1 = -9;
1221 fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1222 Int_t nx = fHistMlem->GetNbinsX();
1223 Int_t ny = fHistMlem->GetNbinsY();
1224 Double_t thresh = fHistMlem->GetMaximum()/10;
1225 Double_t x, y, cont, xq=0, yq=0, qq=0;
1227 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1228 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1229 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1230 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1231 cont = fHistMlem->GetCellContent(j,i);
1232 if (cont < thresh) continue;
1233 if (i != i1) {i1 = i; nsumy++;}
1234 if (j != j1) {j1 = j; nsumx++;}
1235 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1244 Int_t i2 = 0, j2 = 0;
1247 // one bin in Y - add one more (with the largest signal)
1248 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1249 if (i == iymax) continue;
1250 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1251 cont = fHistMlem->GetCellContent(j,i);
1254 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1255 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1264 if (i2 != i1) nsumy++;
1265 if (j2 != j1) nsumx++;
1267 } // if (nsumy == 1)
1270 // one bin in X - add one more (with the largest signal)
1272 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1273 if (j == ixmax) continue;
1274 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1275 cont = fHistMlem->GetCellContent(j,i);
1278 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1279 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1288 if (i2 != i1) nsumy++;
1289 if (j2 != j1) nsumx++;
1291 } // if (nsumx == 1)
1293 xyc[0] = xq/qq; xyc[1] = yq/qq;
1294 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1297 //_____________________________________________________________________________
1298 Int_t AliMUONClusterFinderMLEM::FindNearest(const AliMUONPad *pixPtr0)
1300 /// Find the pixel nearest to the given one
1301 /// (algorithm may be not very efficient)
1303 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1304 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1305 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1308 for (Int_t i = 0; i < nPix; ++i) {
1309 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1310 if (pixPtr == pixPtr0 || pixPtr->Charge() < fLowestPixelCharge) continue;
1311 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1312 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1313 r = dx *dx + dy * dy;
1314 if (r < rmin) { rmin = r; imin = i; }
1319 //_____________________________________________________________________________
1321 AliMUONClusterFinderMLEM::Paint(Option_t*)
1323 /// Paint cluster and pixels
1325 AliMpArea area(fPreCluster->Area());
1327 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1329 gVirtualX->SetFillStyle(1001);
1330 gVirtualX->SetFillColor(3);
1331 gVirtualX->SetLineColor(3);
1335 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1337 AliMUONPad* pixel = Pixel(i);
1339 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1340 pixel->Coord(1)-pixel->Size(1)*s,
1341 pixel->Coord(0)+pixel->Size(0)*s,
1342 pixel->Coord(1)+pixel->Size(1)*s);
1344 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1346 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1347 // pixel->Coord(1) + sign*pixel->Size(1),
1348 // pixel->Coord(0) + pixel->Size(0),
1349 // pixel->Coord(1) - sign*pixel->Size(1)
1355 gVirtualX->SetFillStyle(0);
1357 fPreCluster->Paint();
1359 gVirtualX->SetLineColor(1);
1360 gVirtualX->SetLineWidth(2);
1361 gVirtualX->SetFillStyle(0);
1362 gVirtualX->SetTextColor(1);
1363 gVirtualX->SetTextAlign(22);
1365 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1367 AliMUONPad* pixel = Pixel(i);
1368 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1369 pixel->Coord(1)-pixel->Size(1),
1370 pixel->Coord(0)+pixel->Size(0),
1371 pixel->Coord(1)+pixel->Size(1));
1372 gVirtualX->SetTextSize(pixel->Size(0)*60);
1374 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1378 //_____________________________________________________________________________
1379 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1381 /// Find local maxima in pixel space for large preclusters in order to
1382 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1383 /// or to find additional fitting seeds if clusters were not completely resolved
1385 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1387 Double_t xylim[4] = {999, 999, 999, 999};
1389 Int_t nPix = pixArray->GetEntriesFast();
1390 AliMUONPad *pixPtr = 0;
1391 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1392 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1393 for (Int_t i = 0; i < 4; ++i)
1394 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1396 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
1398 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1399 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1400 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1401 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1402 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1403 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1404 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1406 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1408 Int_t nMax = 0, indx, nxy = ny * nx;
1409 Int_t *isLocalMax = new Int_t[nxy];
1410 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
1412 for (Int_t i = 1; i <= ny; ++i) {
1414 for (Int_t j = 1; j <= nx; ++j) {
1415 if (fHistAnode->GetCellContent(j,i) < fLowestPixelCharge) continue;
1416 //if (isLocalMax[indx+j-1] < 0) continue;
1417 if (isLocalMax[indx+j-1] != 0) continue;
1418 FlagLocalMax(fHistAnode, i, j, isLocalMax);
1422 for (Int_t i = 1; i <= ny; ++i) {
1424 for (Int_t j = 1; j <= nx; ++j) {
1425 if (isLocalMax[indx+j-1] > 0) {
1426 localMax[nMax] = indx + j - 1;
1427 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
1428 ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
1429 if (nMax > 99) break;
1433 AliError(" Too many local maxima !!!");
1437 if (fDebug) cout << " Local max: " << nMax << endl;
1438 delete [] isLocalMax;
1442 //_____________________________________________________________________________
1443 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1445 /// Flag pixels (whether or not local maxima)
1447 Int_t nx = hist->GetNbinsX();
1448 Int_t ny = hist->GetNbinsY();
1449 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1450 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1452 Int_t ie = i + 2, je = j + 2;
1453 for (Int_t i1 = i-1; i1 < ie; ++i1) {
1454 if (i1 < 1 || i1 > ny) continue;
1455 indx1 = (i1 - 1) * nx;
1456 for (Int_t j1 = j-1; j1 < je; ++j1) {
1457 if (j1 < 1 || j1 > nx) continue;
1458 if (i == i1 && j == j1) continue;
1459 indx2 = indx1 + j1 - 1;
1460 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1461 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1462 else if (cont > cont1) isLocalMax[indx2] = -1;
1463 else { // the same charge
1464 isLocalMax[indx] = 1;
1465 if (isLocalMax[indx2] == 0) {
1466 FlagLocalMax(hist, i1, j1, isLocalMax);
1467 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1468 else isLocalMax[indx2] = -1;
1473 isLocalMax[indx] = 1; // local maximum
1476 //_____________________________________________________________________________
1477 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1478 const Int_t *localMax, Int_t iMax)
1480 /// Find pixel cluster around local maximum \a iMax and pick up pads
1481 /// overlapping with it
1484 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1486 hist->Draw("lego1Fb"); // debug
1491 Int_t nx = fHistAnode->GetNbinsX();
1492 Int_t ny = fHistAnode->GetNbinsY();
1493 Int_t ic = localMax[iMax] / nx + 1;
1494 Int_t jc = localMax[iMax] % nx + 1;
1495 Int_t nxy = ny * nx;
1496 Bool_t *used = new Bool_t[nxy];
1497 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
1499 // Drop all pixels from the array - pick up only the ones from the cluster
1500 fPixArray->Delete();
1502 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1503 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1504 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1505 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1506 Double_t cont = fHistAnode->GetCellContent(jc,ic);
1507 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1508 used[(ic-1)*nx+jc-1] = kTRUE;
1509 AddBinSimple(fHistAnode, ic, jc);
1510 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1512 Int_t nPix = fPixArray->GetEntriesFast();
1513 Int_t npad = cluster.Multiplicity();
1515 for (Int_t i = 0; i < nPix; ++i)
1517 AliMUONPad* pixPtr = Pixel(i);
1518 pixPtr->SetSize(0,wx);
1519 pixPtr->SetSize(1,wy);
1522 // Pick up pads which overlap with found pixels
1523 for (Int_t i = 0; i < npad; ++i)
1525 //cluster.Pad(i)->SetStatus(-1);
1526 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1529 for (Int_t i = 0; i < nPix; ++i)
1531 AliMUONPad* pixPtr = Pixel(i);
1532 for (Int_t j = 0; j < npad; ++j)
1534 AliMUONPad* pad = cluster.Pad(j);
1535 //if (pad->Status() == 0) continue;
1536 if (pad->Status() == fgkZero) continue;
1537 if ( Overlap(*pad,*pixPtr) )
1539 //pad->SetStatus(0);
1540 pad->SetStatus(fgkZero);
1541 if (fDebug) { cout << j << " "; pad->Print("full"); }
1549 //_____________________________________________________________________________
1551 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
1553 /// Add adjacent bins (+-1 in X and Y) to the cluster
1555 Int_t nx = hist->GetNbinsX();
1556 Int_t ny = hist->GetNbinsY();
1557 Double_t cont1, cont = hist->GetCellContent(jc,ic);
1558 AliMUONPad *pixPtr = 0;
1560 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1561 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1562 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1563 cont1 = hist->GetCellContent(j,i);
1564 if (cont1 > cont) continue;
1565 if (cont1 < fLowestPixelCharge) continue;
1566 pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j),
1567 hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1568 fPixArray->Add(pixPtr);
1573 //_____________________________________________________________________________
1574 AliMUONClusterFinderMLEM&
1575 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1577 /// Protected assignement operator
1579 if (this == &rhs) return *this;
1581 AliFatal("Not implemented.");
1586 //_____________________________________________________________________________
1587 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1589 /// Add virtual pad (with small charge) to improve fit for clusters
1590 /// with number of pads == 2 per direction
1592 // Find out non-bending and bending planes
1593 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1595 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1596 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1597 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1602 Bool_t same = kFALSE;
1603 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1606 Bool_t check[2] = {kFALSE, kFALSE};
1608 nxy[0] = nxy[1] = 0;
1609 for (Int_t inb = 0; inb < 2; ++inb) {
1610 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1611 if (inb == 0 && AliMp::PairFirst(cn) == 2) check[inb] = kTRUE; // check non-bending plane
1612 else if (inb == 1 && AliMp::PairSecond(cn) == 2) check[inb] = kTRUE; // check bending plane
1614 nxy[0] = TMath::Max (nxy[0], AliMp::PairFirst(cn));
1615 nxy[1] = TMath::Max (nxy[1], AliMp::PairSecond(cn));
1616 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1617 else if (inb == 1 && AliMp::PairSecond(cn) < 2) nonb[inb] = !nonb[inb];
1621 if (nxy[0] > 2) check[0] = kFALSE;
1622 if (nxy[1] > 2) check[1] = kFALSE;
1624 if (!check[0] && !check[1]) return;
1626 for (Int_t inb = 0; inb < 2; ++inb) {
1627 if (!check[inb]) continue;
1629 // Find pads with maximum and next to maximum charges
1630 Int_t maxPads[2] = {-1, -1};
1631 Double_t amax[2] = {0};
1632 Int_t mult = cluster.Multiplicity();
1633 for (Int_t j = 0; j < mult; ++j) {
1634 AliMUONPad *pad = cluster.Pad(j);
1635 if (pad->Cathode() != nonb[inb]) continue;
1636 for (Int_t j2 = 0; j2 < 2; ++j2) {
1637 if (pad->Charge() > amax[j2]) {
1638 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1639 amax[j2] = pad->Charge();
1646 // Find min and max dimensions of the cluster
1647 Double_t limits[2] = {9999, -9999};
1648 for (Int_t j = 0; j < mult; ++j) {
1649 AliMUONPad *pad = cluster.Pad(j);
1650 if (pad->Cathode() != nonb[inb]) continue;
1651 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1652 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1655 // Loop over max and next to max pads
1656 Bool_t add = kFALSE;
1658 for (Int_t j = 0; j < 2; ++j) {
1660 if (maxPads[j] < 0) continue;
1662 if (amax[1] / amax[0] < 0.5) break;
1664 // Check if pad at the cluster limit
1665 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1667 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1668 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1670 //cout << " *** Pad not at the cluster limit: " << j << endl;
1673 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1675 // Add pad (if it exists)
1677 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1678 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1679 //AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1680 AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos.X(), pos.Y(),kFALSE);
1681 if (!mppad.IsValid()) continue; // non-existing pad
1682 AliMUONPad muonPad(fDetElemId, nonb[inb], mppad.GetIx(), mppad.GetIy(),
1683 mppad.GetPositionX(), mppad.GetPositionY(),
1684 mppad.GetDimensionX(), mppad.GetDimensionY(), 0);
1685 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, fLowestPadCharge));
1686 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1687 else muonPad.SetCharge(TMath::Min (amax[j]/15, fLowestPadCharge));
1688 if (muonPad.Charge() < 2.0*fLowestPixelCharge) muonPad.SetCharge(2.0*fLowestPixelCharge);
1689 muonPad.SetReal(kFALSE);
1690 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1691 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1692 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1693 cluster.AddPad(muonPad); // add pad to the cluster
1700 //_____________________________________________________________________________
1701 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1702 Int_t &nInX, Int_t &nInY) const
1704 /// Find number of pads in X and Y-directions (excluding virtual ones and
1707 //Int_t statusToTest = 1;
1708 Int_t statusToTest = fgkUseForFit;
1710 //if ( nInX < 0 ) statusToTest = 0;
1711 if ( nInX < 0 ) statusToTest = fgkZero;
1713 Bool_t mustMatch(kTRUE);
1715 Long_t cn = cluster.NofPads(statusToTest,mustMatch);
1717 nInX = AliMp::PairFirst(cn);
1718 nInY = AliMp::PairSecond(cn);
1721 //_____________________________________________________________________________
1722 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1724 /// Remove pixel at index i
1725 AliMUONPad* pixPtr = Pixel(i);
1726 fPixArray->RemoveAt(i);
1730 //_____________________________________________________________________________
1731 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1733 /// Process simple cluster (small number of pads) without EM-procedure
1735 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1736 Double_t parOk[3] = {0.};
1737 TObjArray *clusters[1];
1738 clusters[0] = fPixArray;
1740 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1742 Int_t mult = cluster.Multiplicity();
1743 for (Int_t i = 0; i < mult; ++i)
1745 AliMUONPad* pad = cluster.Pad(i);
1747 if ( pad->IsSaturated())
1756 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1758 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
1761 //_____________________________________________________________________________
1763 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1765 /// Returns pixel at index i
1766 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1769 //_____________________________________________________________________________
1771 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1774 TString swhat(what);
1776 if ( swhat.Contains("precluster") )
1778 if ( fPreCluster) fPreCluster->Print();
1782 //_____________________________________________________________________________
1784 AliMUONClusterFinderMLEM::SetChargeHints(Double_t lowestPadCharge, Double_t lowestClusterCharge)
1786 /// Set some thresholds we use later on in the algorithm
1787 fLowestPadCharge=lowestPadCharge;
1788 fLowestClusterCharge=lowestClusterCharge;
1789 fLowestPixelCharge=fLowestPadCharge/12.0;