]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONClusterFinderMLEM.cxx
Take into account new class AliAODTracklets.
[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>
b1a19e07 48//#include "AliCodeTimer.h"
a6435fb4 49
c0a16418 50/// \cond CLASSIMP
51ClassImp(AliMUONClusterFinderMLEM)
52/// \endcond
53
54const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
55const Double_t AliMUONClusterFinderMLEM::fgkSaturation = 3000; // average saturation level
56const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
57const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
58const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
59
60 TMinuit* AliMUONClusterFinderMLEM::fgMinuit = 0x0;
61
62//_____________________________________________________________________________
b1a19e07 63AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
c0a16418 64 : AliMUONVClusterFinder(),
b1a19e07 65fPreClusterFinder(clusterFinder),
c0a16418 66fPreCluster(0x0),
67fClusterList(),
68fEventNumber(0),
69fDetElemId(-1),
70fClusterNumber(0),
71fZpad(0.0),
72fReco(1),
73fCathBeg(0),
74fPixArray(new TObjArray(20)),
9bbd7f60 75fDebug(0),
c0a16418 76fPlot(plot),
77fTimers(new TObjArray(kLast)),
78fSplitter(0x0),
79fNClusters(0),
80fNAddVirtualPads(0)
81{
82 /// Constructor
83
84 fSegmentation[1] = fSegmentation[0] = 0x0;
85
86 fPadBeg[0] = fPadBeg[1] = fCathBeg;
87
88 if (!fgMinuit) fgMinuit = new TMinuit(8);
89
90 fTimers->SetOwner(kTRUE);
91
92 for ( Int_t i = 0; i < kLast; ++i )
93 {
94 TStopwatch* t = new TStopwatch;
95 fTimers->AddLast(new TStopwatch);
96 t->Start(kTRUE);
97 t->Stop();
98 }
9bbd7f60 99 if (fPlot) fDebug = 1;
c0a16418 100}
101
102//_____________________________________________________________________________
103AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
104{
105/// Destructor
106 delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
107// delete fDraw;
108 delete fPreClusterFinder;
109 for ( Int_t i = 0; i < kLast; ++i )
110 {
111 AliInfo(Form("Timer %d",i));
112 Timer(i)->Print();
113 }
114 delete fTimers;
115 delete fSplitter;
116 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
117 fNClusters,fNAddVirtualPads));
118}
119
120//_____________________________________________________________________________
121Bool_t
122AliMUONClusterFinderMLEM::Prepare(const AliMpVSegmentation* segmentations[2],
3887d11e 123 const AliMUONVDigitStore& digitStore)
c0a16418 124{
125 /// Prepare for clustering
126
127 for ( Int_t i = 0; i < 2; ++i )
128 {
129 fSegmentation[i] = segmentations[i];
130 }
131
132 // Find out the DetElemId
133 fDetElemId = -1;
134
3887d11e 135 TIter next(digitStore.CreateIterator());
136 AliMUONVDigit* d = static_cast<AliMUONVDigit*>(next());
137
138 if (d)
c0a16418 139 {
3887d11e 140 fDetElemId = d->DetElemId();
c0a16418 141 }
142
143 if ( fDetElemId < 0 )
144 {
145 AliWarning("Could not find DE. Probably no digits at all ?");
146 return kFALSE;
147 }
148
149 delete fSplitter;
150 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
9bbd7f60 151 fSplitter->SetDebug(fDebug);
c0a16418 152
153 // find out current event number, and reset the cluster number
154 fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
155 fClusterNumber = -1;
156 fClusterList.Delete();
157
b1a19e07 158 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
c0a16418 159
3887d11e 160 return fPreClusterFinder->Prepare(segmentations,digitStore);
c0a16418 161}
162
163//_____________________________________________________________________________
164AliMUONCluster*
165AliMUONClusterFinderMLEM::NextCluster()
166{
167 /// Return next cluster
168
c0a16418 169 // if the list of clusters is not void, pick one from there
9bbd7f60 170 TObject* o = fClusterList.At(++fClusterNumber);
171 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
c0a16418 172
173 //FIXME : at this point, must check whether we've used all the digits
174 //from precluster : if not, let the preclustering know about those unused
175 //digits, so it can reuse them
176
177 // if the cluster list is exhausted, we need to go to the next
178 // pre-cluster and treat it
179
180 fPreCluster = fPreClusterFinder->NextCluster();
181
182 if (!fPreCluster)
183 {
184 // we are done
185 return 0x0;
186 }
187
b1a19e07 188// AliCodeTimerAuto("Excluding preclustering")
189
190
c0a16418 191 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
9bbd7f60 192 fClusterNumber = -1; //AZ
c0a16418 193
194 WorkOnPreCluster();
195
196 // WorkOnPreCluster may have used only part of the pads, so we check that
197 // now, and let the unused pads be reused by the preclustering...
198
199 for ( Int_t i = 0; i < fPreCluster->Multiplicity(); ++i )
200 {
201 AliMUONPad* pad = fPreCluster->Pad(i);
202 if ( !pad->IsUsed() )
203 {
204 fPreClusterFinder->UsePad(*pad);
205 }
206 }
207
208 return NextCluster();
209}
210
211//_____________________________________________________________________________
212Bool_t
213AliMUONClusterFinderMLEM::WorkOnPreCluster()
214{
215 /// Starting from a precluster, builds a pixel array, and then
216 /// extract clusters from this array
217
9bbd7f60 218 // Set saturation flag - it is not set if working directly with MC digits (w/out
219 // creating raw data) !!!
220 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
221 AliMUONPad* pad = fPreCluster->Pad(j);
222 if (pad->IsSaturated()) break;
223 if (pad->Charge() > fgkSaturation-1) pad->SetSaturated(kTRUE); //FIXME : remove usage of fgkSaturation
224 }
225
226 if (fDebug) {
227 cout << " *** Event # " << AliRunLoader::GetRunLoader()->GetEventNumber()
228 << " det. elem.: " << fDetElemId << endl;
229 for (Int_t j=0; j<fPreCluster->Multiplicity(); ++j) {
230 AliMUONPad* pad = fPreCluster->Pad(j);
231 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
232 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
233 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
234 }
235 }
236
c0a16418 237 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
c0a16418 238 if (!cluster) return kFALSE;
9bbd7f60 239
c0a16418 240 BuildPixArray(*cluster);
241
242 if ( fPixArray->GetLast() < 0 )
243 {
244 AliDebug(1,"No pixel for the above cluster");
245 delete cluster;
246 return kFALSE;
247 }
248
249 // Use MLEM for cluster finder
250 Int_t nMax = 1, localMax[100], maxPos[100];
251 Double_t maxVal[100];
252
253 if (cluster->Multiplicity() > 50)
254 {
255 nMax = FindLocalMaxima(fPixArray, localMax, maxVal);
256 }
257
258 if (nMax > 1)
259 {
260 TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
261 }
262
263 Int_t iSimple = 0, nInX = -1, nInY;
264
265 PadsInXandY(*cluster,nInX, nInY);
266
267 if (nMax == 1 && nInX < 4 && nInY < 4)
268 {
269 iSimple = 1; //1; // simple cluster
270 }
271
272 for (Int_t i=0; i<nMax; ++i)
273 {
274 if (nMax > 1)
275 {
276 FindCluster(*cluster,localMax, maxPos[i]);
277 }
278 Timer(kMainLoop)->Start(kFALSE);
279 MainLoop(*cluster,iSimple);
280 Timer(kMainLoop)->Stop();
281 if (i < nMax-1)
282 {
283 for (Int_t j=0; j<cluster->Multiplicity(); ++j)
284 {
285 AliMUONPad* pad = cluster->Pad(j);
286 if ( pad->Status() == 0 ) continue; // pad charge was not modified
287 pad->SetStatus(0);
288 pad->RevertCharge(); // use backup charge value
289 }
290 }
291 } // for (Int_t i=0; i<nMax;
292 if (nMax > 1) ((TH2D*) gROOT->FindObject("anode"))->Delete();
293 TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
294 if (mlem) mlem->Delete();
295 delete cluster;
296 return kTRUE;
297}
298
299//_____________________________________________________________________________
300Bool_t
301AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
302{
303 /// Check if the pad and the pixel overlaps
304
305 // make a fake pad from the pixel
306 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
307 pixel.Coord(0),pixel.Coord(1),
308 pixel.Size(0),pixel.Size(1),0);
309
310 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
311}
312
313//_____________________________________________________________________________
314AliMUONCluster*
315AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
316{
317 /// Check precluster in order to attempt to simplify it (mostly for
318 /// two-cathode preclusters)
319
320 if (origCluster.Multiplicity()==1)
321 {
322 // Disregard one-pad clusters (leftovers from splitting)
323 return 0x0;
324 }
325
326 Timer(kCheckPreCluster)->Start(kFALSE);
327
328
329 AliMUONCluster* cluster = static_cast<AliMUONCluster*>(origCluster.Clone());
330
c0a16418 331 AliDebug(2,"Start of CheckPreCluster=");
b1a19e07 332// StdoutToAliDebug(2,cluster->Print("full"));
c0a16418 333
334 // Check if one-cathode precluster
335 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
336 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
337
338 AliMUONCluster* rv(0x0);
339
340 if (i1 != i2)
341 {
342 rv = CheckPreclusterTwoCathodes(cluster);
343 }
344 else
345 {
346 rv = CheckPreclusterOneCathode(cluster);
347 }
348 Timer(kCheckPreCluster)->Stop();
349 return rv;
350}
351
352//_____________________________________________________________________________
353AliMUONCluster*
354AliMUONClusterFinderMLEM::CheckPreclusterOneCathode(AliMUONCluster* cluster)
355{
356 /// Check single-cathode precluster
357 AliWarning("Reimplement me!");
b1a19e07 358 AliDebug(2,"End of CheckPreClusterOneCathode=");
359// StdoutToAliDebug(2,cluster->Print("full"));
c0a16418 360
361 return cluster;
362}
363
364//_____________________________________________________________________________
365AliMUONCluster*
366AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
367{
71a2d3aa 368 /// Check two-cathode cluster
c0a16418 369
370 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
371 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
372
373 Int_t npad = cluster->Multiplicity();
374 Int_t* flags = new Int_t[npad];
375 memset(flags,0,npad*sizeof(Int_t));
376
377 // Check pad overlaps
378 for ( Int_t i=0; i<npad; ++i)
379 {
380 AliMUONPad* padi = cluster->Pad(i);
381 if ( padi->Cathode() != i1 ) continue;
382 for (Int_t j=i+1; j<npad; ++j)
383 {
384 AliMUONPad* padj = cluster->Pad(j);
385 if ( padj->Cathode() != i2 ) continue;
386 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
387 flags[i] = flags[j] = 1; // mark overlapped pads
388 }
389 }
390
391 // Check if all pads overlap
392 Int_t nFlags=0;
393 for (Int_t i=0; i<npad; ++i)
394 {
395 if (flags[i]) continue;
396 ++nFlags;
397 }
398
399 if (nFlags > 0)
400 {
401 // not all pads overlap.
9bbd7f60 402 if (fDebug) cout << " nFlags: " << nFlags << endl;
403 TObjArray toBeRemoved;
c0a16418 404 for (Int_t i=0; i<npad; ++i)
405 {
406 AliMUONPad* pad = cluster->Pad(i);
407 if (flags[i]) continue;
408 Int_t cath = pad->Cathode();
409 Int_t cath1 = TMath::Even(cath);
410 // Check for edge effect (missing pads on the _other_ cathode)
411 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
412 if (!mpPad.IsValid()) continue;
9bbd7f60 413 if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
c0a16418 414 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
415 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
9bbd7f60 416 toBeRemoved.AddLast(pad);
417 //AZ cluster->RemovePad(pad);
c0a16418 418 fPreCluster->Pad(i)->Release();
9bbd7f60 419 //AZ --npad;
420 }
421 for ( Int_t i = 0; i <= toBeRemoved.GetLast(); ++i )
422 {
423 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.At(i)));
c0a16418 424 }
425 }
426
427 // Check correlations of cathode charges
9bbd7f60 428 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
c0a16418 429 {
430 // big difference
431 Int_t cathode = cluster->MaxRawChargeCathode();
432 Int_t imin(0);
433 Int_t imax(0);
434 Double_t cmax(0);
435 Double_t cmin(1E9);
436
437 // get min and max pad charges on the cathode opposite to the
438 // max pad (given by MaxRawChargeCathode())
439 //
440 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
441 {
442 AliMUONPad* pad = cluster->Pad(i);
443 if ( pad->Cathode() != cathode || !pad->IsReal() )
444 {
445 // only consider pads in the opposite cathode, and
446 // onyl consider real pads (i.e. exclude the virtual ones)
447 continue;
448 }
449 if ( pad->Charge() < cmin )
450 {
451 cmin = pad->Charge();
452 imin = i;
453 }
454 if ( pad->Charge() > cmax )
455 {
456 cmax = pad->Charge();
457 imax = i;
458 }
459 }
460 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
461 imin,imax,cmin,cmax));
462 //
463 // arrange pads according to their distance to the max, normalized
464 // to the pad size
465 Double_t* dist = new Double_t[cluster->Multiplicity()];
466 Double_t dxMin(1E9);
467 Double_t dyMin(1E9);
468 Double_t dmin(0);
469
470 AliMUONPad* padmax = cluster->Pad(imax);
471
472 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
473 {
474 dist[i] = 0.0;
475 if ( i == imax) continue;
476 AliMUONPad* pad = cluster->Pad(i);
477 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
478 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
479 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
480 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
481 if ( i == imin )
482 {
483 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
484 dxMin = dx;
485 dyMin = dy;
486 }
487 }
488
489 TMath::Sort(cluster->Multiplicity(),dist,flags,kFALSE); // in ascending order
490 Double_t xmax(-1);
491 TObjArray toBeRemoved;
492
493 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
494 {
495 Int_t indx = flags[i];
496 AliMUONPad* pad = cluster->Pad(indx);
497 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
498 if ( dist[indx] > dmin )
499 {
500 // farther than the minimum pad
501 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
502 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
503 dx *= dxMin;
504 dy *= dyMin;
505 if (dx >= 0 && dy >= 0) continue;
506 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
507 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
508 }
509 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
510 {
511 // release pad
512 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
513 {
514 cmax = TMath::Max(pad->Charge(),cmax);
515 }
516 else
517 {
518 cmax = pad->Charge();
519 }
520 xmax = dist[indx];
521 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
522 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
523 pad->Charge()));
524
525 toBeRemoved.AddLast(pad);
526 fPreCluster->Pad(indx)->Release();
527 }
528 }
529 for ( Int_t i = 0; i <= toBeRemoved.GetLast(); ++i )
530 {
531 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.At(i)));
532 }
533 delete[] dist;
534 }
535
536 delete[] flags;
537
538 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
b1a19e07 539// StdoutToAliDebug(2,cluster->Print("full"));
c0a16418 540
541 return cluster;
542}
543
544//_____________________________________________________________________________
545void
546AliMUONClusterFinderMLEM::CheckOverlaps()
547{
548 /// For debug only : check if some pixels overlap...
549
550 Int_t nPix = fPixArray->GetLast()+1;
551 Int_t dummy(0);
552
553 for ( Int_t i = 0; i < nPix; ++i )
554 {
555 AliMUONPad* pixelI = Pixel(i);
556 AliMUONPad pi(dummy,dummy,dummy,dummy,
557 pixelI->Coord(0),pixelI->Coord(1),
558 pixelI->Size(0),pixelI->Size(1),0.0);
559
560 for ( Int_t j = i+1; j < nPix; ++j )
561 {
562 AliMUONPad* pixelJ = Pixel(j);
563 AliMUONPad pj(dummy,dummy,dummy,dummy,
564 pixelJ->Coord(0),pixelJ->Coord(1),
565 pixelJ->Size(0),pixelJ->Size(1),0.0);
566 AliMpArea area;
567
568 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
569 {
570 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
571 StdoutToAliInfo(pixelI->Print();
572 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
573 pixelJ->Print();
574 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
575 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
576 cout << "-------" << endl;
577 );
578
579 }
580 }
581 }
582}
583
584//_____________________________________________________________________________
585void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
586{
587 /// Build pixel array for MLEM method
588
589 Int_t npad = cluster.Multiplicity();
590 if (npad<=0)
591 {
592 AliWarning("Got no pad at all ?!");
593 }
594
595 fPixArray->Delete();
596
597 if ( !cluster.Multiplicity(0) || !cluster.Multiplicity(1) )
598 {
599 BuildPixArrayOneCathode(cluster);
600 }
601 else
602 {
603 BuildPixArrayTwoCathodes(cluster);
604 }
605
606 fPixArray->Sort(); // FIXME : not really needed, only to compare with ClusterFinderAZ
607
608 Int_t nPix = fPixArray->GetLast()+1;
609
b1a19e07 610 AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
c0a16418 611
612 Double_t xPadMin(1E9);
613 Double_t yPadMin(1E9);
614
615 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i )
616 {
617 AliMUONPad* pad = cluster.Pad(i);
618 xPadMin = TMath::Min (xPadMin, pad->DX());
619 yPadMin = TMath::Min (yPadMin, pad->DY());
620 }
621
622 Double_t wxmin(1E9);
623 Double_t wymin(1E9);
624
625 for ( Int_t i = 0; i < nPix; ++i )
626 {
627 AliMUONPad* pixPtr = Pixel(i);
628 wxmin = TMath::Min(wxmin, pixPtr->Size(0));
629 wymin = TMath::Min(wymin, pixPtr->Size(1));
630 }
631
632 wxmin = TMath::Abs (wxmin - xPadMin/2) > 0.001 ? xPadMin : xPadMin / 2;
633 wymin = TMath::Abs (wymin - yPadMin/2) > 0.001 ? yPadMin : yPadMin / 2;
634
635 // Check if small pixel X-size
636 AdjustPixel(cluster,wxmin, 0);
637 // Check if small pixel Y-size
638 AdjustPixel(cluster,wymin, 1);
639 // Check if large pixel size
640 AdjustPixel(wxmin, wymin);
641
642 // Remove discarded pixels
643 for (Int_t i=0; i<nPix; ++i)
644 {
645 AliMUONPad* pixPtr = Pixel(i);
646 if (pixPtr->Charge() < 1)
647 {
648 AliDebug(2,Form("Removing pixel %d with charge<1 : ",i));
b1a19e07 649// StdoutToAliDebug(2,pixPtr->Print());
c0a16418 650 RemovePixel(i);
651 }
652 }
653
654 fPixArray->Compress();
655 nPix = fPixArray->GetEntriesFast();
656
b1a19e07 657 AliDebug(2,Form("nPix after AdjustPixel=%d",nPix));
c0a16418 658
659 if ( nPix > cluster.Multiplicity() )
660 {
b1a19e07 661 AliDebug(2,Form("Will trim number of pixels to number of pads"));
c0a16418 662
663 // Too many pixels - sort and remove pixels with the lowest signal
664 fPixArray->Sort();
665 for ( Int_t i = cluster.Multiplicity(); i<nPix; ++i )
666 {
667 RemovePixel(i);
668 }
669 fPixArray->Compress();
670 nPix = fPixArray->GetEntriesFast();
671 } // if (nPix > npad)
672
673// StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
674// fPixArray->Print(););
b1a19e07 675// CheckOverlaps();//FIXME : this is for debug only. Remove it.
c0a16418 676}
677
678//_____________________________________________________________________________
679void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
680{
681 /// From a single-cathode cluster, build the pixel array
682
b1a19e07 683 AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
c0a16418 684
685 for ( Int_t j=0; j<cluster.Multiplicity(); ++j)
686 {
687 AliMUONPad* pad = cluster.Pad(j);
688 AliMUONPad* pixPtr = new AliMUONPad(pad->Position(),pad->Dimensions(),
689 pad->Charge());
690 fPixArray->Add(pixPtr);
691 }
692}
693
694//_____________________________________________________________________________
695void AliMUONClusterFinderMLEM::BuildPixArrayTwoCathodes(AliMUONCluster& cluster)
696{
697 /// From a two-cathodes cluster, build the pixel array
698
b1a19e07 699 AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
c0a16418 700
701 Int_t i1 = cluster.Pad(0)->Cathode();
702 Int_t i2 = TMath::Even(i1);
703
704 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i)
705 {
706 AliMUONPad* padi = cluster.Pad(i);
707 if (padi->Cathode() != i1) continue;
708
709 for ( Int_t j = 1; j < cluster.Multiplicity(); ++j)
710 {
711 AliMUONPad* padj = cluster.Pad(j);
712 if (padj->Cathode() != i2) continue;
713
714 AliMpArea overlap;
715
716 if ( AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize,overlap) )
717 {
718 AliMUONPad* pixPtr = new AliMUONPad(overlap.Position(),
719 overlap.Dimensions(),
720 TMath::Min(padi->Charge(),padj->Charge()));
721 if ( ( padi->Charge() <= padj->Charge() && padi->IsSaturated() ) ||
722 ( padj->Charge() < padi->Charge() && padj->IsSaturated() ) )
723 {
724 // if smallest charge of the 2 pads is already above saturation, then
725 // the pixel is saturated...
726 pixPtr->SetSaturated(kTRUE);
727 }
728 pixPtr->SetReal(kFALSE);
729 fPixArray->Add(pixPtr);
730 }
731 }
732 }
733}
734
735//_____________________________________________________________________________
736void AliMUONClusterFinderMLEM::AdjustPixel(AliMUONCluster& cluster,
737 Float_t width, Int_t ixy)
738{
739 /// Check if some pixels have small size (adjust if necessary)
740
741 AliDebug(2,Form("width=%e ixy=%d",width,ixy));
742
743 AliMUONPad *pixPtr, *pixPtr1 = 0;
744 Int_t ixy1 = TMath::Even(ixy);
745 Int_t nPix = fPixArray->GetEntriesFast();
746
747 for (Int_t i=0; i<nPix; i++)
748 {
749 pixPtr = Pixel(i);
750 if (pixPtr->Charge() < 1) continue; // discarded pixel
751 if (pixPtr->Size(ixy)-width < -1.e-4)
752 {
753 // try to merge
754 for (Int_t j=i+1; j<nPix; j++)
755 {
756 pixPtr1 = Pixel(j);
757 if (pixPtr1->Charge() < 1) continue; // discarded pixel
758 if (TMath::Abs(pixPtr1->Size(ixy)-width) < fgkDistancePrecision) continue; // right size
759 if (TMath::Abs(pixPtr1->Coord(ixy1)-pixPtr->Coord(ixy1)) > fgkDistancePrecision) continue; // different rows/columns
760 if (TMath::Abs(pixPtr1->Coord(ixy)-pixPtr->Coord(ixy)) < 2*width)
761 {
762 // merge
763 Double_t tmp = pixPtr->Coord(ixy) + pixPtr1->Size(ixy)*
764 TMath::Sign (1., pixPtr1->Coord(ixy) - pixPtr->Coord(ixy));
765 pixPtr->SetCoord(ixy, tmp);
766 pixPtr->SetSize(ixy, width);
767 pixPtr->SetCharge(TMath::Min (pixPtr->Charge(),pixPtr1->Charge()));
768 pixPtr1->SetCharge(0);
769 pixPtr1 = 0;
770 break;
771 }
772 } // for (Int_t j=i+1;
773 if (pixPtr1 || i == nPix-1)
774 {
775 // edge pixel - just increase its size
776 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
777 {
778 AliMUONPad* pad = cluster.Pad(j);
779 Double_t d = ( ixy == 0 ) ? pad->X() : ( ixy == 1 ) ? pad->Y() : -1E9;
780
781 if (pixPtr->Coord(ixy) < d)
782 {
783 pixPtr->Shift(ixy, pixPtr->Size(ixy)-width);
784 }
785 else
786 {
787 pixPtr->Shift(ixy, -pixPtr->Size(ixy)+width);
788 }
789 pixPtr->SetSize(ixy, width);
790 break;
791 }
792 }
793 } // if (pixPtr->Size(ixy)-width < -1.e-4)
794 } // for (Int_t i=0; i<nPix;
795}
796
797//_____________________________________________________________________________
798void AliMUONClusterFinderMLEM::AdjustPixel(Double_t wxmin, Double_t wymin)
799{
800/// Check if some pixels have large size (adjust if necessary)
801
802 AliDebug(2,Form("wxmin=%e wymin=%e",wxmin,wymin));
803
804 Int_t n1[2], n2[2], iOK = 1, nPix = fPixArray->GetEntriesFast();
805 AliMUONPad *pixPtr, pix;
806 Double_t xy0[2] = {9999, 9999}, wxy[2], dist[2] = {0};
807
808 // Check if large pixel size
809 for (Int_t i = 0; i < nPix; i++) {
810 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
811 if (pixPtr->Charge() < 1) continue; // discarded pixel
812 if (pixPtr->Size(0) - wxmin < 1.e-4) {
813 if (xy0[0] > 9998) xy0[0] = pixPtr->Coord(0); // position of a "normal" pixel
814 if (pixPtr->Size(1) - wymin < 1.e-4) {
815 if (xy0[1] > 9998) xy0[1] = pixPtr->Coord(1); // position of a "normal" pixel
816 continue;
817 } else iOK = 0; // large pixel
818 } else {
819 iOK = 0; // large pixel
820 if (xy0[1] > 9998 && pixPtr->Size(1) - wymin < 1.e-4) xy0[1] = pixPtr->Coord(1); // "normal" pixel
821 }
822 if (xy0[0] < 9998 && xy0[1] < 9998) break;
823 }
824 if (iOK) return;
825
826 wxy[0] = wxmin;
827 wxy[1] = wymin;
828 //cout << xy0[0] << " " << xy0[1] << endl;
829 for (Int_t i = 0; i < nPix; i++) {
830 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
831 if (pixPtr->Charge() < 1) continue; // discarded pixel
832 n1[0] = n1[1] = 999;
833 n2[0] = n2[1] = 1;
834 for (Int_t j = 0; j < 2; j++) {
835 if (pixPtr->Size(j) - wxy[j] < 1.e-4) continue;
836 dist[j] = (pixPtr->Coord(j) - xy0[j]) / wxy[j] / 2; // normalized distance to "normal" pixel
837 n2[j] = TMath::Nint (pixPtr->Size(j) / wxy[j]);
838 n1[j] = n2[j] == 1 ? TMath::Nint(dist[j]) : (Int_t)dist[j];
839 }
840 if (n1[0] > 998 && n1[1] > 998) continue;
841 if (fDebug) cout << " Different " << pixPtr->Size(0) << " " << wxy[0] << " "
842 << pixPtr->Size(1) << " " << wxy[1] <<endl;
843
844 if (n2[0] > 2 || n2[1] > 2) {
845 //cout << n2[0] << " " << n2[1] << endl;
846 if (n2[0] > 2 && n1[0] < 999) n1[0]--;
847 if (n2[1] > 2 && n1[1] < 999) n1[1]--;
848 }
849 //cout << n1[0] << " " << n2[0] << " " << n1[1] << " " << n2[1] << endl;
850 pix = *pixPtr;
851 pix.SetSize(0, wxy[0]); pix.SetSize(1, wxy[1]);
852 //pixPtr->Print();
853 for (Int_t ii = 0; ii < n2[0]; ii++) {
854 if (n1[0] < 999) pix.SetCoord(0, xy0[0] + (n1[0] + TMath::Sign(1.,dist[0]) * ii) * 2 * wxy[0]);
855 for (Int_t jj = 0; jj < n2[1]; jj++) {
856 if (n1[1] < 999) pix.SetCoord(1, xy0[1] + (n1[1] + TMath::Sign(1.,dist[1]) * jj) * 2 * wxy[1]);
857 fPixArray->Add(new AliMUONPad(pix));
858 //pix.Print();
859 }
860 }
861 pixPtr->SetCharge(0);
862 } // for (Int_t i = 0; i < nPix;
863}
864
865//_____________________________________________________________________________
866void
867AliMUONClusterFinderMLEM::Plot(const char* basename)
868{
869 /// Make a plot and save it as png
870
9bbd7f60 871 return; //AZ
c0a16418 872 if (!fPlot) return;
873
874 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
875 c->Draw();
876 Draw();
877 c->Modified();
878 c->Update();
879 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
880 fDetElemId,fClusterNumber));
881}
882
883//_____________________________________________________________________________
884void
885AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
886 Double_t* coef,
887 Double_t* probi)
888{
889 /// Compute coefficients needed for MLEM algorithm
890
891 Int_t nPix = fPixArray->GetLast()+1;
892
893 memset(probi,0,nPix*sizeof(Double_t));
894
895 for ( Int_t j=0; j<cluster.Multiplicity(); ++j )
896 {
897 AliMUONPad* pad = cluster.Pad(j);
898 Int_t indx = j*nPix;
899
900 for ( Int_t ipix=0; ipix<nPix; ++ipix )
901 {
902 Int_t indx1 = indx + ipix;
903 if (pad->Status() < 0)
904 {
905 coef[indx1] = 0;
906 continue;
907 }
908 AliMUONPad* pixPtr = Pixel(ipix);
909 // coef is the charge (given by Mathieson integral) on pad, assuming
910 // the Mathieson is center at pixel.
911 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
b1a19e07 912 AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
913 pad->Ix(),pad->Iy(),
914 pad->X(),pad->Y(),
915 pad->DX(),pad->DY(),
916 pixPtr->Coord(0),pixPtr->Coord(1),
917 pixPtr->Size(0),pixPtr->Size(1),
918 coef[indx1]));
c0a16418 919
920 probi[ipix] += coef[indx1];
921 }
922 }
923}
924
925//_____________________________________________________________________________
926Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
927{
928 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
929
930 Int_t nPix = fPixArray->GetLast()+1;
931
932 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
b1a19e07 933// StdoutToAliDebug(2,cluster.Print("full"););
c0a16418 934
935 if ( nPix < 0 )
936 {
937 AliDebug(1,"No pixels, why am I here then ?");
938 }
939
940 AddVirtualPad(cluster); // add virtual pads if necessary
941
942 Int_t npadTot = cluster.Multiplicity();
943 Int_t npadOK = 0;
944 for (Int_t i = 0; i < npadTot; ++i)
945 {
946 if (cluster.Pad(i)->Status() == 0) ++npadOK;
947 }
948
949 TH2D* mlem(0x0);
950 Double_t* coef(0x0);
951 Double_t* probi(0x0);
952 Int_t lc(0); // loop counter (for debug)
953
b1a19e07 954// Plot("mlem.start");
c0a16418 955
956 while (1)
957 {
958 ++lc;
959 mlem = (TH2D*) gROOT->FindObject("mlem");
960 delete mlem;
961
962 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
963 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
b1a19e07 964// StdoutToAliDebug(2,fPixArray->Print("","full"));
c0a16418 965
966 coef = new Double_t [npadTot*nPix];
967 probi = new Double_t [nPix];
968
969 // Calculate coefficients and pixel visibilities
970 ComputeCoefficients(cluster,coef,probi);
971
972 for (Int_t ipix=0; ipix<nPix; ++ipix)
973 {
974 if (probi[ipix] < 0.01)
975 {
976 AliMUONPad* pixel = Pixel(ipix);
977 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
b1a19e07 978// StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
c0a16418 979 pixel->SetCharge(0); // "invisible" pixel
980 }
981 }
982
983 // MLEM algorithm
984 Mlem(cluster,coef, probi, 15);
985
986 Double_t xylim[4] = {999, 999, 999, 999};
987 AliMUONPad* pixPtr(0x0);
988
989 for ( Int_t ipix=0; ipix<nPix; ++ipix )
990 {
991 pixPtr = Pixel(ipix);
992 for ( Int_t i=0; i<4; ++i )
993 {
994 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
995 }
996 }
997 for (Int_t i=0; i<4; i++)
998 {
999 xylim[i] -= pixPtr->Size(i/2);
1000 }
1001
1002
1003 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1004 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1005
1006// StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
1007 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
1008 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
1009 xylim[0],-xylim[1],xylim[2],-xylim[3]
1010 ));
1011
1012 mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1013
1014 for (Int_t ipix=0; ipix<nPix; ++ipix)
1015 {
1016 AliMUONPad* pixPtr = Pixel(ipix);
1017 mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
1018 }
1019
1020 // Check if the total charge of pixels is too low
1021 Double_t qTot = 0;
1022 for ( Int_t i=0; i<nPix; ++i)
1023 {
1024 qTot += Pixel(i)->Charge();
1025 }
1026
1027 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
1028 {
1029 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
1030 delete [] coef;
1031 delete [] probi;
1032 coef = 0;
1033 probi = 0;
1034 fPixArray->Delete();
1035 for ( Int_t i=0; i<npadTot; ++i)
1036 {
1037 AliMUONPad* pad = cluster.Pad(i);
1038 if ( pad->Status() == 0) pad->SetStatus(-1);
1039 }
1040 return kFALSE;
1041 }
1042
1043 if (iSimple)
1044 {
1045 // Simple cluster - skip further passes thru EM-procedure
1046 Simple(cluster);
1047 delete [] coef;
1048 delete [] probi;
1049 coef = 0;
1050 probi = 0;
1051 fPixArray->Delete();
1052 return kTRUE;
1053 }
1054
1055 // Calculate position of the center-of-gravity around the maximum pixel
1056 Double_t xyCOG[2];
1057 FindCOG(mlem, xyCOG);
1058
1059 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
1060 pixPtr->Size(0) > pixPtr->Size(1)) break;
1061
1062 // Sort pixels according to the charge
1063 fPixArray->Sort();
1064 Double_t pixMin = 0.01*Pixel(0)->Charge();
1065 pixMin = TMath::Min(pixMin,50.);
1066
1067 // Decrease pixel size and shift pixels to make them centered at
1068 // the maximum one
1069 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
1070 Int_t ix(1);
1071 Double_t width = 0;
1072 Double_t shift[2] = { 0.0, 0.0 };
1073 for (Int_t i=0; i<4; i++) xylim[i] = 999;
1074 Int_t nPix1 = nPix;
1075 nPix = 0;
1076 for (Int_t ipix=0; ipix<nPix1; ++ipix)
1077 {
1078 AliMUONPad* pixPtr = Pixel(ipix);
1079 if ( nPix >= npadOK // too many pixels already
1080 ||
1081 pixPtr->Charge() < pixMin // too low charge
1082 )
1083 {
1084 RemovePixel(ipix);
1085 continue;
1086 }
1087 for (Int_t i=0; i<2; ++i)
1088 {
1089 if (!i)
1090 {
1091 pixPtr->SetCharge(10);
1092 pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
1093 width = -pixPtr->Size(indx);
1094 pixPtr->Shift(indx, width);
1095 // Shift pixel position
1096 if (ix)
1097 {
1098 ix = 0;
1099 for (Int_t j=0; j<2; ++j)
1100 {
1101 shift[j] = pixPtr->Coord(j) - xyCOG[j];
1102 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
1103 }
1104 } // if (ix)
1105 pixPtr->Shift(0, -shift[0]);
1106 pixPtr->Shift(1, -shift[1]);
1107 }
1108 else
1109 {
1110 pixPtr = new AliMUONPad(*pixPtr);
1111 pixPtr->Shift(indx, -2*width);
1112 fPixArray->Add(pixPtr);
1113 }
1114 for (Int_t i=0; i<4; ++i)
1115 {
1116 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1117 }
1118 } // for (Int_t i=0; i<2;
1119 nPix += 2;
1120 } // for (Int_t ipix=0;
1121
1122 fPixArray->Compress();
1123 nPix = fPixArray->GetEntriesFast();
1124
1125 AliDebug(2,Form("After shift:"));
b1a19e07 1126// StdoutToAliDebug(2,fPixArray->Print("","full"););
1127// Plot(Form("mlem.lc%d",lc+1));
c0a16418 1128
1129 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1130 xyCOG[0],xyCOG[1],
1131 xylim[0],xylim[1],
1132 xylim[2],xylim[3]));
1133
1134 // Remove excessive pixels
1135 if (nPix > npadOK)
1136 {
1137 for (Int_t ipix=npadOK; ipix<nPix; ++ipix)
1138 {
1139 RemovePixel(ipix);
1140 }
1141 }
1142 else
1143 {
1144 AliMUONPad* pixPtr = Pixel(0);
1145 // add pixels if the maximum is at the limit of pixel area
1146 // start from Y-direction
1147 Int_t j = 0;
1148 for (Int_t i=3; i>-1; --i)
1149 {
1150 if (nPix < npadOK &&
1151 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2))
1152 {
1153 AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1154 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1155 j = TMath::Even (i/2);
1156 p->SetCoord(j, xyCOG[j]);
1157 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
b1a19e07 1158// StdoutToAliDebug(2,cout << " ---- ";
1159// p->Print("corners"););
c0a16418 1160 fPixArray->Add(p);
1161 ++nPix;
1162 }
1163 }
1164 }
1165 fPixArray->Compress();
1166 nPix = fPixArray->GetEntriesFast();
1167 delete [] coef;
1168 delete [] probi;
1169 coef = 0;
1170 probi = 0;
1171 } // while (1)
1172
1173 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
b1a19e07 1174// StdoutToAliDebug(2,fPixArray->Print("","full"););
c0a16418 1175
1176 // remove pixels with low signal or low visibility
1177 // Cuts are empirical !!!
1178 Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1179 thresh = TMath::Min (thresh,50.);
1180 Double_t cmax = -1;
1181 Double_t charge = 0;
1182
1183 for ( Int_t i=0; i<nPix; ++i)
1184 {
1185 cmax = TMath::Max (cmax,probi[i]);
1186 }
1187
1188 // Mark pixels which should be removed
1189 for (Int_t i=0; i<nPix; ++i)
1190 {
1191 AliMUONPad* pixPtr = Pixel(i);
1192 charge = pixPtr->Charge();
1193 if (charge < thresh)
1194 {
1195 pixPtr->SetCharge(-charge);
1196 }
1197 }
1198
1199 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1200 Int_t near = 0;
1201 for (Int_t i=0; i<nPix; ++i)
1202 {
1203 AliMUONPad* pixPtr = Pixel(i);
1204 charge = pixPtr->Charge();
1205 if (charge > 0) continue;
1206 near = FindNearest(pixPtr);
1207 pixPtr->SetCharge(0);
1208 probi[i] = 0; // make it "invisible"
1209 AliMUONPad* pnear = Pixel(near);
1210 pnear->SetCharge(pnear->Charge() + (-charge));
1211 }
1212 Mlem(cluster,coef,probi,2);
1213
1214 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
b1a19e07 1215// StdoutToAliDebug(2,fPixArray->Print("","full"););
1216// Plot("mlem.beforesplit");
c0a16418 1217
1218 // Update histogram
1219 for (Int_t i=0; i<nPix; ++i)
1220 {
1221 AliMUONPad* pixPtr = Pixel(i);
1222 Int_t ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1223 Int_t iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1224 mlem->SetBinContent(ix, iy, pixPtr->Charge());
1225 }
1226
1227 // Try to split into clusters
1228 Bool_t ok = kTRUE;
1229 if (mlem->GetSum() < 1)
1230 {
1231 ok = kFALSE;
1232 }
1233 else
1234 {
1235 fSplitter->Split(cluster,mlem,coef,fClusterList);
1236 }
1237
1238 delete [] coef;
1239 delete [] probi;
1240 coef = 0;
1241 probi = 0;
1242 fPixArray->Delete();
1243
1244 return ok;
1245}
1246
1247//_____________________________________________________________________________
1248void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1249 Double_t* coef, Double_t* probi,
1250 Int_t nIter)
1251{
1252 /// Use MLEM to find pixel charges
1253
1254 Int_t nPix = fPixArray->GetEntriesFast();
1255
1256 Int_t npad = cluster.Multiplicity();
1257
1258 Double_t* probi1 = new Double_t[nPix];
1259 Double_t probMax = 0;
1260 Double_t tmp = TMath::MaxElement(nPix,probi);
1261
1262 for (Int_t ipix=0; ipix<nPix; ++ipix)
1263 {
1264 probMax = TMath::Max(probMax,probi[ipix]);
1265 }
1266
1267 if (probMax!=tmp) { AliWarning(Form("probMax=%e tmp=%e",probMax,tmp)); }
1268
1269 for (Int_t iter=0; iter<nIter; ++iter)
1270 {
1271 // Do iterations
1272 for (Int_t ipix=0; ipix<nPix; ++ipix)
1273 {
9bbd7f60 1274 Pixel(ipix)->SetChargeBackup(0);
c0a16418 1275 // Correct each pixel
1276 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1277 Double_t sum = 0;
1278 probi1[ipix] = probMax;
1279 for (Int_t j=0; j<npad; j++)
1280 {
1281 AliMUONPad* pad = cluster.Pad(j);
1282 if (pad->Status() < 0) continue;
1283 Double_t sum1 = 0;
1284 Int_t indx1 = j*nPix;
1285 Int_t indx = indx1 + ipix;
1286 // Calculate expectation
1287 for (Int_t i=0; i<nPix; i++)
1288 {
1289 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1290 }
1291 if ( pad->Charge() > fgkSaturation-1 && sum1 > pad->Charge() ) //FIXME : remove usage of fgkSaturation
1292 {
1293 if ( !pad->IsSaturated() )
1294 {
1295 AliWarning("Got a pad charge above saturation not backed-up by pad->IsSaturated() function : ");
1296 StdoutToAliWarning(pad->Print("full"));
1297 }
1298 // correct for pad charge overflows
1299 probi1[ipix] -= coef[indx];
1300 continue;
1301 }
1302
1303 if (coef[indx] > 1.e-6)
1304 {
1305 sum += pad->Charge()*coef[indx]/sum1;
1306 }
1307 } // for (Int_t j=0;
1308 AliMUONPad* pixPtr = Pixel(ipix);
1309 if (probi1[ipix] > 1.e-6)
1310 {
9bbd7f60 1311 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1312 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
c0a16418 1313 }
1314 } // for (Int_t ipix=0;
9bbd7f60 1315 for (Int_t i = 0; i < nPix; ++i) {
1316 AliMUONPad* pixPtr = Pixel(i);
1317 pixPtr->SetCharge(pixPtr->ChargeBackup());
1318 }
c0a16418 1319 } // for (Int_t iter=0;
1320 delete [] probi1;
1321}
1322
1323//_____________________________________________________________________________
1324void AliMUONClusterFinderMLEM::FindCOG(TH2D *mlem, Double_t *xyc)
1325{
1326 /// Calculate position of the center-of-gravity around the maximum pixel
1327
1328 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1329 Int_t i1 = -9, j1 = -9;
1330 mlem->GetMaximumBin(ixmax,iymax,ix);
1331 Int_t nx = mlem->GetNbinsX();
1332 Int_t ny = mlem->GetNbinsY();
1333 Double_t thresh = mlem->GetMaximum()/10;
1334 Double_t x, y, cont, xq=0, yq=0, qq=0;
1335
1336 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1337 y = mlem->GetYaxis()->GetBinCenter(i);
1338 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1339 cont = mlem->GetCellContent(j,i);
1340 if (cont < thresh) continue;
1341 if (i != i1) {i1 = i; nsumy++;}
1342 if (j != j1) {j1 = j; nsumx++;}
1343 x = mlem->GetXaxis()->GetBinCenter(j);
1344 xq += x*cont;
1345 yq += y*cont;
1346 qq += cont;
1347 nsum++;
1348 }
1349 }
1350
1351 Double_t cmax = 0;
1352 Int_t i2 = 0, j2 = 0;
1353 x = y = 0;
1354 if (nsumy == 1) {
1355 // one bin in Y - add one more (with the largest signal)
1356 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1357 if (i == iymax) continue;
1358 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1359 cont = mlem->GetCellContent(j,i);
1360 if (cont > cmax) {
1361 cmax = cont;
1362 x = mlem->GetXaxis()->GetBinCenter(j);
1363 y = mlem->GetYaxis()->GetBinCenter(i);
1364 i2 = i;
1365 j2 = j;
1366 }
1367 }
1368 }
1369 xq += x*cmax;
1370 yq += y*cmax;
1371 qq += cmax;
1372 if (i2 != i1) nsumy++;
1373 if (j2 != j1) nsumx++;
1374 nsum++;
1375 } // if (nsumy == 1)
1376
1377 if (nsumx == 1) {
1378 // one bin in X - add one more (with the largest signal)
1379 cmax = x = y = 0;
1380 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1381 if (j == ixmax) continue;
1382 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1383 cont = mlem->GetCellContent(j,i);
1384 if (cont > cmax) {
1385 cmax = cont;
1386 x = mlem->GetXaxis()->GetBinCenter(j);
1387 y = mlem->GetYaxis()->GetBinCenter(i);
1388 i2 = i;
1389 j2 = j;
1390 }
1391 }
1392 }
1393 xq += x*cmax;
1394 yq += y*cmax;
1395 qq += cmax;
1396 if (i2 != i1) nsumy++;
1397 if (j2 != j1) nsumx++;
1398 nsum++;
1399 } // if (nsumx == 1)
1400
1401 xyc[0] = xq/qq; xyc[1] = yq/qq;
1402 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1403}
1404
1405//_____________________________________________________________________________
1406Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1407{
1408/// Find the pixel nearest to the given one
1409/// (algorithm may be not very efficient)
1410
1411 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1412 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1413 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1414 AliMUONPad *pixPtr;
1415
1416 for (Int_t i=0; i<nPix; i++) {
1417 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1418 if (pixPtr->Charge() < 0.5) continue;
1419 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1420 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1421 r = dx *dx + dy * dy;
1422 if (r < rmin) { rmin = r; imin = i; }
1423 }
1424 return imin;
1425}
1426
1427
1428//_____________________________________________________________________________
1429TStopwatch*
1430AliMUONClusterFinderMLEM::Timer(Int_t i) const
1431{
1432 /// Return timer at index i
1433 return static_cast<TStopwatch*>(fTimers->At(i));
1434}
1435
1436//_____________________________________________________________________________
1437void
1438AliMUONClusterFinderMLEM::Paint(Option_t*)
1439{
1440 /// Paint cluster and pixels
1441
1442 AliMpArea area(fPreCluster->Area());
1443
1444 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1445
1446 gVirtualX->SetFillStyle(1001);
1447 gVirtualX->SetFillColor(3);
1448 gVirtualX->SetLineColor(3);
1449
1450 Double_t s(1.0);
1451
1452 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1453 {
1454 AliMUONPad* pixel = Pixel(i);
1455
1456 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1457 pixel->Coord(1)-pixel->Size(1)*s,
1458 pixel->Coord(0)+pixel->Size(0)*s,
1459 pixel->Coord(1)+pixel->Size(1)*s);
1460
1461// for ( Int_t sign = -1; sign < 2; sign +=2 )
1462// {
1463// gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1464// pixel->Coord(1) + sign*pixel->Size(1),
1465// pixel->Coord(0) + pixel->Size(0),
1466// pixel->Coord(1) - sign*pixel->Size(1)
1467// );
1468// }
1469 }
1470
1471
1472 gVirtualX->SetFillStyle(0);
1473
1474 fPreCluster->Paint();
1475
1476 gVirtualX->SetLineColor(1);
1477 gVirtualX->SetLineWidth(2);
1478 gVirtualX->SetFillStyle(0);
1479 gVirtualX->SetTextColor(1);
1480 gVirtualX->SetTextAlign(22);
1481
1482 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1483 {
1484 AliMUONPad* pixel = Pixel(i);
1485 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1486 pixel->Coord(1)-pixel->Size(1),
1487 pixel->Coord(0)+pixel->Size(0),
1488 pixel->Coord(1)+pixel->Size(1));
1489 gVirtualX->SetTextSize(pixel->Size(0)*60);
1490
1491 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1492 }
1493}
1494
1495//_____________________________________________________________________________
1496Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1497{
1498/// Find local maxima in pixel space for large preclusters in order to
1499/// try to split them into smaller pieces (to speed up the MLEM procedure)
1500/// or to find additional fitting seeds if clusters were not completely resolved
1501
1502 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1503
1504 TH2D *hist = NULL;
1505 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
1506 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
1507 //if (hist) hist->Delete();
1508
1509 Double_t xylim[4] = {999, 999, 999, 999};
1510
1511 Int_t nPix = pixArray->GetEntriesFast();
1512 AliMUONPad *pixPtr = 0;
1513 for (Int_t ipix=0; ipix<nPix; ipix++) {
1514 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1515 for (Int_t i=0; i<4; i++)
1516 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1517 }
1518 for (Int_t i=0; i<4; i++) xylim[i] -= pixPtr->Size(i/2);
1519
1520 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1521 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1522 if (pixArray == fPixArray) hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1523 else hist = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1524 for (Int_t ipix=0; ipix<nPix; ipix++) {
1525 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1526 hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1527 }
1528// if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1529
1530 Int_t nMax = 0, indx;
1531 Int_t *isLocalMax = new Int_t[ny*nx];
1532 for (Int_t i=0; i<ny*nx; i++) isLocalMax[i] = 0;
1533
1534 for (Int_t i=1; i<=ny; i++) {
1535 indx = (i-1) * nx;
1536 for (Int_t j=1; j<=nx; j++) {
1537 if (hist->GetCellContent(j,i) < 0.5) continue;
1538 //if (isLocalMax[indx+j-1] < 0) continue;
1539 if (isLocalMax[indx+j-1] != 0) continue;
1540 FlagLocalMax(hist, i, j, isLocalMax);
1541 }
1542 }
1543
1544 for (Int_t i=1; i<=ny; i++) {
1545 indx = (i-1) * nx;
1546 for (Int_t j=1; j<=nx; j++) {
1547 if (isLocalMax[indx+j-1] > 0) {
1548 localMax[nMax] = indx + j - 1;
1549 maxVal[nMax++] = hist->GetCellContent(j,i);
1550 if (nMax > 99) AliFatal(" Too many local maxima !!!");
1551 }
1552 }
1553 }
1554 if (fDebug) cout << " Local max: " << nMax << endl;
1555 delete [] isLocalMax; isLocalMax = 0;
1556 return nMax;
1557}
1558
1559//_____________________________________________________________________________
1560void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1561{
1562/// Flag pixels (whether or not local maxima)
1563
1564 Int_t nx = hist->GetNbinsX();
1565 Int_t ny = hist->GetNbinsY();
1566 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1567 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1568
1569 for (Int_t i1=i-1; i1<i+2; i1++) {
1570 if (i1 < 1 || i1 > ny) continue;
1571 indx1 = (i1 - 1) * nx;
1572 for (Int_t j1=j-1; j1<j+2; j1++) {
1573 if (j1 < 1 || j1 > nx) continue;
1574 if (i == i1 && j == j1) continue;
1575 indx2 = indx1 + j1 - 1;
1576 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1577 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1578 else if (cont > cont1) isLocalMax[indx2] = -1;
1579 else { // the same charge
1580 isLocalMax[indx] = 1;
1581 if (isLocalMax[indx2] == 0) {
1582 FlagLocalMax(hist, i1, j1, isLocalMax);
1583 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1584 else isLocalMax[indx2] = -1;
1585 }
1586 }
1587 }
1588 }
1589 isLocalMax[indx] = 1; // local maximum
1590}
1591
1592//_____________________________________________________________________________
1593void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1594 Int_t *localMax, Int_t iMax)
1595{
1596/// Find pixel cluster around local maximum \a iMax and pick up pads
1597/// overlapping with it
1598
1599 TH2D *hist = (TH2D*) gROOT->FindObject("anode");
1600 Int_t nx = hist->GetNbinsX();
1601 Int_t ny = hist->GetNbinsY();
1602 Int_t ic = localMax[iMax] / nx + 1;
1603 Int_t jc = localMax[iMax] % nx + 1;
1604 Bool_t *used = new Bool_t[ny*nx];
1605 for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE;
1606
1607 // Drop all pixels from the array - pick up only the ones from the cluster
1608 fPixArray->Delete();
1609
1610 Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2;
1611 Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;
1612 Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
1613 Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
1614 Double_t cont = hist->GetCellContent(jc,ic);
1615 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1616 used[(ic-1)*nx+jc-1] = kTRUE;
1617 fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1618
1619 Int_t nPix = fPixArray->GetEntriesFast();
1620 Int_t npad = cluster.Multiplicity();
1621
1622 for (Int_t i=0; i<nPix; ++i)
1623 {
1624 AliMUONPad* pixPtr = Pixel(i);
1625 pixPtr->SetSize(0,wx);
1626 pixPtr->SetSize(1,wy);
1627 }
1628
1629 // Pick up pads which overlap with found pixels
1630 for (Int_t i=0; i<npad; i++)
1631 {
1632 cluster.Pad(i)->SetStatus(-1);
1633 }
1634
1635 for (Int_t i=0; i<nPix; i++)
1636 {
1637 AliMUONPad* pixPtr = Pixel(i);
1638 for (Int_t j=0; j<npad; ++j)
1639 {
1640 AliMUONPad* pad = cluster.Pad(j);
1641 if ( Overlap(*pad,*pixPtr) )
1642 {
1643 pad->SetStatus(0);
1644 }
1645 }
1646 }
1647
1648 delete [] used; used = 0;
1649}
1650
1651//_____________________________________________________________________________
1652AliMUONClusterFinderMLEM&
1653AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1654{
1655/// Protected assignement operator
1656
1657 if (this == &rhs) return *this;
1658
1659 AliFatal("Not implemented.");
1660
1661 return *this;
1662}
1663
1664//_____________________________________________________________________________
1665void
1666AliMUONClusterFinderMLEM::Neighbours(Int_t cathode, Int_t ix, Int_t iy,
1667 Int_t& n, Int_t* xList, Int_t* yList)
1668{
1669 /// Get the list of neighbours of pad at (cathode,ix,iy)
1670 n = 0;
1671
1672 const AliMpVSegmentation* seg = fSegmentation[cathode];
1673
1674 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1675
1676 // Define the region to look into : a region slightly bigger
1677 // than the pad itself (5% bigger), in order to catch first neighbours.
1678
1679 AliMpArea area(pad.Position(),pad.Dimensions()*1.05);
1680
1681 AliMpVPadIterator* it = seg->CreateIterator(area);
1682 it->First();
1683 while ( !it->IsDone() && n < 10 )
1684 {
1685 AliMpPad p = it->CurrentItem();
1686 if ( p != pad ) // skip self
1687 {
1688 xList[n] = p.GetIndices().GetFirst();
1689 yList[n] = p.GetIndices().GetSecond();
1690 ++n;
1691 }
1692 it->Next();
1693 }
1694 delete it;
1695}
1696
1697//_____________________________________________________________________________
1698void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1699{
1700 /// Add virtual pad (with small charge) to improve fit for some
1701 /// clusters (when pad with max charge is at the extreme of the cluster)
1702
1703 // Get number of pads in X and Y-directions
1704 Int_t nInX = -1, nInY;
1705 PadsInXandY(cluster,nInX, nInY);
1706 ++fNClusters;
c0a16418 1707
1708 // Add virtual pad only if number of pads per direction == 2
1709 if (nInX != 2 && nInY != 2) return;
1710
1711 ++fNAddVirtualPads;
1712
1713 // Find pads with max charge
1714 Int_t maxpad[2][2] = {{-1, -1}, {-1, -1}}, cath;
1715 Double_t sigmax[2] = {0}, aamax[2] = {0};
1716 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
1717 {
1718 AliMUONPad* pad = cluster.Pad(j);
1719 if (pad->Status() != 0) continue;
1720 cath = pad->Cathode();
1721 if (pad->Charge() > sigmax[cath])
1722 {
1723 maxpad[cath][1] = maxpad[cath][0];
1724 aamax[cath] = sigmax[cath];
1725 sigmax[cath] = pad->Charge();
1726 maxpad[cath][0] = j;
1727 }
1728 }
1729
1730 if (maxpad[0][0] >= 0 && maxpad[0][1] < 0 || maxpad[1][0] >= 0 && maxpad[1][1] < 0)
1731 {
1732 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
1733 {
1734 AliMUONPad* pad = cluster.Pad(j);
1735 if (pad->Status() != 0) continue;
1736 cath = pad->Cathode();
1737 if (j == maxpad[cath][0] || j == maxpad[cath][1]) continue;
1738 if ( pad->Charge() > aamax[cath])
1739 {
1740 aamax[cath] = pad->Charge();
1741 maxpad[cath][1] = j;
1742 }
1743 }
1744 }
1745
1746 // cout << "-------AddVirtualPad" << endl;
1747// cout << Form("nInX %2d nInY %2d",nInX,nInY) << endl;
1748//
1749// cluster.Print("full");
1750//
1751// for ( Int_t i = 0; i < 2; ++i )
1752// {
1753// for ( Int_t j = 0; j < 2; ++j )
1754// {
1755// cout << Form("maxpad[%d][%d]=%d",i,j,maxpad[i][j]) << endl;
1756// }
1757// }
1758
1759 // Check for mirrors (side X on cathode 0)
1760 Bool_t mirror = kFALSE;
1761 if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1762 {
1763 AliMUONPad* maxPadCath[] = { cluster.Pad(maxpad[0][0]), cluster.Pad(maxpad[1][0]) };
1764 mirror = maxPadCath[0]->DX() < maxPadCath[0]->DY();
1765 if (!mirror && TMath::Abs( maxPadCath[0]->DX() - maxPadCath[1]->DX()) < 0.001)
1766 {
1767 // Special case when pads on both cathodes have the same size
1768 Int_t yud[2] = {0};
1769 for (Int_t j = 0; j < cluster.Multiplicity(); ++j)
1770 {
1771 AliMUONPad* pad = cluster.Pad(j);
1772 cath = pad->Cathode();
1773 if (j == maxpad[cath][0]) continue;
1774 if ( pad->Ix() != maxPadCath[cath]->Ix() ) continue;
1775 if ( TMath::Abs(pad->Iy() - maxPadCath[cath]->Iy()) == 1 )
1776 {
1777 yud[cath]++;
1778 }
1779 }
1780 if (!yud[0]) mirror = kTRUE; // take the other cathode
1781 } // if (!mirror &&...
1782 } // if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1783
1784// // Find neughbours of pads with max charges
1785 Int_t nn, xList[10], yList[10], ix0, iy0, ix, iy, neighb;
1786 for (cath=0; cath<2; cath++)
1787 {
1788 if (!cath && maxpad[0][0] < 0) continue; // one-sided cluster - cathode 1
1789 if (cath && maxpad[1][0] < 0) break; // one-sided cluster - cathode 0
1790 if (maxpad[1][0] >= 0)
1791 {
1792 if (!mirror)
1793 {
1794 if (!cath && nInY != 2) continue;
1795 if (cath && nInX != 2 && (maxpad[0][0] >= 0 || nInY != 2)) continue;
1796 }
1797 else
1798 {
1799 if (!cath && nInX != 2) continue;
1800 if (cath && nInY != 2 && (maxpad[0][0] >= 0 || nInX != 2)) continue;
1801 }
1802 }
1803
1804 Int_t iAddX = 0, iAddY = 0, ix1 = 0, iy1 = 0, iPad = 0;
1805 if (maxpad[0][0] < 0) iPad = 1;
1806
1807 for (iPad=0; iPad<2; iPad++)
1808 {
1809 if (maxpad[cath][iPad] < 0) continue;
1810 if (iPad && !iAddX && !iAddY) break;
1811 if (iPad && cluster.Pad(maxpad[cath][1])->Charge() / sigmax[cath] < 0.5) break;
1812
1813 Int_t neighbx = 0, neighby = 0;
1814 ix0 = cluster.Pad(maxpad[cath][iPad])->Ix();
1815 iy0 = cluster.Pad(maxpad[cath][iPad])->Iy();
1816 Neighbours(cath,ix0,iy0,nn,xList,yList);
1817 //Float_t zpad;
1818 for (Int_t j=0; j<nn; j++) {
1819 if (TMath::Abs(xList[j]-ix0) == 1 || xList[j]*ix0 == -1) neighbx++;
1820 if (TMath::Abs(yList[j]-iy0) == 1 || yList[j]*iy0 == -1) neighby++;
1821 }
1822 if (!mirror) {
1823 if (cath) neighb = neighbx;
1824 else neighb = neighby;
1825 if (maxpad[0][0] < 0) neighb += neighby;
1826 else if (maxpad[1][0] < 0) neighb += neighbx;
1827 } else {
1828 if (!cath) neighb = neighbx;
1829 else neighb = neighby;
1830 if (maxpad[0][0] < 0) neighb += neighbx;
1831 else if (maxpad[1][0] < 0) neighb += neighby;
1832 }
1833
1834 for (Int_t j=0; j< cluster.Multiplicity(); ++j)
1835 {
1836 AliMUONPad* pad = cluster.Pad(j);
1837 if ( pad->Cathode() != cath) continue;
1838 ix = pad->Ix();
1839 iy = pad->Iy();
1840 if (iy == iy0 && ix == ix0) continue;
1841 for (Int_t k=0; k<nn; ++k)
1842 {
1843 if (xList[k] != ix || yList[k] != iy) continue;
1844 if (!mirror)
1845 {
1846 if ((!cath || maxpad[0][0] < 0) &&
1847 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
1848 if (!iPad && TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) ix1 = xList[k]; //19-12-05
1849 xList[k] = yList[k] = 0;
1850 neighb--;
1851 break;
1852 }
1853 if ((cath || maxpad[1][0] < 0) &&
1854 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
1855 if (!iPad) ix1 = xList[k]; //19-12-05
1856 xList[k] = yList[k] = 0;
1857 neighb--;
1858 }
1859 } else {
1860 if ((!cath || maxpad[0][0] < 0) &&
1861 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
1862 if (!iPad) ix1 = xList[k]; //19-12-05
1863 xList[k] = yList[k] = 0;
1864 neighb--;
1865 break;
1866 }
1867 if ((cath || maxpad[1][0] < 0) &&
1868 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
1869 xList[k] = yList[k] = 0;
1870 neighb--;
1871 }
1872 }
1873 break;
1874 } // for (Int_t k=0; k<nn;
1875 if (!neighb) break;
1876 } // for (Int_t j=0; j< cluster.Multiplicity();
1877 if (!neighb) continue;
1878
1879 // Add virtual pad
1880 Int_t npads, isec;
1881 isec = npads = 0;
1882 for (Int_t j=0; j<nn; j++)
1883 {
1884 if (xList[j] == 0 && yList[j] == 0) continue;
1885 // npads = fnPads[0] + fnPads[1];
1886 // fPadIJ[0][npads] = cath;
1887 // fPadIJ[1][npads] = 0;
1888 ix = xList[j];
1889 iy = yList[j];
1890 if (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) {
1891 if (iy != iy0) continue; // new segmentation - check
1892 if (nInX != 2) continue; // new
1893 if (!mirror) {
1894 if (!cath && maxpad[1][0] >= 0) continue;
1895 } else {
1896 if (cath && maxpad[0][0] >= 0) continue;
1897 }
1898 if (iPad && !iAddX) continue;
1899 AliMpPad pad = fSegmentation[cath]->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1900 // fXyq[0][npads] = pad.Position().X();
1901 // fXyq[1][npads] = pad.Position().Y();
9bbd7f60 1902 AliMUONPad muonPad(fDetElemId, cath, ix, iy, pad.Position().X(), pad.Position().Y(),
1903 pad.Dimensions().X(), pad.Dimensions().Y(), 0);
c0a16418 1904 // fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
1905 // if (fXyq[0][npads] > 1.e+5) continue; // temporary fix
1906 if (muonPad.Coord(0) > 1.e+5) continue; // temporary fix
1907 if (ix == ix1) continue; //19-12-05
1908 if (ix1 == ix0) continue;
1909 if (maxpad[1][0] < 0 || mirror && maxpad[0][0] >= 0) {
1910 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/100, 5.);
1911 // else fXyq[2][npads] = TMath::Min (aamax[0]/100, 5.);
1912 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[0]/100, 5.));
1913 else muonPad.SetCharge(TMath::Min (aamax[0]/100, 5.));
1914 }
1915 else {
1916 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/100, 5.);
1917 // else fXyq[2][npads] = TMath::Min (aamax[1]/100, 5.);
1918 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[1]/100, 5.));
1919 else muonPad.SetCharge(TMath::Min (aamax[1]/100, 5.));
1920 }
1921 // fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
1922 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1923 // fXyq[3][npads] = -2; // flag
1924 // fPadIJ[2][npads] = ix;
1925 // fPadIJ[3][npads] = iy;
9bbd7f60 1926 muonPad.SetReal(kFALSE);
c0a16418 1927 // fnPads[1]++;
1928 // iAddX = npads;
1929 iAddX = 1;
9bbd7f60 1930 if (fDebug) printf(" ***** Add virtual pad in X ***** %f %f %f %3d %3d %f %f \n",
1931 muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy,
1932 muonPad.DX(), muonPad.DY());
c0a16418 1933 cluster.AddPad(muonPad); // add pad to the cluster
1934 ix1 = ix0;
1935 continue;
1936 }
1937 if (nInY != 2) continue;
1938 if (!mirror && cath && maxpad[0][0] >= 0) continue;
1939 if (mirror && !cath && maxpad[1][0] >= 0) continue;
1940 if (TMath::Abs(iy-iy0) == 1 || TMath::Abs(iy*iy0) == 1) {
1941 if (ix != ix0) continue; // new segmentation - check
1942 if (iPad && !iAddY) continue;
1943 AliMpPad pad = fSegmentation[cath]->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1944 // fXyq[0][npads] = pad.Position().X();
1945 // fXyq[1][npads] = pad.Position().Y();
1946 // fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
9bbd7f60 1947 AliMUONPad muonPad(fDetElemId, cath, ix, iy, pad.Position().X(), pad.Position().Y(),
1948 pad.Dimensions().X(), pad.Dimensions().Y(), 0);
c0a16418 1949 if (iy1 == iy0) continue;
1950 //if (iPad && iy1 == iy0) continue;
1951 if (maxpad[0][0] < 0 || mirror && maxpad[1][0] >= 0) {
1952 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/15, fgkZeroSuppression);
1953 // else fXyq[2][npads] = TMath::Min (aamax[1]/15, fgkZeroSuppression);
1954 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[1]/15, fgkZeroSuppression));
1955 else muonPad.SetCharge(TMath::Min (aamax[1]/15, fgkZeroSuppression));
1956 }
1957 else {
1958 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/15, fgkZeroSuppression);
1959 // else fXyq[2][npads] = TMath::Min (aamax[0]/15, fgkZeroSuppression);
1960 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[0]/15, fgkZeroSuppression));
1961 else muonPad.SetCharge(TMath::Min (aamax[0]/15, fgkZeroSuppression));
1962 }
1963 // fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
1964 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1965 // fXyq[3][npads] = -2; // flag
1966 // fPadIJ[2][npads] = ix;
1967 // fPadIJ[3][npads] = iy;
9bbd7f60 1968 muonPad.SetReal(kFALSE);
c0a16418 1969 // fnPads[1]++;
1970 // iAddY = npads;
1971 iAddY = 1;
9bbd7f60 1972 if (fDebug) printf(" ***** Add virtual pad in Y ***** %f %f %f %3d %3d %f %f \n",
1973 muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy,
1974 muonPad.DX(), muonPad.DY());
c0a16418 1975 cluster.AddPad(muonPad); // add pad to the cluster
1976 iy1 = iy0;
1977 }
1978 } // for (Int_t j=0; j<nn;
1979 } // for (Int_t iPad=0;
1980 } // for (cath=0; cath<2;
1981}
1982
1983//_____________________________________________________________________________
1984void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1985 Int_t &nInX, Int_t &nInY) const
1986{
1987 /// Find number of pads in X and Y-directions (excluding virtual ones and
1988 /// overflows)
1989
1990 Int_t statusToTest = 1;
1991
1992 if ( nInX < 0 ) statusToTest = 0;
1993
1994 Bool_t mustMatch(kTRUE);
1995
1996 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1997
1998 nInX = cn.GetFirst();
1999 nInY = cn.GetSecond();
2000}
2001
2002//_____________________________________________________________________________
2003void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
2004{
2005 /// Remove pixel at index i
2006 AliMUONPad* pixPtr = Pixel(i);
2007 fPixArray->RemoveAt(i);
2008 delete pixPtr;
2009}
2010
2011//_____________________________________________________________________________
2012void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
2013{
2014/// Process simple cluster (small number of pads) without EM-procedure
2015
2016 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
2017 Double_t parOk[3] = {0.};
2018 TObjArray *clusters[1];
2019 clusters[0] = fPixArray;
2020
2021 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
2022
2023 for (Int_t i = 0; i < cluster.Multiplicity(); ++i)
2024 {
2025 AliMUONPad* pad = cluster.Pad(i);
9bbd7f60 2026 if ( pad->IsSaturated())
c0a16418 2027 {
2028 pad->SetStatus(-9);
2029 }
2030 else
2031 {
2032 pad->SetStatus(1);
2033 }
2034 }
2035 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList);
2036}
2037
2038//_____________________________________________________________________________
2039AliMUONPad*
2040AliMUONClusterFinderMLEM::Pixel(Int_t i) const
2041{
2042 /// Returns pixel at index i
2043 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
2044}
2045
2046//_____________________________________________________________________________
2047void
2048AliMUONClusterFinderMLEM::Print(Option_t* what) const
2049{
2050 /// printout
2051 TString swhat(what);
2052 swhat.ToLower();
2053 if ( swhat.Contains("precluster") )
2054 {
2055 if ( fPreCluster) fPreCluster->Print();
2056 }
2057}
2058
2059