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"
53 ClassImp(AliMUONClusterFinderMLEM)
56 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
57 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
58 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
60 // Status flags for pads
61 const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
62 const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
63 const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
64 const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
65 const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
66 const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
68 //_____________________________________________________________________________
69 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
70 : AliMUONVClusterFinder(),
71 fPreClusterFinder(clusterFinder),
79 fPixArray(new TObjArray(20)),
85 fLowestPixelCharge(0),
87 fLowestClusterCharge(0)
91 fkSegmentation[1] = fkSegmentation[0] = 0x0;
93 if (fPlot) fDebug = 1;
96 //_____________________________________________________________________________
97 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
100 delete fPixArray; fPixArray = 0;
102 delete fPreClusterFinder;
104 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
105 fNClusters,fNAddVirtualPads));
108 //_____________________________________________________________________________
110 AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
112 const AliMpArea& area,
113 const AliMpVSegmentation* seg[2])
115 /// Prepare for clustering
116 // AliCodeTimerAuto("",0)
118 for ( Int_t i = 0; i < 2; ++i )
120 fkSegmentation[i] = seg[i];
123 // Find out the DetElemId
124 fDetElemId = detElemId;
127 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,
131 fLowestClusterCharge);
132 fSplitter->SetDebug(fDebug);
134 // find out current event number, and reset the cluster number
135 AliRunLoader *runLoader = AliRunLoader::Instance();
136 fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
138 fClusterList.Delete();
141 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
143 if ( fPreClusterFinder->NeedSegmentation() )
145 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
149 return fPreClusterFinder->Prepare(detElemId,pads,area);
153 //_____________________________________________________________________________
155 AliMUONClusterFinderMLEM::NextCluster()
157 /// Return next cluster
158 // AliCodeTimerAuto("",0)
160 // if the list of clusters is not void, pick one from there
163 // do we have clusters in our list ?
164 if ( fClusterNumber < fClusterList.GetLast() )
166 o = fClusterList.At(++fClusterNumber);
169 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
171 //FIXME : at this point, must check whether we've used all the digits
172 //from precluster : if not, let the preclustering know about those unused
173 //digits, so it can reuse them
175 // if the cluster list is exhausted, we need to go to the next
176 // pre-cluster and treat it
178 fPreCluster = fPreClusterFinder->NextCluster();
181 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
182 fClusterNumber = -1; //AZ
192 // WorkOnPreCluster may have used only part of the pads, so we check that
193 // now, and let the unused pads be reused by the preclustering...
195 Int_t mult = fPreCluster->Multiplicity();
196 for ( Int_t i = 0; i < mult; ++i )
198 AliMUONPad* pad = fPreCluster->Pad(i);
199 if ( !pad->IsUsed() )
201 fPreClusterFinder->UsePad(*pad);
205 return NextCluster();
208 //_____________________________________________________________________________
210 AliMUONClusterFinderMLEM::WorkOnPreCluster()
212 /// Starting from a precluster, builds a pixel array, and then
213 /// extract clusters from this array
215 // AliCodeTimerAuto("",0)
218 cout << " *** Event # " << fEventNumber
219 << " det. elem.: " << fDetElemId << endl;
220 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
221 AliMUONPad* pad = fPreCluster->Pad(j);
222 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
223 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
224 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
228 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
229 if (!cluster) return kFALSE;
231 BuildPixArray(*cluster);
233 if ( fPixArray->GetLast() < 0 )
235 AliDebug(1,"No pixel for the above cluster");
240 // Use MLEM for cluster finder
241 Int_t nMax = 1, localMax[100], maxPos[100];
242 Double_t maxVal[100];
244 Int_t iSimple = 0, nInX = -1, nInY;
246 PadsInXandY(*cluster,nInX, nInY);
248 if (nInX < 4 && nInY < 4)
250 iSimple = 1; // simple cluster
254 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
256 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
257 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
261 for (Int_t i = 0; i < nMax; ++i)
265 FindCluster(*cluster,localMax, maxPos[i]);
268 MainLoop(*cluster,iSimple);
272 Int_t mult = cluster->Multiplicity();
273 for (Int_t j = 0; j < mult; ++j)
275 AliMUONPad* pad = cluster->Pad(j);
276 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
277 if ( pad->Status() != fgkOver ) continue; // pad was not used
279 pad->SetStatus(fgkZero);
280 pad->RevertCharge(); // use backup charge value
283 } // for (Int_t i=0; i<nMax;
286 fHistMlem = fHistAnode = 0x0;
291 //_____________________________________________________________________________
293 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
295 /// Check if the pad and the pixel overlaps
297 // make a fake pad from the pixel
298 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
299 pixel.Coord(0),pixel.Coord(1),
300 pixel.Size(0),pixel.Size(1),0);
302 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
305 //_____________________________________________________________________________
307 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
309 /// Check precluster in order to attempt to simplify it (mostly for
310 /// two-cathode preclusters)
312 AliCodeTimerAuto("",0)
314 // Disregard small clusters (leftovers from splitting or noise)
315 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
316 origCluster.Charge(0)+origCluster.Charge(1) < fLowestClusterCharge )
321 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
323 AliDebug(2,"Start of CheckPreCluster=");
324 //StdoutToAliDebug(2,cluster->Print("full"));
326 AliMUONCluster* rv(0x0);
328 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
330 rv = CheckPreclusterTwoCathodes(cluster);
339 //_____________________________________________________________________________
341 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
343 /// Check two-cathode cluster
345 Int_t npad = cluster->Multiplicity();
346 Int_t* flags = new Int_t[npad];
347 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
349 // Check pad overlaps
350 for ( Int_t i = 0; i < npad; ++i)
352 AliMUONPad* padi = cluster->Pad(i);
353 if ( padi->Cathode() != 0 ) continue;
354 for (Int_t j = i+1; j < npad; ++j)
356 AliMUONPad* padj = cluster->Pad(j);
357 if ( padj->Cathode() != 1 ) continue;
358 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
359 flags[i] = flags[j] = 1; // mark overlapped pads
363 // Check if all pads overlap
365 for (Int_t i = 0; i < npad; ++i)
367 if (!flags[i]) ++nFlags;
372 // not all pads overlap.
373 if (fDebug) cout << " nFlags: " << nFlags << endl;
374 TObjArray toBeRemoved;
375 for (Int_t i = 0; i < npad; ++i)
377 AliMUONPad* pad = cluster->Pad(i);
378 if (flags[i]) continue;
379 Int_t cath = pad->Cathode();
380 Int_t cath1 = TMath::Even(cath);
381 // Check for edge effect (missing pads on the _other_ cathode)
383 fkSegmentation[cath1]->PadByPosition(pad->Position().X(),
384 pad->Position().Y(),kFALSE);
385 if (!mpPad.IsValid()) continue;
386 if (nFlags == 1 && pad->Charge() < fLowestPadCharge) continue;
387 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
388 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
389 toBeRemoved.AddLast(pad);
390 fPreCluster->Pad(i)->Release();
392 Int_t nRemove = toBeRemoved.GetEntriesFast();
393 for ( Int_t i = 0; i < nRemove; ++i )
395 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
399 // Check correlations of cathode charges
400 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
403 Int_t cathode = cluster->MaxRawChargeCathode();
409 // get min and max pad charges on the cathode opposite to the
410 // max pad (given by MaxRawChargeCathode())
412 Int_t mult = cluster->Multiplicity();
413 for ( Int_t i = 0; i < mult; ++i )
415 AliMUONPad* pad = cluster->Pad(i);
416 if ( pad->Cathode() != cathode || !pad->IsReal() )
418 // only consider pads in the opposite cathode, and
419 // only consider real pads (i.e. exclude the virtual ones)
422 if ( pad->Charge() < cmin )
424 cmin = pad->Charge();
431 else if ( pad->Charge() > cmax )
433 cmax = pad->Charge();
437 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
438 imin,imax,cmin,cmax));
440 // arrange pads according to their distance to the max, normalized
442 Double_t* dist = new Double_t[mult];
447 AliMUONPad* padmax = cluster->Pad(imax);
449 for ( Int_t i = 0; i < mult; ++i )
452 if ( i == imax) continue;
453 AliMUONPad* pad = cluster->Pad(i);
454 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
455 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
456 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
457 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
460 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
466 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
467 Double_t xmax(-1), distPrev(999);
468 TObjArray toBeRemoved;
470 for ( Int_t i = 0; i < mult; ++i )
472 Int_t indx = flags[i];
473 AliMUONPad* pad = cluster->Pad(indx);
474 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
475 if ( dist[indx] > dmin )
477 // farther than the minimum pad
478 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
479 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
482 if (dx >= 0 && dy >= 0) continue;
483 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
484 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
486 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
487 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
490 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
492 cmax = TMath::Max(pad->Charge(),cmax);
496 cmax = pad->Charge();
499 distPrev = dist[indx];
500 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
501 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
504 toBeRemoved.AddLast(pad);
505 fPreCluster->Pad(indx)->Release();
508 Int_t nRemove = toBeRemoved.GetEntriesFast();
509 for ( Int_t i = 0; i < nRemove; ++i )
511 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
518 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
519 //StdoutToAliDebug(2,cluster->Print("full"));
524 //_____________________________________________________________________________
526 AliMUONClusterFinderMLEM::CheckOverlaps()
528 /// For debug only : check if some pixels overlap...
530 Int_t nPix = fPixArray->GetLast()+1;
533 for ( Int_t i = 0; i < nPix; ++i )
535 AliMUONPad* pixelI = Pixel(i);
536 AliMUONPad pi(dummy,dummy,dummy,dummy,
537 pixelI->Coord(0),pixelI->Coord(1),
538 pixelI->Size(0),pixelI->Size(1),0.0);
540 for ( Int_t j = i+1; j < nPix; ++j )
542 AliMUONPad* pixelJ = Pixel(j);
543 AliMUONPad pj(dummy,dummy,dummy,dummy,
544 pixelJ->Coord(0),pixelJ->Coord(1),
545 pixelJ->Size(0),pixelJ->Size(1),0.0);
548 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
550 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
552 StdoutToAliInfo(pixelI->Print();
553 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
555 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
556 cout << " Area surface = " << area.GetDimensionX()*area.GetDimensionY()*4 << endl;
557 cout << "-------" << endl;
565 //_____________________________________________________________________________
566 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
568 /// Build pixel array for MLEM method
570 Int_t npad = cluster.Multiplicity();
573 AliWarning("Got no pad at all ?!");
577 BuildPixArrayOneCathode(cluster);
579 Int_t nPix = fPixArray->GetLast()+1;
581 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
585 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
587 // Too many pixels - sort and remove pixels with the lowest signal
589 for ( Int_t i = npad; i < nPix; ++i )
593 fPixArray->Compress();
594 } // if (nPix > npad)
596 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
597 // fPixArray->Print(););
598 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
601 //_____________________________________________________________________________
602 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
604 /// Build the pixel array
606 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
608 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
609 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2]={99999,99999};
610 Int_t found[2] = {0,0}, mult = cluster.Multiplicity();
612 for ( Int_t i = 0; i < mult; ++i) {
613 AliMUONPad* pad = cluster.Pad(i);
614 for (Int_t j = 0; j < 2; ++j) {
615 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
616 xy0[j] = pad->Coord(j);
620 if (found[0] && found[1]) break;
623 Double_t min[2], max[2];
624 Int_t cath0 = 0, cath1 = 1;
625 if (cluster.Multiplicity(0) == 0) cath0 = 1;
626 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
629 Double_t leftDownX, leftDownY;
630 cluster.Area(cath0).LeftDownCorner(leftDownX, leftDownY);
631 Double_t rightUpX, rightUpY;
632 cluster.Area(cath0).RightUpCorner(rightUpX, rightUpY);
637 if (cath1 != cath0) {
638 cluster.Area(cath1).LeftDownCorner(leftDownX, leftDownY);
639 cluster.Area(cath1).RightUpCorner(rightUpX, rightUpY);
640 min[0] = TMath::Max (min[0], leftDownX);
641 min[1] = TMath::Max (min[1], leftDownY);
642 max[0] = TMath::Min (max[0], rightUpX);
643 max[1] = TMath::Min (max[1], rightUpY);
647 //width[0] /= 2; width[1] /= 2; // just for check
648 Int_t nbins[2]={0,0};
649 for (Int_t i = 0; i < 2; ++i) {
650 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
651 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
652 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
653 + TMath::Sign(0.5,dist)) * width[i] * 2;
654 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
655 if (nbins[i] == 0) ++nbins[i];
656 max[i] = min[i] + nbins[i] * width[i] * 2;
657 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
661 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
662 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
663 TAxis *xaxis = hist1->GetXaxis();
664 TAxis *yaxis = hist1->GetYaxis();
667 for ( Int_t i = 0; i < mult; ++i) {
668 AliMUONPad* pad = cluster.Pad(i);
669 Int_t ix0 = xaxis->FindBin(pad->X());
670 Int_t iy0 = yaxis->FindBin(pad->Y());
671 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
675 for (Int_t i = 1; i <= nbins[0]; ++i) {
676 Double_t x = xaxis->GetBinCenter(i);
677 for (Int_t j = 1; j <= nbins[1]; ++j) {
678 if (hist2->GetBinContent(hist2->GetBin(i,j)) < 0.1) continue;
679 //if (hist2->GetBinContent(hist2->GetBin(i,j)) < 1.1 && cluster.Multiplicity(0) &&
680 // cluster.Multiplicity(1)) continue;
681 if (cath0 != cath1) {
683 Double_t cont = hist2->GetBinContent(hist2->GetBin(i,j));
684 if (cont < 999.) continue;
685 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
687 Double_t y = yaxis->GetBinCenter(j);
688 Double_t charge = hist1->GetBinContent(hist1->GetBin(i,j));
689 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
690 fPixArray->Add(pixPtr);
694 if (fPixArray->GetEntriesFast() == 1) {
695 // Split pixel into 2
696 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
697 pixPtr->SetSize(0,width[0]/2.);
698 pixPtr->Shift(0,-width[0]/4.);
699 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
700 fPixArray->Add(pixPtr);
703 //fPixArray->Print();
708 //_____________________________________________________________________________
709 void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
710 TH2D *hist1, TH2D *hist2)
712 /// "Span" pad over histogram in the direction idir
714 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
715 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
716 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
718 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
720 for (Int_t i = 0; i < nbinPad; ++i) {
721 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
722 if (ixy > nbins) break;
723 Double_t lowEdge = axis->GetBinLowEdge(ixy);
724 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
725 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
728 Double_t cont = pad->Charge();
729 if (hist2->GetBinContent(hist2->GetBin(ix0, ixy)) > 0.1)
730 cont = TMath::Min (hist1->GetBinContent(hist1->GetBin(ix0, ixy)), cont);
731 hist1->SetBinContent(hist1->GetBin(ix0, ixy), cont);
732 //hist2->SetBinContent(hist2->GetBin(ix0, ixy), hist2->GetBinContent(hist2->GetBin(ix0, ixy))+1);
733 hist2->SetBinContent(hist2->GetBin(ix0, ixy), hist2->GetBinContent(hist2->GetBin(ix0, ixy))+amask);
737 for (Int_t i = -1; i > -nbinPad; --i) {
738 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
740 Double_t upEdge = axis->GetBinUpEdge(ixy);
741 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
742 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
745 Double_t cont = pad->Charge();
746 if (hist2->GetBinContent(hist2->GetBin(ix0, ixy)) > 0.1)
747 cont = TMath::Min (hist1->GetBinContent(hist1->GetBin(ix0, ixy)), cont);
748 hist1->SetBinContent(hist1->GetBin(ix0, ixy), cont);
749 //hist2->SetBinContent(hist1->GetBin(ix0, ixy), hist2->GetBinContent(hist2->GetBin(ix0, ixy))+1);
750 hist2->SetBinContent(hist2->GetBin(ix0, ixy), hist2->GetBinContent(hist2->GetBin(ix0, ixy))+amask);
755 //_____________________________________________________________________________
757 AliMUONClusterFinderMLEM::Plot(const char* /*basename*/)
759 /// Make a plot and save it as png
762 // if (!fPlot) return;
764 // TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
769 // c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
770 // fDetElemId,fClusterNumber));
773 //_____________________________________________________________________________
775 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
779 /// Compute coefficients needed for MLEM algorithm
781 Int_t npadTot = cluster.Multiplicity();
782 Int_t nPix = fPixArray->GetLast()+1;
784 //memset(probi,0,nPix*sizeof(Double_t));
785 for (Int_t j = 0; j < npadTot*nPix; ++j) coef[j] = 0.;
786 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
788 Int_t mult = cluster.Multiplicity();
789 for ( Int_t j = 0; j < mult; ++j )
791 AliMUONPad* pad = cluster.Pad(j);
794 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
796 Int_t indx1 = indx + ipix;
797 //if (pad->Status() < 0)
798 if (pad->Status() != fgkZero)
803 AliMUONPad* pixPtr = Pixel(ipix);
804 // coef is the charge (given by Mathieson integral) on pad, assuming
805 // the Mathieson is center at pixel.
806 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
807 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
808 // pad->Ix(),pad->Iy(),
809 // pad->X(),pad->Y(),
810 // pad->DX(),pad->DY(),
811 // pixPtr->Coord(0),pixPtr->Coord(1),
812 // pixPtr->Size(0),pixPtr->Size(1),
815 probi[ipix] += coef[indx1];
820 //_____________________________________________________________________________
821 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
823 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
825 // AliCodeTimerAuto("",0)
827 Int_t nPix = fPixArray->GetLast()+1;
829 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
830 //StdoutToAliDebug(2,cluster.Print("full"););
834 AliDebug(1,"No pixels, why am I here then ?");
837 AddVirtualPad(cluster); // add virtual pads if necessary
839 Int_t npadTot = cluster.Multiplicity();
841 for (Int_t i = 0; i < npadTot; ++i)
843 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
844 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
848 Double_t* probi(0x0);
849 Int_t lc(0); // loop counter
851 //Plot("mlem.start");
852 AliMUONPad* pixPtr = Pixel(0);
853 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
861 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
862 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
863 //StdoutToAliDebug(2,fPixArray->Print("full"));
865 coef = new Double_t [npadTot*nPix];
866 probi = new Double_t [nPix];
868 // Calculate coefficients and pixel visibilities
869 ComputeCoefficients(cluster,coef,probi);
871 for (Int_t ipix = 0; ipix < nPix; ++ipix)
873 if (probi[ipix] < 0.01)
875 AliMUONPad* pixel = Pixel(ipix);
876 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
877 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
878 pixel->SetCharge(0); // "invisible" pixel
883 Mlem(cluster,coef, probi, 15);
885 // Find histogram limits for the 1'st pass only - for others computed below
887 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
889 pixPtr = Pixel(ipix);
890 for ( Int_t i = 0; i < 2; ++i )
893 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
894 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
897 } else pixPtr = Pixel(0);
899 for (Int_t i = 0; i < 4; i++)
901 xylim[i] -= pixPtr->Size(i/2);
904 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
905 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
907 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
908 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
909 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
910 xylim[0],-xylim[1],xylim[2],-xylim[3]
913 AliDebug(2,Form("LowestPadCharge=%e",fLowestPadCharge));
917 fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
919 for (Int_t ipix = 0; ipix < nPix; ++ipix)
921 AliMUONPad* pixPtr2 = Pixel(ipix);
922 fHistMlem->Fill(pixPtr2->Coord(0),pixPtr2->Coord(1),pixPtr2->Charge());
925 // Check if the total charge of pixels is too low
927 for ( Int_t i = 0; i < nPix; ++i)
929 qTot += Pixel(i)->Charge();
932 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < fLowestClusterCharge ) )
934 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
938 for ( Int_t i = 0; i < npadTot; ++i)
940 AliMUONPad* pad = cluster.Pad(i);
941 //if ( pad->Status() == 0) pad->SetStatus(-1);
942 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
949 // Simple cluster - skip further passes thru EM-procedure
957 // Calculate position of the center-of-gravity around the maximum pixel
961 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
962 pixPtr->Size(0) > pixPtr->Size(1)) break;
964 // Sort pixels according to the charge
965 MaskPeaks(1); // mask local maxima
967 MaskPeaks(0); // unmask local maxima
968 Double_t pixMin = 0.01*Pixel(0)->Charge();
969 pixMin = TMath::Min(pixMin,100*fLowestPixelCharge);
971 // Decrease pixel size and shift pixels to make them centered at
973 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
976 Double_t shift[2] = { 0.0, 0.0 };
977 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
980 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
982 AliMUONPad* pixPtr2 = Pixel(ipix);
983 if ( nPix >= npadOK // too many pixels already
985 ((pixPtr2->Charge() < pixMin) && (pixPtr2->Status() != fgkMustKeep)) // too low charge
991 for (Int_t i = 0; i < 2; ++i)
995 pixPtr2->SetCharge(fLowestPadCharge);
996 pixPtr2->SetSize(indx, pixPtr2->Size(indx)/2);
997 width = -pixPtr2->Size(indx);
998 pixPtr2->Shift(indx, width);
999 // Shift pixel position
1003 for (Int_t j = 0; j < 2; ++j)
1005 shift[j] = pixPtr2->Coord(j) - xyCOG[j];
1006 shift[j] -= ((Int_t)(shift[j]/pixPtr2->Size(j)/2))*pixPtr2->Size(j)*2;
1009 pixPtr2->Shift(0, -shift[0]);
1010 pixPtr2->Shift(1, -shift[1]);
1013 else if (nPix < npadOK)
1015 pixPtr2 = new AliMUONPad(*pixPtr2);
1016 pixPtr2->Shift(indx, -2*width);
1017 pixPtr2->SetStatus(fgkZero);
1018 fPixArray->Add(pixPtr2);
1021 else continue; // skip adjustment of histo limits
1022 for (Int_t j = 0; j < 4; ++j)
1024 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr2->Coord(j/2));
1026 } // for (Int_t i=0; i<2;
1027 } // for (Int_t ipix=0;
1029 fPixArray->Compress();
1031 AliDebug(2,Form("After shift:"));
1032 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1033 //Plot(Form("mlem.lc%d",lc+1));
1035 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1038 xylim[2],xylim[3]));
1042 AliMUONPad* pixPtr2 = Pixel(0);
1043 // add pixels if the maximum is at the limit of pixel area:
1044 // start from Y-direction
1046 for (Int_t i = 3; i > -1; --i)
1048 if (nPix < npadOK &&
1049 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr2->Size(i/2))
1051 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1052 AliMUONPad* p = new AliMUONPad(*pixPtr2);
1053 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr2->Size(i/2));
1054 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1055 j = TMath::Even (i/2);
1056 p->SetCoord(j, xyCOG[j]);
1057 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1058 //StdoutToAliDebug(2,cout << " ---- ";
1059 // p->Print("corners"););
1065 nPix = fPixArray->GetEntriesFast();
1070 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1071 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1073 // remove pixels with low signal or low visibility
1074 // Cuts are empirical !!!
1075 Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,2.0*fLowestPixelCharge);
1076 thresh = TMath::Min (thresh,100.0*fLowestPixelCharge);
1077 Double_t charge = 0;
1079 // Mark pixels which should be removed
1080 for (Int_t i = 0; i < nPix; ++i)
1082 AliMUONPad* pixPtr2 = Pixel(i);
1083 charge = pixPtr2->Charge();
1084 if (charge < thresh)
1086 pixPtr2->SetCharge(-charge);
1090 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1092 for (Int_t i = 0; i < nPix; ++i)
1094 AliMUONPad* pixPtr2 = Pixel(i);
1095 charge = pixPtr2->Charge();
1096 if (charge > 0) continue;
1097 near = FindNearest(pixPtr2);
1098 pixPtr2->SetCharge(0);
1099 probi[i] = 0; // make it "invisible"
1100 AliMUONPad* pnear = Pixel(near);
1101 pnear->SetCharge(pnear->Charge() + (-charge));
1103 Mlem(cluster,coef,probi,2);
1105 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1106 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1107 //Plot("mlem.beforesplit");
1110 for (Int_t i = 0; i < nPix; ++i)
1112 AliMUONPad* pixPtr2 = Pixel(i);
1113 Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr2->Coord(0));
1114 Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr2->Coord(1));
1115 fHistMlem->SetBinContent(ix, iy, pixPtr2->Charge());
1118 // Try to split into clusters
1120 if (fHistMlem->GetSum() < 2.0*fLowestPixelCharge)
1126 fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
1131 fPixArray->Delete();
1136 //_____________________________________________________________________________
1137 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1139 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1140 /// - to avoid loosing low charge pixels after sorting)
1142 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1143 AliMUONPad* pix = Pixel(i);
1144 if (pix->Status() == fgkMustKeep) {
1145 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1146 else pix->SetCharge(pix->Charge()-10000.);
1151 //_____________________________________________________________________________
1152 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1153 const Double_t* coef, Double_t* probi,
1156 /// Use MLEM to find pixel charges
1158 Int_t nPix = fPixArray->GetEntriesFast();
1160 Int_t npad = cluster.Multiplicity();
1162 Double_t* probi1 = new Double_t[nPix];
1163 Double_t probMax = TMath::MaxElement(nPix,probi);
1165 for (Int_t iter = 0; iter < nIter; ++iter)
1168 for (Int_t ipix = 0; ipix < nPix; ++ipix)
1170 Pixel(ipix)->SetChargeBackup(0);
1171 // Correct each pixel
1173 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1175 probi1[ipix] = probMax;
1176 for (Int_t j = 0; j < npad; ++j)
1178 AliMUONPad* pad = cluster.Pad(j);
1179 //if (pad->Status() < 0) continue;
1180 if (pad->Status() != fgkZero) continue;
1182 Int_t indx1 = j*nPix;
1183 Int_t indx = indx1 + ipix;
1184 // Calculate expectation
1185 for (Int_t i = 0; i < nPix; ++i)
1187 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1188 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1190 if ( pad->IsSaturated() && sum1 > pad->Charge() )
1192 // correct for pad charge overflows
1193 probi1[ipix] -= coef[indx];
1195 //sum1 = pad->Charge();
1198 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1199 } // for (Int_t j=0;
1200 AliMUONPad* pixPtr = Pixel(ipix);
1201 if (probi1[ipix] > 1.e-6)
1203 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1204 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1206 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1207 } // for (Int_t ipix=0;
1209 for (Int_t i = 0; i < nPix; ++i) {
1210 AliMUONPad* pixPtr = Pixel(i);
1211 pixPtr->RevertCharge();
1212 qTot += pixPtr->Charge();
1215 // Can happen in clusters with large number of overflows - speeding up
1219 } // for (Int_t iter=0;
1223 //_____________________________________________________________________________
1224 void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
1226 /// Calculate position of the center-of-gravity around the maximum pixel
1228 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1229 Int_t i1 = -9, j1 = -9;
1230 fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1231 Int_t nx = fHistMlem->GetNbinsX();
1232 Int_t ny = fHistMlem->GetNbinsY();
1233 Double_t thresh = fHistMlem->GetMaximum()/10;
1234 Double_t x, y, cont, xq=0, yq=0, qq=0;
1236 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1237 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1238 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1239 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1240 cont = fHistMlem->GetBinContent(fHistMlem->GetBin(j,i));
1241 if (cont < thresh) continue;
1242 if (i != i1) {i1 = i; nsumy++;}
1243 if (j != j1) {j1 = j; nsumx++;}
1244 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1253 Int_t i2 = 0, j2 = 0;
1256 // one bin in Y - add one more (with the largest signal)
1257 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1258 if (i == iymax) continue;
1259 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1260 cont = fHistMlem->GetBinContent(fHistMlem->GetBin(j,i));
1263 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1264 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1273 if (i2 != i1) nsumy++;
1274 if (j2 != j1) nsumx++;
1276 } // if (nsumy == 1)
1279 // one bin in X - add one more (with the largest signal)
1281 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1282 if (j == ixmax) continue;
1283 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1284 cont = fHistMlem->GetBinContent(fHistMlem->GetBin(j,i));
1287 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1288 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1297 if (i2 != i1) nsumy++;
1298 if (j2 != j1) nsumx++;
1300 } // if (nsumx == 1)
1302 xyc[0] = xq/qq; xyc[1] = yq/qq;
1303 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1306 //_____________________________________________________________________________
1307 Int_t AliMUONClusterFinderMLEM::FindNearest(const AliMUONPad *pixPtr0)
1309 /// Find the pixel nearest to the given one
1310 /// (algorithm may be not very efficient)
1312 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1313 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1314 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1317 for (Int_t i = 0; i < nPix; ++i) {
1318 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1319 if (pixPtr == pixPtr0 || pixPtr->Charge() < fLowestPixelCharge) continue;
1320 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1321 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1322 r = dx *dx + dy * dy;
1323 if (r < rmin) { rmin = r; imin = i; }
1328 //_____________________________________________________________________________
1330 AliMUONClusterFinderMLEM::Paint(Option_t*)
1332 /// Paint cluster and pixels
1334 AliMpArea area(fPreCluster->Area());
1336 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1338 gVirtualX->SetFillStyle(1001);
1339 gVirtualX->SetFillColor(3);
1340 gVirtualX->SetLineColor(3);
1344 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1346 AliMUONPad* pixel = Pixel(i);
1348 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1349 pixel->Coord(1)-pixel->Size(1)*s,
1350 pixel->Coord(0)+pixel->Size(0)*s,
1351 pixel->Coord(1)+pixel->Size(1)*s);
1353 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1355 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1356 // pixel->Coord(1) + sign*pixel->Size(1),
1357 // pixel->Coord(0) + pixel->Size(0),
1358 // pixel->Coord(1) - sign*pixel->Size(1)
1364 gVirtualX->SetFillStyle(0);
1366 fPreCluster->Paint();
1368 gVirtualX->SetLineColor(1);
1369 gVirtualX->SetLineWidth(2);
1370 gVirtualX->SetFillStyle(0);
1371 gVirtualX->SetTextColor(1);
1372 gVirtualX->SetTextAlign(22);
1374 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1376 AliMUONPad* pixel = Pixel(i);
1377 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1378 pixel->Coord(1)-pixel->Size(1),
1379 pixel->Coord(0)+pixel->Size(0),
1380 pixel->Coord(1)+pixel->Size(1));
1381 gVirtualX->SetTextSize(pixel->Size(0)*60);
1383 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1387 //_____________________________________________________________________________
1388 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1390 /// Find local maxima in pixel space for large preclusters in order to
1391 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1392 /// or to find additional fitting seeds if clusters were not completely resolved
1394 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1396 Double_t xylim[4] = {999, 999, 999, 999};
1398 Int_t nPix = pixArray->GetEntriesFast();
1400 if ( nPix <= 0 ) return 0;
1402 AliMUONPad *pixPtr = 0;
1403 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1404 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1405 for (Int_t i = 0; i < 4; ++i)
1406 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1408 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
1410 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1411 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1412 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1413 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1414 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1415 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1416 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1418 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1420 Int_t nMax = 0, indx, nxy = ny * nx;
1421 Int_t *isLocalMax = new Int_t[nxy];
1422 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
1424 for (Int_t i = 1; i <= ny; ++i) {
1426 for (Int_t j = 1; j <= nx; ++j) {
1427 if (fHistAnode->GetBinContent(fHistAnode->GetBin(j,i)) < fLowestPixelCharge) continue;
1428 //if (isLocalMax[indx+j-1] < 0) continue;
1429 if (isLocalMax[indx+j-1] != 0) continue;
1430 FlagLocalMax(fHistAnode, i, j, isLocalMax);
1434 for (Int_t i = 1; i <= ny; ++i) {
1436 for (Int_t j = 1; j <= nx; ++j) {
1437 if (isLocalMax[indx+j-1] > 0) {
1438 localMax[nMax] = indx + j - 1;
1439 maxVal[nMax++] = fHistAnode->GetBinContent(fHistAnode->GetBin(j,i));
1440 ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
1441 if (nMax > 99) break;
1445 AliError(" Too many local maxima !!!");
1449 if (fDebug) cout << " Local max: " << nMax << endl;
1450 delete [] isLocalMax;
1454 //_____________________________________________________________________________
1455 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1457 /// Flag pixels (whether or not local maxima)
1459 Int_t nx = hist->GetNbinsX();
1460 Int_t ny = hist->GetNbinsY();
1461 Int_t cont = TMath::Nint (hist->GetBinContent(hist->GetBin(j,i)));
1462 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1464 Int_t ie = i + 2, je = j + 2;
1465 for (Int_t i1 = i-1; i1 < ie; ++i1) {
1466 if (i1 < 1 || i1 > ny) continue;
1467 indx1 = (i1 - 1) * nx;
1468 for (Int_t j1 = j-1; j1 < je; ++j1) {
1469 if (j1 < 1 || j1 > nx) continue;
1470 if (i == i1 && j == j1) continue;
1471 indx2 = indx1 + j1 - 1;
1472 cont1 = TMath::Nint (hist->GetBinContent(hist->GetBin(j1,i1)));
1473 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1474 else if (cont > cont1) isLocalMax[indx2] = -1;
1475 else { // the same charge
1476 isLocalMax[indx] = 1;
1477 if (isLocalMax[indx2] == 0) {
1478 FlagLocalMax(hist, i1, j1, isLocalMax);
1479 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1480 else isLocalMax[indx2] = -1;
1485 isLocalMax[indx] = 1; // local maximum
1488 //_____________________________________________________________________________
1489 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1490 const Int_t *localMax, Int_t iMax)
1492 /// Find pixel cluster around local maximum \a iMax and pick up pads
1493 /// overlapping with it
1496 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1498 hist->Draw("lego1Fb"); // debug
1503 Int_t nx = fHistAnode->GetNbinsX();
1504 Int_t ny = fHistAnode->GetNbinsY();
1505 Int_t ic = localMax[iMax] / nx + 1;
1506 Int_t jc = localMax[iMax] % nx + 1;
1507 Int_t nxy = ny * nx;
1508 Bool_t *used = new Bool_t[nxy];
1509 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
1511 // Drop all pixels from the array - pick up only the ones from the cluster
1512 fPixArray->Delete();
1514 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1515 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1516 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1517 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1518 Double_t cont = fHistAnode->GetBinContent( fHistAnode->GetBin(jc,ic));
1519 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1520 used[(ic-1)*nx+jc-1] = kTRUE;
1521 AddBinSimple(fHistAnode, ic, jc);
1522 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1524 Int_t nPix = fPixArray->GetEntriesFast();
1525 Int_t npad = cluster.Multiplicity();
1527 for (Int_t i = 0; i < nPix; ++i)
1529 AliMUONPad* pixPtr = Pixel(i);
1530 pixPtr->SetSize(0,wx);
1531 pixPtr->SetSize(1,wy);
1534 // Pick up pads which overlap with found pixels
1535 for (Int_t i = 0; i < npad; ++i)
1537 //cluster.Pad(i)->SetStatus(-1);
1538 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1541 for (Int_t i = 0; i < nPix; ++i)
1543 AliMUONPad* pixPtr = Pixel(i);
1544 for (Int_t j = 0; j < npad; ++j)
1546 AliMUONPad* pad = cluster.Pad(j);
1547 //if (pad->Status() == 0) continue;
1548 if (pad->Status() == fgkZero) continue;
1549 if ( Overlap(*pad,*pixPtr) )
1551 //pad->SetStatus(0);
1552 pad->SetStatus(fgkZero);
1553 if (fDebug) { cout << j << " "; pad->Print("full"); }
1561 //_____________________________________________________________________________
1563 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
1565 /// Add adjacent bins (+-1 in X and Y) to the cluster
1567 Int_t nx = hist->GetNbinsX();
1568 Int_t ny = hist->GetNbinsY();
1569 Double_t cont1, cont = hist->GetBinContent(hist->GetBin(jc,ic));
1570 AliMUONPad *pixPtr = 0;
1572 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1573 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1574 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1575 cont1 = hist->GetBinContent(hist->GetBin(j,i));
1576 if (cont1 > cont) continue;
1577 if (cont1 < fLowestPixelCharge) continue;
1578 pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j),
1579 hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1580 fPixArray->Add(pixPtr);
1585 //_____________________________________________________________________________
1586 AliMUONClusterFinderMLEM&
1587 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1589 /// Protected assignement operator
1591 if (this == &rhs) return *this;
1593 AliFatal("Not implemented.");
1598 //_____________________________________________________________________________
1599 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1601 /// Add virtual pad (with small charge) to improve fit for clusters
1602 /// with number of pads == 2 per direction
1604 // Find out non-bending and bending planes
1605 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1607 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1608 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1609 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1614 Bool_t same = kFALSE;
1615 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1618 Bool_t check[2] = {kFALSE, kFALSE};
1620 nxy[0] = nxy[1] = 0;
1621 for (Int_t inb = 0; inb < 2; ++inb) {
1622 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1623 if (inb == 0 && AliMp::PairFirst(cn) == 2) check[inb] = kTRUE; // check non-bending plane
1624 else if (inb == 1 && AliMp::PairSecond(cn) == 2) check[inb] = kTRUE; // check bending plane
1626 nxy[0] = TMath::Max (nxy[0], AliMp::PairFirst(cn));
1627 nxy[1] = TMath::Max (nxy[1], AliMp::PairSecond(cn));
1628 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1629 else if (inb == 1 && AliMp::PairSecond(cn) < 2) nonb[inb] = !nonb[inb];
1633 if (nxy[0] > 2) check[0] = kFALSE;
1634 if (nxy[1] > 2) check[1] = kFALSE;
1636 if (!check[0] && !check[1]) return;
1638 for (Int_t inb = 0; inb < 2; ++inb) {
1639 if (!check[inb]) continue;
1641 // Find pads with maximum and next to maximum charges
1642 Int_t maxPads[2] = {-1, -1};
1643 Double_t amax[2] = {0};
1644 Int_t mult = cluster.Multiplicity();
1645 for (Int_t j = 0; j < mult; ++j) {
1646 AliMUONPad *pad = cluster.Pad(j);
1647 if (pad->Cathode() != nonb[inb]) continue;
1648 for (Int_t j2 = 0; j2 < 2; ++j2) {
1649 if (pad->Charge() > amax[j2]) {
1650 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1651 amax[j2] = pad->Charge();
1658 // Find min and max dimensions of the cluster
1659 Double_t limits[2] = {9999, -9999};
1660 for (Int_t j = 0; j < mult; ++j) {
1661 AliMUONPad *pad = cluster.Pad(j);
1662 if (pad->Cathode() != nonb[inb]) continue;
1663 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1664 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1667 // Loop over max and next to max pads
1668 Bool_t add = kFALSE;
1670 for (Int_t j = 0; j < 2; ++j) {
1672 if (maxPads[j] < 0) continue;
1674 if (amax[1] / amax[0] < 0.5) break;
1676 // Check if pad at the cluster limit
1677 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1679 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1680 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1682 //cout << " *** Pad not at the cluster limit: " << j << endl;
1685 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1687 // Add pad (if it exists)
1689 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1690 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1691 //AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1692 AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos.X(), pos.Y(),kFALSE);
1693 if (!mppad.IsValid()) continue; // non-existing pad
1694 AliMUONPad muonPad(fDetElemId, nonb[inb], mppad.GetIx(), mppad.GetIy(),
1695 mppad.GetPositionX(), mppad.GetPositionY(),
1696 mppad.GetDimensionX(), mppad.GetDimensionY(), 0);
1697 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, fLowestPadCharge));
1698 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1699 else muonPad.SetCharge(TMath::Min (amax[j]/15, fLowestPadCharge));
1700 if (muonPad.Charge() < 2.0*fLowestPixelCharge) muonPad.SetCharge(2.0*fLowestPixelCharge);
1701 muonPad.SetReal(kFALSE);
1702 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1703 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1704 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1705 cluster.AddPad(muonPad); // add pad to the cluster
1712 //_____________________________________________________________________________
1713 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1714 Int_t &nInX, Int_t &nInY) const
1716 /// Find number of pads in X and Y-directions (excluding virtual ones and
1719 //Int_t statusToTest = 1;
1720 Int_t statusToTest = fgkUseForFit;
1722 //if ( nInX < 0 ) statusToTest = 0;
1723 if ( nInX < 0 ) statusToTest = fgkZero;
1725 Bool_t mustMatch(kTRUE);
1727 Long_t cn = cluster.NofPads(statusToTest,mustMatch);
1729 nInX = AliMp::PairFirst(cn);
1730 nInY = AliMp::PairSecond(cn);
1733 //_____________________________________________________________________________
1734 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1736 /// Remove pixel at index i
1737 AliMUONPad* pixPtr = Pixel(i);
1738 fPixArray->RemoveAt(i);
1742 //_____________________________________________________________________________
1743 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1745 /// Process simple cluster (small number of pads) without EM-procedure
1747 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1748 Double_t parOk[3] = {0.};
1749 TObjArray *clusters[1];
1750 clusters[0] = fPixArray;
1752 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1754 Int_t mult = cluster.Multiplicity();
1755 for (Int_t i = 0; i < mult; ++i)
1757 AliMUONPad* pad = cluster.Pad(i);
1759 if ( pad->IsSaturated())
1768 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1770 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
1773 //_____________________________________________________________________________
1775 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1777 /// Returns pixel at index i
1778 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1781 //_____________________________________________________________________________
1783 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1786 TString swhat(what);
1788 if ( swhat.Contains("precluster") )
1790 if ( fPreCluster) fPreCluster->Print();
1794 //_____________________________________________________________________________
1796 AliMUONClusterFinderMLEM::SetChargeHints(Double_t lowestPadCharge, Double_t lowestClusterCharge)
1798 /// Set some thresholds we use later on in the algorithm
1799 fLowestPadCharge=lowestPadCharge;
1800 fLowestClusterCharge=lowestClusterCharge;
1801 fLowestPixelCharge=fLowestPadCharge/12.0;