]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONClusterFinderMLEM.cxx
Disabling cluster sorting - breaks the MLEM algorithm.
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderMLEM.cxx
CommitLineData
c0a16418 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/* $Id$ */
17
18/// \class AliMUONClusterFinderMLEM
19///
20/// Clusterizer class based on the Expectation-Maximization algorithm
21///
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).
25///
26/// \author Laurent Aphecetche (for the "new" C++ structure) and
27/// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
28
9edefa04 29#include "AliMUONClusterFinderMLEM.h"
c0a16418 30#include "AliLog.h"
31#include "AliMUONCluster.h"
32#include "AliMUONClusterSplitterMLEM.h"
3887d11e 33#include "AliMUONVDigit.h"
c0a16418 34#include "AliMUONPad.h"
35#include "AliMUONPreClusterFinder.h"
36#include "AliMpPad.h"
37#include "AliMpVPadIterator.h"
38#include "AliMpVSegmentation.h"
39#include "AliRunLoader.h"
3887d11e 40#include "AliMUONVDigitStore.h"
a6435fb4 41#include <Riostream.h>
42#include <TH2.h>
43#include <TMinuit.h>
44#include <TCanvas.h>
45#include <TStopwatch.h>
46#include <TMath.h>
a1e17193 47#include <TROOT.h>
a6435fb4 48
c0a16418 49/// \cond CLASSIMP
50ClassImp(AliMUONClusterFinderMLEM)
51/// \endcond
52
53const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
54const Double_t AliMUONClusterFinderMLEM::fgkSaturation = 3000; // average saturation level
55const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
56const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
57const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
58
59 TMinuit* AliMUONClusterFinderMLEM::fgMinuit = 0x0;
60
61//_____________________________________________________________________________
62AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot)
63 : AliMUONVClusterFinder(),
64fPreClusterFinder(new AliMUONPreClusterFinder),
65fPreCluster(0x0),
66fClusterList(),
67fEventNumber(0),
68fDetElemId(-1),
69fClusterNumber(0),
70fZpad(0.0),
71fReco(1),
72fCathBeg(0),
73fPixArray(new TObjArray(20)),
74fDebug(1),
75fPlot(plot),
76fTimers(new TObjArray(kLast)),
77fSplitter(0x0),
78fNClusters(0),
79fNAddVirtualPads(0)
80{
81 /// Constructor
82
83 fSegmentation[1] = fSegmentation[0] = 0x0;
84
85 fPadBeg[0] = fPadBeg[1] = fCathBeg;
86
87 if (!fgMinuit) fgMinuit = new TMinuit(8);
88
89 fTimers->SetOwner(kTRUE);
90
91 for ( Int_t i = 0; i < kLast; ++i )
92 {
93 TStopwatch* t = new TStopwatch;
94 fTimers->AddLast(new TStopwatch);
95 t->Start(kTRUE);
96 t->Stop();
97 }
98}
99
100//_____________________________________________________________________________
101AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
102{
103/// Destructor
104 delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
105// delete fDraw;
106 delete fPreClusterFinder;
107 for ( Int_t i = 0; i < kLast; ++i )
108 {
109 AliInfo(Form("Timer %d",i));
110 Timer(i)->Print();
111 }
112 delete fTimers;
113 delete fSplitter;
114 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
115 fNClusters,fNAddVirtualPads));
116}
117
118//_____________________________________________________________________________
119Bool_t
120AliMUONClusterFinderMLEM::Prepare(const AliMpVSegmentation* segmentations[2],
3887d11e 121 const AliMUONVDigitStore& digitStore)
c0a16418 122{
123 /// Prepare for clustering
124
125 for ( Int_t i = 0; i < 2; ++i )
126 {
127 fSegmentation[i] = segmentations[i];
128 }
129
130 // Find out the DetElemId
131 fDetElemId = -1;
132
3887d11e 133 TIter next(digitStore.CreateIterator());
134 AliMUONVDigit* d = static_cast<AliMUONVDigit*>(next());
135
136 if (d)
c0a16418 137 {
3887d11e 138 fDetElemId = d->DetElemId();
c0a16418 139 }
140
141 if ( fDetElemId < 0 )
142 {
143 AliWarning("Could not find DE. Probably no digits at all ?");
144 return kFALSE;
145 }
146
147 delete fSplitter;
148 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
149
150 // find out current event number, and reset the cluster number
151 fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
152 fClusterNumber = -1;
153 fClusterList.Delete();
154
155// AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
156
3887d11e 157 return fPreClusterFinder->Prepare(segmentations,digitStore);
c0a16418 158}
159
160//_____________________________________________________________________________
161AliMUONCluster*
162AliMUONClusterFinderMLEM::NextCluster()
163{
164 /// Return next cluster
165
166 ++fClusterNumber;
167
168 // if the list of clusters is not void, pick one from there
169 if ( fClusterList.GetLast() >= 0 )
170 {
171 TObject* o = fClusterList.At(0);
172 fClusterList.RemoveAt(0);
173 return static_cast<AliMUONCluster*>(o);
174 }
175
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
179
180 // if the cluster list is exhausted, we need to go to the next
181 // pre-cluster and treat it
182
183 fPreCluster = fPreClusterFinder->NextCluster();
184
185 if (!fPreCluster)
186 {
187 // we are done
188 return 0x0;
189 }
190
191 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
192
193 WorkOnPreCluster();
194
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...
197
198 for ( Int_t i = 0; i < fPreCluster->Multiplicity(); ++i )
199 {
200 AliMUONPad* pad = fPreCluster->Pad(i);
201 if ( !pad->IsUsed() )
202 {
203 fPreClusterFinder->UsePad(*pad);
204 }
205 }
206
207 return NextCluster();
208}
209
210//_____________________________________________________________________________
211Bool_t
212AliMUONClusterFinderMLEM::WorkOnPreCluster()
213{
214 /// Starting from a precluster, builds a pixel array, and then
215 /// extract clusters from this array
216
217 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
218
219 if (!cluster) return kFALSE;
220
221 BuildPixArray(*cluster);
222
223 if ( fPixArray->GetLast() < 0 )
224 {
225 AliDebug(1,"No pixel for the above cluster");
226 delete cluster;
227 return kFALSE;
228 }
229
230 // Use MLEM for cluster finder
231 Int_t nMax = 1, localMax[100], maxPos[100];
232 Double_t maxVal[100];
233
234 if (cluster->Multiplicity() > 50)
235 {
236 nMax = FindLocalMaxima(fPixArray, localMax, maxVal);
237 }
238
239 if (nMax > 1)
240 {
241 TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
242 }
243
244 Int_t iSimple = 0, nInX = -1, nInY;
245
246 PadsInXandY(*cluster,nInX, nInY);
247
248 if (nMax == 1 && nInX < 4 && nInY < 4)
249 {
250 iSimple = 1; //1; // simple cluster
251 }
252
253 for (Int_t i=0; i<nMax; ++i)
254 {
255 if (nMax > 1)
256 {
257 FindCluster(*cluster,localMax, maxPos[i]);
258 }
259 Timer(kMainLoop)->Start(kFALSE);
260 MainLoop(*cluster,iSimple);
261 Timer(kMainLoop)->Stop();
262 if (i < nMax-1)
263 {
264 for (Int_t j=0; j<cluster->Multiplicity(); ++j)
265 {
266 AliMUONPad* pad = cluster->Pad(j);
267 if ( pad->Status() == 0 ) continue; // pad charge was not modified
268 pad->SetStatus(0);
269 pad->RevertCharge(); // use backup charge value
270 }
271 }
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();
276 delete cluster;
277 return kTRUE;
278}
279
280//_____________________________________________________________________________
281Bool_t
282AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
283{
284 /// Check if the pad and the pixel overlaps
285
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);
290
291 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
292}
293
294//_____________________________________________________________________________
295AliMUONCluster*
296AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
297{
298 /// Check precluster in order to attempt to simplify it (mostly for
299 /// two-cathode preclusters)
300
301 if (origCluster.Multiplicity()==1)
302 {
303 // Disregard one-pad clusters (leftovers from splitting)
304 return 0x0;
305 }
306
307 Timer(kCheckPreCluster)->Start(kFALSE);
308
309
310 AliMUONCluster* cluster = static_cast<AliMUONCluster*>(origCluster.Clone());
311
312 cluster->Sort();
313
314 AliDebug(2,"Start of CheckPreCluster=");
315 StdoutToAliDebug(2,cluster->Print("full"));
316
317 // Check if one-cathode precluster
318 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
319 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
320
321 AliMUONCluster* rv(0x0);
322
323 if (i1 != i2)
324 {
325 rv = CheckPreclusterTwoCathodes(cluster);
326 }
327 else
328 {
329 rv = CheckPreclusterOneCathode(cluster);
330 }
331 Timer(kCheckPreCluster)->Stop();
332 return rv;
333}
334
335//_____________________________________________________________________________
336AliMUONCluster*
337AliMUONClusterFinderMLEM::CheckPreclusterOneCathode(AliMUONCluster* cluster)
338{
339 /// Check single-cathode precluster
340 AliWarning("Reimplement me!");
341 AliDebug(2,"End of CheckPreClusterOneCathode=");
342 StdoutToAliDebug(2,cluster->Print("full"));
343
344 return cluster;
345}
346
347//_____________________________________________________________________________
348AliMUONCluster*
349AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
350{
71a2d3aa 351 /// Check two-cathode cluster
c0a16418 352
353 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
354 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
355
356 Int_t npad = cluster->Multiplicity();
357 Int_t* flags = new Int_t[npad];
358 memset(flags,0,npad*sizeof(Int_t));
359
360 // Check pad overlaps
361 for ( Int_t i=0; i<npad; ++i)
362 {
363 AliMUONPad* padi = cluster->Pad(i);
364 if ( padi->Cathode() != i1 ) continue;
365 for (Int_t j=i+1; j<npad; ++j)
366 {
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
371 }
372 }
373
374 // Check if all pads overlap
375 Int_t nFlags=0;
376 for (Int_t i=0; i<npad; ++i)
377 {
378 if (flags[i]) continue;
379 ++nFlags;
380 }
381
382 if (nFlags > 0)
383 {
384 // not all pads overlap.
385 for (Int_t i=0; i<npad; ++i)
386 {
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();
398 --npad;
399 }
400 }
401
402 // Check correlations of cathode charges
403 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry()*2 > 1 )
404 {
405 // big difference
406 Int_t cathode = cluster->MaxRawChargeCathode();
407 Int_t imin(0);
408 Int_t imax(0);
409 Double_t cmax(0);
410 Double_t cmin(1E9);
411
412 // get min and max pad charges on the cathode opposite to the
413 // max pad (given by MaxRawChargeCathode())
414 //
415 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
416 {
417 AliMUONPad* pad = cluster->Pad(i);
418 if ( pad->Cathode() != cathode || !pad->IsReal() )
419 {
420 // only consider pads in the opposite cathode, and
421 // onyl consider real pads (i.e. exclude the virtual ones)
422 continue;
423 }
424 if ( pad->Charge() < cmin )
425 {
426 cmin = pad->Charge();
427 imin = i;
428 }
429 if ( pad->Charge() > cmax )
430 {
431 cmax = pad->Charge();
432 imax = i;
433 }
434 }
435 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
436 imin,imax,cmin,cmax));
437 //
438 // arrange pads according to their distance to the max, normalized
439 // to the pad size
440 Double_t* dist = new Double_t[cluster->Multiplicity()];
441 Double_t dxMin(1E9);
442 Double_t dyMin(1E9);
443 Double_t dmin(0);
444
445 AliMUONPad* padmax = cluster->Pad(imax);
446
447 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
448 {
449 dist[i] = 0.0;
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);
456 if ( i == imin )
457 {
458 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
459 dxMin = dx;
460 dyMin = dy;
461 }
462 }
463
464 TMath::Sort(cluster->Multiplicity(),dist,flags,kFALSE); // in ascending order
465 Double_t xmax(-1);
466 TObjArray toBeRemoved;
467
468 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
469 {
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 )
474 {
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;
478 dx *= dxMin;
479 dy *= dyMin;
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;
483 }
484 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
485 {
486 // release pad
487 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
488 {
489 cmax = TMath::Max(pad->Charge(),cmax);
490 }
491 else
492 {
493 cmax = pad->Charge();
494 }
495 xmax = dist[indx];
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(),
498 pad->Charge()));
499
500 toBeRemoved.AddLast(pad);
501 fPreCluster->Pad(indx)->Release();
502 }
503 }
504 for ( Int_t i = 0; i <= toBeRemoved.GetLast(); ++i )
505 {
506 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.At(i)));
507 }
508 delete[] dist;
509 }
510
511 delete[] flags;
512
513 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
514 StdoutToAliDebug(2,cluster->Print("full"));
515
516 return cluster;
517}
518
519//_____________________________________________________________________________
520void
521AliMUONClusterFinderMLEM::CheckOverlaps()
522{
523 /// For debug only : check if some pixels overlap...
524
525 Int_t nPix = fPixArray->GetLast()+1;
526 Int_t dummy(0);
527
528 for ( Int_t i = 0; i < nPix; ++i )
529 {
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);
534
535 for ( Int_t j = i+1; j < nPix; ++j )
536 {
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);
541 AliMpArea area;
542
543 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
544 {
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;
548 pixelJ->Print();
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;
552 );
553
554 }
555 }
556 }
557}
558
559//_____________________________________________________________________________
560void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
561{
562 /// Build pixel array for MLEM method
563
564 Int_t npad = cluster.Multiplicity();
565 if (npad<=0)
566 {
567 AliWarning("Got no pad at all ?!");
568 }
569
570 fPixArray->Delete();
571
572 if ( !cluster.Multiplicity(0) || !cluster.Multiplicity(1) )
573 {
574 BuildPixArrayOneCathode(cluster);
575 }
576 else
577 {
578 BuildPixArrayTwoCathodes(cluster);
579 }
580
581 fPixArray->Sort(); // FIXME : not really needed, only to compare with ClusterFinderAZ
582
583 Int_t nPix = fPixArray->GetLast()+1;
584
585// AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
586
587 Double_t xPadMin(1E9);
588 Double_t yPadMin(1E9);
589
590 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i )
591 {
592 AliMUONPad* pad = cluster.Pad(i);
593 xPadMin = TMath::Min (xPadMin, pad->DX());
594 yPadMin = TMath::Min (yPadMin, pad->DY());
595 }
596
597 Double_t wxmin(1E9);
598 Double_t wymin(1E9);
599
600 for ( Int_t i = 0; i < nPix; ++i )
601 {
602 AliMUONPad* pixPtr = Pixel(i);
603 wxmin = TMath::Min(wxmin, pixPtr->Size(0));
604 wymin = TMath::Min(wymin, pixPtr->Size(1));
605 }
606
607 wxmin = TMath::Abs (wxmin - xPadMin/2) > 0.001 ? xPadMin : xPadMin / 2;
608 wymin = TMath::Abs (wymin - yPadMin/2) > 0.001 ? yPadMin : yPadMin / 2;
609
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);
616
617 // Remove discarded pixels
618 for (Int_t i=0; i<nPix; ++i)
619 {
620 AliMUONPad* pixPtr = Pixel(i);
621 if (pixPtr->Charge() < 1)
622 {
623 AliDebug(2,Form("Removing pixel %d with charge<1 : ",i));
624 StdoutToAliDebug(2,pixPtr->Print());
625 RemovePixel(i);
626 }
627 }
628
629 fPixArray->Compress();
630 nPix = fPixArray->GetEntriesFast();
631
632// AliDebug(2,Form("nPix after AdjustPixel=%d",nPix));
633
634 if ( nPix > cluster.Multiplicity() )
635 {
636// AliDebug(2,Form("Will trim number of pixels to number of pads"));
637
638 // Too many pixels - sort and remove pixels with the lowest signal
639 fPixArray->Sort();
640 for ( Int_t i = cluster.Multiplicity(); i<nPix; ++i )
641 {
642 RemovePixel(i);
643 }
644 fPixArray->Compress();
645 nPix = fPixArray->GetEntriesFast();
646 } // if (nPix > npad)
647
648// StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
649// fPixArray->Print(););
650 CheckOverlaps();//FIXME : this is for debug only. Remove it.
651}
652
653//_____________________________________________________________________________
654void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
655{
656 /// From a single-cathode cluster, build the pixel array
657
658// AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
659
660 for ( Int_t j=0; j<cluster.Multiplicity(); ++j)
661 {
662 AliMUONPad* pad = cluster.Pad(j);
663 AliMUONPad* pixPtr = new AliMUONPad(pad->Position(),pad->Dimensions(),
664 pad->Charge());
665 fPixArray->Add(pixPtr);
666 }
667}
668
669//_____________________________________________________________________________
670void AliMUONClusterFinderMLEM::BuildPixArrayTwoCathodes(AliMUONCluster& cluster)
671{
672 /// From a two-cathodes cluster, build the pixel array
673
674// AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
675
676 Int_t i1 = cluster.Pad(0)->Cathode();
677 Int_t i2 = TMath::Even(i1);
678
679 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i)
680 {
681 AliMUONPad* padi = cluster.Pad(i);
682 if (padi->Cathode() != i1) continue;
683
684 for ( Int_t j = 1; j < cluster.Multiplicity(); ++j)
685 {
686 AliMUONPad* padj = cluster.Pad(j);
687 if (padj->Cathode() != i2) continue;
688
689 AliMpArea overlap;
690
691 if ( AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize,overlap) )
692 {
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() ) )
698 {
699 // if smallest charge of the 2 pads is already above saturation, then
700 // the pixel is saturated...
701 pixPtr->SetSaturated(kTRUE);
702 }
703 pixPtr->SetReal(kFALSE);
704 fPixArray->Add(pixPtr);
705 }
706 }
707 }
708}
709
710//_____________________________________________________________________________
711void AliMUONClusterFinderMLEM::AdjustPixel(AliMUONCluster& cluster,
712 Float_t width, Int_t ixy)
713{
714 /// Check if some pixels have small size (adjust if necessary)
715
716 AliDebug(2,Form("width=%e ixy=%d",width,ixy));
717
718 AliMUONPad *pixPtr, *pixPtr1 = 0;
719 Int_t ixy1 = TMath::Even(ixy);
720 Int_t nPix = fPixArray->GetEntriesFast();
721
722 for (Int_t i=0; i<nPix; i++)
723 {
724 pixPtr = Pixel(i);
725 if (pixPtr->Charge() < 1) continue; // discarded pixel
726 if (pixPtr->Size(ixy)-width < -1.e-4)
727 {
728 // try to merge
729 for (Int_t j=i+1; j<nPix; j++)
730 {
731 pixPtr1 = Pixel(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)
736 {
737 // merge
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);
744 pixPtr1 = 0;
745 break;
746 }
747 } // for (Int_t j=i+1;
748 if (pixPtr1 || i == nPix-1)
749 {
750 // edge pixel - just increase its size
751 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
752 {
753 AliMUONPad* pad = cluster.Pad(j);
754 Double_t d = ( ixy == 0 ) ? pad->X() : ( ixy == 1 ) ? pad->Y() : -1E9;
755
756 if (pixPtr->Coord(ixy) < d)
757 {
758 pixPtr->Shift(ixy, pixPtr->Size(ixy)-width);
759 }
760 else
761 {
762 pixPtr->Shift(ixy, -pixPtr->Size(ixy)+width);
763 }
764 pixPtr->SetSize(ixy, width);
765 break;
766 }
767 }
768 } // if (pixPtr->Size(ixy)-width < -1.e-4)
769 } // for (Int_t i=0; i<nPix;
770}
771
772//_____________________________________________________________________________
773void AliMUONClusterFinderMLEM::AdjustPixel(Double_t wxmin, Double_t wymin)
774{
775/// Check if some pixels have large size (adjust if necessary)
776
777 AliDebug(2,Form("wxmin=%e wymin=%e",wxmin,wymin));
778
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};
782
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
791 continue;
792 } else iOK = 0; // large pixel
793 } else {
794 iOK = 0; // large pixel
795 if (xy0[1] > 9998 && pixPtr->Size(1) - wymin < 1.e-4) xy0[1] = pixPtr->Coord(1); // "normal" pixel
796 }
797 if (xy0[0] < 9998 && xy0[1] < 9998) break;
798 }
799 if (iOK) return;
800
801 wxy[0] = wxmin;
802 wxy[1] = wymin;
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
807 n1[0] = n1[1] = 999;
808 n2[0] = n2[1] = 1;
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];
814 }
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;
818
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]--;
823 }
824 //cout << n1[0] << " " << n2[0] << " " << n1[1] << " " << n2[1] << endl;
825 pix = *pixPtr;
826 pix.SetSize(0, wxy[0]); pix.SetSize(1, wxy[1]);
827 //pixPtr->Print();
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));
833 //pix.Print();
834 }
835 }
836 pixPtr->SetCharge(0);
837 } // for (Int_t i = 0; i < nPix;
838}
839
840//_____________________________________________________________________________
841void
842AliMUONClusterFinderMLEM::Plot(const char* basename)
843{
844 /// Make a plot and save it as png
845
846 if (!fPlot) return;
847
848 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
849 c->Draw();
850 Draw();
851 c->Modified();
852 c->Update();
853 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
854 fDetElemId,fClusterNumber));
855}
856
857//_____________________________________________________________________________
858void
859AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
860 Double_t* coef,
861 Double_t* probi)
862{
863 /// Compute coefficients needed for MLEM algorithm
864
865 Int_t nPix = fPixArray->GetLast()+1;
866
867 memset(probi,0,nPix*sizeof(Double_t));
868
869 for ( Int_t j=0; j<cluster.Multiplicity(); ++j )
870 {
871 AliMUONPad* pad = cluster.Pad(j);
872 Int_t indx = j*nPix;
873
874 for ( Int_t ipix=0; ipix<nPix; ++ipix )
875 {
876 Int_t indx1 = indx + ipix;
877 if (pad->Status() < 0)
878 {
879 coef[indx1] = 0;
880 continue;
881 }
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),
892// coef[indx1]));
893
894 probi[ipix] += coef[indx1];
895 }
896 }
897}
898
899//_____________________________________________________________________________
900Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
901{
902 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
903
904 Int_t nPix = fPixArray->GetLast()+1;
905
906 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
907 StdoutToAliDebug(2,cluster.Print("full"););
908
909 if ( nPix < 0 )
910 {
911 AliDebug(1,"No pixels, why am I here then ?");
912 }
913
914 AddVirtualPad(cluster); // add virtual pads if necessary
915
916 Int_t npadTot = cluster.Multiplicity();
917 Int_t npadOK = 0;
918 for (Int_t i = 0; i < npadTot; ++i)
919 {
920 if (cluster.Pad(i)->Status() == 0) ++npadOK;
921 }
922
923 TH2D* mlem(0x0);
924 Double_t* coef(0x0);
925 Double_t* probi(0x0);
926 Int_t lc(0); // loop counter (for debug)
927
928 Plot("mlem.start");
929
930 while (1)
931 {
932 ++lc;
933 mlem = (TH2D*) gROOT->FindObject("mlem");
934 delete mlem;
935
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"));
939
940 coef = new Double_t [npadTot*nPix];
941 probi = new Double_t [nPix];
942
943 // Calculate coefficients and pixel visibilities
944 ComputeCoefficients(cluster,coef,probi);
945
946 for (Int_t ipix=0; ipix<nPix; ++ipix)
947 {
948 if (probi[ipix] < 0.01)
949 {
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
954 }
955 }
956
957 // MLEM algorithm
958 Mlem(cluster,coef, probi, 15);
959
960 Double_t xylim[4] = {999, 999, 999, 999};
961 AliMUONPad* pixPtr(0x0);
962
963 for ( Int_t ipix=0; ipix<nPix; ++ipix )
964 {
965 pixPtr = Pixel(ipix);
966 for ( Int_t i=0; i<4; ++i )
967 {
968 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
969 }
970 }
971 for (Int_t i=0; i<4; i++)
972 {
973 xylim[i] -= pixPtr->Size(i/2);
974 }
975
976
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);
979
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]
984 ));
985
986 mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
987
988 for (Int_t ipix=0; ipix<nPix; ++ipix)
989 {
990 AliMUONPad* pixPtr = Pixel(ipix);
991 mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
992 }
993
994 // Check if the total charge of pixels is too low
995 Double_t qTot = 0;
996 for ( Int_t i=0; i<nPix; ++i)
997 {
998 qTot += Pixel(i)->Charge();
999 }
1000
1001 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
1002 {
1003 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
1004 delete [] coef;
1005 delete [] probi;
1006 coef = 0;
1007 probi = 0;
1008 fPixArray->Delete();
1009 for ( Int_t i=0; i<npadTot; ++i)
1010 {
1011 AliMUONPad* pad = cluster.Pad(i);
1012 if ( pad->Status() == 0) pad->SetStatus(-1);
1013 }
1014 return kFALSE;
1015 }
1016
1017 if (iSimple)
1018 {
1019 // Simple cluster - skip further passes thru EM-procedure
1020 Simple(cluster);
1021 delete [] coef;
1022 delete [] probi;
1023 coef = 0;
1024 probi = 0;
1025 fPixArray->Delete();
1026 return kTRUE;
1027 }
1028
1029 // Calculate position of the center-of-gravity around the maximum pixel
1030 Double_t xyCOG[2];
1031 FindCOG(mlem, xyCOG);
1032
1033 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
1034 pixPtr->Size(0) > pixPtr->Size(1)) break;
1035
1036 // Sort pixels according to the charge
1037 fPixArray->Sort();
1038 Double_t pixMin = 0.01*Pixel(0)->Charge();
1039 pixMin = TMath::Min(pixMin,50.);
1040
1041 // Decrease pixel size and shift pixels to make them centered at
1042 // the maximum one
1043 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
1044 Int_t ix(1);
1045 Double_t width = 0;
1046 Double_t shift[2] = { 0.0, 0.0 };
1047 for (Int_t i=0; i<4; i++) xylim[i] = 999;
1048 Int_t nPix1 = nPix;
1049 nPix = 0;
1050 for (Int_t ipix=0; ipix<nPix1; ++ipix)
1051 {
1052 AliMUONPad* pixPtr = Pixel(ipix);
1053 if ( nPix >= npadOK // too many pixels already
1054 ||
1055 pixPtr->Charge() < pixMin // too low charge
1056 )
1057 {
1058 RemovePixel(ipix);
1059 continue;
1060 }
1061 for (Int_t i=0; i<2; ++i)
1062 {
1063 if (!i)
1064 {
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
1070 if (ix)
1071 {
1072 ix = 0;
1073 for (Int_t j=0; j<2; ++j)
1074 {
1075 shift[j] = pixPtr->Coord(j) - xyCOG[j];
1076 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
1077 }
1078 } // if (ix)
1079 pixPtr->Shift(0, -shift[0]);
1080 pixPtr->Shift(1, -shift[1]);
1081 }
1082 else
1083 {
1084 pixPtr = new AliMUONPad(*pixPtr);
1085 pixPtr->Shift(indx, -2*width);
1086 fPixArray->Add(pixPtr);
1087 }
1088 for (Int_t i=0; i<4; ++i)
1089 {
1090 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1091 }
1092 } // for (Int_t i=0; i<2;
1093 nPix += 2;
1094 } // for (Int_t ipix=0;
1095
1096 fPixArray->Compress();
1097 nPix = fPixArray->GetEntriesFast();
1098
1099 AliDebug(2,Form("After shift:"));
1100 StdoutToAliDebug(2,fPixArray->Print("","full"););
1101 Plot(Form("mlem.lc%d",lc+1));
1102
1103 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1104 xyCOG[0],xyCOG[1],
1105 xylim[0],xylim[1],
1106 xylim[2],xylim[3]));
1107
1108 // Remove excessive pixels
1109 if (nPix > npadOK)
1110 {
1111 for (Int_t ipix=npadOK; ipix<nPix; ++ipix)
1112 {
1113 RemovePixel(ipix);
1114 }
1115 }
1116 else
1117 {
1118 AliMUONPad* pixPtr = Pixel(0);
1119 // add pixels if the maximum is at the limit of pixel area
1120 // start from Y-direction
1121 Int_t j = 0;
1122 for (Int_t i=3; i>-1; --i)
1123 {
1124 if (nPix < npadOK &&
1125 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2))
1126 {
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"););
1134 fPixArray->Add(p);
1135 ++nPix;
1136 }
1137 }
1138 }
1139 fPixArray->Compress();
1140 nPix = fPixArray->GetEntriesFast();
1141 delete [] coef;
1142 delete [] probi;
1143 coef = 0;
1144 probi = 0;
1145 } // while (1)
1146
1147 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1148 StdoutToAliDebug(2,fPixArray->Print("","full"););
1149
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.);
1154 Double_t cmax = -1;
1155 Double_t charge = 0;
1156
1157 for ( Int_t i=0; i<nPix; ++i)
1158 {
1159 cmax = TMath::Max (cmax,probi[i]);
1160 }
1161
1162 // Mark pixels which should be removed
1163 for (Int_t i=0; i<nPix; ++i)
1164 {
1165 AliMUONPad* pixPtr = Pixel(i);
1166 charge = pixPtr->Charge();
1167 if (charge < thresh)
1168 {
1169 pixPtr->SetCharge(-charge);
1170 }
1171 }
1172
1173 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1174 Int_t near = 0;
1175 for (Int_t i=0; i<nPix; ++i)
1176 {
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));
1185 }
1186 Mlem(cluster,coef,probi,2);
1187
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");
1191
1192 // Update histogram
1193 for (Int_t i=0; i<nPix; ++i)
1194 {
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());
1199 }
1200
1201 // Try to split into clusters
1202 Bool_t ok = kTRUE;
1203 if (mlem->GetSum() < 1)
1204 {
1205 ok = kFALSE;
1206 }
1207 else
1208 {
1209 fSplitter->Split(cluster,mlem,coef,fClusterList);
1210 }
1211
1212 delete [] coef;
1213 delete [] probi;
1214 coef = 0;
1215 probi = 0;
1216 fPixArray->Delete();
1217
1218 return ok;
1219}
1220
1221//_____________________________________________________________________________
1222void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1223 Double_t* coef, Double_t* probi,
1224 Int_t nIter)
1225{
1226 /// Use MLEM to find pixel charges
1227
1228 Int_t nPix = fPixArray->GetEntriesFast();
1229
1230 Int_t npad = cluster.Multiplicity();
1231
1232 Double_t* probi1 = new Double_t[nPix];
1233 Double_t probMax = 0;
1234 Double_t tmp = TMath::MaxElement(nPix,probi);
1235
1236 for (Int_t ipix=0; ipix<nPix; ++ipix)
1237 {
1238 probMax = TMath::Max(probMax,probi[ipix]);
1239 }
1240
1241 if (probMax!=tmp) { AliWarning(Form("probMax=%e tmp=%e",probMax,tmp)); }
1242
1243 for (Int_t iter=0; iter<nIter; ++iter)
1244 {
1245 // Do iterations
1246 for (Int_t ipix=0; ipix<nPix; ++ipix)
1247 {
1248 // Correct each pixel
1249 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1250 Double_t sum = 0;
1251 probi1[ipix] = probMax;
1252 for (Int_t j=0; j<npad; j++)
1253 {
1254 AliMUONPad* pad = cluster.Pad(j);
1255 if (pad->Status() < 0) continue;
1256 Double_t sum1 = 0;
1257 Int_t indx1 = j*nPix;
1258 Int_t indx = indx1 + ipix;
1259 // Calculate expectation
1260 for (Int_t i=0; i<nPix; i++)
1261 {
1262 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1263 }
1264 if ( pad->Charge() > fgkSaturation-1 && sum1 > pad->Charge() ) //FIXME : remove usage of fgkSaturation
1265 {
1266 if ( !pad->IsSaturated() )
1267 {
1268 AliWarning("Got a pad charge above saturation not backed-up by pad->IsSaturated() function : ");
1269 StdoutToAliWarning(pad->Print("full"));
1270 }
1271 // correct for pad charge overflows
1272 probi1[ipix] -= coef[indx];
1273 continue;
1274 }
1275
1276 if (coef[indx] > 1.e-6)
1277 {
1278 sum += pad->Charge()*coef[indx]/sum1;
1279 }
1280 } // for (Int_t j=0;
1281 AliMUONPad* pixPtr = Pixel(ipix);
1282 if (probi1[ipix] > 1.e-6)
1283 {
1284 pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1285 }
1286 } // for (Int_t ipix=0;
1287 } // for (Int_t iter=0;
1288 delete [] probi1;
1289}
1290
1291//_____________________________________________________________________________
1292void AliMUONClusterFinderMLEM::FindCOG(TH2D *mlem, Double_t *xyc)
1293{
1294 /// Calculate position of the center-of-gravity around the maximum pixel
1295
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;
1303
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);
1312 xq += x*cont;
1313 yq += y*cont;
1314 qq += cont;
1315 nsum++;
1316 }
1317 }
1318
1319 Double_t cmax = 0;
1320 Int_t i2 = 0, j2 = 0;
1321 x = y = 0;
1322 if (nsumy == 1) {
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);
1328 if (cont > cmax) {
1329 cmax = cont;
1330 x = mlem->GetXaxis()->GetBinCenter(j);
1331 y = mlem->GetYaxis()->GetBinCenter(i);
1332 i2 = i;
1333 j2 = j;
1334 }
1335 }
1336 }
1337 xq += x*cmax;
1338 yq += y*cmax;
1339 qq += cmax;
1340 if (i2 != i1) nsumy++;
1341 if (j2 != j1) nsumx++;
1342 nsum++;
1343 } // if (nsumy == 1)
1344
1345 if (nsumx == 1) {
1346 // one bin in X - add one more (with the largest signal)
1347 cmax = x = y = 0;
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);
1352 if (cont > cmax) {
1353 cmax = cont;
1354 x = mlem->GetXaxis()->GetBinCenter(j);
1355 y = mlem->GetYaxis()->GetBinCenter(i);
1356 i2 = i;
1357 j2 = j;
1358 }
1359 }
1360 }
1361 xq += x*cmax;
1362 yq += y*cmax;
1363 qq += cmax;
1364 if (i2 != i1) nsumy++;
1365 if (j2 != j1) nsumx++;
1366 nsum++;
1367 } // if (nsumx == 1)
1368
1369 xyc[0] = xq/qq; xyc[1] = yq/qq;
1370 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1371}
1372
1373//_____________________________________________________________________________
1374Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1375{
1376/// Find the pixel nearest to the given one
1377/// (algorithm may be not very efficient)
1378
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);
1382 AliMUONPad *pixPtr;
1383
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; }
1391 }
1392 return imin;
1393}
1394
1395
1396//_____________________________________________________________________________
1397TStopwatch*
1398AliMUONClusterFinderMLEM::Timer(Int_t i) const
1399{
1400 /// Return timer at index i
1401 return static_cast<TStopwatch*>(fTimers->At(i));
1402}
1403
1404//_____________________________________________________________________________
1405void
1406AliMUONClusterFinderMLEM::Paint(Option_t*)
1407{
1408 /// Paint cluster and pixels
1409
1410 AliMpArea area(fPreCluster->Area());
1411
1412 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1413
1414 gVirtualX->SetFillStyle(1001);
1415 gVirtualX->SetFillColor(3);
1416 gVirtualX->SetLineColor(3);
1417
1418 Double_t s(1.0);
1419
1420 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1421 {
1422 AliMUONPad* pixel = Pixel(i);
1423
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);
1428
1429// for ( Int_t sign = -1; sign < 2; sign +=2 )
1430// {
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)
1435// );
1436// }
1437 }
1438
1439
1440 gVirtualX->SetFillStyle(0);
1441
1442 fPreCluster->Paint();
1443
1444 gVirtualX->SetLineColor(1);
1445 gVirtualX->SetLineWidth(2);
1446 gVirtualX->SetFillStyle(0);
1447 gVirtualX->SetTextColor(1);
1448 gVirtualX->SetTextAlign(22);
1449
1450 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1451 {
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);
1458
1459 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1460 }
1461}
1462
1463//_____________________________________________________________________________
1464Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1465{
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
1469
1470 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1471
1472 TH2D *hist = NULL;
1473 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
1474 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
1475 //if (hist) hist->Delete();
1476
1477 Double_t xylim[4] = {999, 999, 999, 999};
1478
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));
1485 }
1486 for (Int_t i=0; i<4; i++) xylim[i] -= pixPtr->Size(i/2);
1487
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());
1495 }
1496// if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1497
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;
1501
1502 for (Int_t i=1; i<=ny; i++) {
1503 indx = (i-1) * nx;
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);
1509 }
1510 }
1511
1512 for (Int_t i=1; i<=ny; i++) {
1513 indx = (i-1) * nx;
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 !!!");
1519 }
1520 }
1521 }
1522 if (fDebug) cout << " Local max: " << nMax << endl;
1523 delete [] isLocalMax; isLocalMax = 0;
1524 return nMax;
1525}
1526
1527//_____________________________________________________________________________
1528void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1529{
1530/// Flag pixels (whether or not local maxima)
1531
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;
1536
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;
1553 }
1554 }
1555 }
1556 }
1557 isLocalMax[indx] = 1; // local maximum
1558}
1559
1560//_____________________________________________________________________________
1561void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1562 Int_t *localMax, Int_t iMax)
1563{
1564/// Find pixel cluster around local maximum \a iMax and pick up pads
1565/// overlapping with it
1566
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;
1574
1575 // Drop all pixels from the array - pick up only the ones from the cluster
1576 fPixArray->Delete();
1577
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
1586
1587 Int_t nPix = fPixArray->GetEntriesFast();
1588 Int_t npad = cluster.Multiplicity();
1589
1590 for (Int_t i=0; i<nPix; ++i)
1591 {
1592 AliMUONPad* pixPtr = Pixel(i);
1593 pixPtr->SetSize(0,wx);
1594 pixPtr->SetSize(1,wy);
1595 }
1596
1597 // Pick up pads which overlap with found pixels
1598 for (Int_t i=0; i<npad; i++)
1599 {
1600 cluster.Pad(i)->SetStatus(-1);
1601 }
1602
1603 for (Int_t i=0; i<nPix; i++)
1604 {
1605 AliMUONPad* pixPtr = Pixel(i);
1606 for (Int_t j=0; j<npad; ++j)
1607 {
1608 AliMUONPad* pad = cluster.Pad(j);
1609 if ( Overlap(*pad,*pixPtr) )
1610 {
1611 pad->SetStatus(0);
1612 }
1613 }
1614 }
1615
1616 delete [] used; used = 0;
1617}
1618
1619//_____________________________________________________________________________
1620AliMUONClusterFinderMLEM&
1621AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1622{
1623/// Protected assignement operator
1624
1625 if (this == &rhs) return *this;
1626
1627 AliFatal("Not implemented.");
1628
1629 return *this;
1630}
1631
1632//_____________________________________________________________________________
1633void
1634AliMUONClusterFinderMLEM::Neighbours(Int_t cathode, Int_t ix, Int_t iy,
1635 Int_t& n, Int_t* xList, Int_t* yList)
1636{
1637 /// Get the list of neighbours of pad at (cathode,ix,iy)
1638 n = 0;
1639
1640 const AliMpVSegmentation* seg = fSegmentation[cathode];
1641
1642 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1643
1644 // Define the region to look into : a region slightly bigger
1645 // than the pad itself (5% bigger), in order to catch first neighbours.
1646
1647 AliMpArea area(pad.Position(),pad.Dimensions()*1.05);
1648
1649 AliMpVPadIterator* it = seg->CreateIterator(area);
1650 it->First();
1651 while ( !it->IsDone() && n < 10 )
1652 {
1653 AliMpPad p = it->CurrentItem();
1654 if ( p != pad ) // skip self
1655 {
1656 xList[n] = p.GetIndices().GetFirst();
1657 yList[n] = p.GetIndices().GetSecond();
1658 ++n;
1659 }
1660 it->Next();
1661 }
1662 delete it;
1663}
1664
1665//_____________________________________________________________________________
1666void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1667{
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)
1670
1671 // Get number of pads in X and Y-directions
1672 Int_t nInX = -1, nInY;
1673 PadsInXandY(cluster,nInX, nInY);
1674 ++fNClusters;
1675 if (fDebug) cout << " Chamber: " << fDetElemId / 100 - 1 << " " << nInX << " " << nInY << endl;
1676
1677 // Add virtual pad only if number of pads per direction == 2
1678 if (nInX != 2 && nInY != 2) return;
1679
1680 ++fNAddVirtualPads;
1681
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)
1686 {
1687 AliMUONPad* pad = cluster.Pad(j);
1688 if (pad->Status() != 0) continue;
1689 cath = pad->Cathode();
1690 if (pad->Charge() > sigmax[cath])
1691 {
1692 maxpad[cath][1] = maxpad[cath][0];
1693 aamax[cath] = sigmax[cath];
1694 sigmax[cath] = pad->Charge();
1695 maxpad[cath][0] = j;
1696 }
1697 }
1698
1699 if (maxpad[0][0] >= 0 && maxpad[0][1] < 0 || maxpad[1][0] >= 0 && maxpad[1][1] < 0)
1700 {
1701 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
1702 {
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])
1708 {
1709 aamax[cath] = pad->Charge();
1710 maxpad[cath][1] = j;
1711 }
1712 }
1713 }
1714
1715 // cout << "-------AddVirtualPad" << endl;
1716// cout << Form("nInX %2d nInY %2d",nInX,nInY) << endl;
1717//
1718// cluster.Print("full");
1719//
1720// for ( Int_t i = 0; i < 2; ++i )
1721// {
1722// for ( Int_t j = 0; j < 2; ++j )
1723// {
1724// cout << Form("maxpad[%d][%d]=%d",i,j,maxpad[i][j]) << endl;
1725// }
1726// }
1727
1728 // Check for mirrors (side X on cathode 0)
1729 Bool_t mirror = kFALSE;
1730 if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1731 {
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)
1735 {
1736 // Special case when pads on both cathodes have the same size
1737 Int_t yud[2] = {0};
1738 for (Int_t j = 0; j < cluster.Multiplicity(); ++j)
1739 {
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 )
1745 {
1746 yud[cath]++;
1747 }
1748 }
1749 if (!yud[0]) mirror = kTRUE; // take the other cathode
1750 } // if (!mirror &&...
1751 } // if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1752
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++)
1756 {
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)
1760 {
1761 if (!mirror)
1762 {
1763 if (!cath && nInY != 2) continue;
1764 if (cath && nInX != 2 && (maxpad[0][0] >= 0 || nInY != 2)) continue;
1765 }
1766 else
1767 {
1768 if (!cath && nInX != 2) continue;
1769 if (cath && nInY != 2 && (maxpad[0][0] >= 0 || nInX != 2)) continue;
1770 }
1771 }
1772
1773 Int_t iAddX = 0, iAddY = 0, ix1 = 0, iy1 = 0, iPad = 0;
1774 if (maxpad[0][0] < 0) iPad = 1;
1775
1776 for (iPad=0; iPad<2; iPad++)
1777 {
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;
1781
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);
1786 //Float_t zpad;
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++;
1790 }
1791 if (!mirror) {
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;
1796 } else {
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;
1801 }
1802
1803 for (Int_t j=0; j< cluster.Multiplicity(); ++j)
1804 {
1805 AliMUONPad* pad = cluster.Pad(j);
1806 if ( pad->Cathode() != cath) continue;
1807 ix = pad->Ix();
1808 iy = pad->Iy();
1809 if (iy == iy0 && ix == ix0) continue;
1810 for (Int_t k=0; k<nn; ++k)
1811 {
1812 if (xList[k] != ix || yList[k] != iy) continue;
1813 if (!mirror)
1814 {
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;
1819 neighb--;
1820 break;
1821 }
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;
1826 neighb--;
1827 }
1828 } else {
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;
1833 neighb--;
1834 break;
1835 }
1836 if ((cath || maxpad[1][0] < 0) &&
1837 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
1838 xList[k] = yList[k] = 0;
1839 neighb--;
1840 }
1841 }
1842 break;
1843 } // for (Int_t k=0; k<nn;
1844 if (!neighb) break;
1845 } // for (Int_t j=0; j< cluster.Multiplicity();
1846 if (!neighb) continue;
1847
1848 // Add virtual pad
1849 Int_t npads, isec;
1850 isec = npads = 0;
1851 for (Int_t j=0; j<nn; j++)
1852 {
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;
1857 ix = xList[j];
1858 iy = yList[j];
1859 if (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) {
1860 if (iy != iy0) continue; // new segmentation - check
1861 if (nInX != 2) continue; // new
1862 if (!mirror) {
1863 if (!cath && maxpad[1][0] >= 0) continue;
1864 } else {
1865 if (cath && maxpad[0][0] >= 0) continue;
1866 }
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.));
1882 }
1883 else {
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.));
1888 }
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
1895 // fnPads[1]++;
1896 // iAddX = npads;
1897 iAddX = 1;
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
1904 ix1 = ix0;
1905 continue;
1906 }
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));
1925 }
1926 else {
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));
1931 }
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
1938 // fnPads[1]++;
1939 // iAddY = npads;
1940 iAddY = 1;
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
1944 iy1 = iy0;
1945 }
1946 } // for (Int_t j=0; j<nn;
1947 } // for (Int_t iPad=0;
1948 } // for (cath=0; cath<2;
1949}
1950
1951//_____________________________________________________________________________
1952void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1953 Int_t &nInX, Int_t &nInY) const
1954{
1955 /// Find number of pads in X and Y-directions (excluding virtual ones and
1956 /// overflows)
1957
1958 Int_t statusToTest = 1;
1959
1960 if ( nInX < 0 ) statusToTest = 0;
1961
1962 Bool_t mustMatch(kTRUE);
1963
1964 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1965
1966 nInX = cn.GetFirst();
1967 nInY = cn.GetSecond();
1968}
1969
1970//_____________________________________________________________________________
1971void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1972{
1973 /// Remove pixel at index i
1974 AliMUONPad* pixPtr = Pixel(i);
1975 fPixArray->RemoveAt(i);
1976 delete pixPtr;
1977}
1978
1979//_____________________________________________________________________________
1980void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1981{
1982/// Process simple cluster (small number of pads) without EM-procedure
1983
1984 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1985 Double_t parOk[3] = {0.};
1986 TObjArray *clusters[1];
1987 clusters[0] = fPixArray;
1988
1989 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1990
1991 for (Int_t i = 0; i < cluster.Multiplicity(); ++i)
1992 {
1993 AliMUONPad* pad = cluster.Pad(i);
1994 if ( pad->Charge() > fgkSaturation-1) //FIXME : remove usage of fgkSaturation
1995 {
1996 pad->SetStatus(-9);
1997 }
1998 else
1999 {
2000 pad->SetStatus(1);
2001 }
2002 }
2003 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList);
2004}
2005
2006//_____________________________________________________________________________
2007AliMUONPad*
2008AliMUONClusterFinderMLEM::Pixel(Int_t i) const
2009{
2010 /// Returns pixel at index i
2011 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
2012}
2013
2014//_____________________________________________________________________________
2015void
2016AliMUONClusterFinderMLEM::Print(Option_t* what) const
2017{
2018 /// printout
2019 TString swhat(what);
2020 swhat.ToLower();
2021 if ( swhat.Contains("precluster") )
2022 {
2023 if ( fPreCluster) fPreCluster->Print();
2024 }
2025}
2026
2027