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