1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONClusterFinderMLEM
21 /// Clusterizer class based on the Expectation-Maximization algorithm
23 /// Pre-clustering is handled by AliMUONPreClusterFinder
24 /// From a precluster a pixel array is built, and from this pixel array
25 /// a list of clusters is output (using AliMUONClusterSplitterMLEM).
27 /// \author Laurent Aphecetche (for the "new" C++ structure) and
28 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
29 //-----------------------------------------------------------------------------
31 #include "AliMUONClusterFinderMLEM.h"
33 #include "AliMUONCluster.h"
34 #include "AliMUONClusterSplitterMLEM.h"
35 #include "AliMUONVDigit.h"
36 #include "AliMUONPad.h"
37 #include "AliMUONPreClusterFinder.h"
39 #include "AliMpVPadIterator.h"
40 #include "AliMpVSegmentation.h"
41 #include "AliRunLoader.h"
42 #include "AliMUONVDigitStore.h"
43 #include <Riostream.h>
48 #include "AliCodeTimer.h"
51 ClassImp(AliMUONClusterFinderMLEM)
54 const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
55 //const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
56 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
57 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
58 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
60 // Status flags for pads
61 const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
62 const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
63 const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
64 const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
65 const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
66 const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
68 //_____________________________________________________________________________
69 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
70 : AliMUONVClusterFinder(),
71 fPreClusterFinder(clusterFinder),
79 fPixArray(new TObjArray(20)),
88 fSegmentation[1] = fSegmentation[0] = 0x0;
90 if (fPlot) fDebug = 1;
93 //_____________________________________________________________________________
94 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
97 delete fPixArray; fPixArray = 0;
99 delete fPreClusterFinder;
101 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
102 fNClusters,fNAddVirtualPads));
105 //_____________________________________________________________________________
107 AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
108 TClonesArray* pads[2],
109 const AliMpArea& area,
110 const AliMpVSegmentation* seg[2])
112 /// Prepare for clustering
113 // AliCodeTimerAuto("")
115 for ( Int_t i = 0; i < 2; ++i )
117 fSegmentation[i] = seg[i];
120 // Find out the DetElemId
121 fDetElemId = detElemId;
124 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
125 fSplitter->SetDebug(fDebug);
127 // find out current event number, and reset the cluster number
128 AliRunLoader *runLoader = AliRunLoader::Instance();
129 fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
131 fClusterList.Delete();
133 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
135 if ( fPreClusterFinder->NeedSegmentation() )
137 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
141 return fPreClusterFinder->Prepare(detElemId,pads,area);
145 //_____________________________________________________________________________
147 AliMUONClusterFinderMLEM::NextCluster()
149 /// Return next cluster
150 // AliCodeTimerAuto("")
152 // if the list of clusters is not void, pick one from there
153 TObject* o = fClusterList.At(++fClusterNumber);
154 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
156 //FIXME : at this point, must check whether we've used all the digits
157 //from precluster : if not, let the preclustering know about those unused
158 //digits, so it can reuse them
160 // if the cluster list is exhausted, we need to go to the next
161 // pre-cluster and treat it
163 fPreCluster = fPreClusterFinder->NextCluster();
165 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
166 fClusterNumber = -1; //AZ
176 // WorkOnPreCluster may have used only part of the pads, so we check that
177 // now, and let the unused pads be reused by the preclustering...
179 Int_t mult = fPreCluster->Multiplicity();
180 for ( Int_t i = 0; i < mult; ++i )
182 AliMUONPad* pad = fPreCluster->Pad(i);
183 if ( !pad->IsUsed() )
185 fPreClusterFinder->UsePad(*pad);
189 return NextCluster();
192 //_____________________________________________________________________________
194 AliMUONClusterFinderMLEM::WorkOnPreCluster()
196 /// Starting from a precluster, builds a pixel array, and then
197 /// extract clusters from this array
199 // AliCodeTimerAuto("")
202 cout << " *** Event # " << fEventNumber
203 << " det. elem.: " << fDetElemId << endl;
204 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
205 AliMUONPad* pad = fPreCluster->Pad(j);
206 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
207 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
208 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
212 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
213 if (!cluster) return kFALSE;
215 BuildPixArray(*cluster);
217 if ( fPixArray->GetLast() < 0 )
219 AliDebug(1,"No pixel for the above cluster");
224 // Use MLEM for cluster finder
225 Int_t nMax = 1, localMax[100], maxPos[100];
226 Double_t maxVal[100];
228 Int_t iSimple = 0, nInX = -1, nInY;
230 PadsInXandY(*cluster,nInX, nInY);
232 if (nInX < 4 && nInY < 4)
234 iSimple = 1; // simple cluster
238 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
240 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
241 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
245 for (Int_t i = 0; i < nMax; ++i)
249 FindCluster(*cluster,localMax, maxPos[i]);
252 MainLoop(*cluster,iSimple);
256 Int_t mult = cluster->Multiplicity();
257 for (Int_t j = 0; j < mult; ++j)
259 AliMUONPad* pad = cluster->Pad(j);
260 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
261 if ( pad->Status() != fgkOver ) continue; // pad was not used
263 pad->SetStatus(fgkZero);
264 pad->RevertCharge(); // use backup charge value
267 } // for (Int_t i=0; i<nMax;
270 fHistMlem = fHistAnode = 0x0;
275 //_____________________________________________________________________________
277 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
279 /// Check if the pad and the pixel overlaps
281 // make a fake pad from the pixel
282 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
283 pixel.Coord(0),pixel.Coord(1),
284 pixel.Size(0),pixel.Size(1),0);
286 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
289 //_____________________________________________________________________________
291 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
293 /// Check precluster in order to attempt to simplify it (mostly for
294 /// two-cathode preclusters)
298 // Disregard small clusters (leftovers from splitting or noise)
299 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
300 origCluster.Charge(0)+origCluster.Charge(1) < 10)
305 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
307 AliDebug(2,"Start of CheckPreCluster=");
308 //StdoutToAliDebug(2,cluster->Print("full"));
310 AliMUONCluster* rv(0x0);
312 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
314 rv = CheckPreclusterTwoCathodes(cluster);
323 //_____________________________________________________________________________
325 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
327 /// Check two-cathode cluster
329 Int_t npad = cluster->Multiplicity();
330 Int_t* flags = new Int_t[npad];
331 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
333 // Check pad overlaps
334 for ( Int_t i = 0; i < npad; ++i)
336 AliMUONPad* padi = cluster->Pad(i);
337 if ( padi->Cathode() != 0 ) continue;
338 for (Int_t j = i+1; j < npad; ++j)
340 AliMUONPad* padj = cluster->Pad(j);
341 if ( padj->Cathode() != 1 ) continue;
342 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
343 flags[i] = flags[j] = 1; // mark overlapped pads
347 // Check if all pads overlap
349 for (Int_t i = 0; i < npad; ++i)
351 if (!flags[i]) ++nFlags;
356 // not all pads overlap.
357 if (fDebug) cout << " nFlags: " << nFlags << endl;
358 TObjArray toBeRemoved;
359 for (Int_t i = 0; i < npad; ++i)
361 AliMUONPad* pad = cluster->Pad(i);
362 if (flags[i]) continue;
363 Int_t cath = pad->Cathode();
364 Int_t cath1 = TMath::Even(cath);
365 // Check for edge effect (missing pads on the _other_ cathode)
366 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
367 if (!mpPad.IsValid()) continue;
368 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
369 if (nFlags == 1 && pad->Charge() < 20) continue;
370 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
371 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
372 toBeRemoved.AddLast(pad);
373 fPreCluster->Pad(i)->Release();
375 Int_t nRemove = toBeRemoved.GetEntriesFast();
376 for ( Int_t i = 0; i < nRemove; ++i )
378 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
382 // Check correlations of cathode charges
383 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
386 Int_t cathode = cluster->MaxRawChargeCathode();
392 // get min and max pad charges on the cathode opposite to the
393 // max pad (given by MaxRawChargeCathode())
395 Int_t mult = cluster->Multiplicity();
396 for ( Int_t i = 0; i < mult; ++i )
398 AliMUONPad* pad = cluster->Pad(i);
399 if ( pad->Cathode() != cathode || !pad->IsReal() )
401 // only consider pads in the opposite cathode, and
402 // only consider real pads (i.e. exclude the virtual ones)
405 if ( pad->Charge() < cmin )
407 cmin = pad->Charge();
414 else if ( pad->Charge() > cmax )
416 cmax = pad->Charge();
420 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
421 imin,imax,cmin,cmax));
423 // arrange pads according to their distance to the max, normalized
425 Double_t* dist = new Double_t[mult];
430 AliMUONPad* padmax = cluster->Pad(imax);
432 for ( Int_t i = 0; i < mult; ++i )
435 if ( i == imax) continue;
436 AliMUONPad* pad = cluster->Pad(i);
437 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
438 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
439 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
440 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
443 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
449 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
450 Double_t xmax(-1), distPrev(999);
451 TObjArray toBeRemoved;
453 for ( Int_t i = 0; i < mult; ++i )
455 Int_t indx = flags[i];
456 AliMUONPad* pad = cluster->Pad(indx);
457 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
458 if ( dist[indx] > dmin )
460 // farther than the minimum pad
461 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
462 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
465 if (dx >= 0 && dy >= 0) continue;
466 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
467 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
469 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
470 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
473 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
475 cmax = TMath::Max(pad->Charge(),cmax);
479 cmax = pad->Charge();
482 distPrev = dist[indx];
483 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
484 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
487 toBeRemoved.AddLast(pad);
488 fPreCluster->Pad(indx)->Release();
491 Int_t nRemove = toBeRemoved.GetEntriesFast();
492 for ( Int_t i = 0; i < nRemove; ++i )
494 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
501 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
502 //StdoutToAliDebug(2,cluster->Print("full"));
507 //_____________________________________________________________________________
509 AliMUONClusterFinderMLEM::CheckOverlaps()
511 /// For debug only : check if some pixels overlap...
513 Int_t nPix = fPixArray->GetLast()+1;
516 for ( Int_t i = 0; i < nPix; ++i )
518 AliMUONPad* pixelI = Pixel(i);
519 AliMUONPad pi(dummy,dummy,dummy,dummy,
520 pixelI->Coord(0),pixelI->Coord(1),
521 pixelI->Size(0),pixelI->Size(1),0.0);
523 for ( Int_t j = i+1; j < nPix; ++j )
525 AliMUONPad* pixelJ = Pixel(j);
526 AliMUONPad pj(dummy,dummy,dummy,dummy,
527 pixelJ->Coord(0),pixelJ->Coord(1),
528 pixelJ->Size(0),pixelJ->Size(1),0.0);
531 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
533 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
535 StdoutToAliInfo(pixelI->Print();
536 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
538 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
539 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
540 cout << "-------" << endl;
548 //_____________________________________________________________________________
549 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
551 /// Build pixel array for MLEM method
553 Int_t npad = cluster.Multiplicity();
556 AliWarning("Got no pad at all ?!");
560 BuildPixArrayOneCathode(cluster);
562 Int_t nPix = fPixArray->GetLast()+1;
564 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
568 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
570 // Too many pixels - sort and remove pixels with the lowest signal
572 for ( Int_t i = npad; i < nPix; ++i )
576 fPixArray->Compress();
577 } // if (nPix > npad)
579 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
580 // fPixArray->Print(););
581 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
584 //_____________________________________________________________________________
585 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
587 /// Build the pixel array
589 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
591 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
592 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2]={99999,99999};
593 Int_t found[2] = {0,0}, mult = cluster.Multiplicity();
595 for ( Int_t i = 0; i < mult; ++i) {
596 AliMUONPad* pad = cluster.Pad(i);
597 for (Int_t j = 0; j < 2; ++j) {
598 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
599 xy0[j] = pad->Coord(j);
603 if (found[0] && found[1]) break;
606 Double_t min[2], max[2];
607 Int_t cath0 = 0, cath1 = 1;
608 if (cluster.Multiplicity(0) == 0) cath0 = 1;
609 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
611 TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
612 TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
613 min[0] = leftDown.X();
614 min[1] = leftDown.Y();
615 max[0] = rightUp.X();
616 max[1] = rightUp.Y();
617 if (cath1 != cath0) {
618 leftDown = cluster.Area(cath1).LeftDownCorner();
619 rightUp = cluster.Area(cath1).RightUpCorner();
620 min[0] = TMath::Max (min[0], leftDown.X());
621 min[1] = TMath::Max (min[1], leftDown.Y());
622 max[0] = TMath::Min (max[0], rightUp.X());
623 max[1] = TMath::Min (max[1], rightUp.Y());
627 //width[0] /= 2; width[1] /= 2; // just for check
628 Int_t nbins[2]={0,0};
629 for (Int_t i = 0; i < 2; ++i) {
630 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
631 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
632 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
633 + TMath::Sign(0.5,dist)) * width[i] * 2;
634 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
635 if (nbins[i] == 0) ++nbins[i];
636 max[i] = min[i] + nbins[i] * width[i] * 2;
637 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
641 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
642 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
643 TAxis *xaxis = hist1->GetXaxis();
644 TAxis *yaxis = hist1->GetYaxis();
647 for ( Int_t i = 0; i < mult; ++i) {
648 AliMUONPad* pad = cluster.Pad(i);
649 Int_t ix0 = xaxis->FindBin(pad->X());
650 Int_t iy0 = yaxis->FindBin(pad->Y());
651 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
655 for (Int_t i = 1; i <= nbins[0]; ++i) {
656 Double_t x = xaxis->GetBinCenter(i);
657 for (Int_t j = 1; j <= nbins[1]; ++j) {
658 if (hist2->GetCellContent(i,j) < 0.1) continue;
659 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
660 // cluster.Multiplicity(1)) continue;
661 if (cath0 != cath1) {
663 Double_t cont = hist2->GetCellContent(i,j);
664 if (cont < 999.) continue;
665 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
667 Double_t y = yaxis->GetBinCenter(j);
668 Double_t charge = hist1->GetCellContent(i,j);
669 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
670 fPixArray->Add(pixPtr);
674 if (fPixArray->GetEntriesFast() == 1) {
675 // Split pixel into 2
676 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
677 pixPtr->SetSize(0,width[0]/2.);
678 pixPtr->Shift(0,-width[0]/4.);
679 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
680 fPixArray->Add(pixPtr);
683 //fPixArray->Print();
688 //_____________________________________________________________________________
689 void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
690 TH2D *hist1, TH2D *hist2)
692 /// "Span" pad over histogram in the direction idir
694 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
695 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
696 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
698 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
700 for (Int_t i = 0; i < nbinPad; ++i) {
701 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
702 if (ixy > nbins) break;
703 Double_t lowEdge = axis->GetBinLowEdge(ixy);
704 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
705 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
708 Double_t cont = pad->Charge();
709 if (hist2->GetCellContent(ix0, ixy) > 0.1)
710 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
711 hist1->SetCellContent(ix0, ixy, cont);
712 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
713 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
717 for (Int_t i = -1; i > -nbinPad; --i) {
718 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
720 Double_t upEdge = axis->GetBinUpEdge(ixy);
721 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
722 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
725 Double_t cont = pad->Charge();
726 if (hist2->GetCellContent(ix0, ixy) > 0.1)
727 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
728 hist1->SetCellContent(ix0, ixy, cont);
729 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
730 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
735 //_____________________________________________________________________________
737 AliMUONClusterFinderMLEM::Plot(const char* basename)
739 /// Make a plot and save it as png
744 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
749 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
750 fDetElemId,fClusterNumber));
753 //_____________________________________________________________________________
755 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
759 /// Compute coefficients needed for MLEM algorithm
761 Int_t npadTot = cluster.Multiplicity();
762 Int_t nPix = fPixArray->GetLast()+1;
764 //memset(probi,0,nPix*sizeof(Double_t));
765 for (Int_t j = 0; j < npadTot*nPix; ++j) coef[j] = 0.;
766 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
768 Int_t mult = cluster.Multiplicity();
769 for ( Int_t j = 0; j < mult; ++j )
771 AliMUONPad* pad = cluster.Pad(j);
774 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
776 Int_t indx1 = indx + ipix;
777 //if (pad->Status() < 0)
778 if (pad->Status() != fgkZero)
783 AliMUONPad* pixPtr = Pixel(ipix);
784 // coef is the charge (given by Mathieson integral) on pad, assuming
785 // the Mathieson is center at pixel.
786 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
787 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
788 // pad->Ix(),pad->Iy(),
789 // pad->X(),pad->Y(),
790 // pad->DX(),pad->DY(),
791 // pixPtr->Coord(0),pixPtr->Coord(1),
792 // pixPtr->Size(0),pixPtr->Size(1),
795 probi[ipix] += coef[indx1];
800 //_____________________________________________________________________________
801 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
803 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
805 // AliCodeTimerAuto("")
807 Int_t nPix = fPixArray->GetLast()+1;
809 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
810 //StdoutToAliDebug(2,cluster.Print("full"););
814 AliDebug(1,"No pixels, why am I here then ?");
817 AddVirtualPad(cluster); // add virtual pads if necessary
819 Int_t npadTot = cluster.Multiplicity();
821 for (Int_t i = 0; i < npadTot; ++i)
823 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
824 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
828 Double_t* probi(0x0);
829 Int_t lc(0); // loop counter
831 //Plot("mlem.start");
832 AliMUONPad* pixPtr = Pixel(0);
833 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
840 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
841 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
842 //StdoutToAliDebug(2,fPixArray->Print("","full"));
844 coef = new Double_t [npadTot*nPix];
845 probi = new Double_t [nPix];
847 // Calculate coefficients and pixel visibilities
848 ComputeCoefficients(cluster,coef,probi);
850 for (Int_t ipix = 0; ipix < nPix; ++ipix)
852 if (probi[ipix] < 0.01)
854 AliMUONPad* pixel = Pixel(ipix);
855 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
856 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
857 pixel->SetCharge(0); // "invisible" pixel
862 Mlem(cluster,coef, probi, 15);
864 // Find histogram limits for the 1'st pass only - for others computed below
866 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
868 pixPtr = Pixel(ipix);
869 for ( Int_t i = 0; i < 2; ++i )
872 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
873 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
876 } else pixPtr = Pixel(0);
878 for (Int_t i = 0; i < 4; i++)
880 xylim[i] -= pixPtr->Size(i/2);
883 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
884 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
886 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
887 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
888 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
889 xylim[0],-xylim[1],xylim[2],-xylim[3]
892 fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
894 for (Int_t ipix = 0; ipix < nPix; ++ipix)
896 AliMUONPad* pixPtr2 = Pixel(ipix);
897 fHistMlem->Fill(pixPtr2->Coord(0),pixPtr2->Coord(1),pixPtr2->Charge());
900 // Check if the total charge of pixels is too low
902 for ( Int_t i = 0; i < nPix; ++i)
904 qTot += Pixel(i)->Charge();
907 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
909 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
913 for ( Int_t i = 0; i < npadTot; ++i)
915 AliMUONPad* pad = cluster.Pad(i);
916 //if ( pad->Status() == 0) pad->SetStatus(-1);
917 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
924 // Simple cluster - skip further passes thru EM-procedure
932 // Calculate position of the center-of-gravity around the maximum pixel
936 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
937 pixPtr->Size(0) > pixPtr->Size(1)) break;
939 // Sort pixels according to the charge
940 MaskPeaks(1); // mask local maxima
942 MaskPeaks(0); // unmask local maxima
943 Double_t pixMin = 0.01*Pixel(0)->Charge();
944 pixMin = TMath::Min(pixMin,50.);
946 // Decrease pixel size and shift pixels to make them centered at
948 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
951 Double_t shift[2] = { 0.0, 0.0 };
952 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
955 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
957 AliMUONPad* pixPtr2 = Pixel(ipix);
958 if ( nPix >= npadOK // too many pixels already
960 ((pixPtr2->Charge() < pixMin) && (pixPtr2->Status() != fgkMustKeep)) // too low charge
966 for (Int_t i = 0; i < 2; ++i)
970 pixPtr2->SetCharge(10);
971 pixPtr2->SetSize(indx, pixPtr2->Size(indx)/2);
972 width = -pixPtr2->Size(indx);
973 pixPtr2->Shift(indx, width);
974 // Shift pixel position
978 for (Int_t j = 0; j < 2; ++j)
980 shift[j] = pixPtr2->Coord(j) - xyCOG[j];
981 shift[j] -= ((Int_t)(shift[j]/pixPtr2->Size(j)/2))*pixPtr2->Size(j)*2;
984 pixPtr2->Shift(0, -shift[0]);
985 pixPtr2->Shift(1, -shift[1]);
988 else if (nPix < npadOK)
990 pixPtr2 = new AliMUONPad(*pixPtr2);
991 pixPtr2->Shift(indx, -2*width);
992 pixPtr2->SetStatus(fgkZero);
993 fPixArray->Add(pixPtr2);
996 else continue; // skip adjustment of histo limits
997 for (Int_t j = 0; j < 4; ++j)
999 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr2->Coord(j/2));
1001 } // for (Int_t i=0; i<2;
1002 } // for (Int_t ipix=0;
1004 fPixArray->Compress();
1006 AliDebug(2,Form("After shift:"));
1007 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1008 //Plot(Form("mlem.lc%d",lc+1));
1010 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1013 xylim[2],xylim[3]));
1017 AliMUONPad* pixPtr2 = Pixel(0);
1018 // add pixels if the maximum is at the limit of pixel area:
1019 // start from Y-direction
1021 for (Int_t i = 3; i > -1; --i)
1023 if (nPix < npadOK &&
1024 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr2->Size(i/2))
1026 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1027 AliMUONPad* p = new AliMUONPad(*pixPtr2);
1028 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr2->Size(i/2));
1029 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1030 j = TMath::Even (i/2);
1031 p->SetCoord(j, xyCOG[j]);
1032 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1033 //StdoutToAliDebug(2,cout << " ---- ";
1034 // p->Print("corners"););
1040 nPix = fPixArray->GetEntriesFast();
1045 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1046 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1048 // remove pixels with low signal or low visibility
1049 // Cuts are empirical !!!
1050 Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,1.);
1051 thresh = TMath::Min (thresh,50.);
1052 Double_t charge = 0;
1054 // Mark pixels which should be removed
1055 for (Int_t i = 0; i < nPix; ++i)
1057 AliMUONPad* pixPtr2 = Pixel(i);
1058 charge = pixPtr2->Charge();
1059 if (charge < thresh)
1061 pixPtr2->SetCharge(-charge);
1065 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1067 for (Int_t i = 0; i < nPix; ++i)
1069 AliMUONPad* pixPtr2 = Pixel(i);
1070 charge = pixPtr2->Charge();
1071 if (charge > 0) continue;
1072 near = FindNearest(pixPtr2);
1073 pixPtr2->SetCharge(0);
1074 probi[i] = 0; // make it "invisible"
1075 AliMUONPad* pnear = Pixel(near);
1076 pnear->SetCharge(pnear->Charge() + (-charge));
1078 Mlem(cluster,coef,probi,2);
1080 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1081 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1082 //Plot("mlem.beforesplit");
1085 for (Int_t i = 0; i < nPix; ++i)
1087 AliMUONPad* pixPtr2 = Pixel(i);
1088 Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr2->Coord(0));
1089 Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr2->Coord(1));
1090 fHistMlem->SetBinContent(ix, iy, pixPtr2->Charge());
1093 // Try to split into clusters
1095 if (fHistMlem->GetSum() < 1)
1101 fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
1106 fPixArray->Delete();
1111 //_____________________________________________________________________________
1112 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1114 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1115 /// - to avoid loosing low charge pixels after sorting)
1117 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1118 AliMUONPad* pix = Pixel(i);
1119 if (pix->Status() == fgkMustKeep) {
1120 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1121 else pix->SetCharge(pix->Charge()-10000.);
1126 //_____________________________________________________________________________
1127 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1128 Double_t* coef, Double_t* probi,
1131 /// Use MLEM to find pixel charges
1133 Int_t nPix = fPixArray->GetEntriesFast();
1135 Int_t npad = cluster.Multiplicity();
1137 Double_t* probi1 = new Double_t[nPix];
1138 Double_t probMax = TMath::MaxElement(nPix,probi);
1140 for (Int_t iter = 0; iter < nIter; ++iter)
1143 for (Int_t ipix = 0; ipix < nPix; ++ipix)
1145 Pixel(ipix)->SetChargeBackup(0);
1146 // Correct each pixel
1148 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1150 probi1[ipix] = probMax;
1151 for (Int_t j = 0; j < npad; ++j)
1153 AliMUONPad* pad = cluster.Pad(j);
1154 //if (pad->Status() < 0) continue;
1155 if (pad->Status() != fgkZero) continue;
1157 Int_t indx1 = j*nPix;
1158 Int_t indx = indx1 + ipix;
1159 // Calculate expectation
1160 for (Int_t i = 0; i < nPix; ++i)
1162 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1163 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1165 if ( pad->IsSaturated() && sum1 > pad->Charge() )
1167 // correct for pad charge overflows
1168 probi1[ipix] -= coef[indx];
1170 //sum1 = pad->Charge();
1173 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1174 } // for (Int_t j=0;
1175 AliMUONPad* pixPtr = Pixel(ipix);
1176 if (probi1[ipix] > 1.e-6)
1178 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1179 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1181 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1182 } // for (Int_t ipix=0;
1184 for (Int_t i = 0; i < nPix; ++i) {
1185 AliMUONPad* pixPtr = Pixel(i);
1186 pixPtr->RevertCharge();
1187 qTot += pixPtr->Charge();
1190 // Can happen in clusters with large number of overflows - speeding up
1194 } // for (Int_t iter=0;
1198 //_____________________________________________________________________________
1199 void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
1201 /// Calculate position of the center-of-gravity around the maximum pixel
1203 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1204 Int_t i1 = -9, j1 = -9;
1205 fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1206 Int_t nx = fHistMlem->GetNbinsX();
1207 Int_t ny = fHistMlem->GetNbinsY();
1208 Double_t thresh = fHistMlem->GetMaximum()/10;
1209 Double_t x, y, cont, xq=0, yq=0, qq=0;
1211 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1212 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1213 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1214 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1215 cont = fHistMlem->GetCellContent(j,i);
1216 if (cont < thresh) continue;
1217 if (i != i1) {i1 = i; nsumy++;}
1218 if (j != j1) {j1 = j; nsumx++;}
1219 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1228 Int_t i2 = 0, j2 = 0;
1231 // one bin in Y - add one more (with the largest signal)
1232 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1233 if (i == iymax) continue;
1234 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1235 cont = fHistMlem->GetCellContent(j,i);
1238 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1239 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1248 if (i2 != i1) nsumy++;
1249 if (j2 != j1) nsumx++;
1251 } // if (nsumy == 1)
1254 // one bin in X - add one more (with the largest signal)
1256 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1257 if (j == ixmax) continue;
1258 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1259 cont = fHistMlem->GetCellContent(j,i);
1262 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1263 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1272 if (i2 != i1) nsumy++;
1273 if (j2 != j1) nsumx++;
1275 } // if (nsumx == 1)
1277 xyc[0] = xq/qq; xyc[1] = yq/qq;
1278 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1281 //_____________________________________________________________________________
1282 Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1284 /// Find the pixel nearest to the given one
1285 /// (algorithm may be not very efficient)
1287 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1288 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1289 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1292 for (Int_t i = 0; i < nPix; ++i) {
1293 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1294 if (pixPtr == pixPtr0 || pixPtr->Charge() < 0.5) continue;
1295 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1296 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1297 r = dx *dx + dy * dy;
1298 if (r < rmin) { rmin = r; imin = i; }
1303 //_____________________________________________________________________________
1305 AliMUONClusterFinderMLEM::Paint(Option_t*)
1307 /// Paint cluster and pixels
1309 AliMpArea area(fPreCluster->Area());
1311 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1313 gVirtualX->SetFillStyle(1001);
1314 gVirtualX->SetFillColor(3);
1315 gVirtualX->SetLineColor(3);
1319 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1321 AliMUONPad* pixel = Pixel(i);
1323 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1324 pixel->Coord(1)-pixel->Size(1)*s,
1325 pixel->Coord(0)+pixel->Size(0)*s,
1326 pixel->Coord(1)+pixel->Size(1)*s);
1328 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1330 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1331 // pixel->Coord(1) + sign*pixel->Size(1),
1332 // pixel->Coord(0) + pixel->Size(0),
1333 // pixel->Coord(1) - sign*pixel->Size(1)
1339 gVirtualX->SetFillStyle(0);
1341 fPreCluster->Paint();
1343 gVirtualX->SetLineColor(1);
1344 gVirtualX->SetLineWidth(2);
1345 gVirtualX->SetFillStyle(0);
1346 gVirtualX->SetTextColor(1);
1347 gVirtualX->SetTextAlign(22);
1349 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1351 AliMUONPad* pixel = Pixel(i);
1352 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1353 pixel->Coord(1)-pixel->Size(1),
1354 pixel->Coord(0)+pixel->Size(0),
1355 pixel->Coord(1)+pixel->Size(1));
1356 gVirtualX->SetTextSize(pixel->Size(0)*60);
1358 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1362 //_____________________________________________________________________________
1363 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1365 /// Find local maxima in pixel space for large preclusters in order to
1366 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1367 /// or to find additional fitting seeds if clusters were not completely resolved
1369 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1371 Double_t xylim[4] = {999, 999, 999, 999};
1373 Int_t nPix = pixArray->GetEntriesFast();
1374 AliMUONPad *pixPtr = 0;
1375 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1376 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1377 for (Int_t i = 0; i < 4; ++i)
1378 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1380 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
1382 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1383 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1384 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1385 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1386 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1387 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1388 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1390 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1392 Int_t nMax = 0, indx, nxy = ny * nx;
1393 Int_t *isLocalMax = new Int_t[nxy];
1394 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
1396 for (Int_t i = 1; i <= ny; ++i) {
1398 for (Int_t j = 1; j <= nx; ++j) {
1399 if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
1400 //if (isLocalMax[indx+j-1] < 0) continue;
1401 if (isLocalMax[indx+j-1] != 0) continue;
1402 FlagLocalMax(fHistAnode, i, j, isLocalMax);
1406 for (Int_t i = 1; i <= ny; ++i) {
1408 for (Int_t j = 1; j <= nx; ++j) {
1409 if (isLocalMax[indx+j-1] > 0) {
1410 localMax[nMax] = indx + j - 1;
1411 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
1412 ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
1413 if (nMax > 99) break;
1417 AliError(" Too many local maxima !!!");
1421 if (fDebug) cout << " Local max: " << nMax << endl;
1422 delete [] isLocalMax;
1426 //_____________________________________________________________________________
1427 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1429 /// Flag pixels (whether or not local maxima)
1431 Int_t nx = hist->GetNbinsX();
1432 Int_t ny = hist->GetNbinsY();
1433 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1434 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1436 Int_t ie = i + 2, je = j + 2;
1437 for (Int_t i1 = i-1; i1 < ie; ++i1) {
1438 if (i1 < 1 || i1 > ny) continue;
1439 indx1 = (i1 - 1) * nx;
1440 for (Int_t j1 = j-1; j1 < je; ++j1) {
1441 if (j1 < 1 || j1 > nx) continue;
1442 if (i == i1 && j == j1) continue;
1443 indx2 = indx1 + j1 - 1;
1444 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1445 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1446 else if (cont > cont1) isLocalMax[indx2] = -1;
1447 else { // the same charge
1448 isLocalMax[indx] = 1;
1449 if (isLocalMax[indx2] == 0) {
1450 FlagLocalMax(hist, i1, j1, isLocalMax);
1451 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1452 else isLocalMax[indx2] = -1;
1457 isLocalMax[indx] = 1; // local maximum
1460 //_____________________________________________________________________________
1461 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1462 Int_t *localMax, Int_t iMax)
1464 /// Find pixel cluster around local maximum \a iMax and pick up pads
1465 /// overlapping with it
1468 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1470 hist->Draw("lego1Fb"); // debug
1475 Int_t nx = fHistAnode->GetNbinsX();
1476 Int_t ny = fHistAnode->GetNbinsY();
1477 Int_t ic = localMax[iMax] / nx + 1;
1478 Int_t jc = localMax[iMax] % nx + 1;
1479 Int_t nxy = ny * nx;
1480 Bool_t *used = new Bool_t[nxy];
1481 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
1483 // Drop all pixels from the array - pick up only the ones from the cluster
1484 fPixArray->Delete();
1486 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1487 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1488 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1489 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1490 Double_t cont = fHistAnode->GetCellContent(jc,ic);
1491 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1492 used[(ic-1)*nx+jc-1] = kTRUE;
1493 AddBinSimple(fHistAnode, ic, jc);
1494 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1496 Int_t nPix = fPixArray->GetEntriesFast();
1497 Int_t npad = cluster.Multiplicity();
1499 for (Int_t i = 0; i < nPix; ++i)
1501 AliMUONPad* pixPtr = Pixel(i);
1502 pixPtr->SetSize(0,wx);
1503 pixPtr->SetSize(1,wy);
1506 // Pick up pads which overlap with found pixels
1507 for (Int_t i = 0; i < npad; ++i)
1509 //cluster.Pad(i)->SetStatus(-1);
1510 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1513 for (Int_t i = 0; i < nPix; ++i)
1515 AliMUONPad* pixPtr = Pixel(i);
1516 for (Int_t j = 0; j < npad; ++j)
1518 AliMUONPad* pad = cluster.Pad(j);
1519 //if (pad->Status() == 0) continue;
1520 if (pad->Status() == fgkZero) continue;
1521 if ( Overlap(*pad,*pixPtr) )
1523 //pad->SetStatus(0);
1524 pad->SetStatus(fgkZero);
1525 if (fDebug) { cout << j << " "; pad->Print("full"); }
1533 //_____________________________________________________________________________
1535 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
1537 /// Add adjacent bins (+-1 in X and Y) to the cluster
1539 Int_t nx = hist->GetNbinsX();
1540 Int_t ny = hist->GetNbinsY();
1541 Double_t cont1, cont = hist->GetCellContent(jc,ic);
1542 AliMUONPad *pixPtr = 0;
1544 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1545 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1546 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1547 cont1 = hist->GetCellContent(j,i);
1548 if (cont1 > cont) continue;
1549 if (cont1 < 0.5) continue;
1550 pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j),
1551 hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1552 fPixArray->Add(pixPtr);
1557 //_____________________________________________________________________________
1558 AliMUONClusterFinderMLEM&
1559 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1561 /// Protected assignement operator
1563 if (this == &rhs) return *this;
1565 AliFatal("Not implemented.");
1570 //_____________________________________________________________________________
1571 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1573 /// Add virtual pad (with small charge) to improve fit for clusters
1574 /// with number of pads == 2 per direction
1576 // Find out non-bending and bending planes
1577 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1579 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1580 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1581 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1586 Bool_t same = kFALSE;
1587 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1590 Bool_t check[2] = {kFALSE, kFALSE};
1592 nxy[0] = nxy[1] = 0;
1593 for (Int_t inb = 0; inb < 2; ++inb) {
1594 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1595 if (inb == 0 && cn.GetFirst() == 2) check[inb] = kTRUE; // check non-bending plane
1596 else if (inb == 1 && cn.GetSecond() == 2) check[inb] = kTRUE; // check bending plane
1598 nxy[0] = TMath::Max (nxy[0], cn.GetFirst());
1599 nxy[1] = TMath::Max (nxy[1], cn.GetSecond());
1600 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1601 else if (inb == 1 && cn.GetSecond() < 2) nonb[inb] = !nonb[inb];
1605 if (nxy[0] > 2) check[0] = kFALSE;
1606 if (nxy[1] > 2) check[1] = kFALSE;
1608 if (!check[0] && !check[1]) return;
1610 for (Int_t inb = 0; inb < 2; ++inb) {
1611 if (!check[inb]) continue;
1613 // Find pads with maximum and next to maximum charges
1614 Int_t maxPads[2] = {-1, -1};
1615 Double_t amax[2] = {0};
1616 Int_t mult = cluster.Multiplicity();
1617 for (Int_t j = 0; j < mult; ++j) {
1618 AliMUONPad *pad = cluster.Pad(j);
1619 if (pad->Cathode() != nonb[inb]) continue;
1620 for (Int_t j2 = 0; j2 < 2; ++j2) {
1621 if (pad->Charge() > amax[j2]) {
1622 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1623 amax[j2] = pad->Charge();
1630 // Find min and max dimensions of the cluster
1631 Double_t limits[2] = {9999, -9999};
1632 for (Int_t j = 0; j < mult; ++j) {
1633 AliMUONPad *pad = cluster.Pad(j);
1634 if (pad->Cathode() != nonb[inb]) continue;
1635 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1636 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1639 // Loop over max and next to max pads
1640 Bool_t add = kFALSE;
1642 for (Int_t j = 0; j < 2; ++j) {
1644 if (maxPads[j] < 0) continue;
1646 if (amax[1] / amax[0] < 0.5) break;
1648 // Check if pad at the cluster limit
1649 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1651 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1652 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1654 //cout << " *** Pad not at the cluster limit: " << j << endl;
1657 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1659 // Add pad (if it exists)
1661 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1662 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1663 //AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1664 AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kFALSE);
1665 if (!mppad.IsValid()) continue; // non-existing pad
1666 cn = mppad.GetIndices();
1667 AliMUONPad muonPad(fDetElemId, nonb[inb], cn.GetFirst(), cn.GetSecond(),
1668 mppad.Position().X(), mppad.Position().Y(),
1669 mppad.Dimensions().X(), mppad.Dimensions().Y(), 0);
1670 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, 5.));
1671 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1672 else muonPad.SetCharge(TMath::Min (amax[j]/15, 6.));
1673 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1674 muonPad.SetReal(kFALSE);
1675 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1676 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1677 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1678 cluster.AddPad(muonPad); // add pad to the cluster
1685 //_____________________________________________________________________________
1686 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1687 Int_t &nInX, Int_t &nInY) const
1689 /// Find number of pads in X and Y-directions (excluding virtual ones and
1692 //Int_t statusToTest = 1;
1693 Int_t statusToTest = fgkUseForFit;
1695 //if ( nInX < 0 ) statusToTest = 0;
1696 if ( nInX < 0 ) statusToTest = fgkZero;
1698 Bool_t mustMatch(kTRUE);
1700 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1702 nInX = cn.GetFirst();
1703 nInY = cn.GetSecond();
1706 //_____________________________________________________________________________
1707 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1709 /// Remove pixel at index i
1710 AliMUONPad* pixPtr = Pixel(i);
1711 fPixArray->RemoveAt(i);
1715 //_____________________________________________________________________________
1716 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1718 /// Process simple cluster (small number of pads) without EM-procedure
1720 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1721 Double_t parOk[3] = {0.};
1722 TObjArray *clusters[1];
1723 clusters[0] = fPixArray;
1725 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1727 Int_t mult = cluster.Multiplicity();
1728 for (Int_t i = 0; i < mult; ++i)
1730 AliMUONPad* pad = cluster.Pad(i);
1732 if ( pad->IsSaturated())
1741 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1743 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
1746 //_____________________________________________________________________________
1748 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1750 /// Returns pixel at index i
1751 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1754 //_____________________________________________________________________________
1756 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1759 TString swhat(what);
1761 if ( swhat.Contains("precluster") )
1763 if ( fPreCluster) fPreCluster->Print();