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 /// \class AliMUONClusterFinderMLEM
20 /// Clusterizer class based on the Expectation-Maximization algorithm
22 /// Pre-clustering is handled by AliMUONPreClusterFinder
23 /// From a precluster a pixel array is built, and from this pixel array
24 /// a list of clusters is output (using AliMUONClusterSplitterMLEM).
26 /// \author Laurent Aphecetche (for the "new" C++ structure) and
27 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
29 #include "AliMUONClusterFinderMLEM.h"
32 #include "AliMUONCluster.h"
33 #include "AliMUONClusterSplitterMLEM.h"
34 #include "AliMUONDigit.h"
35 #include "AliMUONPad.h"
36 #include "AliMUONPreClusterFinder.h"
38 #include "AliMpVPadIterator.h"
39 #include "AliMpVSegmentation.h"
40 #include "AliRunLoader.h"
41 #include <Riostream.h>
45 #include "TStopwatch.h"
48 ClassImp(AliMUONClusterFinderMLEM)
51 const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
52 const Double_t AliMUONClusterFinderMLEM::fgkSaturation = 3000; // average saturation level
53 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
54 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
55 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
57 TMinuit* AliMUONClusterFinderMLEM::fgMinuit = 0x0;
59 //_____________________________________________________________________________
60 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot)
61 : AliMUONVClusterFinder(),
62 fPreClusterFinder(new AliMUONPreClusterFinder),
71 fPixArray(new TObjArray(20)),
74 fTimers(new TObjArray(kLast)),
81 fSegmentation[1] = fSegmentation[0] = 0x0;
83 fPadBeg[0] = fPadBeg[1] = fCathBeg;
85 if (!fgMinuit) fgMinuit = new TMinuit(8);
87 fTimers->SetOwner(kTRUE);
89 for ( Int_t i = 0; i < kLast; ++i )
91 TStopwatch* t = new TStopwatch;
92 fTimers->AddLast(new TStopwatch);
98 //_____________________________________________________________________________
99 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
102 delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
104 delete fPreClusterFinder;
105 for ( Int_t i = 0; i < kLast; ++i )
107 AliInfo(Form("Timer %d",i));
112 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
113 fNClusters,fNAddVirtualPads));
116 //_____________________________________________________________________________
118 AliMUONClusterFinderMLEM::Prepare(const AliMpVSegmentation* segmentations[2],
119 TClonesArray* digits[2])
121 /// Prepare for clustering
123 for ( Int_t i = 0; i < 2; ++i )
125 fSegmentation[i] = segmentations[i];
128 // Find out the DetElemId
131 for ( Int_t i = 0; i < 2; ++i )
133 AliMUONDigit* d = static_cast<AliMUONDigit*>(digits[i]->First());
136 fDetElemId = d->DetElemId();
141 if ( fDetElemId < 0 )
143 AliWarning("Could not find DE. Probably no digits at all ?");
148 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
150 // find out current event number, and reset the cluster number
151 fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
153 fClusterList.Delete();
155 // AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
157 return fPreClusterFinder->Prepare(segmentations,digits);
160 //_____________________________________________________________________________
162 AliMUONClusterFinderMLEM::NextCluster()
164 /// Return next cluster
168 // if the list of clusters is not void, pick one from there
169 if ( fClusterList.GetLast() >= 0 )
171 TObject* o = fClusterList.At(0);
172 fClusterList.RemoveAt(0);
173 return static_cast<AliMUONCluster*>(o);
176 //FIXME : at this point, must check whether we've used all the digits
177 //from precluster : if not, let the preclustering know about those unused
178 //digits, so it can reuse them
180 // if the cluster list is exhausted, we need to go to the next
181 // pre-cluster and treat it
183 fPreCluster = fPreClusterFinder->NextCluster();
191 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
195 // WorkOnPreCluster may have used only part of the pads, so we check that
196 // now, and let the unused pads be reused by the preclustering...
198 for ( Int_t i = 0; i < fPreCluster->Multiplicity(); ++i )
200 AliMUONPad* pad = fPreCluster->Pad(i);
201 if ( !pad->IsUsed() )
203 fPreClusterFinder->UsePad(*pad);
207 return NextCluster();
210 //_____________________________________________________________________________
212 AliMUONClusterFinderMLEM::WorkOnPreCluster()
214 /// Starting from a precluster, builds a pixel array, and then
215 /// extract clusters from this array
217 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
219 if (!cluster) return kFALSE;
221 BuildPixArray(*cluster);
223 if ( fPixArray->GetLast() < 0 )
225 AliDebug(1,"No pixel for the above cluster");
230 // Use MLEM for cluster finder
231 Int_t nMax = 1, localMax[100], maxPos[100];
232 Double_t maxVal[100];
234 if (cluster->Multiplicity() > 50)
236 nMax = FindLocalMaxima(fPixArray, localMax, maxVal);
241 TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
244 Int_t iSimple = 0, nInX = -1, nInY;
246 PadsInXandY(*cluster,nInX, nInY);
248 if (nMax == 1 && nInX < 4 && nInY < 4)
250 iSimple = 1; //1; // simple cluster
253 for (Int_t i=0; i<nMax; ++i)
257 FindCluster(*cluster,localMax, maxPos[i]);
259 Timer(kMainLoop)->Start(kFALSE);
260 MainLoop(*cluster,iSimple);
261 Timer(kMainLoop)->Stop();
264 for (Int_t j=0; j<cluster->Multiplicity(); ++j)
266 AliMUONPad* pad = cluster->Pad(j);
267 if ( pad->Status() == 0 ) continue; // pad charge was not modified
269 pad->RevertCharge(); // use backup charge value
272 } // for (Int_t i=0; i<nMax;
273 if (nMax > 1) ((TH2D*) gROOT->FindObject("anode"))->Delete();
274 TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
275 if (mlem) mlem->Delete();
280 //_____________________________________________________________________________
282 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
284 /// Check if the pad and the pixel overlaps
286 // make a fake pad from the pixel
287 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
288 pixel.Coord(0),pixel.Coord(1),
289 pixel.Size(0),pixel.Size(1),0);
291 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
294 //_____________________________________________________________________________
296 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
298 /// Check precluster in order to attempt to simplify it (mostly for
299 /// two-cathode preclusters)
301 if (origCluster.Multiplicity()==1)
303 // Disregard one-pad clusters (leftovers from splitting)
307 Timer(kCheckPreCluster)->Start(kFALSE);
310 AliMUONCluster* cluster = static_cast<AliMUONCluster*>(origCluster.Clone());
314 AliDebug(2,"Start of CheckPreCluster=");
315 StdoutToAliDebug(2,cluster->Print("full"));
317 // Check if one-cathode precluster
318 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
319 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
321 AliMUONCluster* rv(0x0);
325 rv = CheckPreclusterTwoCathodes(cluster);
329 rv = CheckPreclusterOneCathode(cluster);
331 Timer(kCheckPreCluster)->Stop();
335 //_____________________________________________________________________________
337 AliMUONClusterFinderMLEM::CheckPreclusterOneCathode(AliMUONCluster* cluster)
339 /// Check single-cathode precluster
340 AliWarning("Reimplement me!");
341 AliDebug(2,"End of CheckPreClusterOneCathode=");
342 StdoutToAliDebug(2,cluster->Print("full"));
347 //_____________________________________________________________________________
349 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
351 // Check two-cathode cluster
353 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
354 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
356 Int_t npad = cluster->Multiplicity();
357 Int_t* flags = new Int_t[npad];
358 memset(flags,0,npad*sizeof(Int_t));
360 // Check pad overlaps
361 for ( Int_t i=0; i<npad; ++i)
363 AliMUONPad* padi = cluster->Pad(i);
364 if ( padi->Cathode() != i1 ) continue;
365 for (Int_t j=i+1; j<npad; ++j)
367 AliMUONPad* padj = cluster->Pad(j);
368 if ( padj->Cathode() != i2 ) continue;
369 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
370 flags[i] = flags[j] = 1; // mark overlapped pads
374 // Check if all pads overlap
376 for (Int_t i=0; i<npad; ++i)
378 if (flags[i]) continue;
384 // not all pads overlap.
385 for (Int_t i=0; i<npad; ++i)
387 AliMUONPad* pad = cluster->Pad(i);
388 if (flags[i]) continue;
389 Int_t cath = pad->Cathode();
390 Int_t cath1 = TMath::Even(cath);
391 // Check for edge effect (missing pads on the _other_ cathode)
392 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
393 if (!mpPad.IsValid()) continue;
394 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
395 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
396 cluster->RemovePad(pad);
397 fPreCluster->Pad(i)->Release();
402 // Check correlations of cathode charges
403 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry()*2 > 1 )
406 Int_t cathode = cluster->MaxRawChargeCathode();
412 // get min and max pad charges on the cathode opposite to the
413 // max pad (given by MaxRawChargeCathode())
415 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
417 AliMUONPad* pad = cluster->Pad(i);
418 if ( pad->Cathode() != cathode || !pad->IsReal() )
420 // only consider pads in the opposite cathode, and
421 // onyl consider real pads (i.e. exclude the virtual ones)
424 if ( pad->Charge() < cmin )
426 cmin = pad->Charge();
429 if ( pad->Charge() > cmax )
431 cmax = pad->Charge();
435 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
436 imin,imax,cmin,cmax));
438 // arrange pads according to their distance to the max, normalized
440 Double_t* dist = new Double_t[cluster->Multiplicity()];
445 AliMUONPad* padmax = cluster->Pad(imax);
447 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
450 if ( i == imax) continue;
451 AliMUONPad* pad = cluster->Pad(i);
452 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
453 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
454 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
455 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
458 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
464 TMath::Sort(cluster->Multiplicity(),dist,flags,kFALSE); // in ascending order
466 TObjArray toBeRemoved;
468 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
470 Int_t indx = flags[i];
471 AliMUONPad* pad = cluster->Pad(indx);
472 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
473 if ( dist[indx] > dmin )
475 // farther than the minimum pad
476 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
477 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
480 if (dx >= 0 && dy >= 0) continue;
481 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
482 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
484 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
487 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
489 cmax = TMath::Max(pad->Charge(),cmax);
493 cmax = pad->Charge();
496 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
497 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
500 toBeRemoved.AddLast(pad);
501 fPreCluster->Pad(indx)->Release();
504 for ( Int_t i = 0; i <= toBeRemoved.GetLast(); ++i )
506 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.At(i)));
513 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
514 StdoutToAliDebug(2,cluster->Print("full"));
519 //_____________________________________________________________________________
521 AliMUONClusterFinderMLEM::CheckOverlaps()
523 /// For debug only : check if some pixels overlap...
525 Int_t nPix = fPixArray->GetLast()+1;
528 for ( Int_t i = 0; i < nPix; ++i )
530 AliMUONPad* pixelI = Pixel(i);
531 AliMUONPad pi(dummy,dummy,dummy,dummy,
532 pixelI->Coord(0),pixelI->Coord(1),
533 pixelI->Size(0),pixelI->Size(1),0.0);
535 for ( Int_t j = i+1; j < nPix; ++j )
537 AliMUONPad* pixelJ = Pixel(j);
538 AliMUONPad pj(dummy,dummy,dummy,dummy,
539 pixelJ->Coord(0),pixelJ->Coord(1),
540 pixelJ->Size(0),pixelJ->Size(1),0.0);
543 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
545 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
546 StdoutToAliInfo(pixelI->Print();
547 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
549 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
550 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
551 cout << "-------" << endl;
559 //_____________________________________________________________________________
560 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
562 /// Build pixel array for MLEM method
564 Int_t npad = cluster.Multiplicity();
567 AliWarning("Got no pad at all ?!");
572 if ( !cluster.Multiplicity(0) || !cluster.Multiplicity(1) )
574 BuildPixArrayOneCathode(cluster);
578 BuildPixArrayTwoCathodes(cluster);
581 fPixArray->Sort(); // FIXME : not really needed, only to compare with ClusterFinderAZ
583 Int_t nPix = fPixArray->GetLast()+1;
585 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
587 Double_t xPadMin(1E9);
588 Double_t yPadMin(1E9);
590 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i )
592 AliMUONPad* pad = cluster.Pad(i);
593 xPadMin = TMath::Min (xPadMin, pad->DX());
594 yPadMin = TMath::Min (yPadMin, pad->DY());
600 for ( Int_t i = 0; i < nPix; ++i )
602 AliMUONPad* pixPtr = Pixel(i);
603 wxmin = TMath::Min(wxmin, pixPtr->Size(0));
604 wymin = TMath::Min(wymin, pixPtr->Size(1));
607 wxmin = TMath::Abs (wxmin - xPadMin/2) > 0.001 ? xPadMin : xPadMin / 2;
608 wymin = TMath::Abs (wymin - yPadMin/2) > 0.001 ? yPadMin : yPadMin / 2;
610 // Check if small pixel X-size
611 AdjustPixel(cluster,wxmin, 0);
612 // Check if small pixel Y-size
613 AdjustPixel(cluster,wymin, 1);
614 // Check if large pixel size
615 AdjustPixel(wxmin, wymin);
617 // Remove discarded pixels
618 for (Int_t i=0; i<nPix; ++i)
620 AliMUONPad* pixPtr = Pixel(i);
621 if (pixPtr->Charge() < 1)
623 AliDebug(2,Form("Removing pixel %d with charge<1 : ",i));
624 StdoutToAliDebug(2,pixPtr->Print());
629 fPixArray->Compress();
630 nPix = fPixArray->GetEntriesFast();
632 // AliDebug(2,Form("nPix after AdjustPixel=%d",nPix));
634 if ( nPix > cluster.Multiplicity() )
636 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
638 // Too many pixels - sort and remove pixels with the lowest signal
640 for ( Int_t i = cluster.Multiplicity(); i<nPix; ++i )
644 fPixArray->Compress();
645 nPix = fPixArray->GetEntriesFast();
646 } // if (nPix > npad)
648 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
649 // fPixArray->Print(););
650 CheckOverlaps();//FIXME : this is for debug only. Remove it.
653 //_____________________________________________________________________________
654 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
656 /// From a single-cathode cluster, build the pixel array
658 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
660 for ( Int_t j=0; j<cluster.Multiplicity(); ++j)
662 AliMUONPad* pad = cluster.Pad(j);
663 AliMUONPad* pixPtr = new AliMUONPad(pad->Position(),pad->Dimensions(),
665 fPixArray->Add(pixPtr);
669 //_____________________________________________________________________________
670 void AliMUONClusterFinderMLEM::BuildPixArrayTwoCathodes(AliMUONCluster& cluster)
672 /// From a two-cathodes cluster, build the pixel array
674 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
676 Int_t i1 = cluster.Pad(0)->Cathode();
677 Int_t i2 = TMath::Even(i1);
679 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i)
681 AliMUONPad* padi = cluster.Pad(i);
682 if (padi->Cathode() != i1) continue;
684 for ( Int_t j = 1; j < cluster.Multiplicity(); ++j)
686 AliMUONPad* padj = cluster.Pad(j);
687 if (padj->Cathode() != i2) continue;
691 if ( AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize,overlap) )
693 AliMUONPad* pixPtr = new AliMUONPad(overlap.Position(),
694 overlap.Dimensions(),
695 TMath::Min(padi->Charge(),padj->Charge()));
696 if ( ( padi->Charge() <= padj->Charge() && padi->IsSaturated() ) ||
697 ( padj->Charge() < padi->Charge() && padj->IsSaturated() ) )
699 // if smallest charge of the 2 pads is already above saturation, then
700 // the pixel is saturated...
701 pixPtr->SetSaturated(kTRUE);
703 pixPtr->SetReal(kFALSE);
704 fPixArray->Add(pixPtr);
710 //_____________________________________________________________________________
711 void AliMUONClusterFinderMLEM::AdjustPixel(AliMUONCluster& cluster,
712 Float_t width, Int_t ixy)
714 /// Check if some pixels have small size (adjust if necessary)
716 AliDebug(2,Form("width=%e ixy=%d",width,ixy));
718 AliMUONPad *pixPtr, *pixPtr1 = 0;
719 Int_t ixy1 = TMath::Even(ixy);
720 Int_t nPix = fPixArray->GetEntriesFast();
722 for (Int_t i=0; i<nPix; i++)
725 if (pixPtr->Charge() < 1) continue; // discarded pixel
726 if (pixPtr->Size(ixy)-width < -1.e-4)
729 for (Int_t j=i+1; j<nPix; j++)
732 if (pixPtr1->Charge() < 1) continue; // discarded pixel
733 if (TMath::Abs(pixPtr1->Size(ixy)-width) < fgkDistancePrecision) continue; // right size
734 if (TMath::Abs(pixPtr1->Coord(ixy1)-pixPtr->Coord(ixy1)) > fgkDistancePrecision) continue; // different rows/columns
735 if (TMath::Abs(pixPtr1->Coord(ixy)-pixPtr->Coord(ixy)) < 2*width)
738 Double_t tmp = pixPtr->Coord(ixy) + pixPtr1->Size(ixy)*
739 TMath::Sign (1., pixPtr1->Coord(ixy) - pixPtr->Coord(ixy));
740 pixPtr->SetCoord(ixy, tmp);
741 pixPtr->SetSize(ixy, width);
742 pixPtr->SetCharge(TMath::Min (pixPtr->Charge(),pixPtr1->Charge()));
743 pixPtr1->SetCharge(0);
747 } // for (Int_t j=i+1;
748 if (pixPtr1 || i == nPix-1)
750 // edge pixel - just increase its size
751 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
753 AliMUONPad* pad = cluster.Pad(j);
754 Double_t d = ( ixy == 0 ) ? pad->X() : ( ixy == 1 ) ? pad->Y() : -1E9;
756 if (pixPtr->Coord(ixy) < d)
758 pixPtr->Shift(ixy, pixPtr->Size(ixy)-width);
762 pixPtr->Shift(ixy, -pixPtr->Size(ixy)+width);
764 pixPtr->SetSize(ixy, width);
768 } // if (pixPtr->Size(ixy)-width < -1.e-4)
769 } // for (Int_t i=0; i<nPix;
772 //_____________________________________________________________________________
773 void AliMUONClusterFinderMLEM::AdjustPixel(Double_t wxmin, Double_t wymin)
775 /// Check if some pixels have large size (adjust if necessary)
777 AliDebug(2,Form("wxmin=%e wymin=%e",wxmin,wymin));
779 Int_t n1[2], n2[2], iOK = 1, nPix = fPixArray->GetEntriesFast();
780 AliMUONPad *pixPtr, pix;
781 Double_t xy0[2] = {9999, 9999}, wxy[2], dist[2] = {0};
783 // Check if large pixel size
784 for (Int_t i = 0; i < nPix; i++) {
785 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
786 if (pixPtr->Charge() < 1) continue; // discarded pixel
787 if (pixPtr->Size(0) - wxmin < 1.e-4) {
788 if (xy0[0] > 9998) xy0[0] = pixPtr->Coord(0); // position of a "normal" pixel
789 if (pixPtr->Size(1) - wymin < 1.e-4) {
790 if (xy0[1] > 9998) xy0[1] = pixPtr->Coord(1); // position of a "normal" pixel
792 } else iOK = 0; // large pixel
794 iOK = 0; // large pixel
795 if (xy0[1] > 9998 && pixPtr->Size(1) - wymin < 1.e-4) xy0[1] = pixPtr->Coord(1); // "normal" pixel
797 if (xy0[0] < 9998 && xy0[1] < 9998) break;
803 //cout << xy0[0] << " " << xy0[1] << endl;
804 for (Int_t i = 0; i < nPix; i++) {
805 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
806 if (pixPtr->Charge() < 1) continue; // discarded pixel
809 for (Int_t j = 0; j < 2; j++) {
810 if (pixPtr->Size(j) - wxy[j] < 1.e-4) continue;
811 dist[j] = (pixPtr->Coord(j) - xy0[j]) / wxy[j] / 2; // normalized distance to "normal" pixel
812 n2[j] = TMath::Nint (pixPtr->Size(j) / wxy[j]);
813 n1[j] = n2[j] == 1 ? TMath::Nint(dist[j]) : (Int_t)dist[j];
815 if (n1[0] > 998 && n1[1] > 998) continue;
816 if (fDebug) cout << " Different " << pixPtr->Size(0) << " " << wxy[0] << " "
817 << pixPtr->Size(1) << " " << wxy[1] <<endl;
819 if (n2[0] > 2 || n2[1] > 2) {
820 //cout << n2[0] << " " << n2[1] << endl;
821 if (n2[0] > 2 && n1[0] < 999) n1[0]--;
822 if (n2[1] > 2 && n1[1] < 999) n1[1]--;
824 //cout << n1[0] << " " << n2[0] << " " << n1[1] << " " << n2[1] << endl;
826 pix.SetSize(0, wxy[0]); pix.SetSize(1, wxy[1]);
828 for (Int_t ii = 0; ii < n2[0]; ii++) {
829 if (n1[0] < 999) pix.SetCoord(0, xy0[0] + (n1[0] + TMath::Sign(1.,dist[0]) * ii) * 2 * wxy[0]);
830 for (Int_t jj = 0; jj < n2[1]; jj++) {
831 if (n1[1] < 999) pix.SetCoord(1, xy0[1] + (n1[1] + TMath::Sign(1.,dist[1]) * jj) * 2 * wxy[1]);
832 fPixArray->Add(new AliMUONPad(pix));
836 pixPtr->SetCharge(0);
837 } // for (Int_t i = 0; i < nPix;
840 //_____________________________________________________________________________
842 AliMUONClusterFinderMLEM::Plot(const char* basename)
844 /// Make a plot and save it as png
848 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
853 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
854 fDetElemId,fClusterNumber));
857 //_____________________________________________________________________________
859 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
863 /// Compute coefficients needed for MLEM algorithm
865 Int_t nPix = fPixArray->GetLast()+1;
867 memset(probi,0,nPix*sizeof(Double_t));
869 for ( Int_t j=0; j<cluster.Multiplicity(); ++j )
871 AliMUONPad* pad = cluster.Pad(j);
874 for ( Int_t ipix=0; ipix<nPix; ++ipix )
876 Int_t indx1 = indx + ipix;
877 if (pad->Status() < 0)
882 AliMUONPad* pixPtr = Pixel(ipix);
883 // coef is the charge (given by Mathieson integral) on pad, assuming
884 // the Mathieson is center at pixel.
885 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
886 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
887 // pad->Ix(),pad->Iy(),
888 // pad->X(),pad->Y(),
889 // pad->DX(),pad->DY(),
890 // pixPtr->Coord(0),pixPtr->Coord(1),
891 // pixPtr->Size(0),pixPtr->Size(1),
894 probi[ipix] += coef[indx1];
899 //_____________________________________________________________________________
900 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
902 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
904 Int_t nPix = fPixArray->GetLast()+1;
906 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
907 StdoutToAliDebug(2,cluster.Print("full"););
911 AliDebug(1,"No pixels, why am I here then ?");
914 AddVirtualPad(cluster); // add virtual pads if necessary
916 Int_t npadTot = cluster.Multiplicity();
918 for (Int_t i = 0; i < npadTot; ++i)
920 if (cluster.Pad(i)->Status() == 0) ++npadOK;
925 Double_t* probi(0x0);
926 Int_t lc(0); // loop counter (for debug)
933 mlem = (TH2D*) gROOT->FindObject("mlem");
936 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
937 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
938 StdoutToAliDebug(2,fPixArray->Print("","full"));
940 coef = new Double_t [npadTot*nPix];
941 probi = new Double_t [nPix];
943 // Calculate coefficients and pixel visibilities
944 ComputeCoefficients(cluster,coef,probi);
946 for (Int_t ipix=0; ipix<nPix; ++ipix)
948 if (probi[ipix] < 0.01)
950 AliMUONPad* pixel = Pixel(ipix);
951 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
952 StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
953 pixel->SetCharge(0); // "invisible" pixel
958 Mlem(cluster,coef, probi, 15);
960 Double_t xylim[4] = {999, 999, 999, 999};
961 AliMUONPad* pixPtr(0x0);
963 for ( Int_t ipix=0; ipix<nPix; ++ipix )
965 pixPtr = Pixel(ipix);
966 for ( Int_t i=0; i<4; ++i )
968 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
971 for (Int_t i=0; i<4; i++)
973 xylim[i] -= pixPtr->Size(i/2);
977 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
978 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
980 // StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
981 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
982 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
983 xylim[0],-xylim[1],xylim[2],-xylim[3]
986 mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
988 for (Int_t ipix=0; ipix<nPix; ++ipix)
990 AliMUONPad* pixPtr = Pixel(ipix);
991 mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
994 // Check if the total charge of pixels is too low
996 for ( Int_t i=0; i<nPix; ++i)
998 qTot += Pixel(i)->Charge();
1001 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
1003 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
1008 fPixArray->Delete();
1009 for ( Int_t i=0; i<npadTot; ++i)
1011 AliMUONPad* pad = cluster.Pad(i);
1012 if ( pad->Status() == 0) pad->SetStatus(-1);
1019 // Simple cluster - skip further passes thru EM-procedure
1025 fPixArray->Delete();
1029 // Calculate position of the center-of-gravity around the maximum pixel
1031 FindCOG(mlem, xyCOG);
1033 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
1034 pixPtr->Size(0) > pixPtr->Size(1)) break;
1036 // Sort pixels according to the charge
1038 Double_t pixMin = 0.01*Pixel(0)->Charge();
1039 pixMin = TMath::Min(pixMin,50.);
1041 // Decrease pixel size and shift pixels to make them centered at
1043 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
1046 Double_t shift[2] = { 0.0, 0.0 };
1047 for (Int_t i=0; i<4; i++) xylim[i] = 999;
1050 for (Int_t ipix=0; ipix<nPix1; ++ipix)
1052 AliMUONPad* pixPtr = Pixel(ipix);
1053 if ( nPix >= npadOK // too many pixels already
1055 pixPtr->Charge() < pixMin // too low charge
1061 for (Int_t i=0; i<2; ++i)
1065 pixPtr->SetCharge(10);
1066 pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
1067 width = -pixPtr->Size(indx);
1068 pixPtr->Shift(indx, width);
1069 // Shift pixel position
1073 for (Int_t j=0; j<2; ++j)
1075 shift[j] = pixPtr->Coord(j) - xyCOG[j];
1076 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
1079 pixPtr->Shift(0, -shift[0]);
1080 pixPtr->Shift(1, -shift[1]);
1084 pixPtr = new AliMUONPad(*pixPtr);
1085 pixPtr->Shift(indx, -2*width);
1086 fPixArray->Add(pixPtr);
1088 for (Int_t i=0; i<4; ++i)
1090 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1092 } // for (Int_t i=0; i<2;
1094 } // for (Int_t ipix=0;
1096 fPixArray->Compress();
1097 nPix = fPixArray->GetEntriesFast();
1099 AliDebug(2,Form("After shift:"));
1100 StdoutToAliDebug(2,fPixArray->Print("","full"););
1101 Plot(Form("mlem.lc%d",lc+1));
1103 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1106 xylim[2],xylim[3]));
1108 // Remove excessive pixels
1111 for (Int_t ipix=npadOK; ipix<nPix; ++ipix)
1118 AliMUONPad* pixPtr = Pixel(0);
1119 // add pixels if the maximum is at the limit of pixel area
1120 // start from Y-direction
1122 for (Int_t i=3; i>-1; --i)
1124 if (nPix < npadOK &&
1125 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2))
1127 AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1128 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1129 j = TMath::Even (i/2);
1130 p->SetCoord(j, xyCOG[j]);
1131 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1132 StdoutToAliDebug(2,cout << " ---- ";
1133 p->Print("corners"););
1139 fPixArray->Compress();
1140 nPix = fPixArray->GetEntriesFast();
1147 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1148 StdoutToAliDebug(2,fPixArray->Print("","full"););
1150 // remove pixels with low signal or low visibility
1151 // Cuts are empirical !!!
1152 Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1153 thresh = TMath::Min (thresh,50.);
1155 Double_t charge = 0;
1157 for ( Int_t i=0; i<nPix; ++i)
1159 cmax = TMath::Max (cmax,probi[i]);
1162 // Mark pixels which should be removed
1163 for (Int_t i=0; i<nPix; ++i)
1165 AliMUONPad* pixPtr = Pixel(i);
1166 charge = pixPtr->Charge();
1167 if (charge < thresh)
1169 pixPtr->SetCharge(-charge);
1173 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1175 for (Int_t i=0; i<nPix; ++i)
1177 AliMUONPad* pixPtr = Pixel(i);
1178 charge = pixPtr->Charge();
1179 if (charge > 0) continue;
1180 near = FindNearest(pixPtr);
1181 pixPtr->SetCharge(0);
1182 probi[i] = 0; // make it "invisible"
1183 AliMUONPad* pnear = Pixel(near);
1184 pnear->SetCharge(pnear->Charge() + (-charge));
1186 Mlem(cluster,coef,probi,2);
1188 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1189 StdoutToAliDebug(2,fPixArray->Print("","full"););
1190 Plot("mlem.beforesplit");
1193 for (Int_t i=0; i<nPix; ++i)
1195 AliMUONPad* pixPtr = Pixel(i);
1196 Int_t ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1197 Int_t iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1198 mlem->SetBinContent(ix, iy, pixPtr->Charge());
1201 // Try to split into clusters
1203 if (mlem->GetSum() < 1)
1209 fSplitter->Split(cluster,mlem,coef,fClusterList);
1216 fPixArray->Delete();
1221 //_____________________________________________________________________________
1222 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1223 Double_t* coef, Double_t* probi,
1226 /// Use MLEM to find pixel charges
1228 Int_t nPix = fPixArray->GetEntriesFast();
1230 Int_t npad = cluster.Multiplicity();
1232 Double_t* probi1 = new Double_t[nPix];
1233 Double_t probMax = 0;
1234 Double_t tmp = TMath::MaxElement(nPix,probi);
1236 for (Int_t ipix=0; ipix<nPix; ++ipix)
1238 probMax = TMath::Max(probMax,probi[ipix]);
1241 if (probMax!=tmp) { AliWarning(Form("probMax=%e tmp=%e",probMax,tmp)); }
1243 for (Int_t iter=0; iter<nIter; ++iter)
1246 for (Int_t ipix=0; ipix<nPix; ++ipix)
1248 // Correct each pixel
1249 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1251 probi1[ipix] = probMax;
1252 for (Int_t j=0; j<npad; j++)
1254 AliMUONPad* pad = cluster.Pad(j);
1255 if (pad->Status() < 0) continue;
1257 Int_t indx1 = j*nPix;
1258 Int_t indx = indx1 + ipix;
1259 // Calculate expectation
1260 for (Int_t i=0; i<nPix; i++)
1262 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1264 if ( pad->Charge() > fgkSaturation-1 && sum1 > pad->Charge() ) //FIXME : remove usage of fgkSaturation
1266 if ( !pad->IsSaturated() )
1268 AliWarning("Got a pad charge above saturation not backed-up by pad->IsSaturated() function : ");
1269 StdoutToAliWarning(pad->Print("full"));
1271 // correct for pad charge overflows
1272 probi1[ipix] -= coef[indx];
1276 if (coef[indx] > 1.e-6)
1278 sum += pad->Charge()*coef[indx]/sum1;
1280 } // for (Int_t j=0;
1281 AliMUONPad* pixPtr = Pixel(ipix);
1282 if (probi1[ipix] > 1.e-6)
1284 pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1286 } // for (Int_t ipix=0;
1287 } // for (Int_t iter=0;
1291 //_____________________________________________________________________________
1292 void AliMUONClusterFinderMLEM::FindCOG(TH2D *mlem, Double_t *xyc)
1294 /// Calculate position of the center-of-gravity around the maximum pixel
1296 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1297 Int_t i1 = -9, j1 = -9;
1298 mlem->GetMaximumBin(ixmax,iymax,ix);
1299 Int_t nx = mlem->GetNbinsX();
1300 Int_t ny = mlem->GetNbinsY();
1301 Double_t thresh = mlem->GetMaximum()/10;
1302 Double_t x, y, cont, xq=0, yq=0, qq=0;
1304 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1305 y = mlem->GetYaxis()->GetBinCenter(i);
1306 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1307 cont = mlem->GetCellContent(j,i);
1308 if (cont < thresh) continue;
1309 if (i != i1) {i1 = i; nsumy++;}
1310 if (j != j1) {j1 = j; nsumx++;}
1311 x = mlem->GetXaxis()->GetBinCenter(j);
1320 Int_t i2 = 0, j2 = 0;
1323 // one bin in Y - add one more (with the largest signal)
1324 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1325 if (i == iymax) continue;
1326 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1327 cont = mlem->GetCellContent(j,i);
1330 x = mlem->GetXaxis()->GetBinCenter(j);
1331 y = mlem->GetYaxis()->GetBinCenter(i);
1340 if (i2 != i1) nsumy++;
1341 if (j2 != j1) nsumx++;
1343 } // if (nsumy == 1)
1346 // one bin in X - add one more (with the largest signal)
1348 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1349 if (j == ixmax) continue;
1350 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1351 cont = mlem->GetCellContent(j,i);
1354 x = mlem->GetXaxis()->GetBinCenter(j);
1355 y = mlem->GetYaxis()->GetBinCenter(i);
1364 if (i2 != i1) nsumy++;
1365 if (j2 != j1) nsumx++;
1367 } // if (nsumx == 1)
1369 xyc[0] = xq/qq; xyc[1] = yq/qq;
1370 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1373 //_____________________________________________________________________________
1374 Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1376 /// Find the pixel nearest to the given one
1377 /// (algorithm may be not very efficient)
1379 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1380 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1381 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1384 for (Int_t i=0; i<nPix; i++) {
1385 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1386 if (pixPtr->Charge() < 0.5) continue;
1387 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1388 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1389 r = dx *dx + dy * dy;
1390 if (r < rmin) { rmin = r; imin = i; }
1396 //_____________________________________________________________________________
1398 AliMUONClusterFinderMLEM::Timer(Int_t i) const
1400 /// Return timer at index i
1401 return static_cast<TStopwatch*>(fTimers->At(i));
1404 //_____________________________________________________________________________
1406 AliMUONClusterFinderMLEM::Paint(Option_t*)
1408 /// Paint cluster and pixels
1410 AliMpArea area(fPreCluster->Area());
1412 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1414 gVirtualX->SetFillStyle(1001);
1415 gVirtualX->SetFillColor(3);
1416 gVirtualX->SetLineColor(3);
1420 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1422 AliMUONPad* pixel = Pixel(i);
1424 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1425 pixel->Coord(1)-pixel->Size(1)*s,
1426 pixel->Coord(0)+pixel->Size(0)*s,
1427 pixel->Coord(1)+pixel->Size(1)*s);
1429 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1431 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1432 // pixel->Coord(1) + sign*pixel->Size(1),
1433 // pixel->Coord(0) + pixel->Size(0),
1434 // pixel->Coord(1) - sign*pixel->Size(1)
1440 gVirtualX->SetFillStyle(0);
1442 fPreCluster->Paint();
1444 gVirtualX->SetLineColor(1);
1445 gVirtualX->SetLineWidth(2);
1446 gVirtualX->SetFillStyle(0);
1447 gVirtualX->SetTextColor(1);
1448 gVirtualX->SetTextAlign(22);
1450 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1452 AliMUONPad* pixel = Pixel(i);
1453 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1454 pixel->Coord(1)-pixel->Size(1),
1455 pixel->Coord(0)+pixel->Size(0),
1456 pixel->Coord(1)+pixel->Size(1));
1457 gVirtualX->SetTextSize(pixel->Size(0)*60);
1459 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1463 //_____________________________________________________________________________
1464 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1466 /// Find local maxima in pixel space for large preclusters in order to
1467 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1468 /// or to find additional fitting seeds if clusters were not completely resolved
1470 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1473 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
1474 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
1475 //if (hist) hist->Delete();
1477 Double_t xylim[4] = {999, 999, 999, 999};
1479 Int_t nPix = pixArray->GetEntriesFast();
1480 AliMUONPad *pixPtr = 0;
1481 for (Int_t ipix=0; ipix<nPix; ipix++) {
1482 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1483 for (Int_t i=0; i<4; i++)
1484 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1486 for (Int_t i=0; i<4; i++) xylim[i] -= pixPtr->Size(i/2);
1488 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1489 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1490 if (pixArray == fPixArray) hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1491 else hist = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1492 for (Int_t ipix=0; ipix<nPix; ipix++) {
1493 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1494 hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1496 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1498 Int_t nMax = 0, indx;
1499 Int_t *isLocalMax = new Int_t[ny*nx];
1500 for (Int_t i=0; i<ny*nx; i++) isLocalMax[i] = 0;
1502 for (Int_t i=1; i<=ny; i++) {
1504 for (Int_t j=1; j<=nx; j++) {
1505 if (hist->GetCellContent(j,i) < 0.5) continue;
1506 //if (isLocalMax[indx+j-1] < 0) continue;
1507 if (isLocalMax[indx+j-1] != 0) continue;
1508 FlagLocalMax(hist, i, j, isLocalMax);
1512 for (Int_t i=1; i<=ny; i++) {
1514 for (Int_t j=1; j<=nx; j++) {
1515 if (isLocalMax[indx+j-1] > 0) {
1516 localMax[nMax] = indx + j - 1;
1517 maxVal[nMax++] = hist->GetCellContent(j,i);
1518 if (nMax > 99) AliFatal(" Too many local maxima !!!");
1522 if (fDebug) cout << " Local max: " << nMax << endl;
1523 delete [] isLocalMax; isLocalMax = 0;
1527 //_____________________________________________________________________________
1528 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1530 /// Flag pixels (whether or not local maxima)
1532 Int_t nx = hist->GetNbinsX();
1533 Int_t ny = hist->GetNbinsY();
1534 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1535 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1537 for (Int_t i1=i-1; i1<i+2; i1++) {
1538 if (i1 < 1 || i1 > ny) continue;
1539 indx1 = (i1 - 1) * nx;
1540 for (Int_t j1=j-1; j1<j+2; j1++) {
1541 if (j1 < 1 || j1 > nx) continue;
1542 if (i == i1 && j == j1) continue;
1543 indx2 = indx1 + j1 - 1;
1544 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1545 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1546 else if (cont > cont1) isLocalMax[indx2] = -1;
1547 else { // the same charge
1548 isLocalMax[indx] = 1;
1549 if (isLocalMax[indx2] == 0) {
1550 FlagLocalMax(hist, i1, j1, isLocalMax);
1551 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1552 else isLocalMax[indx2] = -1;
1557 isLocalMax[indx] = 1; // local maximum
1560 //_____________________________________________________________________________
1561 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1562 Int_t *localMax, Int_t iMax)
1564 /// Find pixel cluster around local maximum \a iMax and pick up pads
1565 /// overlapping with it
1567 TH2D *hist = (TH2D*) gROOT->FindObject("anode");
1568 Int_t nx = hist->GetNbinsX();
1569 Int_t ny = hist->GetNbinsY();
1570 Int_t ic = localMax[iMax] / nx + 1;
1571 Int_t jc = localMax[iMax] % nx + 1;
1572 Bool_t *used = new Bool_t[ny*nx];
1573 for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE;
1575 // Drop all pixels from the array - pick up only the ones from the cluster
1576 fPixArray->Delete();
1578 Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2;
1579 Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;
1580 Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
1581 Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
1582 Double_t cont = hist->GetCellContent(jc,ic);
1583 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1584 used[(ic-1)*nx+jc-1] = kTRUE;
1585 fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1587 Int_t nPix = fPixArray->GetEntriesFast();
1588 Int_t npad = cluster.Multiplicity();
1590 for (Int_t i=0; i<nPix; ++i)
1592 AliMUONPad* pixPtr = Pixel(i);
1593 pixPtr->SetSize(0,wx);
1594 pixPtr->SetSize(1,wy);
1597 // Pick up pads which overlap with found pixels
1598 for (Int_t i=0; i<npad; i++)
1600 cluster.Pad(i)->SetStatus(-1);
1603 for (Int_t i=0; i<nPix; i++)
1605 AliMUONPad* pixPtr = Pixel(i);
1606 for (Int_t j=0; j<npad; ++j)
1608 AliMUONPad* pad = cluster.Pad(j);
1609 if ( Overlap(*pad,*pixPtr) )
1616 delete [] used; used = 0;
1619 //_____________________________________________________________________________
1620 AliMUONClusterFinderMLEM&
1621 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1623 /// Protected assignement operator
1625 if (this == &rhs) return *this;
1627 AliFatal("Not implemented.");
1632 //_____________________________________________________________________________
1634 AliMUONClusterFinderMLEM::Neighbours(Int_t cathode, Int_t ix, Int_t iy,
1635 Int_t& n, Int_t* xList, Int_t* yList)
1637 /// Get the list of neighbours of pad at (cathode,ix,iy)
1640 const AliMpVSegmentation* seg = fSegmentation[cathode];
1642 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1644 // Define the region to look into : a region slightly bigger
1645 // than the pad itself (5% bigger), in order to catch first neighbours.
1647 AliMpArea area(pad.Position(),pad.Dimensions()*1.05);
1649 AliMpVPadIterator* it = seg->CreateIterator(area);
1651 while ( !it->IsDone() && n < 10 )
1653 AliMpPad p = it->CurrentItem();
1654 if ( p != pad ) // skip self
1656 xList[n] = p.GetIndices().GetFirst();
1657 yList[n] = p.GetIndices().GetSecond();
1665 //_____________________________________________________________________________
1666 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1668 /// Add virtual pad (with small charge) to improve fit for some
1669 /// clusters (when pad with max charge is at the extreme of the cluster)
1671 // Get number of pads in X and Y-directions
1672 Int_t nInX = -1, nInY;
1673 PadsInXandY(cluster,nInX, nInY);
1675 if (fDebug) cout << " Chamber: " << fDetElemId / 100 - 1 << " " << nInX << " " << nInY << endl;
1677 // Add virtual pad only if number of pads per direction == 2
1678 if (nInX != 2 && nInY != 2) return;
1682 // Find pads with max charge
1683 Int_t maxpad[2][2] = {{-1, -1}, {-1, -1}}, cath;
1684 Double_t sigmax[2] = {0}, aamax[2] = {0};
1685 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
1687 AliMUONPad* pad = cluster.Pad(j);
1688 if (pad->Status() != 0) continue;
1689 cath = pad->Cathode();
1690 if (pad->Charge() > sigmax[cath])
1692 maxpad[cath][1] = maxpad[cath][0];
1693 aamax[cath] = sigmax[cath];
1694 sigmax[cath] = pad->Charge();
1695 maxpad[cath][0] = j;
1699 if (maxpad[0][0] >= 0 && maxpad[0][1] < 0 || maxpad[1][0] >= 0 && maxpad[1][1] < 0)
1701 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
1703 AliMUONPad* pad = cluster.Pad(j);
1704 if (pad->Status() != 0) continue;
1705 cath = pad->Cathode();
1706 if (j == maxpad[cath][0] || j == maxpad[cath][1]) continue;
1707 if ( pad->Charge() > aamax[cath])
1709 aamax[cath] = pad->Charge();
1710 maxpad[cath][1] = j;
1715 // cout << "-------AddVirtualPad" << endl;
1716 // cout << Form("nInX %2d nInY %2d",nInX,nInY) << endl;
1718 // cluster.Print("full");
1720 // for ( Int_t i = 0; i < 2; ++i )
1722 // for ( Int_t j = 0; j < 2; ++j )
1724 // cout << Form("maxpad[%d][%d]=%d",i,j,maxpad[i][j]) << endl;
1728 // Check for mirrors (side X on cathode 0)
1729 Bool_t mirror = kFALSE;
1730 if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1732 AliMUONPad* maxPadCath[] = { cluster.Pad(maxpad[0][0]), cluster.Pad(maxpad[1][0]) };
1733 mirror = maxPadCath[0]->DX() < maxPadCath[0]->DY();
1734 if (!mirror && TMath::Abs( maxPadCath[0]->DX() - maxPadCath[1]->DX()) < 0.001)
1736 // Special case when pads on both cathodes have the same size
1738 for (Int_t j = 0; j < cluster.Multiplicity(); ++j)
1740 AliMUONPad* pad = cluster.Pad(j);
1741 cath = pad->Cathode();
1742 if (j == maxpad[cath][0]) continue;
1743 if ( pad->Ix() != maxPadCath[cath]->Ix() ) continue;
1744 if ( TMath::Abs(pad->Iy() - maxPadCath[cath]->Iy()) == 1 )
1749 if (!yud[0]) mirror = kTRUE; // take the other cathode
1750 } // if (!mirror &&...
1751 } // if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1753 // // Find neughbours of pads with max charges
1754 Int_t nn, xList[10], yList[10], ix0, iy0, ix, iy, neighb;
1755 for (cath=0; cath<2; cath++)
1757 if (!cath && maxpad[0][0] < 0) continue; // one-sided cluster - cathode 1
1758 if (cath && maxpad[1][0] < 0) break; // one-sided cluster - cathode 0
1759 if (maxpad[1][0] >= 0)
1763 if (!cath && nInY != 2) continue;
1764 if (cath && nInX != 2 && (maxpad[0][0] >= 0 || nInY != 2)) continue;
1768 if (!cath && nInX != 2) continue;
1769 if (cath && nInY != 2 && (maxpad[0][0] >= 0 || nInX != 2)) continue;
1773 Int_t iAddX = 0, iAddY = 0, ix1 = 0, iy1 = 0, iPad = 0;
1774 if (maxpad[0][0] < 0) iPad = 1;
1776 for (iPad=0; iPad<2; iPad++)
1778 if (maxpad[cath][iPad] < 0) continue;
1779 if (iPad && !iAddX && !iAddY) break;
1780 if (iPad && cluster.Pad(maxpad[cath][1])->Charge() / sigmax[cath] < 0.5) break;
1782 Int_t neighbx = 0, neighby = 0;
1783 ix0 = cluster.Pad(maxpad[cath][iPad])->Ix();
1784 iy0 = cluster.Pad(maxpad[cath][iPad])->Iy();
1785 Neighbours(cath,ix0,iy0,nn,xList,yList);
1787 for (Int_t j=0; j<nn; j++) {
1788 if (TMath::Abs(xList[j]-ix0) == 1 || xList[j]*ix0 == -1) neighbx++;
1789 if (TMath::Abs(yList[j]-iy0) == 1 || yList[j]*iy0 == -1) neighby++;
1792 if (cath) neighb = neighbx;
1793 else neighb = neighby;
1794 if (maxpad[0][0] < 0) neighb += neighby;
1795 else if (maxpad[1][0] < 0) neighb += neighbx;
1797 if (!cath) neighb = neighbx;
1798 else neighb = neighby;
1799 if (maxpad[0][0] < 0) neighb += neighbx;
1800 else if (maxpad[1][0] < 0) neighb += neighby;
1803 for (Int_t j=0; j< cluster.Multiplicity(); ++j)
1805 AliMUONPad* pad = cluster.Pad(j);
1806 if ( pad->Cathode() != cath) continue;
1809 if (iy == iy0 && ix == ix0) continue;
1810 for (Int_t k=0; k<nn; ++k)
1812 if (xList[k] != ix || yList[k] != iy) continue;
1815 if ((!cath || maxpad[0][0] < 0) &&
1816 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
1817 if (!iPad && TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) ix1 = xList[k]; //19-12-05
1818 xList[k] = yList[k] = 0;
1822 if ((cath || maxpad[1][0] < 0) &&
1823 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
1824 if (!iPad) ix1 = xList[k]; //19-12-05
1825 xList[k] = yList[k] = 0;
1829 if ((!cath || maxpad[0][0] < 0) &&
1830 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
1831 if (!iPad) ix1 = xList[k]; //19-12-05
1832 xList[k] = yList[k] = 0;
1836 if ((cath || maxpad[1][0] < 0) &&
1837 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
1838 xList[k] = yList[k] = 0;
1843 } // for (Int_t k=0; k<nn;
1845 } // for (Int_t j=0; j< cluster.Multiplicity();
1846 if (!neighb) continue;
1851 for (Int_t j=0; j<nn; j++)
1853 if (xList[j] == 0 && yList[j] == 0) continue;
1854 // npads = fnPads[0] + fnPads[1];
1855 // fPadIJ[0][npads] = cath;
1856 // fPadIJ[1][npads] = 0;
1859 if (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) {
1860 if (iy != iy0) continue; // new segmentation - check
1861 if (nInX != 2) continue; // new
1863 if (!cath && maxpad[1][0] >= 0) continue;
1865 if (cath && maxpad[0][0] >= 0) continue;
1867 if (iPad && !iAddX) continue;
1868 AliMpPad pad = fSegmentation[cath]->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1869 // fXyq[0][npads] = pad.Position().X();
1870 // fXyq[1][npads] = pad.Position().Y();
1871 AliMUONPad muonPad(fDetElemId, cath, ix, iy, pad.Position().X(), pad.Position().Y(), 0, 0, 0);
1872 // fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
1873 // if (fXyq[0][npads] > 1.e+5) continue; // temporary fix
1874 if (muonPad.Coord(0) > 1.e+5) continue; // temporary fix
1875 if (ix == ix1) continue; //19-12-05
1876 if (ix1 == ix0) continue;
1877 if (maxpad[1][0] < 0 || mirror && maxpad[0][0] >= 0) {
1878 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/100, 5.);
1879 // else fXyq[2][npads] = TMath::Min (aamax[0]/100, 5.);
1880 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[0]/100, 5.));
1881 else muonPad.SetCharge(TMath::Min (aamax[0]/100, 5.));
1884 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/100, 5.);
1885 // else fXyq[2][npads] = TMath::Min (aamax[1]/100, 5.);
1886 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[1]/100, 5.));
1887 else muonPad.SetCharge(TMath::Min (aamax[1]/100, 5.));
1889 // fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
1890 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1891 // fXyq[3][npads] = -2; // flag
1892 // fPadIJ[2][npads] = ix;
1893 // fPadIJ[3][npads] = iy;
1894 muonPad.SetSize(0,-2.); //flag
1898 //AliDebug(1,Form("Add virtual pad in X %f %f %f %3d %3d \n",
1899 // fXyq[2][npads], fXyq[0][npads], fXyq[1][npads], ix, iy));
1900 //muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy));
1901 if (fDebug) printf(" ***** Add virtual pad in X ***** %f %f %f %3d %3d \n",
1902 muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy);
1903 cluster.AddPad(muonPad); // add pad to the cluster
1907 if (nInY != 2) continue;
1908 if (!mirror && cath && maxpad[0][0] >= 0) continue;
1909 if (mirror && !cath && maxpad[1][0] >= 0) continue;
1910 if (TMath::Abs(iy-iy0) == 1 || TMath::Abs(iy*iy0) == 1) {
1911 if (ix != ix0) continue; // new segmentation - check
1912 if (iPad && !iAddY) continue;
1913 AliMpPad pad = fSegmentation[cath]->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1914 // fXyq[0][npads] = pad.Position().X();
1915 // fXyq[1][npads] = pad.Position().Y();
1916 // fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
1917 AliMUONPad muonPad(fDetElemId, cath, ix, iy, pad.Position().X(), pad.Position().Y(), 0, 0, 0);
1918 if (iy1 == iy0) continue;
1919 //if (iPad && iy1 == iy0) continue;
1920 if (maxpad[0][0] < 0 || mirror && maxpad[1][0] >= 0) {
1921 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/15, fgkZeroSuppression);
1922 // else fXyq[2][npads] = TMath::Min (aamax[1]/15, fgkZeroSuppression);
1923 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[1]/15, fgkZeroSuppression));
1924 else muonPad.SetCharge(TMath::Min (aamax[1]/15, fgkZeroSuppression));
1927 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/15, fgkZeroSuppression);
1928 // else fXyq[2][npads] = TMath::Min (aamax[0]/15, fgkZeroSuppression);
1929 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[0]/15, fgkZeroSuppression));
1930 else muonPad.SetCharge(TMath::Min (aamax[0]/15, fgkZeroSuppression));
1932 // fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
1933 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1934 // fXyq[3][npads] = -2; // flag
1935 // fPadIJ[2][npads] = ix;
1936 // fPadIJ[3][npads] = iy;
1937 muonPad.SetSize(0,-2.); //flag
1941 if (fDebug) printf(" ***** Add virtual pad in Y ***** %f %f %f %3d %3d \n",
1942 muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy);
1943 cluster.AddPad(muonPad); // add pad to the cluster
1946 } // for (Int_t j=0; j<nn;
1947 } // for (Int_t iPad=0;
1948 } // for (cath=0; cath<2;
1951 //_____________________________________________________________________________
1952 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1953 Int_t &nInX, Int_t &nInY) const
1955 /// Find number of pads in X and Y-directions (excluding virtual ones and
1958 Int_t statusToTest = 1;
1960 if ( nInX < 0 ) statusToTest = 0;
1962 Bool_t mustMatch(kTRUE);
1964 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1966 nInX = cn.GetFirst();
1967 nInY = cn.GetSecond();
1970 //_____________________________________________________________________________
1971 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1973 /// Remove pixel at index i
1974 AliMUONPad* pixPtr = Pixel(i);
1975 fPixArray->RemoveAt(i);
1979 //_____________________________________________________________________________
1980 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1982 /// Process simple cluster (small number of pads) without EM-procedure
1984 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1985 Double_t parOk[3] = {0.};
1986 TObjArray *clusters[1];
1987 clusters[0] = fPixArray;
1989 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1991 for (Int_t i = 0; i < cluster.Multiplicity(); ++i)
1993 AliMUONPad* pad = cluster.Pad(i);
1994 if ( pad->Charge() > fgkSaturation-1) //FIXME : remove usage of fgkSaturation
2003 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList);
2006 //_____________________________________________________________________________
2008 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
2010 /// Returns pixel at index i
2011 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
2014 //_____________________________________________________________________________
2016 AliMUONClusterFinderMLEM::Print(Option_t* what) const
2019 TString swhat(what);
2021 if ( swhat.Contains("precluster") )
2023 if ( fPreCluster) fPreCluster->Print();