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,
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
161 // do we have clusters in our list ?
162 if ( fClusterNumber < fClusterList.GetLast() )
164 o = fClusterList.At(++fClusterNumber);
167 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
169 //FIXME : at this point, must check whether we've used all the digits
170 //from precluster : if not, let the preclustering know about those unused
171 //digits, so it can reuse them
173 // if the cluster list is exhausted, we need to go to the next
174 // pre-cluster and treat it
176 fPreCluster = fPreClusterFinder->NextCluster();
179 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
180 fClusterNumber = -1; //AZ
190 // WorkOnPreCluster may have used only part of the pads, so we check that
191 // now, and let the unused pads be reused by the preclustering...
193 Int_t mult = fPreCluster->Multiplicity();
194 for ( Int_t i = 0; i < mult; ++i )
196 AliMUONPad* pad = fPreCluster->Pad(i);
197 if ( !pad->IsUsed() )
199 fPreClusterFinder->UsePad(*pad);
203 return NextCluster();
206 //_____________________________________________________________________________
208 AliMUONClusterFinderMLEM::WorkOnPreCluster()
210 /// Starting from a precluster, builds a pixel array, and then
211 /// extract clusters from this array
213 // AliCodeTimerAuto("",0)
216 cout << " *** Event # " << fEventNumber
217 << " det. elem.: " << fDetElemId << endl;
218 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
219 AliMUONPad* pad = fPreCluster->Pad(j);
220 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
221 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
222 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
226 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
227 if (!cluster) return kFALSE;
229 BuildPixArray(*cluster);
231 if ( fPixArray->GetLast() < 0 )
233 AliDebug(1,"No pixel for the above cluster");
238 // Use MLEM for cluster finder
239 Int_t nMax = 1, localMax[100], maxPos[100];
240 Double_t maxVal[100];
242 Int_t iSimple = 0, nInX = -1, nInY;
244 PadsInXandY(*cluster,nInX, nInY);
246 if (nInX < 4 && nInY < 4)
248 iSimple = 1; // simple cluster
252 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
254 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
255 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
259 for (Int_t i = 0; i < nMax; ++i)
263 FindCluster(*cluster,localMax, maxPos[i]);
266 MainLoop(*cluster,iSimple);
270 Int_t mult = cluster->Multiplicity();
271 for (Int_t j = 0; j < mult; ++j)
273 AliMUONPad* pad = cluster->Pad(j);
274 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
275 if ( pad->Status() != fgkOver ) continue; // pad was not used
277 pad->SetStatus(fgkZero);
278 pad->RevertCharge(); // use backup charge value
281 } // for (Int_t i=0; i<nMax;
284 fHistMlem = fHistAnode = 0x0;
289 //_____________________________________________________________________________
291 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
293 /// Check if the pad and the pixel overlaps
295 // make a fake pad from the pixel
296 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
297 pixel.Coord(0),pixel.Coord(1),
298 pixel.Size(0),pixel.Size(1),0);
300 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
303 //_____________________________________________________________________________
305 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
307 /// Check precluster in order to attempt to simplify it (mostly for
308 /// two-cathode preclusters)
310 AliCodeTimerAuto("",0)
312 // Disregard small clusters (leftovers from splitting or noise)
313 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
314 origCluster.Charge(0)+origCluster.Charge(1) < fLowestClusterCharge )
319 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
321 AliDebug(2,"Start of CheckPreCluster=");
322 //StdoutToAliDebug(2,cluster->Print("full"));
324 AliMUONCluster* rv(0x0);
326 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
328 rv = CheckPreclusterTwoCathodes(cluster);
337 //_____________________________________________________________________________
339 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
341 /// Check two-cathode cluster
343 Int_t npad = cluster->Multiplicity();
344 Int_t* flags = new Int_t[npad];
345 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
347 // Check pad overlaps
348 for ( Int_t i = 0; i < npad; ++i)
350 AliMUONPad* padi = cluster->Pad(i);
351 if ( padi->Cathode() != 0 ) continue;
352 for (Int_t j = i+1; j < npad; ++j)
354 AliMUONPad* padj = cluster->Pad(j);
355 if ( padj->Cathode() != 1 ) continue;
356 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
357 flags[i] = flags[j] = 1; // mark overlapped pads
361 // Check if all pads overlap
363 for (Int_t i = 0; i < npad; ++i)
365 if (!flags[i]) ++nFlags;
370 // not all pads overlap.
371 if (fDebug) cout << " nFlags: " << nFlags << endl;
372 TObjArray toBeRemoved;
373 for (Int_t i = 0; i < npad; ++i)
375 AliMUONPad* pad = cluster->Pad(i);
376 if (flags[i]) continue;
377 Int_t cath = pad->Cathode();
378 Int_t cath1 = TMath::Even(cath);
379 // Check for edge effect (missing pads on the _other_ cathode)
381 fkSegmentation[cath1]->PadByPosition(pad->Position().X(),
382 pad->Position().Y(),kFALSE);
383 if (!mpPad.IsValid()) continue;
384 if (nFlags == 1 && pad->Charge() < fLowestPadCharge) continue;
385 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
386 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
387 toBeRemoved.AddLast(pad);
388 fPreCluster->Pad(i)->Release();
390 Int_t nRemove = toBeRemoved.GetEntriesFast();
391 for ( Int_t i = 0; i < nRemove; ++i )
393 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
397 // Check correlations of cathode charges
398 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
401 Int_t cathode = cluster->MaxRawChargeCathode();
407 // get min and max pad charges on the cathode opposite to the
408 // max pad (given by MaxRawChargeCathode())
410 Int_t mult = cluster->Multiplicity();
411 for ( Int_t i = 0; i < mult; ++i )
413 AliMUONPad* pad = cluster->Pad(i);
414 if ( pad->Cathode() != cathode || !pad->IsReal() )
416 // only consider pads in the opposite cathode, and
417 // only consider real pads (i.e. exclude the virtual ones)
420 if ( pad->Charge() < cmin )
422 cmin = pad->Charge();
429 else if ( pad->Charge() > cmax )
431 cmax = pad->Charge();
435 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
436 imin,imax,cmin,cmax));
438 // arrange pads according to their distance to the max, normalized
440 Double_t* dist = new Double_t[mult];
445 AliMUONPad* padmax = cluster->Pad(imax);
447 for ( Int_t i = 0; i < mult; ++i )
450 if ( i == imax) continue;
451 AliMUONPad* pad = cluster->Pad(i);
452 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
453 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
454 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
455 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
458 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
464 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
465 Double_t xmax(-1), distPrev(999);
466 TObjArray toBeRemoved;
468 for ( Int_t i = 0; i < mult; ++i )
470 Int_t indx = flags[i];
471 AliMUONPad* pad = cluster->Pad(indx);
472 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
473 if ( dist[indx] > dmin )
475 // farther than the minimum pad
476 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
477 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
480 if (dx >= 0 && dy >= 0) continue;
481 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
482 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
484 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
485 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
488 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
490 cmax = TMath::Max(pad->Charge(),cmax);
494 cmax = pad->Charge();
497 distPrev = dist[indx];
498 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
499 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
502 toBeRemoved.AddLast(pad);
503 fPreCluster->Pad(indx)->Release();
506 Int_t nRemove = toBeRemoved.GetEntriesFast();
507 for ( Int_t i = 0; i < nRemove; ++i )
509 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
516 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
517 //StdoutToAliDebug(2,cluster->Print("full"));
522 //_____________________________________________________________________________
524 AliMUONClusterFinderMLEM::CheckOverlaps()
526 /// For debug only : check if some pixels overlap...
528 Int_t nPix = fPixArray->GetLast()+1;
531 for ( Int_t i = 0; i < nPix; ++i )
533 AliMUONPad* pixelI = Pixel(i);
534 AliMUONPad pi(dummy,dummy,dummy,dummy,
535 pixelI->Coord(0),pixelI->Coord(1),
536 pixelI->Size(0),pixelI->Size(1),0.0);
538 for ( Int_t j = i+1; j < nPix; ++j )
540 AliMUONPad* pixelJ = Pixel(j);
541 AliMUONPad pj(dummy,dummy,dummy,dummy,
542 pixelJ->Coord(0),pixelJ->Coord(1),
543 pixelJ->Size(0),pixelJ->Size(1),0.0);
546 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
548 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
550 StdoutToAliInfo(pixelI->Print();
551 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
553 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
554 cout << " Area surface = " << area.GetDimensionX()*area.GetDimensionY()*4 << endl;
555 cout << "-------" << endl;
563 //_____________________________________________________________________________
564 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
566 /// Build pixel array for MLEM method
568 Int_t npad = cluster.Multiplicity();
571 AliWarning("Got no pad at all ?!");
575 BuildPixArrayOneCathode(cluster);
577 Int_t nPix = fPixArray->GetLast()+1;
579 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
583 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
585 // Too many pixels - sort and remove pixels with the lowest signal
587 for ( Int_t i = npad; i < nPix; ++i )
591 fPixArray->Compress();
592 } // if (nPix > npad)
594 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
595 // fPixArray->Print(););
596 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
599 //_____________________________________________________________________________
600 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
602 /// Build the pixel array
604 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
606 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
607 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2]={99999,99999};
608 Int_t found[2] = {0,0}, mult = cluster.Multiplicity();
610 for ( Int_t i = 0; i < mult; ++i) {
611 AliMUONPad* pad = cluster.Pad(i);
612 for (Int_t j = 0; j < 2; ++j) {
613 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
614 xy0[j] = pad->Coord(j);
618 if (found[0] && found[1]) break;
621 Double_t min[2], max[2];
622 Int_t cath0 = 0, cath1 = 1;
623 if (cluster.Multiplicity(0) == 0) cath0 = 1;
624 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
627 Double_t leftDownX, leftDownY;
628 cluster.Area(cath0).LeftDownCorner(leftDownX, leftDownY);
629 Double_t rightUpX, rightUpY;
630 cluster.Area(cath0).RightUpCorner(rightUpX, rightUpY);
635 if (cath1 != cath0) {
636 cluster.Area(cath1).LeftDownCorner(leftDownX, leftDownY);
637 cluster.Area(cath1).RightUpCorner(rightUpX, rightUpY);
638 min[0] = TMath::Max (min[0], leftDownX);
639 min[1] = TMath::Max (min[1], leftDownY);
640 max[0] = TMath::Min (max[0], rightUpX);
641 max[1] = TMath::Min (max[1], rightUpY);
645 //width[0] /= 2; width[1] /= 2; // just for check
646 Int_t nbins[2]={0,0};
647 for (Int_t i = 0; i < 2; ++i) {
648 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
649 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
650 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
651 + TMath::Sign(0.5,dist)) * width[i] * 2;
652 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
653 if (nbins[i] == 0) ++nbins[i];
654 max[i] = min[i] + nbins[i] * width[i] * 2;
655 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
659 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
660 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
661 TAxis *xaxis = hist1->GetXaxis();
662 TAxis *yaxis = hist1->GetYaxis();
665 for ( Int_t i = 0; i < mult; ++i) {
666 AliMUONPad* pad = cluster.Pad(i);
667 Int_t ix0 = xaxis->FindBin(pad->X());
668 Int_t iy0 = yaxis->FindBin(pad->Y());
669 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
673 for (Int_t i = 1; i <= nbins[0]; ++i) {
674 Double_t x = xaxis->GetBinCenter(i);
675 for (Int_t j = 1; j <= nbins[1]; ++j) {
676 if (hist2->GetCellContent(i,j) < 0.1) continue;
677 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
678 // cluster.Multiplicity(1)) continue;
679 if (cath0 != cath1) {
681 Double_t cont = hist2->GetCellContent(i,j);
682 if (cont < 999.) continue;
683 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
685 Double_t y = yaxis->GetBinCenter(j);
686 Double_t charge = hist1->GetCellContent(i,j);
687 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
688 fPixArray->Add(pixPtr);
692 if (fPixArray->GetEntriesFast() == 1) {
693 // Split pixel into 2
694 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
695 pixPtr->SetSize(0,width[0]/2.);
696 pixPtr->Shift(0,-width[0]/4.);
697 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
698 fPixArray->Add(pixPtr);
701 //fPixArray->Print();
706 //_____________________________________________________________________________
707 void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
708 TH2D *hist1, TH2D *hist2)
710 /// "Span" pad over histogram in the direction idir
712 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
713 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
714 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
716 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
718 for (Int_t i = 0; i < nbinPad; ++i) {
719 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
720 if (ixy > nbins) break;
721 Double_t lowEdge = axis->GetBinLowEdge(ixy);
722 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
723 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
726 Double_t cont = pad->Charge();
727 if (hist2->GetCellContent(ix0, ixy) > 0.1)
728 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
729 hist1->SetCellContent(ix0, ixy, cont);
730 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
731 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
735 for (Int_t i = -1; i > -nbinPad; --i) {
736 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
738 Double_t upEdge = axis->GetBinUpEdge(ixy);
739 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
740 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
743 Double_t cont = pad->Charge();
744 if (hist2->GetCellContent(ix0, ixy) > 0.1)
745 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
746 hist1->SetCellContent(ix0, ixy, cont);
747 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
748 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
753 //_____________________________________________________________________________
755 AliMUONClusterFinderMLEM::Plot(const char* /*basename*/)
757 /// Make a plot and save it as png
760 // if (!fPlot) return;
762 // TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
767 // c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
768 // fDetElemId,fClusterNumber));
771 //_____________________________________________________________________________
773 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
777 /// Compute coefficients needed for MLEM algorithm
779 Int_t npadTot = cluster.Multiplicity();
780 Int_t nPix = fPixArray->GetLast()+1;
782 //memset(probi,0,nPix*sizeof(Double_t));
783 for (Int_t j = 0; j < npadTot*nPix; ++j) coef[j] = 0.;
784 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
786 Int_t mult = cluster.Multiplicity();
787 for ( Int_t j = 0; j < mult; ++j )
789 AliMUONPad* pad = cluster.Pad(j);
792 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
794 Int_t indx1 = indx + ipix;
795 //if (pad->Status() < 0)
796 if (pad->Status() != fgkZero)
801 AliMUONPad* pixPtr = Pixel(ipix);
802 // coef is the charge (given by Mathieson integral) on pad, assuming
803 // the Mathieson is center at pixel.
804 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
805 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
806 // pad->Ix(),pad->Iy(),
807 // pad->X(),pad->Y(),
808 // pad->DX(),pad->DY(),
809 // pixPtr->Coord(0),pixPtr->Coord(1),
810 // pixPtr->Size(0),pixPtr->Size(1),
813 probi[ipix] += coef[indx1];
818 //_____________________________________________________________________________
819 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
821 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
823 // AliCodeTimerAuto("",0)
825 Int_t nPix = fPixArray->GetLast()+1;
827 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
828 //StdoutToAliDebug(2,cluster.Print("full"););
832 AliDebug(1,"No pixels, why am I here then ?");
835 AddVirtualPad(cluster); // add virtual pads if necessary
837 Int_t npadTot = cluster.Multiplicity();
839 for (Int_t i = 0; i < npadTot; ++i)
841 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
842 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
846 Double_t* probi(0x0);
847 Int_t lc(0); // loop counter
849 //Plot("mlem.start");
850 AliMUONPad* pixPtr = Pixel(0);
851 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
859 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
860 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
861 //StdoutToAliDebug(2,fPixArray->Print("full"));
863 coef = new Double_t [npadTot*nPix];
864 probi = new Double_t [nPix];
866 // Calculate coefficients and pixel visibilities
867 ComputeCoefficients(cluster,coef,probi);
869 for (Int_t ipix = 0; ipix < nPix; ++ipix)
871 if (probi[ipix] < 0.01)
873 AliMUONPad* pixel = Pixel(ipix);
874 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
875 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
876 pixel->SetCharge(0); // "invisible" pixel
881 Mlem(cluster,coef, probi, 15);
883 // Find histogram limits for the 1'st pass only - for others computed below
885 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
887 pixPtr = Pixel(ipix);
888 for ( Int_t i = 0; i < 2; ++i )
891 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
892 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
895 } else pixPtr = Pixel(0);
897 for (Int_t i = 0; i < 4; i++)
899 xylim[i] -= pixPtr->Size(i/2);
902 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
903 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
905 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
906 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
907 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
908 xylim[0],-xylim[1],xylim[2],-xylim[3]
911 AliDebug(2,Form("LowestPadCharge=%e",fLowestPadCharge));
915 fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
917 for (Int_t ipix = 0; ipix < nPix; ++ipix)
919 AliMUONPad* pixPtr2 = Pixel(ipix);
920 fHistMlem->Fill(pixPtr2->Coord(0),pixPtr2->Coord(1),pixPtr2->Charge());
923 // Check if the total charge of pixels is too low
925 for ( Int_t i = 0; i < nPix; ++i)
927 qTot += Pixel(i)->Charge();
930 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < fLowestClusterCharge ) )
932 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
936 for ( Int_t i = 0; i < npadTot; ++i)
938 AliMUONPad* pad = cluster.Pad(i);
939 //if ( pad->Status() == 0) pad->SetStatus(-1);
940 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
947 // Simple cluster - skip further passes thru EM-procedure
955 // Calculate position of the center-of-gravity around the maximum pixel
959 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
960 pixPtr->Size(0) > pixPtr->Size(1)) break;
962 // Sort pixels according to the charge
963 MaskPeaks(1); // mask local maxima
965 MaskPeaks(0); // unmask local maxima
966 Double_t pixMin = 0.01*Pixel(0)->Charge();
967 pixMin = TMath::Min(pixMin,100*fLowestPixelCharge);
969 // Decrease pixel size and shift pixels to make them centered at
971 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
974 Double_t shift[2] = { 0.0, 0.0 };
975 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
978 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
980 AliMUONPad* pixPtr2 = Pixel(ipix);
981 if ( nPix >= npadOK // too many pixels already
983 ((pixPtr2->Charge() < pixMin) && (pixPtr2->Status() != fgkMustKeep)) // too low charge
989 for (Int_t i = 0; i < 2; ++i)
993 pixPtr2->SetCharge(fLowestPadCharge);
994 pixPtr2->SetSize(indx, pixPtr2->Size(indx)/2);
995 width = -pixPtr2->Size(indx);
996 pixPtr2->Shift(indx, width);
997 // Shift pixel position
1001 for (Int_t j = 0; j < 2; ++j)
1003 shift[j] = pixPtr2->Coord(j) - xyCOG[j];
1004 shift[j] -= ((Int_t)(shift[j]/pixPtr2->Size(j)/2))*pixPtr2->Size(j)*2;
1007 pixPtr2->Shift(0, -shift[0]);
1008 pixPtr2->Shift(1, -shift[1]);
1011 else if (nPix < npadOK)
1013 pixPtr2 = new AliMUONPad(*pixPtr2);
1014 pixPtr2->Shift(indx, -2*width);
1015 pixPtr2->SetStatus(fgkZero);
1016 fPixArray->Add(pixPtr2);
1019 else continue; // skip adjustment of histo limits
1020 for (Int_t j = 0; j < 4; ++j)
1022 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr2->Coord(j/2));
1024 } // for (Int_t i=0; i<2;
1025 } // for (Int_t ipix=0;
1027 fPixArray->Compress();
1029 AliDebug(2,Form("After shift:"));
1030 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1031 //Plot(Form("mlem.lc%d",lc+1));
1033 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1036 xylim[2],xylim[3]));
1040 AliMUONPad* pixPtr2 = Pixel(0);
1041 // add pixels if the maximum is at the limit of pixel area:
1042 // start from Y-direction
1044 for (Int_t i = 3; i > -1; --i)
1046 if (nPix < npadOK &&
1047 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr2->Size(i/2))
1049 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1050 AliMUONPad* p = new AliMUONPad(*pixPtr2);
1051 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr2->Size(i/2));
1052 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1053 j = TMath::Even (i/2);
1054 p->SetCoord(j, xyCOG[j]);
1055 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1056 //StdoutToAliDebug(2,cout << " ---- ";
1057 // p->Print("corners"););
1063 nPix = fPixArray->GetEntriesFast();
1068 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1069 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1071 // remove pixels with low signal or low visibility
1072 // Cuts are empirical !!!
1073 Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,2.0*fLowestPixelCharge);
1074 thresh = TMath::Min (thresh,100.0*fLowestPixelCharge);
1075 Double_t charge = 0;
1077 // Mark pixels which should be removed
1078 for (Int_t i = 0; i < nPix; ++i)
1080 AliMUONPad* pixPtr2 = Pixel(i);
1081 charge = pixPtr2->Charge();
1082 if (charge < thresh)
1084 pixPtr2->SetCharge(-charge);
1088 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1090 for (Int_t i = 0; i < nPix; ++i)
1092 AliMUONPad* pixPtr2 = Pixel(i);
1093 charge = pixPtr2->Charge();
1094 if (charge > 0) continue;
1095 near = FindNearest(pixPtr2);
1096 pixPtr2->SetCharge(0);
1097 probi[i] = 0; // make it "invisible"
1098 AliMUONPad* pnear = Pixel(near);
1099 pnear->SetCharge(pnear->Charge() + (-charge));
1101 Mlem(cluster,coef,probi,2);
1103 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1104 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1105 //Plot("mlem.beforesplit");
1108 for (Int_t i = 0; i < nPix; ++i)
1110 AliMUONPad* pixPtr2 = Pixel(i);
1111 Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr2->Coord(0));
1112 Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr2->Coord(1));
1113 fHistMlem->SetBinContent(ix, iy, pixPtr2->Charge());
1116 // Try to split into clusters
1118 if (fHistMlem->GetSum() < 2.0*fLowestPixelCharge)
1124 fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
1129 fPixArray->Delete();
1134 //_____________________________________________________________________________
1135 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1137 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1138 /// - to avoid loosing low charge pixels after sorting)
1140 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1141 AliMUONPad* pix = Pixel(i);
1142 if (pix->Status() == fgkMustKeep) {
1143 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1144 else pix->SetCharge(pix->Charge()-10000.);
1149 //_____________________________________________________________________________
1150 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1151 const Double_t* coef, Double_t* probi,
1154 /// Use MLEM to find pixel charges
1156 Int_t nPix = fPixArray->GetEntriesFast();
1158 Int_t npad = cluster.Multiplicity();
1160 Double_t* probi1 = new Double_t[nPix];
1161 Double_t probMax = TMath::MaxElement(nPix,probi);
1163 for (Int_t iter = 0; iter < nIter; ++iter)
1166 for (Int_t ipix = 0; ipix < nPix; ++ipix)
1168 Pixel(ipix)->SetChargeBackup(0);
1169 // Correct each pixel
1171 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1173 probi1[ipix] = probMax;
1174 for (Int_t j = 0; j < npad; ++j)
1176 AliMUONPad* pad = cluster.Pad(j);
1177 //if (pad->Status() < 0) continue;
1178 if (pad->Status() != fgkZero) continue;
1180 Int_t indx1 = j*nPix;
1181 Int_t indx = indx1 + ipix;
1182 // Calculate expectation
1183 for (Int_t i = 0; i < nPix; ++i)
1185 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1186 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1188 if ( pad->IsSaturated() && sum1 > pad->Charge() )
1190 // correct for pad charge overflows
1191 probi1[ipix] -= coef[indx];
1193 //sum1 = pad->Charge();
1196 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1197 } // for (Int_t j=0;
1198 AliMUONPad* pixPtr = Pixel(ipix);
1199 if (probi1[ipix] > 1.e-6)
1201 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1202 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1204 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1205 } // for (Int_t ipix=0;
1207 for (Int_t i = 0; i < nPix; ++i) {
1208 AliMUONPad* pixPtr = Pixel(i);
1209 pixPtr->RevertCharge();
1210 qTot += pixPtr->Charge();
1213 // Can happen in clusters with large number of overflows - speeding up
1217 } // for (Int_t iter=0;
1221 //_____________________________________________________________________________
1222 void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
1224 /// Calculate position of the center-of-gravity around the maximum pixel
1226 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1227 Int_t i1 = -9, j1 = -9;
1228 fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1229 Int_t nx = fHistMlem->GetNbinsX();
1230 Int_t ny = fHistMlem->GetNbinsY();
1231 Double_t thresh = fHistMlem->GetMaximum()/10;
1232 Double_t x, y, cont, xq=0, yq=0, qq=0;
1234 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1235 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1236 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1237 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1238 cont = fHistMlem->GetCellContent(j,i);
1239 if (cont < thresh) continue;
1240 if (i != i1) {i1 = i; nsumy++;}
1241 if (j != j1) {j1 = j; nsumx++;}
1242 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1251 Int_t i2 = 0, j2 = 0;
1254 // one bin in Y - add one more (with the largest signal)
1255 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1256 if (i == iymax) continue;
1257 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1258 cont = fHistMlem->GetCellContent(j,i);
1261 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1262 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1271 if (i2 != i1) nsumy++;
1272 if (j2 != j1) nsumx++;
1274 } // if (nsumy == 1)
1277 // one bin in X - add one more (with the largest signal)
1279 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1280 if (j == ixmax) continue;
1281 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1282 cont = fHistMlem->GetCellContent(j,i);
1285 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1286 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1295 if (i2 != i1) nsumy++;
1296 if (j2 != j1) nsumx++;
1298 } // if (nsumx == 1)
1300 xyc[0] = xq/qq; xyc[1] = yq/qq;
1301 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1304 //_____________________________________________________________________________
1305 Int_t AliMUONClusterFinderMLEM::FindNearest(const AliMUONPad *pixPtr0)
1307 /// Find the pixel nearest to the given one
1308 /// (algorithm may be not very efficient)
1310 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1311 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1312 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1315 for (Int_t i = 0; i < nPix; ++i) {
1316 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1317 if (pixPtr == pixPtr0 || pixPtr->Charge() < fLowestPixelCharge) continue;
1318 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1319 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1320 r = dx *dx + dy * dy;
1321 if (r < rmin) { rmin = r; imin = i; }
1326 //_____________________________________________________________________________
1328 AliMUONClusterFinderMLEM::Paint(Option_t*)
1330 /// Paint cluster and pixels
1332 AliMpArea area(fPreCluster->Area());
1334 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1336 gVirtualX->SetFillStyle(1001);
1337 gVirtualX->SetFillColor(3);
1338 gVirtualX->SetLineColor(3);
1342 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1344 AliMUONPad* pixel = Pixel(i);
1346 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1347 pixel->Coord(1)-pixel->Size(1)*s,
1348 pixel->Coord(0)+pixel->Size(0)*s,
1349 pixel->Coord(1)+pixel->Size(1)*s);
1351 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1353 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1354 // pixel->Coord(1) + sign*pixel->Size(1),
1355 // pixel->Coord(0) + pixel->Size(0),
1356 // pixel->Coord(1) - sign*pixel->Size(1)
1362 gVirtualX->SetFillStyle(0);
1364 fPreCluster->Paint();
1366 gVirtualX->SetLineColor(1);
1367 gVirtualX->SetLineWidth(2);
1368 gVirtualX->SetFillStyle(0);
1369 gVirtualX->SetTextColor(1);
1370 gVirtualX->SetTextAlign(22);
1372 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1374 AliMUONPad* pixel = Pixel(i);
1375 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1376 pixel->Coord(1)-pixel->Size(1),
1377 pixel->Coord(0)+pixel->Size(0),
1378 pixel->Coord(1)+pixel->Size(1));
1379 gVirtualX->SetTextSize(pixel->Size(0)*60);
1381 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1385 //_____________________________________________________________________________
1386 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1388 /// Find local maxima in pixel space for large preclusters in order to
1389 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1390 /// or to find additional fitting seeds if clusters were not completely resolved
1392 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1394 Double_t xylim[4] = {999, 999, 999, 999};
1396 Int_t nPix = pixArray->GetEntriesFast();
1398 if ( nPix <= 0 ) return 0;
1400 AliMUONPad *pixPtr = 0;
1401 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1402 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1403 for (Int_t i = 0; i < 4; ++i)
1404 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1406 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
1408 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1409 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1410 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1411 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1412 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1413 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1414 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1416 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1418 Int_t nMax = 0, indx, nxy = ny * nx;
1419 Int_t *isLocalMax = new Int_t[nxy];
1420 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
1422 for (Int_t i = 1; i <= ny; ++i) {
1424 for (Int_t j = 1; j <= nx; ++j) {
1425 if (fHistAnode->GetCellContent(j,i) < fLowestPixelCharge) continue;
1426 //if (isLocalMax[indx+j-1] < 0) continue;
1427 if (isLocalMax[indx+j-1] != 0) continue;
1428 FlagLocalMax(fHistAnode, i, j, isLocalMax);
1432 for (Int_t i = 1; i <= ny; ++i) {
1434 for (Int_t j = 1; j <= nx; ++j) {
1435 if (isLocalMax[indx+j-1] > 0) {
1436 localMax[nMax] = indx + j - 1;
1437 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
1438 ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
1439 if (nMax > 99) break;
1443 AliError(" Too many local maxima !!!");
1447 if (fDebug) cout << " Local max: " << nMax << endl;
1448 delete [] isLocalMax;
1452 //_____________________________________________________________________________
1453 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1455 /// Flag pixels (whether or not local maxima)
1457 Int_t nx = hist->GetNbinsX();
1458 Int_t ny = hist->GetNbinsY();
1459 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1460 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1462 Int_t ie = i + 2, je = j + 2;
1463 for (Int_t i1 = i-1; i1 < ie; ++i1) {
1464 if (i1 < 1 || i1 > ny) continue;
1465 indx1 = (i1 - 1) * nx;
1466 for (Int_t j1 = j-1; j1 < je; ++j1) {
1467 if (j1 < 1 || j1 > nx) continue;
1468 if (i == i1 && j == j1) continue;
1469 indx2 = indx1 + j1 - 1;
1470 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1471 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1472 else if (cont > cont1) isLocalMax[indx2] = -1;
1473 else { // the same charge
1474 isLocalMax[indx] = 1;
1475 if (isLocalMax[indx2] == 0) {
1476 FlagLocalMax(hist, i1, j1, isLocalMax);
1477 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1478 else isLocalMax[indx2] = -1;
1483 isLocalMax[indx] = 1; // local maximum
1486 //_____________________________________________________________________________
1487 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1488 const Int_t *localMax, Int_t iMax)
1490 /// Find pixel cluster around local maximum \a iMax and pick up pads
1491 /// overlapping with it
1494 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1496 hist->Draw("lego1Fb"); // debug
1501 Int_t nx = fHistAnode->GetNbinsX();
1502 Int_t ny = fHistAnode->GetNbinsY();
1503 Int_t ic = localMax[iMax] / nx + 1;
1504 Int_t jc = localMax[iMax] % nx + 1;
1505 Int_t nxy = ny * nx;
1506 Bool_t *used = new Bool_t[nxy];
1507 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
1509 // Drop all pixels from the array - pick up only the ones from the cluster
1510 fPixArray->Delete();
1512 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1513 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1514 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1515 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1516 Double_t cont = fHistAnode->GetCellContent(jc,ic);
1517 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1518 used[(ic-1)*nx+jc-1] = kTRUE;
1519 AddBinSimple(fHistAnode, ic, jc);
1520 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1522 Int_t nPix = fPixArray->GetEntriesFast();
1523 Int_t npad = cluster.Multiplicity();
1525 for (Int_t i = 0; i < nPix; ++i)
1527 AliMUONPad* pixPtr = Pixel(i);
1528 pixPtr->SetSize(0,wx);
1529 pixPtr->SetSize(1,wy);
1532 // Pick up pads which overlap with found pixels
1533 for (Int_t i = 0; i < npad; ++i)
1535 //cluster.Pad(i)->SetStatus(-1);
1536 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1539 for (Int_t i = 0; i < nPix; ++i)
1541 AliMUONPad* pixPtr = Pixel(i);
1542 for (Int_t j = 0; j < npad; ++j)
1544 AliMUONPad* pad = cluster.Pad(j);
1545 //if (pad->Status() == 0) continue;
1546 if (pad->Status() == fgkZero) continue;
1547 if ( Overlap(*pad,*pixPtr) )
1549 //pad->SetStatus(0);
1550 pad->SetStatus(fgkZero);
1551 if (fDebug) { cout << j << " "; pad->Print("full"); }
1559 //_____________________________________________________________________________
1561 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
1563 /// Add adjacent bins (+-1 in X and Y) to the cluster
1565 Int_t nx = hist->GetNbinsX();
1566 Int_t ny = hist->GetNbinsY();
1567 Double_t cont1, cont = hist->GetCellContent(jc,ic);
1568 AliMUONPad *pixPtr = 0;
1570 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1571 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1572 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1573 cont1 = hist->GetCellContent(j,i);
1574 if (cont1 > cont) continue;
1575 if (cont1 < fLowestPixelCharge) continue;
1576 pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j),
1577 hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1578 fPixArray->Add(pixPtr);
1583 //_____________________________________________________________________________
1584 AliMUONClusterFinderMLEM&
1585 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1587 /// Protected assignement operator
1589 if (this == &rhs) return *this;
1591 AliFatal("Not implemented.");
1596 //_____________________________________________________________________________
1597 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1599 /// Add virtual pad (with small charge) to improve fit for clusters
1600 /// with number of pads == 2 per direction
1602 // Find out non-bending and bending planes
1603 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1605 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1606 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1607 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1612 Bool_t same = kFALSE;
1613 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1616 Bool_t check[2] = {kFALSE, kFALSE};
1618 nxy[0] = nxy[1] = 0;
1619 for (Int_t inb = 0; inb < 2; ++inb) {
1620 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1621 if (inb == 0 && AliMp::PairFirst(cn) == 2) check[inb] = kTRUE; // check non-bending plane
1622 else if (inb == 1 && AliMp::PairSecond(cn) == 2) check[inb] = kTRUE; // check bending plane
1624 nxy[0] = TMath::Max (nxy[0], AliMp::PairFirst(cn));
1625 nxy[1] = TMath::Max (nxy[1], AliMp::PairSecond(cn));
1626 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1627 else if (inb == 1 && AliMp::PairSecond(cn) < 2) nonb[inb] = !nonb[inb];
1631 if (nxy[0] > 2) check[0] = kFALSE;
1632 if (nxy[1] > 2) check[1] = kFALSE;
1634 if (!check[0] && !check[1]) return;
1636 for (Int_t inb = 0; inb < 2; ++inb) {
1637 if (!check[inb]) continue;
1639 // Find pads with maximum and next to maximum charges
1640 Int_t maxPads[2] = {-1, -1};
1641 Double_t amax[2] = {0};
1642 Int_t mult = cluster.Multiplicity();
1643 for (Int_t j = 0; j < mult; ++j) {
1644 AliMUONPad *pad = cluster.Pad(j);
1645 if (pad->Cathode() != nonb[inb]) continue;
1646 for (Int_t j2 = 0; j2 < 2; ++j2) {
1647 if (pad->Charge() > amax[j2]) {
1648 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1649 amax[j2] = pad->Charge();
1656 // Find min and max dimensions of the cluster
1657 Double_t limits[2] = {9999, -9999};
1658 for (Int_t j = 0; j < mult; ++j) {
1659 AliMUONPad *pad = cluster.Pad(j);
1660 if (pad->Cathode() != nonb[inb]) continue;
1661 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1662 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1665 // Loop over max and next to max pads
1666 Bool_t add = kFALSE;
1668 for (Int_t j = 0; j < 2; ++j) {
1670 if (maxPads[j] < 0) continue;
1672 if (amax[1] / amax[0] < 0.5) break;
1674 // Check if pad at the cluster limit
1675 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1677 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1678 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1680 //cout << " *** Pad not at the cluster limit: " << j << endl;
1683 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1685 // Add pad (if it exists)
1687 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1688 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1689 //AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1690 AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos.X(), pos.Y(),kFALSE);
1691 if (!mppad.IsValid()) continue; // non-existing pad
1692 AliMUONPad muonPad(fDetElemId, nonb[inb], mppad.GetIx(), mppad.GetIy(),
1693 mppad.GetPositionX(), mppad.GetPositionY(),
1694 mppad.GetDimensionX(), mppad.GetDimensionY(), 0);
1695 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, fLowestPadCharge));
1696 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1697 else muonPad.SetCharge(TMath::Min (amax[j]/15, fLowestPadCharge));
1698 if (muonPad.Charge() < 2.0*fLowestPixelCharge) muonPad.SetCharge(2.0*fLowestPixelCharge);
1699 muonPad.SetReal(kFALSE);
1700 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1701 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1702 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1703 cluster.AddPad(muonPad); // add pad to the cluster
1710 //_____________________________________________________________________________
1711 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1712 Int_t &nInX, Int_t &nInY) const
1714 /// Find number of pads in X and Y-directions (excluding virtual ones and
1717 //Int_t statusToTest = 1;
1718 Int_t statusToTest = fgkUseForFit;
1720 //if ( nInX < 0 ) statusToTest = 0;
1721 if ( nInX < 0 ) statusToTest = fgkZero;
1723 Bool_t mustMatch(kTRUE);
1725 Long_t cn = cluster.NofPads(statusToTest,mustMatch);
1727 nInX = AliMp::PairFirst(cn);
1728 nInY = AliMp::PairSecond(cn);
1731 //_____________________________________________________________________________
1732 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1734 /// Remove pixel at index i
1735 AliMUONPad* pixPtr = Pixel(i);
1736 fPixArray->RemoveAt(i);
1740 //_____________________________________________________________________________
1741 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1743 /// Process simple cluster (small number of pads) without EM-procedure
1745 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1746 Double_t parOk[3] = {0.};
1747 TObjArray *clusters[1];
1748 clusters[0] = fPixArray;
1750 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1752 Int_t mult = cluster.Multiplicity();
1753 for (Int_t i = 0; i < mult; ++i)
1755 AliMUONPad* pad = cluster.Pad(i);
1757 if ( pad->IsSaturated())
1766 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1768 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
1771 //_____________________________________________________________________________
1773 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1775 /// Returns pixel at index i
1776 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1779 //_____________________________________________________________________________
1781 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1784 TString swhat(what);
1786 if ( swhat.Contains("precluster") )
1788 if ( fPreCluster) fPreCluster->Print();
1792 //_____________________________________________________________________________
1794 AliMUONClusterFinderMLEM::SetChargeHints(Double_t lowestPadCharge, Double_t lowestClusterCharge)
1796 /// Set some thresholds we use later on in the algorithm
1797 fLowestPadCharge=lowestPadCharge;
1798 fLowestClusterCharge=lowestClusterCharge;
1799 fLowestPixelCharge=fLowestPadCharge/12.0;