]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONClusterFinderMLEM.cxx
Define default values of cluster resolution per chamber in the
[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
3d1463c8 18//-----------------------------------------------------------------------------
c0a16418 19/// \class AliMUONClusterFinderMLEM
20///
21/// Clusterizer class based on the Expectation-Maximization algorithm
22///
23/// Pre-clustering is handled by AliMUONPreClusterFinder
24/// From a precluster a pixel array is built, and from this pixel array
25/// a list of clusters is output (using AliMUONClusterSplitterMLEM).
26///
27/// \author Laurent Aphecetche (for the "new" C++ structure) and
28/// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
3d1463c8 29//-----------------------------------------------------------------------------
c0a16418 30
9edefa04 31#include "AliMUONClusterFinderMLEM.h"
c0a16418 32#include "AliLog.h"
33#include "AliMUONCluster.h"
34#include "AliMUONClusterSplitterMLEM.h"
3887d11e 35#include "AliMUONVDigit.h"
c0a16418 36#include "AliMUONPad.h"
37#include "AliMUONPreClusterFinder.h"
38#include "AliMpPad.h"
39#include "AliMpVPadIterator.h"
40#include "AliMpVSegmentation.h"
41#include "AliRunLoader.h"
3887d11e 42#include "AliMUONVDigitStore.h"
a6435fb4 43#include <Riostream.h>
44#include <TH2.h>
45#include <TMinuit.h>
46#include <TCanvas.h>
a6435fb4 47#include <TMath.h>
f6c291ef 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
2abdae6e 55//const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
56const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
c0a16418 57const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
58const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
59
8c718ff8 60// Status flags for pads
61const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
62const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
63const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
64const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
65const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
66const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
c0a16418 67
68//_____________________________________________________________________________
b1a19e07 69AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
c0a16418 70 : AliMUONVClusterFinder(),
b1a19e07 71fPreClusterFinder(clusterFinder),
c0a16418 72fPreCluster(0x0),
73fClusterList(),
74fEventNumber(0),
75fDetElemId(-1),
76fClusterNumber(0),
c80b29a4 77fHistMlem(0x0),
78fHistAnode(0x0),
c0a16418 79fPixArray(new TObjArray(20)),
9bbd7f60 80fDebug(0),
c0a16418 81fPlot(plot),
c0a16418 82fSplitter(0x0),
83fNClusters(0),
84fNAddVirtualPads(0)
85{
86 /// Constructor
87
88 fSegmentation[1] = fSegmentation[0] = 0x0;
89
9bbd7f60 90 if (fPlot) fDebug = 1;
c0a16418 91}
92
93//_____________________________________________________________________________
94AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
95{
96/// Destructor
c80b29a4 97 delete fPixArray; fPixArray = 0;
c0a16418 98// delete fDraw;
99 delete fPreClusterFinder;
c0a16418 100 delete fSplitter;
101 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
102 fNClusters,fNAddVirtualPads));
103}
104
105//_____________________________________________________________________________
106Bool_t
24935e58 107AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
108 TClonesArray* pads[2],
109 const AliMpArea& area,
110 const AliMpVSegmentation* seg[2])
c0a16418 111{
112 /// Prepare for clustering
2060b217 113// AliCodeTimerAuto("")
c0a16418 114
115 for ( Int_t i = 0; i < 2; ++i )
116 {
24935e58 117 fSegmentation[i] = seg[i];
c0a16418 118 }
119
120 // Find out the DetElemId
24935e58 121 fDetElemId = detElemId;
c0a16418 122
123 delete fSplitter;
124 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
9bbd7f60 125 fSplitter->SetDebug(fDebug);
c0a16418 126
127 // find out current event number, and reset the cluster number
7deb8eb0 128 AliRunLoader *runLoader = AliRunLoader::GetRunLoader();
129 fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
c0a16418 130 fClusterNumber = -1;
131 fClusterList.Delete();
132
b1a19e07 133 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
c0a16418 134
24935e58 135 if ( fPreClusterFinder->NeedSegmentation() )
136 {
137 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
138 }
139 else
140 {
141 return fPreClusterFinder->Prepare(detElemId,pads,area);
142 }
c0a16418 143}
144
145//_____________________________________________________________________________
146AliMUONCluster*
147AliMUONClusterFinderMLEM::NextCluster()
148{
149 /// Return next cluster
f6c291ef 150// AliCodeTimerAuto("")
c0a16418 151
c0a16418 152 // if the list of clusters is not void, pick one from there
9bbd7f60 153 TObject* o = fClusterList.At(++fClusterNumber);
154 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
c0a16418 155
156 //FIXME : at this point, must check whether we've used all the digits
157 //from precluster : if not, let the preclustering know about those unused
158 //digits, so it can reuse them
159
160 // if the cluster list is exhausted, we need to go to the next
161 // pre-cluster and treat it
162
163 fPreCluster = fPreClusterFinder->NextCluster();
9bf6860b 164
165 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
166 fClusterNumber = -1; //AZ
167
c0a16418 168 if (!fPreCluster)
169 {
170 // we are done
171 return 0x0;
172 }
173
c0a16418 174 WorkOnPreCluster();
175
176 // WorkOnPreCluster may have used only part of the pads, so we check that
177 // now, and let the unused pads be reused by the preclustering...
178
2abdae6e 179 Int_t mult = fPreCluster->Multiplicity();
180 for ( Int_t i = 0; i < mult; ++i )
c0a16418 181 {
182 AliMUONPad* pad = fPreCluster->Pad(i);
183 if ( !pad->IsUsed() )
184 {
185 fPreClusterFinder->UsePad(*pad);
186 }
187 }
188
189 return NextCluster();
190}
191
192//_____________________________________________________________________________
193Bool_t
194AliMUONClusterFinderMLEM::WorkOnPreCluster()
195{
196 /// Starting from a precluster, builds a pixel array, and then
197 /// extract clusters from this array
198
2abdae6e 199 // AliCodeTimerAuto("")
9bbd7f60 200
201 if (fDebug) {
7deb8eb0 202 cout << " *** Event # " << fEventNumber
9bbd7f60 203 << " det. elem.: " << fDetElemId << endl;
2abdae6e 204 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
9bbd7f60 205 AliMUONPad* pad = fPreCluster->Pad(j);
206 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
207 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
208 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
209 }
210 }
211
c0a16418 212 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
c0a16418 213 if (!cluster) return kFALSE;
9bbd7f60 214
c0a16418 215 BuildPixArray(*cluster);
216
217 if ( fPixArray->GetLast() < 0 )
218 {
219 AliDebug(1,"No pixel for the above cluster");
220 delete cluster;
221 return kFALSE;
222 }
223
224 // Use MLEM for cluster finder
225 Int_t nMax = 1, localMax[100], maxPos[100];
226 Double_t maxVal[100];
227
c0a16418 228 Int_t iSimple = 0, nInX = -1, nInY;
229
230 PadsInXandY(*cluster,nInX, nInY);
231
8c718ff8 232 if (nInX < 4 && nInY < 4)
233 {
234 iSimple = 1; // simple cluster
235 }
236 else
c0a16418 237 {
8c718ff8 238 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
239 if (nMax > 1) {
240 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
241 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
242 }
c0a16418 243 }
244
2abdae6e 245 for (Int_t i = 0; i < nMax; ++i)
c0a16418 246 {
247 if (nMax > 1)
248 {
249 FindCluster(*cluster,localMax, maxPos[i]);
250 }
2abdae6e 251
c0a16418 252 MainLoop(*cluster,iSimple);
2abdae6e 253
c0a16418 254 if (i < nMax-1)
255 {
2abdae6e 256 Int_t mult = cluster->Multiplicity();
257 for (Int_t j = 0; j < mult; ++j)
c0a16418 258 {
259 AliMUONPad* pad = cluster->Pad(j);
8c718ff8 260 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
261 if ( pad->Status() != fgkOver ) continue; // pad was not used
262 //pad->SetStatus(0);
263 pad->SetStatus(fgkZero);
c0a16418 264 pad->RevertCharge(); // use backup charge value
265 }
266 }
267 } // for (Int_t i=0; i<nMax;
c80b29a4 268 delete fHistMlem;
269 delete fHistAnode;
270 fHistMlem = fHistAnode = 0x0;
c0a16418 271 delete cluster;
272 return kTRUE;
273}
274
275//_____________________________________________________________________________
276Bool_t
277AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
278{
279 /// Check if the pad and the pixel overlaps
280
281 // make a fake pad from the pixel
282 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
283 pixel.Coord(0),pixel.Coord(1),
284 pixel.Size(0),pixel.Size(1),0);
285
286 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
287}
288
289//_____________________________________________________________________________
290AliMUONCluster*
291AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
292{
293 /// Check precluster in order to attempt to simplify it (mostly for
294 /// two-cathode preclusters)
2abdae6e 295
24935e58 296 AliCodeTimerAuto("")
f6c291ef 297
2abdae6e 298 // Disregard small clusters (leftovers from splitting or noise)
299 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
300 origCluster.Charge(0)+origCluster.Charge(1) < 10)
c0a16418 301 {
c0a16418 302 return 0x0;
303 }
304
8c718ff8 305 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
c0a16418 306
c0a16418 307 AliDebug(2,"Start of CheckPreCluster=");
2abdae6e 308 //StdoutToAliDebug(2,cluster->Print("full"));
c0a16418 309
c0a16418 310 AliMUONCluster* rv(0x0);
311
2abdae6e 312 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
c0a16418 313 {
314 rv = CheckPreclusterTwoCathodes(cluster);
315 }
316 else
317 {
f6c291ef 318 rv = cluster;
c0a16418 319 }
c0a16418 320 return rv;
321}
322
c0a16418 323//_____________________________________________________________________________
324AliMUONCluster*
325AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
326{
71a2d3aa 327 /// Check two-cathode cluster
c0a16418 328
c0a16418 329 Int_t npad = cluster->Multiplicity();
330 Int_t* flags = new Int_t[npad];
2abdae6e 331 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
c0a16418 332
333 // Check pad overlaps
2abdae6e 334 for ( Int_t i = 0; i < npad; ++i)
c0a16418 335 {
336 AliMUONPad* padi = cluster->Pad(i);
8c718ff8 337 if ( padi->Cathode() != 0 ) continue;
2abdae6e 338 for (Int_t j = i+1; j < npad; ++j)
c0a16418 339 {
340 AliMUONPad* padj = cluster->Pad(j);
8c718ff8 341 if ( padj->Cathode() != 1 ) continue;
c0a16418 342 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
343 flags[i] = flags[j] = 1; // mark overlapped pads
344 }
345 }
346
347 // Check if all pads overlap
348 Int_t nFlags=0;
2abdae6e 349 for (Int_t i = 0; i < npad; ++i)
c0a16418 350 {
2abdae6e 351 if (!flags[i]) ++nFlags;
c0a16418 352 }
353
354 if (nFlags > 0)
355 {
356 // not all pads overlap.
9bbd7f60 357 if (fDebug) cout << " nFlags: " << nFlags << endl;
358 TObjArray toBeRemoved;
2abdae6e 359 for (Int_t i = 0; i < npad; ++i)
c0a16418 360 {
361 AliMUONPad* pad = cluster->Pad(i);
362 if (flags[i]) continue;
363 Int_t cath = pad->Cathode();
364 Int_t cath1 = TMath::Even(cath);
365 // Check for edge effect (missing pads on the _other_ cathode)
366 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
367 if (!mpPad.IsValid()) continue;
8c718ff8 368 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
369 if (nFlags == 1 && pad->Charge() < 20) continue;
c0a16418 370 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
371 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
9bbd7f60 372 toBeRemoved.AddLast(pad);
c0a16418 373 fPreCluster->Pad(i)->Release();
9bbd7f60 374 }
2abdae6e 375 Int_t nRemove = toBeRemoved.GetEntriesFast();
376 for ( Int_t i = 0; i < nRemove; ++i )
9bbd7f60 377 {
2abdae6e 378 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
c0a16418 379 }
380 }
381
382 // Check correlations of cathode charges
9bbd7f60 383 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
c0a16418 384 {
385 // big difference
386 Int_t cathode = cluster->MaxRawChargeCathode();
8c718ff8 387 Int_t imin(-1);
388 Int_t imax(-1);
c0a16418 389 Double_t cmax(0);
390 Double_t cmin(1E9);
391
392 // get min and max pad charges on the cathode opposite to the
393 // max pad (given by MaxRawChargeCathode())
394 //
2abdae6e 395 Int_t mult = cluster->Multiplicity();
396 for ( Int_t i = 0; i < mult; ++i )
c0a16418 397 {
398 AliMUONPad* pad = cluster->Pad(i);
399 if ( pad->Cathode() != cathode || !pad->IsReal() )
400 {
401 // only consider pads in the opposite cathode, and
2abdae6e 402 // only consider real pads (i.e. exclude the virtual ones)
c0a16418 403 continue;
404 }
405 if ( pad->Charge() < cmin )
406 {
407 cmin = pad->Charge();
408 imin = i;
8c718ff8 409 if (imax < 0) {
410 imax = imin;
411 cmax = cmin;
412 }
c0a16418 413 }
8c718ff8 414 else if ( pad->Charge() > cmax )
c0a16418 415 {
416 cmax = pad->Charge();
417 imax = i;
418 }
419 }
420 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
421 imin,imax,cmin,cmax));
422 //
423 // arrange pads according to their distance to the max, normalized
424 // to the pad size
2abdae6e 425 Double_t* dist = new Double_t[mult];
c0a16418 426 Double_t dxMin(1E9);
427 Double_t dyMin(1E9);
428 Double_t dmin(0);
429
430 AliMUONPad* padmax = cluster->Pad(imax);
431
2abdae6e 432 for ( Int_t i = 0; i < mult; ++i )
c0a16418 433 {
434 dist[i] = 0.0;
435 if ( i == imax) continue;
436 AliMUONPad* pad = cluster->Pad(i);
437 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
438 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
439 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
440 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
441 if ( i == imin )
442 {
443 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
444 dxMin = dx;
445 dyMin = dy;
446 }
447 }
448
2abdae6e 449 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
450 Double_t xmax(-1), distPrev(999);
c0a16418 451 TObjArray toBeRemoved;
452
2abdae6e 453 for ( Int_t i = 0; i < mult; ++i )
c0a16418 454 {
455 Int_t indx = flags[i];
456 AliMUONPad* pad = cluster->Pad(indx);
457 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
458 if ( dist[indx] > dmin )
459 {
460 // farther than the minimum pad
461 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
462 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
463 dx *= dxMin;
464 dy *= dyMin;
465 if (dx >= 0 && dy >= 0) continue;
466 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
467 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
468 }
2abdae6e 469 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
c0a16418 470 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
471 {
472 // release pad
473 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
474 {
475 cmax = TMath::Max(pad->Charge(),cmax);
476 }
477 else
478 {
479 cmax = pad->Charge();
480 }
481 xmax = dist[indx];
2abdae6e 482 distPrev = dist[indx];
c0a16418 483 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
484 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
485 pad->Charge()));
486
487 toBeRemoved.AddLast(pad);
488 fPreCluster->Pad(indx)->Release();
489 }
490 }
2abdae6e 491 Int_t nRemove = toBeRemoved.GetEntriesFast();
492 for ( Int_t i = 0; i < nRemove; ++i )
c0a16418 493 {
2abdae6e 494 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
c0a16418 495 }
496 delete[] dist;
497 }
498
499 delete[] flags;
500
501 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
2abdae6e 502 //StdoutToAliDebug(2,cluster->Print("full"));
c0a16418 503
504 return cluster;
505}
506
507//_____________________________________________________________________________
508void
509AliMUONClusterFinderMLEM::CheckOverlaps()
510{
511 /// For debug only : check if some pixels overlap...
512
513 Int_t nPix = fPixArray->GetLast()+1;
514 Int_t dummy(0);
515
516 for ( Int_t i = 0; i < nPix; ++i )
517 {
518 AliMUONPad* pixelI = Pixel(i);
519 AliMUONPad pi(dummy,dummy,dummy,dummy,
520 pixelI->Coord(0),pixelI->Coord(1),
521 pixelI->Size(0),pixelI->Size(1),0.0);
522
523 for ( Int_t j = i+1; j < nPix; ++j )
524 {
525 AliMUONPad* pixelJ = Pixel(j);
526 AliMUONPad pj(dummy,dummy,dummy,dummy,
527 pixelJ->Coord(0),pixelJ->Coord(1),
528 pixelJ->Size(0),pixelJ->Size(1),0.0);
529 AliMpArea area;
530
531 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
532 {
533 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
2abdae6e 534 /*
c0a16418 535 StdoutToAliInfo(pixelI->Print();
536 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
537 pixelJ->Print();
538 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
539 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
540 cout << "-------" << endl;
541 );
2abdae6e 542 */
c0a16418 543 }
544 }
545 }
546}
547
548//_____________________________________________________________________________
549void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
550{
551 /// Build pixel array for MLEM method
552
553 Int_t npad = cluster.Multiplicity();
554 if (npad<=0)
555 {
556 AliWarning("Got no pad at all ?!");
557 }
558
559 fPixArray->Delete();
8c718ff8 560 BuildPixArrayOneCathode(cluster);
c0a16418 561
562 Int_t nPix = fPixArray->GetLast()+1;
563
2abdae6e 564// AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
c0a16418 565
2abdae6e 566 if ( nPix > npad )
c0a16418 567 {
2abdae6e 568// AliDebug(2,Form("Will trim number of pixels to number of pads"));
c0a16418 569
570 // Too many pixels - sort and remove pixels with the lowest signal
571 fPixArray->Sort();
8c718ff8 572 for ( Int_t i = npad; i < nPix; ++i )
c0a16418 573 {
574 RemovePixel(i);
575 }
576 fPixArray->Compress();
c0a16418 577 } // if (nPix > npad)
578
579// StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
580// fPixArray->Print(););
2abdae6e 581 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
c0a16418 582}
583
2abdae6e 584//_____________________________________________________________________________
585void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
586{
587 /// Build the pixel array
588
589// AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
590
2abdae6e 591 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
592 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
8c718ff8 593 Int_t found[2] = {0}, mult = cluster.Multiplicity();
2abdae6e 594
8c718ff8 595 for ( Int_t i = 0; i < mult; ++i) {
2abdae6e 596 AliMUONPad* pad = cluster.Pad(i);
2abdae6e 597 for (Int_t j = 0; j < 2; ++j) {
598 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
599 xy0[j] = pad->Coord(j);
600 found[j] = 1;
601 }
602 }
8c718ff8 603 if (found[0] && found[1]) break;
2abdae6e 604 }
2abdae6e 605
2abdae6e 606 Double_t min[2], max[2];
607 Int_t cath0 = 0, cath1 = 1;
608 if (cluster.Multiplicity(0) == 0) cath0 = 1;
609 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
8c718ff8 610
611 TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
612 TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
613 min[0] = leftDown.X();
614 min[1] = leftDown.Y();
615 max[0] = rightUp.X();
616 max[1] = rightUp.Y();
617 if (cath1 != cath0) {
618 leftDown = cluster.Area(cath1).LeftDownCorner();
619 rightUp = cluster.Area(cath1).RightUpCorner();
620 min[0] = TMath::Max (min[0], leftDown.X());
621 min[1] = TMath::Max (min[1], leftDown.Y());
622 max[0] = TMath::Min (max[0], rightUp.X());
623 max[1] = TMath::Min (max[1], rightUp.Y());
624 }
2abdae6e 625
626 // Adjust limits
627 //width[0] /= 2; width[1] /= 2; // just for check
628 Int_t nbins[2];
629 for (Int_t i = 0; i < 2; ++i) {
630 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
8c718ff8 631 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
2abdae6e 632 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
633 + TMath::Sign(0.5,dist)) * width[i] * 2;
634 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
635 if (nbins[i] == 0) ++nbins[i];
636 max[i] = min[i] + nbins[i] * width[i] * 2;
637 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
638 }
639
640 // Book histogram
641 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
642 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
643 TAxis *xaxis = hist1->GetXaxis();
644 TAxis *yaxis = hist1->GetYaxis();
645
646 // Fill histogram
2abdae6e 647 for ( Int_t i = 0; i < mult; ++i) {
648 AliMUONPad* pad = cluster.Pad(i);
649 Int_t ix0 = xaxis->FindBin(pad->X());
650 Int_t iy0 = yaxis->FindBin(pad->Y());
c80b29a4 651 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
2abdae6e 652 }
653
654 // Store pixels
655 for (Int_t i = 1; i <= nbins[0]; ++i) {
656 Double_t x = xaxis->GetBinCenter(i);
657 for (Int_t j = 1; j <= nbins[1]; ++j) {
658 if (hist2->GetCellContent(i,j) < 0.1) continue;
8c718ff8 659 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
660 // cluster.Multiplicity(1)) continue;
661 if (cath0 != cath1) {
662 // Two-sided cluster
663 Double_t cont = hist2->GetCellContent(i,j);
664 if (cont < 999.) continue;
665 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
666 }
2abdae6e 667 Double_t y = yaxis->GetBinCenter(j);
668 Double_t charge = hist1->GetCellContent(i,j);
669 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
670 fPixArray->Add(pixPtr);
671 }
672 }
8c718ff8 673 //*
674 if (fPixArray->GetEntriesFast() == 1) {
675 // Split pixel into 2
676 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
677 pixPtr->SetSize(0,width[0]/2.);
678 pixPtr->Shift(0,-width[0]/4.);
679 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
680 fPixArray->Add(pixPtr);
681 }
682 //*/
2abdae6e 683 //fPixArray->Print();
684 delete hist1;
685 delete hist2;
686}
687
688//_____________________________________________________________________________
c80b29a4 689void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
690 TH2D *hist1, TH2D *hist2)
2abdae6e 691{
692 /// "Span" pad over histogram in the direction idir
693
2abdae6e 694 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
8c718ff8 695 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
696 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
2abdae6e 697
698 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
699
700 for (Int_t i = 0; i < nbinPad; ++i) {
701 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
702 if (ixy > nbins) break;
703 Double_t lowEdge = axis->GetBinLowEdge(ixy);
704 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
c80b29a4 705 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
2abdae6e 706 else {
707 // Fill histogram
708 Double_t cont = pad->Charge();
709 if (hist2->GetCellContent(ix0, ixy) > 0.1)
710 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
711 hist1->SetCellContent(ix0, ixy, cont);
8c718ff8 712 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
713 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
2abdae6e 714 }
715 }
716
717 for (Int_t i = -1; i > -nbinPad; --i) {
718 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
719 if (ixy < 1) break;
720 Double_t upEdge = axis->GetBinUpEdge(ixy);
721 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
c80b29a4 722 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
2abdae6e 723 else {
724 // Fill histogram
725 Double_t cont = pad->Charge();
726 if (hist2->GetCellContent(ix0, ixy) > 0.1)
727 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
728 hist1->SetCellContent(ix0, ixy, cont);
8c718ff8 729 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
730 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
2abdae6e 731 }
732 }
733}
734
c0a16418 735//_____________________________________________________________________________
736void
737AliMUONClusterFinderMLEM::Plot(const char* basename)
738{
739 /// Make a plot and save it as png
740
9bbd7f60 741 return; //AZ
c0a16418 742 if (!fPlot) return;
743
744 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
745 c->Draw();
746 Draw();
747 c->Modified();
748 c->Update();
749 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
750 fDetElemId,fClusterNumber));
751}
752
753//_____________________________________________________________________________
754void
755AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
756 Double_t* coef,
757 Double_t* probi)
758{
759 /// Compute coefficients needed for MLEM algorithm
760
761 Int_t nPix = fPixArray->GetLast()+1;
762
2abdae6e 763 //memset(probi,0,nPix*sizeof(Double_t));
764 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
c0a16418 765
2abdae6e 766 Int_t mult = cluster.Multiplicity();
767 for ( Int_t j = 0; j < mult; ++j )
c0a16418 768 {
769 AliMUONPad* pad = cluster.Pad(j);
770 Int_t indx = j*nPix;
771
2abdae6e 772 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
c0a16418 773 {
774 Int_t indx1 = indx + ipix;
8c718ff8 775 //if (pad->Status() < 0)
776 if (pad->Status() != fgkZero)
c0a16418 777 {
778 coef[indx1] = 0;
779 continue;
780 }
781 AliMUONPad* pixPtr = Pixel(ipix);
782 // coef is the charge (given by Mathieson integral) on pad, assuming
783 // the Mathieson is center at pixel.
784 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
2abdae6e 785// AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
786// pad->Ix(),pad->Iy(),
787// pad->X(),pad->Y(),
788// pad->DX(),pad->DY(),
789// pixPtr->Coord(0),pixPtr->Coord(1),
790// pixPtr->Size(0),pixPtr->Size(1),
791// coef[indx1]));
c0a16418 792
793 probi[ipix] += coef[indx1];
794 }
795 }
796}
797
798//_____________________________________________________________________________
799Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
800{
801 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
802
2abdae6e 803 // AliCodeTimerAuto("")
f6c291ef 804
c0a16418 805 Int_t nPix = fPixArray->GetLast()+1;
806
807 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
2abdae6e 808 //StdoutToAliDebug(2,cluster.Print("full"););
c0a16418 809
810 if ( nPix < 0 )
811 {
812 AliDebug(1,"No pixels, why am I here then ?");
813 }
814
815 AddVirtualPad(cluster); // add virtual pads if necessary
816
817 Int_t npadTot = cluster.Multiplicity();
818 Int_t npadOK = 0;
819 for (Int_t i = 0; i < npadTot; ++i)
820 {
8c718ff8 821 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
822 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
c0a16418 823 }
824
c0a16418 825 Double_t* coef(0x0);
826 Double_t* probi(0x0);
8c718ff8 827 Int_t lc(0); // loop counter
c0a16418 828
2abdae6e 829 //Plot("mlem.start");
8c718ff8 830 AliMUONPad* pixPtr = Pixel(0);
831 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
832
c0a16418 833 while (1)
834 {
835 ++lc;
c80b29a4 836 delete fHistMlem;
c0a16418 837
838 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
839 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
2abdae6e 840 //StdoutToAliDebug(2,fPixArray->Print("","full"));
c0a16418 841
842 coef = new Double_t [npadTot*nPix];
843 probi = new Double_t [nPix];
844
845 // Calculate coefficients and pixel visibilities
846 ComputeCoefficients(cluster,coef,probi);
847
2abdae6e 848 for (Int_t ipix = 0; ipix < nPix; ++ipix)
c0a16418 849 {
850 if (probi[ipix] < 0.01)
851 {
852 AliMUONPad* pixel = Pixel(ipix);
853 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
2abdae6e 854 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
c0a16418 855 pixel->SetCharge(0); // "invisible" pixel
856 }
857 }
858
859 // MLEM algorithm
860 Mlem(cluster,coef, probi, 15);
861
8c718ff8 862 // Find histogram limits for the 1'st pass only - for others computed below
863 if (lc == 1) {
864 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
865 {
866 pixPtr = Pixel(ipix);
867 for ( Int_t i = 0; i < 2; ++i )
868 {
869 Int_t indx = i * 2;
870 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
871 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
872 }
873 }
874 } else pixPtr = Pixel(0);
875
2abdae6e 876 for (Int_t i = 0; i < 4; i++)
c0a16418 877 {
878 xylim[i] -= pixPtr->Size(i/2);
879 }
c0a16418 880
881 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
882 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
883
2abdae6e 884 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
c0a16418 885 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
886 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
887 xylim[0],-xylim[1],xylim[2],-xylim[3]
888 ));
889
c80b29a4 890 fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
c0a16418 891
2abdae6e 892 for (Int_t ipix = 0; ipix < nPix; ++ipix)
c0a16418 893 {
bf0d3528 894 AliMUONPad* pixPtr2 = Pixel(ipix);
895 fHistMlem->Fill(pixPtr2->Coord(0),pixPtr2->Coord(1),pixPtr2->Charge());
c0a16418 896 }
897
898 // Check if the total charge of pixels is too low
899 Double_t qTot = 0;
2abdae6e 900 for ( Int_t i = 0; i < nPix; ++i)
c0a16418 901 {
902 qTot += Pixel(i)->Charge();
903 }
904
905 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
906 {
907 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
908 delete [] coef;
909 delete [] probi;
c0a16418 910 fPixArray->Delete();
2abdae6e 911 for ( Int_t i = 0; i < npadTot; ++i)
c0a16418 912 {
913 AliMUONPad* pad = cluster.Pad(i);
8c718ff8 914 //if ( pad->Status() == 0) pad->SetStatus(-1);
915 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
c0a16418 916 }
917 return kFALSE;
918 }
919
920 if (iSimple)
921 {
922 // Simple cluster - skip further passes thru EM-procedure
923 Simple(cluster);
924 delete [] coef;
925 delete [] probi;
c0a16418 926 fPixArray->Delete();
927 return kTRUE;
928 }
929
930 // Calculate position of the center-of-gravity around the maximum pixel
931 Double_t xyCOG[2];
c80b29a4 932 FindCOG(xyCOG);
c0a16418 933
934 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
935 pixPtr->Size(0) > pixPtr->Size(1)) break;
936
937 // Sort pixels according to the charge
8c718ff8 938 MaskPeaks(1); // mask local maxima
c0a16418 939 fPixArray->Sort();
8c718ff8 940 MaskPeaks(0); // unmask local maxima
c0a16418 941 Double_t pixMin = 0.01*Pixel(0)->Charge();
942 pixMin = TMath::Min(pixMin,50.);
943
944 // Decrease pixel size and shift pixels to make them centered at
945 // the maximum one
946 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
947 Int_t ix(1);
948 Double_t width = 0;
949 Double_t shift[2] = { 0.0, 0.0 };
2abdae6e 950 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
c0a16418 951 Int_t nPix1 = nPix;
952 nPix = 0;
2abdae6e 953 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
c0a16418 954 {
bf0d3528 955 AliMUONPad* pixPtr2 = Pixel(ipix);
c0a16418 956 if ( nPix >= npadOK // too many pixels already
957 ||
bf0d3528 958 pixPtr2->Charge() < pixMin && pixPtr2->Status() != fgkMustKeep // too low charge
c0a16418 959 )
960 {
961 RemovePixel(ipix);
962 continue;
963 }
2abdae6e 964 for (Int_t i = 0; i < 2; ++i)
c0a16418 965 {
966 if (!i)
967 {
bf0d3528 968 pixPtr2->SetCharge(10);
969 pixPtr2->SetSize(indx, pixPtr2->Size(indx)/2);
970 width = -pixPtr2->Size(indx);
971 pixPtr2->Shift(indx, width);
c0a16418 972 // Shift pixel position
973 if (ix)
974 {
975 ix = 0;
2abdae6e 976 for (Int_t j = 0; j < 2; ++j)
c0a16418 977 {
bf0d3528 978 shift[j] = pixPtr2->Coord(j) - xyCOG[j];
979 shift[j] -= ((Int_t)(shift[j]/pixPtr2->Size(j)/2))*pixPtr2->Size(j)*2;
c0a16418 980 }
981 } // if (ix)
bf0d3528 982 pixPtr2->Shift(0, -shift[0]);
983 pixPtr2->Shift(1, -shift[1]);
8c718ff8 984 ++nPix;
c0a16418 985 }
8c718ff8 986 else if (nPix < npadOK)
c0a16418 987 {
bf0d3528 988 pixPtr2 = new AliMUONPad(*pixPtr2);
989 pixPtr2->Shift(indx, -2*width);
990 pixPtr2->SetStatus(fgkZero);
991 fPixArray->Add(pixPtr2);
8c718ff8 992 ++nPix;
c0a16418 993 }
8c718ff8 994 else continue; // skip adjustment of histo limits
995 for (Int_t j = 0; j < 4; ++j)
c0a16418 996 {
bf0d3528 997 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr2->Coord(j/2));
c0a16418 998 }
999 } // for (Int_t i=0; i<2;
c0a16418 1000 } // for (Int_t ipix=0;
1001
1002 fPixArray->Compress();
c0a16418 1003
1004 AliDebug(2,Form("After shift:"));
2abdae6e 1005 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1006 //Plot(Form("mlem.lc%d",lc+1));
c0a16418 1007
1008 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1009 xyCOG[0],xyCOG[1],
1010 xylim[0],xylim[1],
1011 xylim[2],xylim[3]));
1012
8c718ff8 1013 if (nPix < npadOK)
c0a16418 1014 {
bf0d3528 1015 AliMUONPad* pixPtr2 = Pixel(0);
8c718ff8 1016 // add pixels if the maximum is at the limit of pixel area:
c0a16418 1017 // start from Y-direction
1018 Int_t j = 0;
2abdae6e 1019 for (Int_t i = 3; i > -1; --i)
c0a16418 1020 {
1021 if (nPix < npadOK &&
bf0d3528 1022 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr2->Size(i/2))
c0a16418 1023 {
8c718ff8 1024 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
bf0d3528 1025 AliMUONPad* p = new AliMUONPad(*pixPtr2);
1026 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr2->Size(i/2));
8c718ff8 1027 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
c0a16418 1028 j = TMath::Even (i/2);
1029 p->SetCoord(j, xyCOG[j]);
1030 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
2abdae6e 1031 //StdoutToAliDebug(2,cout << " ---- ";
1032 // p->Print("corners"););
c0a16418 1033 fPixArray->Add(p);
1034 ++nPix;
1035 }
1036 }
1037 }
c0a16418 1038 nPix = fPixArray->GetEntriesFast();
1039 delete [] coef;
1040 delete [] probi;
c0a16418 1041 } // while (1)
1042
1043 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
2abdae6e 1044 //StdoutToAliDebug(2,fPixArray->Print("","full"););
c0a16418 1045
1046 // remove pixels with low signal or low visibility
1047 // Cuts are empirical !!!
c80b29a4 1048 Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,1.);
c0a16418 1049 thresh = TMath::Min (thresh,50.);
c0a16418 1050 Double_t charge = 0;
1051
c0a16418 1052 // Mark pixels which should be removed
2abdae6e 1053 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1054 {
bf0d3528 1055 AliMUONPad* pixPtr2 = Pixel(i);
1056 charge = pixPtr2->Charge();
c0a16418 1057 if (charge < thresh)
1058 {
bf0d3528 1059 pixPtr2->SetCharge(-charge);
c0a16418 1060 }
1061 }
1062
1063 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1064 Int_t near = 0;
2abdae6e 1065 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1066 {
bf0d3528 1067 AliMUONPad* pixPtr2 = Pixel(i);
1068 charge = pixPtr2->Charge();
c0a16418 1069 if (charge > 0) continue;
bf0d3528 1070 near = FindNearest(pixPtr2);
1071 pixPtr2->SetCharge(0);
c0a16418 1072 probi[i] = 0; // make it "invisible"
1073 AliMUONPad* pnear = Pixel(near);
1074 pnear->SetCharge(pnear->Charge() + (-charge));
1075 }
1076 Mlem(cluster,coef,probi,2);
1077
1078 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
2abdae6e 1079 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1080 //Plot("mlem.beforesplit");
c0a16418 1081
8c718ff8 1082 // Update histogram
2abdae6e 1083 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1084 {
bf0d3528 1085 AliMUONPad* pixPtr2 = Pixel(i);
1086 Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr2->Coord(0));
1087 Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr2->Coord(1));
1088 fHistMlem->SetBinContent(ix, iy, pixPtr2->Charge());
c0a16418 1089 }
1090
1091 // Try to split into clusters
1092 Bool_t ok = kTRUE;
c80b29a4 1093 if (fHistMlem->GetSum() < 1)
c0a16418 1094 {
1095 ok = kFALSE;
1096 }
1097 else
1098 {
c80b29a4 1099 fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
c0a16418 1100 }
1101
1102 delete [] coef;
1103 delete [] probi;
c0a16418 1104 fPixArray->Delete();
1105
1106 return ok;
1107}
1108
8c718ff8 1109//_____________________________________________________________________________
1110void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1111{
1112 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1113 /// - to avoid loosing low charge pixels after sorting)
1114
1115 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1116 AliMUONPad* pix = Pixel(i);
1117 if (pix->Status() == fgkMustKeep) {
1118 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1119 else pix->SetCharge(pix->Charge()-10000.);
1120 }
1121 }
1122}
1123
c0a16418 1124//_____________________________________________________________________________
1125void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1126 Double_t* coef, Double_t* probi,
1127 Int_t nIter)
1128{
1129 /// Use MLEM to find pixel charges
1130
1131 Int_t nPix = fPixArray->GetEntriesFast();
1132
1133 Int_t npad = cluster.Multiplicity();
1134
1135 Double_t* probi1 = new Double_t[nPix];
2abdae6e 1136 Double_t probMax = TMath::MaxElement(nPix,probi);
c0a16418 1137
2abdae6e 1138 for (Int_t iter = 0; iter < nIter; ++iter)
c0a16418 1139 {
1140 // Do iterations
2abdae6e 1141 for (Int_t ipix = 0; ipix < nPix; ++ipix)
c0a16418 1142 {
9bbd7f60 1143 Pixel(ipix)->SetChargeBackup(0);
c0a16418 1144 // Correct each pixel
1145 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1146 Double_t sum = 0;
1147 probi1[ipix] = probMax;
2abdae6e 1148 for (Int_t j = 0; j < npad; ++j)
c0a16418 1149 {
1150 AliMUONPad* pad = cluster.Pad(j);
8c718ff8 1151 //if (pad->Status() < 0) continue;
1152 if (pad->Status() != fgkZero) continue;
c0a16418 1153 Double_t sum1 = 0;
1154 Int_t indx1 = j*nPix;
1155 Int_t indx = indx1 + ipix;
1156 // Calculate expectation
2abdae6e 1157 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1158 {
1159 sum1 += Pixel(i)->Charge()*coef[indx1+i];
2abdae6e 1160 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
c0a16418 1161 }
2abdae6e 1162 if ( pad->IsSaturated() && sum1 > pad->Charge() )
c0a16418 1163 {
c0a16418 1164 // correct for pad charge overflows
1165 probi1[ipix] -= coef[indx];
1166 continue;
2abdae6e 1167 //sum1 = pad->Charge();
c0a16418 1168 }
1169
2abdae6e 1170 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
c0a16418 1171 } // for (Int_t j=0;
1172 AliMUONPad* pixPtr = Pixel(ipix);
1173 if (probi1[ipix] > 1.e-6)
1174 {
9bbd7f60 1175 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1176 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
c0a16418 1177 }
2abdae6e 1178 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
c0a16418 1179 } // for (Int_t ipix=0;
2abdae6e 1180 Double_t qTot = 0;
9bbd7f60 1181 for (Int_t i = 0; i < nPix; ++i) {
1182 AliMUONPad* pixPtr = Pixel(i);
2abdae6e 1183 pixPtr->RevertCharge();
1184 qTot += pixPtr->Charge();
1185 }
1186 if (qTot < 1.e-6) {
1187 // Can happen in clusters with large number of overflows - speeding up
1188 delete [] probi1;
1189 return;
9bbd7f60 1190 }
c0a16418 1191 } // for (Int_t iter=0;
1192 delete [] probi1;
1193}
1194
1195//_____________________________________________________________________________
c80b29a4 1196void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
c0a16418 1197{
1198 /// Calculate position of the center-of-gravity around the maximum pixel
1199
1200 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1201 Int_t i1 = -9, j1 = -9;
c80b29a4 1202 fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1203 Int_t nx = fHistMlem->GetNbinsX();
1204 Int_t ny = fHistMlem->GetNbinsY();
1205 Double_t thresh = fHistMlem->GetMaximum()/10;
c0a16418 1206 Double_t x, y, cont, xq=0, yq=0, qq=0;
1207
2abdae6e 1208 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1209 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
c80b29a4 1210 y = fHistMlem->GetYaxis()->GetBinCenter(i);
2abdae6e 1211 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
c80b29a4 1212 cont = fHistMlem->GetCellContent(j,i);
c0a16418 1213 if (cont < thresh) continue;
1214 if (i != i1) {i1 = i; nsumy++;}
1215 if (j != j1) {j1 = j; nsumx++;}
c80b29a4 1216 x = fHistMlem->GetXaxis()->GetBinCenter(j);
c0a16418 1217 xq += x*cont;
1218 yq += y*cont;
1219 qq += cont;
1220 nsum++;
1221 }
1222 }
1223
1224 Double_t cmax = 0;
1225 Int_t i2 = 0, j2 = 0;
1226 x = y = 0;
1227 if (nsumy == 1) {
1228 // one bin in Y - add one more (with the largest signal)
2abdae6e 1229 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
c0a16418 1230 if (i == iymax) continue;
2abdae6e 1231 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
c80b29a4 1232 cont = fHistMlem->GetCellContent(j,i);
c0a16418 1233 if (cont > cmax) {
1234 cmax = cont;
c80b29a4 1235 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1236 y = fHistMlem->GetYaxis()->GetBinCenter(i);
c0a16418 1237 i2 = i;
1238 j2 = j;
1239 }
1240 }
1241 }
1242 xq += x*cmax;
1243 yq += y*cmax;
1244 qq += cmax;
1245 if (i2 != i1) nsumy++;
1246 if (j2 != j1) nsumx++;
1247 nsum++;
1248 } // if (nsumy == 1)
1249
1250 if (nsumx == 1) {
1251 // one bin in X - add one more (with the largest signal)
1252 cmax = x = y = 0;
2abdae6e 1253 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
c0a16418 1254 if (j == ixmax) continue;
2abdae6e 1255 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
c80b29a4 1256 cont = fHistMlem->GetCellContent(j,i);
c0a16418 1257 if (cont > cmax) {
1258 cmax = cont;
c80b29a4 1259 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1260 y = fHistMlem->GetYaxis()->GetBinCenter(i);
c0a16418 1261 i2 = i;
1262 j2 = j;
1263 }
1264 }
1265 }
1266 xq += x*cmax;
1267 yq += y*cmax;
1268 qq += cmax;
1269 if (i2 != i1) nsumy++;
1270 if (j2 != j1) nsumx++;
1271 nsum++;
1272 } // if (nsumx == 1)
1273
1274 xyc[0] = xq/qq; xyc[1] = yq/qq;
1275 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1276}
1277
1278//_____________________________________________________________________________
1279Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1280{
1281/// Find the pixel nearest to the given one
1282/// (algorithm may be not very efficient)
1283
1284 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1285 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1286 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1287 AliMUONPad *pixPtr;
1288
2abdae6e 1289 for (Int_t i = 0; i < nPix; ++i) {
c0a16418 1290 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
8c718ff8 1291 if (pixPtr == pixPtr0 || pixPtr->Charge() < 0.5) continue;
c0a16418 1292 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1293 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1294 r = dx *dx + dy * dy;
1295 if (r < rmin) { rmin = r; imin = i; }
1296 }
1297 return imin;
1298}
1299
c0a16418 1300//_____________________________________________________________________________
1301void
1302AliMUONClusterFinderMLEM::Paint(Option_t*)
1303{
1304 /// Paint cluster and pixels
1305
1306 AliMpArea area(fPreCluster->Area());
1307
1308 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1309
1310 gVirtualX->SetFillStyle(1001);
1311 gVirtualX->SetFillColor(3);
1312 gVirtualX->SetLineColor(3);
1313
1314 Double_t s(1.0);
1315
1316 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1317 {
1318 AliMUONPad* pixel = Pixel(i);
1319
1320 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1321 pixel->Coord(1)-pixel->Size(1)*s,
1322 pixel->Coord(0)+pixel->Size(0)*s,
1323 pixel->Coord(1)+pixel->Size(1)*s);
1324
1325// for ( Int_t sign = -1; sign < 2; sign +=2 )
1326// {
1327// gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1328// pixel->Coord(1) + sign*pixel->Size(1),
1329// pixel->Coord(0) + pixel->Size(0),
1330// pixel->Coord(1) - sign*pixel->Size(1)
1331// );
1332// }
1333 }
1334
1335
1336 gVirtualX->SetFillStyle(0);
1337
1338 fPreCluster->Paint();
1339
1340 gVirtualX->SetLineColor(1);
1341 gVirtualX->SetLineWidth(2);
1342 gVirtualX->SetFillStyle(0);
1343 gVirtualX->SetTextColor(1);
1344 gVirtualX->SetTextAlign(22);
1345
1346 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1347 {
1348 AliMUONPad* pixel = Pixel(i);
1349 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1350 pixel->Coord(1)-pixel->Size(1),
1351 pixel->Coord(0)+pixel->Size(0),
1352 pixel->Coord(1)+pixel->Size(1));
1353 gVirtualX->SetTextSize(pixel->Size(0)*60);
1354
1355 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1356 }
1357}
1358
1359//_____________________________________________________________________________
1360Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1361{
1362/// Find local maxima in pixel space for large preclusters in order to
1363/// try to split them into smaller pieces (to speed up the MLEM procedure)
1364/// or to find additional fitting seeds if clusters were not completely resolved
1365
1366 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1367
c0a16418 1368 Double_t xylim[4] = {999, 999, 999, 999};
1369
1370 Int_t nPix = pixArray->GetEntriesFast();
1371 AliMUONPad *pixPtr = 0;
2abdae6e 1372 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
c0a16418 1373 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
2abdae6e 1374 for (Int_t i = 0; i < 4; ++i)
c0a16418 1375 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1376 }
2abdae6e 1377 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
c0a16418 1378
1379 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1380 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
c80b29a4 1381 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1382 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
2abdae6e 1383 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
c0a16418 1384 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
c80b29a4 1385 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
c0a16418 1386 }
1387// if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1388
2abdae6e 1389 Int_t nMax = 0, indx, nxy = ny * nx;
1390 Int_t *isLocalMax = new Int_t[nxy];
1391 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
c0a16418 1392
2abdae6e 1393 for (Int_t i = 1; i <= ny; ++i) {
c0a16418 1394 indx = (i-1) * nx;
2abdae6e 1395 for (Int_t j = 1; j <= nx; ++j) {
c80b29a4 1396 if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
c0a16418 1397 //if (isLocalMax[indx+j-1] < 0) continue;
1398 if (isLocalMax[indx+j-1] != 0) continue;
c80b29a4 1399 FlagLocalMax(fHistAnode, i, j, isLocalMax);
c0a16418 1400 }
1401 }
1402
2abdae6e 1403 for (Int_t i = 1; i <= ny; ++i) {
c0a16418 1404 indx = (i-1) * nx;
2abdae6e 1405 for (Int_t j = 1; j <= nx; ++j) {
c0a16418 1406 if (isLocalMax[indx+j-1] > 0) {
1407 localMax[nMax] = indx + j - 1;
c80b29a4 1408 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
1409 ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
85a72902 1410 if (nMax > 99) break;
c0a16418 1411 }
1412 }
85a72902 1413 if (nMax > 99) {
1414 AliError(" Too many local maxima !!!");
1415 break;
1416 }
c0a16418 1417 }
1418 if (fDebug) cout << " Local max: " << nMax << endl;
2abdae6e 1419 delete [] isLocalMax;
c0a16418 1420 return nMax;
1421}
1422
1423//_____________________________________________________________________________
1424void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1425{
1426/// Flag pixels (whether or not local maxima)
1427
1428 Int_t nx = hist->GetNbinsX();
1429 Int_t ny = hist->GetNbinsY();
1430 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1431 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1432
2abdae6e 1433 Int_t ie = i + 2, je = j + 2;
1434 for (Int_t i1 = i-1; i1 < ie; ++i1) {
c0a16418 1435 if (i1 < 1 || i1 > ny) continue;
1436 indx1 = (i1 - 1) * nx;
2abdae6e 1437 for (Int_t j1 = j-1; j1 < je; ++j1) {
c0a16418 1438 if (j1 < 1 || j1 > nx) continue;
1439 if (i == i1 && j == j1) continue;
1440 indx2 = indx1 + j1 - 1;
1441 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1442 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1443 else if (cont > cont1) isLocalMax[indx2] = -1;
1444 else { // the same charge
1445 isLocalMax[indx] = 1;
1446 if (isLocalMax[indx2] == 0) {
1447 FlagLocalMax(hist, i1, j1, isLocalMax);
1448 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1449 else isLocalMax[indx2] = -1;
1450 }
1451 }
1452 }
1453 }
1454 isLocalMax[indx] = 1; // local maximum
1455}
1456
1457//_____________________________________________________________________________
1458void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1459 Int_t *localMax, Int_t iMax)
1460{
1461/// Find pixel cluster around local maximum \a iMax and pick up pads
1462/// overlapping with it
1463
2abdae6e 1464 /* Just for check
1465 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1466 c->cd();
1467 hist->Draw("lego1Fb"); // debug
1468 c->Update();
1469 Int_t tmp;
1470 cin >> tmp;
1471 */
c80b29a4 1472 Int_t nx = fHistAnode->GetNbinsX();
1473 Int_t ny = fHistAnode->GetNbinsY();
c0a16418 1474 Int_t ic = localMax[iMax] / nx + 1;
1475 Int_t jc = localMax[iMax] % nx + 1;
2abdae6e 1476 Int_t nxy = ny * nx;
1477 Bool_t *used = new Bool_t[nxy];
1478 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
c0a16418 1479
1480 // Drop all pixels from the array - pick up only the ones from the cluster
1481 fPixArray->Delete();
1482
c80b29a4 1483 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1484 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1485 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1486 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1487 Double_t cont = fHistAnode->GetCellContent(jc,ic);
c0a16418 1488 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1489 used[(ic-1)*nx+jc-1] = kTRUE;
c80b29a4 1490 AddBinSimple(fHistAnode, ic, jc);
2abdae6e 1491 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
c0a16418 1492
1493 Int_t nPix = fPixArray->GetEntriesFast();
1494 Int_t npad = cluster.Multiplicity();
1495
2abdae6e 1496 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1497 {
1498 AliMUONPad* pixPtr = Pixel(i);
1499 pixPtr->SetSize(0,wx);
1500 pixPtr->SetSize(1,wy);
1501 }
1502
1503 // Pick up pads which overlap with found pixels
2abdae6e 1504 for (Int_t i = 0; i < npad; ++i)
c0a16418 1505 {
8c718ff8 1506 //cluster.Pad(i)->SetStatus(-1);
1507 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
c0a16418 1508 }
1509
2abdae6e 1510 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1511 {
1512 AliMUONPad* pixPtr = Pixel(i);
2abdae6e 1513 for (Int_t j = 0; j < npad; ++j)
c0a16418 1514 {
1515 AliMUONPad* pad = cluster.Pad(j);
8c718ff8 1516 //if (pad->Status() == 0) continue;
1517 if (pad->Status() == fgkZero) continue;
c0a16418 1518 if ( Overlap(*pad,*pixPtr) )
1519 {
8c718ff8 1520 //pad->SetStatus(0);
1521 pad->SetStatus(fgkZero);
2abdae6e 1522 if (fDebug) { cout << j << " "; pad->Print("full"); }
c0a16418 1523 }
1524 }
1525 }
1526
2abdae6e 1527 delete [] used;
1528}
1529
1530//_____________________________________________________________________________
1531void
c80b29a4 1532AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
2abdae6e 1533{
1534 /// Add adjacent bins (+-1 in X and Y) to the cluster
1535
c80b29a4 1536 Int_t nx = hist->GetNbinsX();
1537 Int_t ny = hist->GetNbinsY();
1538 Double_t cont1, cont = hist->GetCellContent(jc,ic);
2abdae6e 1539 AliMUONPad *pixPtr = 0;
1540
1541 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1542 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1543 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
c80b29a4 1544 cont1 = hist->GetCellContent(j,i);
2abdae6e 1545 if (cont1 > cont) continue;
1546 if (cont1 < 0.5) continue;
c80b29a4 1547 pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j),
1548 hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
2abdae6e 1549 fPixArray->Add(pixPtr);
1550 }
1551 }
c0a16418 1552}
1553
1554//_____________________________________________________________________________
1555AliMUONClusterFinderMLEM&
1556AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1557{
1558/// Protected assignement operator
1559
1560 if (this == &rhs) return *this;
1561
1562 AliFatal("Not implemented.");
1563
1564 return *this;
1565}
1566
c0a16418 1567//_____________________________________________________________________________
1568void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1569{
2abdae6e 1570 /// Add virtual pad (with small charge) to improve fit for clusters
1571 /// with number of pads == 2 per direction
c0a16418 1572
2abdae6e 1573 // Find out non-bending and bending planes
1574 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
c0a16418 1575
2abdae6e 1576 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1577 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1578 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1579 nonb[0] = 0;
1580 nonb[1] = 1;
1581 }
1582
1583 Bool_t same = kFALSE;
1584 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1585
1586 AliMpIntPair cn;
1587 Bool_t check[2] = {kFALSE, kFALSE};
1588 Int_t nxy[2];
1589 nxy[0] = nxy[1] = 0;
1590 for (Int_t inb = 0; inb < 2; ++inb) {
1591 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1592 if (inb == 0 && cn.GetFirst() == 2) check[inb] = kTRUE; // check non-bending plane
1593 else if (inb == 1 && cn.GetSecond() == 2) check[inb] = kTRUE; // check bending plane
1594 if (same) {
1595 nxy[0] = TMath::Max (nxy[0], cn.GetFirst());
1596 nxy[1] = TMath::Max (nxy[1], cn.GetSecond());
1597 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1598 else if (inb == 1 && cn.GetSecond() < 2) nonb[inb] = !nonb[inb];
c0a16418 1599 }
1600 }
2abdae6e 1601 if (same) {
1602 if (nxy[0] > 2) check[0] = kFALSE;
1603 if (nxy[1] > 2) check[1] = kFALSE;
1604 }
1605 if (!check[0] && !check[1]) return;
1606
1607 for (Int_t inb = 0; inb < 2; ++inb) {
1608 if (!check[inb]) continue;
1609
1610 // Find pads with maximum and next to maximum charges
1611 Int_t maxPads[2] = {-1, -1};
1612 Double_t amax[2] = {0};
1613 Int_t mult = cluster.Multiplicity();
1614 for (Int_t j = 0; j < mult; ++j) {
1615 AliMUONPad *pad = cluster.Pad(j);
1616 if (pad->Cathode() != nonb[inb]) continue;
1617 for (Int_t j2 = 0; j2 < 2; ++j2) {
1618 if (pad->Charge() > amax[j2]) {
1619 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1620 amax[j2] = pad->Charge();
1621 maxPads[j2] = j;
1622 break;
1623 }
c0a16418 1624 }
1625 }
c0a16418 1626
2abdae6e 1627 // Find min and max dimensions of the cluster
1628 Double_t limits[2] = {9999, -9999};
1629 for (Int_t j = 0; j < mult; ++j) {
1630 AliMUONPad *pad = cluster.Pad(j);
1631 if (pad->Cathode() != nonb[inb]) continue;
1632 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1633 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
c0a16418 1634 }
2abdae6e 1635
1636 // Loop over max and next to max pads
1637 Bool_t add = kFALSE;
1638 Int_t idirAdd = 0;
1639 for (Int_t j = 0; j < 2; ++j) {
1640 if (j == 1) {
1641 if (maxPads[j] < 0) continue;
1642 if (!add) break;
1643 if (amax[1] / amax[0] < 0.5) break;
c0a16418 1644 }
2abdae6e 1645 // Check if pad at the cluster limit
1646 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1647 Int_t idir = 0;
1648 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1649 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1650 else {
1651 //cout << " *** Pad not at the cluster limit: " << j << endl;
1652 break;
c0a16418 1653 }
2abdae6e 1654 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1655
1656 // Add pad (if it exists)
1657 TVector2 pos;
1658 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1659 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1660 //AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1661 AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kFALSE);
1662 if (!mppad.IsValid()) continue; // non-existing pad
1663 cn = mppad.GetIndices();
1664 AliMUONPad muonPad(fDetElemId, nonb[inb], cn.GetFirst(), cn.GetSecond(),
1665 mppad.Position().X(), mppad.Position().Y(),
1666 mppad.Dimensions().X(), mppad.Dimensions().Y(), 0);
1667 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, 5.));
8c718ff8 1668 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1669 else muonPad.SetCharge(TMath::Min (amax[j]/15, 6.));
2abdae6e 1670 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1671 muonPad.SetReal(kFALSE);
1672 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1673 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1674 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1675 cluster.AddPad(muonPad); // add pad to the cluster
1676 add = kTRUE;
1677 idirAdd = idir;
1678 }
1679 }
c0a16418 1680}
1681
1682//_____________________________________________________________________________
1683void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1684 Int_t &nInX, Int_t &nInY) const
1685{
1686 /// Find number of pads in X and Y-directions (excluding virtual ones and
1687 /// overflows)
1688
8c718ff8 1689 //Int_t statusToTest = 1;
1690 Int_t statusToTest = fgkUseForFit;
c0a16418 1691
8c718ff8 1692 //if ( nInX < 0 ) statusToTest = 0;
1693 if ( nInX < 0 ) statusToTest = fgkZero;
c0a16418 1694
1695 Bool_t mustMatch(kTRUE);
1696
1697 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1698
1699 nInX = cn.GetFirst();
1700 nInY = cn.GetSecond();
1701}
1702
1703//_____________________________________________________________________________
1704void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1705{
1706 /// Remove pixel at index i
1707 AliMUONPad* pixPtr = Pixel(i);
1708 fPixArray->RemoveAt(i);
1709 delete pixPtr;
1710}
1711
1712//_____________________________________________________________________________
1713void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1714{
1715/// Process simple cluster (small number of pads) without EM-procedure
1716
1717 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1718 Double_t parOk[3] = {0.};
1719 TObjArray *clusters[1];
1720 clusters[0] = fPixArray;
1721
1722 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1723
2abdae6e 1724 Int_t mult = cluster.Multiplicity();
1725 for (Int_t i = 0; i < mult; ++i)
c0a16418 1726 {
1727 AliMUONPad* pad = cluster.Pad(i);
8c718ff8 1728 /*
9bbd7f60 1729 if ( pad->IsSaturated())
c0a16418 1730 {
1731 pad->SetStatus(-9);
1732 }
1733 else
1734 {
1735 pad->SetStatus(1);
1736 }
8c718ff8 1737 */
1738 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
c0a16418 1739 }
c80b29a4 1740 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
c0a16418 1741}
1742
1743//_____________________________________________________________________________
1744AliMUONPad*
1745AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1746{
1747 /// Returns pixel at index i
1748 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1749}
1750
1751//_____________________________________________________________________________
1752void
1753AliMUONClusterFinderMLEM::Print(Option_t* what) const
1754{
1755 /// printout
1756 TString swhat(what);
1757 swhat.ToLower();
1758 if ( swhat.Contains("precluster") )
1759 {
1760 if ( fPreCluster) fPreCluster->Print();
1761 }
1762}
1763
1764