1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 //-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
19 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
20 // of new IO (à la PHOS)
21 // Mar 2007, Aleksei Pavlinov - new algoritmh of pseudo clusters
22 //////////////////////////////////////////////////////////////////////////////
23 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
24 // unfolds the clusters having several local maxima.
25 // Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
26 // EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all
27 // parameters including input digits branch title, thresholds etc.)
28 // This TTask is normally called from Reconstructioner, but can as well be used in
31 // root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")
32 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
33 // //reads gAlice from header file "..."
34 // root [1] cl->ExecuteTask()
35 // //finds RecPoints in all events stored in galice.root
36 // root [2] cl->SetDigitsBranch("digits2")
37 // //sets another title for Digitis (input) branch
38 // root [3] cl->SetRecPointsBranch("recp2")
39 // //sets another title four output branches
40 // root [4] cl->SetTowerLocalMaxCut(0.03)
41 // //set clusterization parameters
42 // root [5] cl->ExecuteTask("deb all time")
43 // //once more finds RecPoints options are
44 // // deb - print number of found rec points
45 // // deb all - print number of found RecPoints and some their characteristics
46 // // time - print benchmarking results
48 // --- ROOT system ---
59 #include <TBenchmark.h>
63 // --- Standard library ---
66 // --- AliRoot header files ---
67 #include "AliRunLoader.h"
70 #include "AliEMCALClusterizerv1.h"
71 #include "AliEMCALRecPoint.h"
72 #include "AliEMCALDigit.h"
73 #include "AliEMCALDigitizer.h"
75 #include "AliEMCALGeometry.h"
76 #include "AliEMCALHistoUtilities.h"
77 #include "AliEMCALRecParam.h"
78 #include "AliEMCALReconstructor.h"
79 #include "AliCDBManager.h"
82 #include "AliCDBEntry.h"
84 ClassImp(AliEMCALClusterizerv1)
86 //____________________________________________________________________________
87 AliEMCALClusterizerv1::AliEMCALClusterizerv1()
88 : AliEMCALClusterizer(),
89 fHists(0),fPointE(0),fPointL1(0),fPointL2(0),
90 fPointDis(0),fPointMult(0),fDigitAmp(0),fMaxE(0),
91 fMaxL1(0),fMaxL2(0),fMaxDis(0),fGeom(0),
94 fNumberOfECAClusters(0),fNTowerInGroup(0),fCalibData(0),
95 fADCchannelECA(0.),fADCpedestalECA(0.),fECAClusteringThreshold(0.),fECALocMaxCut(0.),
96 fECAW0(0.),fTimeCut(0.),fMinECut(0.)
98 // ctor with the indication of the file where header Tree and digits Tree are stored
104 //____________________________________________________________________________
105 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
110 //____________________________________________________________________________
111 Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp, Int_t AbsId)
114 // Convert digitized amplitude into energy.
115 // Calibration parameters are taken from calibration data base for raw data,
116 // or from digitizer parameters for simulated data.
121 AliFatal("Did not get geometry from EMCALLoader") ;
130 Bool_t bCell = fGeom->GetCellIndex(AbsId, iSupMod, nModule, nIphi, nIeta) ;
132 fGeom->PrintGeometry();
133 Error("Calibrate()"," Wrong cell id number : %i", AbsId);
137 fGeom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
139 fADCchannelECA = fCalibData->GetADCchannel (iSupMod,ieta,iphi);
140 fADCpedestalECA = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
142 return -fADCpedestalECA + amp * fADCchannelECA ;
145 else //Return energy with default parameters if calibration is not available
146 return -fADCpedestalECA + amp * fADCchannelECA ;
150 //____________________________________________________________________________
151 void AliEMCALClusterizerv1::Digits2Clusters(Option_t * option)
153 // Steering method to perform clusterization for the current event
156 if(strstr(option,"tim"))
157 gBenchmark->Start("EMCALClusterizer");
159 if(strstr(option,"print"))
162 //Get calibration parameters from file or digitizer default values.
163 GetCalibrationParameters() ;
166 fNumberOfECAClusters = 0;
168 if(strstr(option,"pseudo"))
169 MakeClusters("pseudo") ; //both types
171 MakeClusters("") ; //only the real clusters
178 if(strstr(option,"deb") || strstr(option,"all"))
179 PrintRecPoints(option) ;
181 AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
183 if(strstr(option,"tim")){
184 gBenchmark->Stop("EMCALClusterizer");
185 printf("Exec took %f seconds for Clusterizing",
186 gBenchmark->GetCpuTime("EMCALClusterizer"));
190 //____________________________________________________________________________
191 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
192 Int_t nPar, Float_t * fitparameters) const
194 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
195 // The initial values for fitting procedure are set equal to the positions of local maxima.
196 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
198 gMinuit->mncler(); // Reset Minuit's list of paramters
199 gMinuit->SetPrintLevel(-1) ; // No Printout
200 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
201 // To set the address of the minimization function
202 TList * toMinuit = new TList();
203 toMinuit->AddAt(emcRP,0) ;
204 toMinuit->AddAt(fDigitsArr,1) ;
206 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
208 // filling initial values for fit parameters
209 AliEMCALDigit * digit ;
213 Int_t nDigits = (Int_t) nPar / 3 ;
217 for(iDigit = 0; iDigit < nDigits; iDigit++){
218 digit = maxAt[iDigit];
222 // have to be tune for TRD1; May 31,06
224 // fGeom->AbsToRelNumbering(digit->GetId(), relid) ; // obsolete method
225 // fGeom->PosInAlice(relid, x, z) ; // obsolete method
227 Float_t energy = maxAtEnergy[iDigit] ;
229 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
232 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : x = %f", x ) ;
235 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
238 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : z = %f", z) ;
241 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
244 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : energy = %f", energy) ;
249 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
254 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
255 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
256 gMinuit->SetMaxIterations(5);
257 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
259 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
261 if(ierflg == 4){ // Minimum not found
262 Error("FindFit", "EMCAL Unfolding Fit not converged, cluster abandoned " ) ;
265 for(index = 0; index < nPar; index++){
268 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
269 fitparameters[index] = val ;
277 //____________________________________________________________________________
278 void AliEMCALClusterizerv1::GetCalibrationParameters()
280 // Set calibration parameters:
281 // if calibration database exists, they are read from database,
282 // otherwise, they are taken from digitizer.
284 // It is a user responsilibity to open CDB before reconstruction,
286 // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
288 //Check if calibration is stored in data base
290 if(!fCalibData && (AliCDBManager::Instance()->IsDefaultStorageSet()))
292 AliCDBEntry *entry = (AliCDBEntry*)
293 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
294 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
298 AliFatal("Calibration parameters not found in CDB!");
300 // Please fix it!! Or better just remove it...
303 // //If calibration is not available use default parameters
305 // if ( !emcalLoader->Digitizer() )
306 // emcalLoader->LoadDigitizer();
307 // AliEMCALDigitizer * dig = dynamic_cast<AliEMCALDigitizer*>(emcalLoader->Digitizer());
309 // fADCchannelECA = dig->GetECAchannel() ;
310 // fADCpedestalECA = dig->GetECApedestal();
314 //____________________________________________________________________________
315 void AliEMCALClusterizerv1::Init()
317 // Make all memory allocations which can not be done in default constructor.
318 // Attach the Clusterizer task to the list of EMCAL tasks
320 AliRunLoader *rl = AliRunLoader::GetRunLoader();
321 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL"))
322 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
324 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaulGeometryName());
326 fGeom->GetTransformationForSM(); // Global <-> Local
327 AliInfo(Form("geom 0x%x",fGeom));
330 gMinuit = new TMinuit(100) ;
332 fHists = BookHists();
335 //____________________________________________________________________________
336 void AliEMCALClusterizerv1::InitParameters()
338 // Initializes the parameters for the Clusterizer
339 fNumberOfECAClusters = 0;
341 fNTowerInGroup = 36; //Produces maximum of 80 pseudoclusters per event
343 fECALocMaxCut = 0.03; // ??
345 fTimeCut = 300e-9 ; // 300 ns time cut (to be tuned)
350 const AliEMCALRecParam* recParam = AliEMCALReconstructor::GetRecParam();
352 AliFatal("Reconstruction parameters for EMCAL not set!");
355 fECAClusteringThreshold = recParam->GetClusteringThreshold();
356 fECAW0 = recParam->GetW0();
357 fMinECut = recParam->GetMinECut();
358 AliDebug(1,Form("Reconstruction parameters: fECAClusteringThreshold=%.3f, fECAW=%.3f, fMinECut=%.3f",
359 fECAClusteringThreshold,fECAW0,fMinECut));
364 //____________________________________________________________________________
365 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2) const
367 // Gives the neighbourness of two digits = 0 are not neighbour ; continue searching
369 // = 2 is in different SM; continue searching
370 // neighbours are defined as digits having at least a common vertex
371 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
372 // which is compared to a digit (d2) not yet in a cluster
375 static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
376 static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
377 static Int_t rowdiff, coldiff;
380 fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
381 fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
382 if(nSupMod1 != nSupMod2) return 2; // different SM
384 fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
385 fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
387 rowdiff = TMath::Abs(iphi1 - iphi2);
388 coldiff = TMath::Abs(ieta1 - ieta2) ;
390 // neighbours with at least commom side; May 11, 2007
391 if ((coldiff==0 && abs(rowdiff)==1) || (rowdiff==0 && abs(coldiff)==1)) rv = 1;
393 if (gDebug == 2 && rv==1)
394 printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n",
395 rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);
400 //____________________________________________________________________________
401 Int_t AliEMCALClusterizerv1::AreInGroup(AliEMCALDigit * d1, AliEMCALDigit * d2) const
403 // Tells whether two digits fall within the same supermodule and
404 // tower grouping. The number of towers in a group is controlled by
405 // the parameter nTowersInGroup
406 // = 0 are not in same group but continue searching
408 // = 2 is in different SM, quit from searching
409 // = 3 different tower group, quit from searching
411 // The order of d1 and d2 is important: first (d1) should be a digit
412 // already in a cluster which is compared to a digit (d2) not yet in a cluster
414 //JLK Question: does the quit from searching assume that the digits
415 //are ordered, so that once you are in a different SM, you'll not
416 //find another in the list that will match? How about my TowerGroup search?
419 static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
420 static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
421 static Int_t towerGroup1 = -1, towerGroup2 = -1;
424 fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
425 fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
426 if(nSupMod1 != nSupMod2) return 2; // different SM
428 static Int_t nTowerInSM = fGeom->GetNCellsInSupMod()/fGeom->GetNCellsInModule();
430 //figure out which tower grouping each digit belongs to
431 for(int it = 0; it < nTowerInSM/fNTowerInGroup; it++) {
432 if(nModule1 <= nTowerInSM - it*fNTowerInGroup) towerGroup1 = it;
433 if(nModule2 <= nTowerInSM - it*fNTowerInGroup) towerGroup2 = it;
435 if(towerGroup1 != towerGroup2) return 3; //different Towergroup
437 //same SM, same towergroup, we're happy
438 if(towerGroup1 == towerGroup2 && towerGroup2 >= 0)
441 if (gDebug == 2 && rv==1)
442 printf("AreInGroup: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n",
443 rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);
448 //____________________________________________________________________________
449 void AliEMCALClusterizerv1::MakeClusters(char* option)
451 // Steering method to construct the clusters stored in a list of Reconstructed Points
452 // A cluster is defined as a list of neighbour digits
453 // Mar 03, 2007 by PAI
455 if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
459 // Set up TObjArray with pointers to digits to work on
460 TObjArray *digitsC = new TObjArray();
461 TIter nextdigit(fDigitsArr);
462 AliEMCALDigit *digit;
463 while ( (digit = dynamic_cast<AliEMCALDigit*>(nextdigit())) ) {
464 digitsC->AddLast(digit);
467 //Start with pseudoclusters, if option
468 if(strstr(option,"pseudo")) {
470 // New algorithm : will be created one pseudo cluster per module
472 AliDebug(1,Form("Pseudo clustering #digits : %i ",fDigitsArr->GetEntries()));
474 AliEMCALRecPoint *recPoints[12]; // max size is 12 : see fGeom->GetNumberOfSuperModules();
475 for(int i=0; i<12; i++) recPoints[i] = 0;
476 TIter nextdigitC(digitsC) ;
478 // PseudoClusterization starts
479 int nSM = 0; // # of SM
480 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
481 if(fGeom->CheckAbsCellId(digit->GetId()) ) { //Is this an EMCAL digit? Just maing sure...
482 nSM = fGeom->GetSuperModuleNumber(digit->GetId());
483 if(recPoints[nSM] == 0) {
484 recPoints[nSM] = new AliEMCALRecPoint(Form("PC%2.2i", nSM));
485 recPoints[nSM]->SetClusterType(AliESDCaloCluster::kPseudoCluster);
487 recPoints[nSM]->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId()));
490 fNumberOfECAClusters = 0;
491 for(int i=0; i<fGeom->GetNumberOfSuperModules(); i++) { // put non empty rec.points to container
492 if(recPoints[i]) fRecPoints->AddAt(recPoints[i], fNumberOfECAClusters++);
494 AliDebug(1,Form(" Number of PC %d ", fNumberOfECAClusters));
498 // Now do real clusters
501 double e = 0.0, ehs = 0.0;
502 TIter nextdigitC(digitsC);
504 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // clean up digits
505 e = Calibrate(digit->GetAmp(), digit->GetId());
506 AliEMCALHistoUtilities::FillH1(fHists, 10, digit->GetAmp());
507 AliEMCALHistoUtilities::FillH1(fHists, 11, e);
508 if ( e < fMinECut || digit->GetTimeR() > fTimeCut )
509 digitsC->Remove(digit);
513 AliDebug(1,Form("MakeClusters: Number of digits %d -> (e %f), ehs %d\n",
514 fDigitsArr->GetEntries(),fMinECut,ehs));
518 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
519 TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
521 if(fGeom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmp(), digit->GetId()) > fECAClusteringThreshold ) ){
522 // start a new Tower RecPoint
523 if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1) ;
524 AliEMCALRecPoint *recPoint = new AliEMCALRecPoint("") ;
525 fRecPoints->AddAt(recPoint, fNumberOfECAClusters) ;
526 recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters)) ;
527 fNumberOfECAClusters++ ;
529 recPoint->SetClusterType(AliESDCaloCluster::kClusterv1);
531 recPoint->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId())) ;
532 TObjArray clusterDigits;
533 clusterDigits.AddLast(digit);
534 digitsC->Remove(digit) ;
536 AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(),
537 Calibrate(digit->GetAmp(),digit->GetId()), fECAClusteringThreshold));
539 // Grow cluster by finding neighbours
540 TIter nextClusterDigit(&clusterDigits);
541 while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) ) { // scan over digits in cluster
542 TIter nextdigitN(digitsC);
543 AliEMCALDigit *digitN = 0; // digi neighbor
544 while ( (digitN = (AliEMCALDigit *)nextdigitN()) ) { // scan over all digits to look for neighbours
545 if (AreNeighbours(digit, digitN)==1) { // call (digit,digitN) in THAT oder !!!!!
546 recPoint->AddDigit(*digitN, Calibrate(digitN->GetAmp(),digitN->GetId()) ) ;
547 clusterDigits.AddLast(digitN) ;
548 digitsC->Remove(digitN) ;
550 } // scan over digits
551 } // scan over digits already in cluster
553 AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy()));
559 AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast()));
562 void AliEMCALClusterizerv1::MakeUnfolding() const
564 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
567 //____________________________________________________________________________
568 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
570 // Shape of the shower (see EMCAL TDR)
571 // If you change this function, change also the gradient evaluation in ChiSquare()
573 Double_t r4 = r*r*r*r ;
574 Double_t r295 = TMath::Power(r, 2.95) ;
575 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
579 //____________________________________________________________________________
580 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * /*iniTower*/,
582 AliEMCALDigit ** /*maxAt*/,
583 Float_t * /*maxAtEnergy*/) const
585 // Performs the unfolding of a cluster with nMax overlapping showers
587 Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
591 //_____________________________________________________________________________
592 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
594 Double_t * /*x*/, Int_t /*iflag*/)
596 // Calculates the Chi square for the cluster unfolding minimization
597 // Number of parameters, Gradient, Chi squared, parameters, what to do
599 ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
601 //____________________________________________________________________________
602 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
604 // Print clusterizer parameters
606 TString message("\n") ;
608 if( strcmp(GetName(), "") !=0 ){
612 TString taskName(Version()) ;
614 printf("--------------- ");
615 printf(taskName.Data()) ;
617 printf("Clusterizing digits: ");
618 printf("\n ECA Local Maximum cut = %f", fECALocMaxCut);
619 printf("\n ECA Logarithmic weight = %f", fECAW0);
621 printf("\nUnfolding on\n");
623 printf("\nUnfolding off\n");
625 printf("------------------------------------------------------------------");
628 printf("AliEMCALClusterizerv1 not initialized ") ;
631 //____________________________________________________________________________
632 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
634 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
635 if(strstr(option,"deb")) {
636 printf("PrintRecPoints: Clusterization result:") ;
638 printf(" Found %d ECA Rec Points\n ",
639 fRecPoints->GetEntriesFast()) ;
642 if(strstr(option,"all")) {
643 if(strstr(option,"deb")) {
644 printf("\n-----------------------------------------------------------------------\n") ;
645 printf("Clusters in ECAL section\n") ;
646 printf("Index Ene(GeV) Multi Module GX GY GZ lX lY lZ Dispersion Lambda 1 Lambda 2 # of prim Primaries list\n") ;
654 AliEMCALHistoUtilities::FillH1(fHists, 12, double(fRecPoints->GetEntries()));
656 for (index = 0 ; index < fRecPoints->GetEntries() ; index++) {
657 AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(fRecPoints->At(index)) ;
659 //rp->GetGlobalPosition(globalpos);
661 rp->GetLocalPosition(localpos);
663 rp->GetElipsAxis(lambda);
666 primaries = rp->GetPrimaries(nprimaries);
667 if(strstr(option,"deb"))
668 printf("\n%6d %8.4f %3d %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4f %4f %2d : ",
669 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
670 globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(),
671 rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ;
673 if(rp->GetEnergy()>maxE){
674 maxE=rp->GetEnergy();
677 maxDis=rp->GetDispersion();
679 fPointE->Fill(rp->GetEnergy());
680 fPointL1->Fill(lambda[0]);
681 fPointL2->Fill(lambda[1]);
682 fPointDis->Fill(rp->GetDispersion());
683 fPointMult->Fill(rp->GetMultiplicity());
685 if(strstr(option,"deb")){
686 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
687 printf("%d ", primaries[iprimary] ) ;
695 fMaxDis->Fill(maxDis);
697 if(strstr(option,"deb"))
698 printf("\n-----------------------------------------------------------------------\n");
701 TList* AliEMCALClusterizerv1::BookHists()
703 //set up histograms for monitoring clusterizer performance
707 fPointE = new TH1F("00_pointE","point energy", 2000, 0.0, 150.);
708 fPointL1 = new TH1F("01_pointL1","point L1", 1000, 0.0, 3.);
709 fPointL2 = new TH1F("02_pointL2","point L2", 1000, 0.0, 3.);
710 fPointDis = new TH1F("03_pointDisp","point dispersion", 1000, 0.0, 10.);
711 fPointMult = new TH1F("04_pointMult","#cell in point(cluster)", 101, -0.5, 100.5);
712 fDigitAmp = new TH1F("05_digitAmp","Digit Amplitude", 2000, 0.0, 5000.);
713 fMaxE = new TH1F("06_maxE","Max point energy", 2000, 0.0, 150.);
714 fMaxL1 = new TH1F("07_maxL1","Largest (first) of eigenvalue of covariance matrix", 1000, 0.0, 3.);
715 fMaxL2 = new TH1F("08_maxL2","Smalest (second) of eigenvalue of covariace matrix", 1000, 0.0, 3.);
716 fMaxDis = new TH1F("09_maxDis","Point dispersion", 1000, 0.0, 10.); // 9
718 new TH1F("10_adcOfDigits","adc of digits(threshold control)", 1001, -0.5, 1000.5); // 10
719 new TH1F("11_energyOfDigits","energy of digits(threshold control)", 1000, 0.0, 1.); // 11
720 new TH1F("12_numberOfPoints","number of points(clusters)", 101, -0.5, 100.5); // 12
722 return AliEMCALHistoUtilities::MoveHistsToList("EmcalClusterizerv1ControlHists", kFALSE);
725 void AliEMCALClusterizerv1::SaveHists(const char *fn)
727 AliEMCALHistoUtilities::SaveListOfHists(fHists, fn, kTRUE);
730 void AliEMCALClusterizerv1::PrintRecoInfo()
732 printf(" AliEMCALClusterizerv1::PrintRecoInfo() : version %s \n", Version() );
733 TH1F *h = (TH1F*)fHists->At(12);
735 printf(" ## Multiplicity of RecPoints ## \n");
736 for(int i=1; i<=h->GetNbinsX(); i++) {
737 int nbin = int((*h)[i]);
738 int mult = int(h->GetBinCenter(i));
739 if(nbin > 0) printf(" %i : %5.5i %6.3f %% \n", mult, nbin, 100.*nbin/h->GetEntries());
744 void AliEMCALClusterizerv1::DrawLambdasHists()
748 if(fMaxL2) fMaxL2->Draw("same");
750 fMaxDis->Draw("same");