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 //Evaluate position, dispersion and other RecPoint properties for EC section
179 for(index = 0; index < fRecPoints->GetEntries(); index++) {
180 if (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->GetClusterType() != AliESDCaloCluster::kPseudoCluster)
181 dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->EvalAll(fECAW0,fDigitsArr) ;
186 for(index = 0; index < fRecPoints->GetEntries(); index++) {
187 (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->SetIndexInList(index) ;
188 (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->Print();
193 if(strstr(option,"deb") || strstr(option,"all"))
194 PrintRecPoints(option) ;
196 AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
198 if(strstr(option,"tim")){
199 gBenchmark->Stop("EMCALClusterizer");
200 printf("Exec took %f seconds for Clusterizing",
201 gBenchmark->GetCpuTime("EMCALClusterizer"));
205 //____________________________________________________________________________
206 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
207 Int_t nPar, Float_t * fitparameters) const
209 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
210 // The initial values for fitting procedure are set equal to the positions of local maxima.
211 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
213 gMinuit->mncler(); // Reset Minuit's list of paramters
214 gMinuit->SetPrintLevel(-1) ; // No Printout
215 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
216 // To set the address of the minimization function
217 TList * toMinuit = new TList();
218 toMinuit->AddAt(emcRP,0) ;
219 toMinuit->AddAt(fDigitsArr,1) ;
221 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
223 // filling initial values for fit parameters
224 AliEMCALDigit * digit ;
228 Int_t nDigits = (Int_t) nPar / 3 ;
232 for(iDigit = 0; iDigit < nDigits; iDigit++){
233 digit = maxAt[iDigit];
237 // have to be tune for TRD1; May 31,06
239 // fGeom->AbsToRelNumbering(digit->GetId(), relid) ; // obsolete method
240 // fGeom->PosInAlice(relid, x, z) ; // obsolete method
242 Float_t energy = maxAtEnergy[iDigit] ;
244 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
247 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : x = %f", x ) ;
250 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
253 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : z = %f", z) ;
256 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
259 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : energy = %f", energy) ;
264 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
269 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
270 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
271 gMinuit->SetMaxIterations(5);
272 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
274 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
276 if(ierflg == 4){ // Minimum not found
277 Error("FindFit", "EMCAL Unfolding Fit not converged, cluster abandoned " ) ;
280 for(index = 0; index < nPar; index++){
283 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
284 fitparameters[index] = val ;
292 //____________________________________________________________________________
293 void AliEMCALClusterizerv1::GetCalibrationParameters()
295 // Set calibration parameters:
296 // if calibration database exists, they are read from database,
297 // otherwise, they are taken from digitizer.
299 // It is a user responsilibity to open CDB before reconstruction,
301 // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
303 //Check if calibration is stored in data base
305 if(!fCalibData && (AliCDBManager::Instance()->IsDefaultStorageSet()))
307 AliCDBEntry *entry = (AliCDBEntry*)
308 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
309 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
313 AliFatal("Calibration parameters not found in CDB!");
315 // Please fix it!! Or better just remove it...
318 // //If calibration is not available use default parameters
320 // if ( !emcalLoader->Digitizer() )
321 // emcalLoader->LoadDigitizer();
322 // AliEMCALDigitizer * dig = dynamic_cast<AliEMCALDigitizer*>(emcalLoader->Digitizer());
324 // fADCchannelECA = dig->GetECAchannel() ;
325 // fADCpedestalECA = dig->GetECApedestal();
329 //____________________________________________________________________________
330 void AliEMCALClusterizerv1::Init()
332 // Make all memory allocations which can not be done in default constructor.
333 // Attach the Clusterizer task to the list of EMCAL tasks
335 AliRunLoader *rl = AliRunLoader::GetRunLoader();
336 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL"))
337 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
339 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaulGeometryName());
341 AliInfo(Form("geom 0x%x",fGeom));
344 gMinuit = new TMinuit(100) ;
346 fHists = BookHists();
349 //____________________________________________________________________________
350 void AliEMCALClusterizerv1::InitParameters()
352 // Initializes the parameters for the Clusterizer
353 fNumberOfECAClusters = 0;
355 fNTowerInGroup = 36; //Produces maximum of 80 pseudoclusters per event
357 fECALocMaxCut = 0.03; // ??
359 fTimeCut = 300e-9 ; // 300 ns time cut (to be tuned)
364 const AliEMCALRecParam* recParam = AliEMCALReconstructor::GetRecParam();
366 AliFatal("Reconstruction parameters for EMCAL not set!");
369 fECAClusteringThreshold = recParam->GetClusteringThreshold();
370 fECAW0 = recParam->GetW0();
371 fMinECut = recParam->GetMinECut();
372 AliDebug(1,Form("Reconstruction parameters: fECAClusteringThreshold=%.3f, fECAW=%.3f, fMinECut=%.3f",
373 fECAClusteringThreshold,fECAW0,fMinECut));
378 //____________________________________________________________________________
379 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2) const
381 // Gives the neighbourness of two digits = 0 are not neighbour ; continue searching
383 // = 2 is in different SM; continue searching
384 // neighbours are defined as digits having at least a common vertex
385 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
386 // which is compared to a digit (d2) not yet in a cluster
389 static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
390 static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
391 static Int_t rowdiff, coldiff;
394 fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
395 fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
396 if(nSupMod1 != nSupMod2) return 2; // different SM
398 fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
399 fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
401 rowdiff = TMath::Abs(iphi1 - iphi2);
402 coldiff = TMath::Abs(ieta1 - ieta2) ;
404 // neighbours with at least commom side; May 11, 2007
405 if ((coldiff==0 && abs(rowdiff)==1) || (rowdiff==0 && abs(coldiff)==1)) rv = 1;
407 if (gDebug == 2 && rv==1)
408 printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n",
409 rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);
414 //____________________________________________________________________________
415 Int_t AliEMCALClusterizerv1::AreInGroup(AliEMCALDigit * d1, AliEMCALDigit * d2) const
417 // Tells whether two digits fall within the same supermodule and
418 // tower grouping. The number of towers in a group is controlled by
419 // the parameter nTowersInGroup
420 // = 0 are not in same group but continue searching
422 // = 2 is in different SM, quit from searching
423 // = 3 different tower group, quit from searching
425 // The order of d1 and d2 is important: first (d1) should be a digit
426 // already in a cluster which is compared to a digit (d2) not yet in a cluster
428 //JLK Question: does the quit from searching assume that the digits
429 //are ordered, so that once you are in a different SM, you'll not
430 //find another in the list that will match? How about my TowerGroup search?
433 static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
434 static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
435 static Int_t towerGroup1 = -1, towerGroup2 = -1;
438 fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
439 fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
440 if(nSupMod1 != nSupMod2) return 2; // different SM
442 static Int_t nTowerInSM = fGeom->GetNCellsInSupMod()/fGeom->GetNCellsInModule();
444 //figure out which tower grouping each digit belongs to
445 for(int it = 0; it < nTowerInSM/fNTowerInGroup; it++) {
446 if(nModule1 <= nTowerInSM - it*fNTowerInGroup) towerGroup1 = it;
447 if(nModule2 <= nTowerInSM - it*fNTowerInGroup) towerGroup2 = it;
449 if(towerGroup1 != towerGroup2) return 3; //different Towergroup
451 //same SM, same towergroup, we're happy
452 if(towerGroup1 == towerGroup2 && towerGroup2 >= 0)
455 if (gDebug == 2 && rv==1)
456 printf("AreInGroup: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n",
457 rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);
462 //____________________________________________________________________________
463 void AliEMCALClusterizerv1::MakeClusters(char* option)
465 // Steering method to construct the clusters stored in a list of Reconstructed Points
466 // A cluster is defined as a list of neighbour digits
467 // Mar 03, 2007 by PAI
469 if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
473 // Set up TObjArray with pointers to digits to work on
474 TObjArray *digitsC = new TObjArray();
475 TIter nextdigit(fDigitsArr);
476 AliEMCALDigit *digit;
477 while ( (digit = dynamic_cast<AliEMCALDigit*>(nextdigit())) ) {
478 digitsC->AddLast(digit);
481 //Start with pseudoclusters, if option
482 if(strstr(option,"pseudo")) {
484 // New algorithm : will be created one pseudo cluster per module
486 AliDebug(1,Form("Pseudo clustering #digits : %i ",fDigitsArr->GetEntries()));
488 AliEMCALRecPoint *recPoints[12]; // max size is 12 : see fGeom->GetNumberOfSuperModules();
489 for(int i=0; i<12; i++) recPoints[i] = 0;
490 TIter nextdigitC(digitsC) ;
492 // PseudoClusterization starts
493 int nSM = 0; // # of SM
494 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
495 if(fGeom->CheckAbsCellId(digit->GetId()) ) { //Is this an EMCAL digit? Just maing sure...
496 nSM = fGeom->GetSuperModuleNumber(digit->GetId());
497 if(recPoints[nSM] == 0) {
498 recPoints[nSM] = new AliEMCALRecPoint(Form("PC%2.2i", nSM));
499 recPoints[nSM]->SetClusterType(AliESDCaloCluster::kPseudoCluster);
501 recPoints[nSM]->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId()));
504 fNumberOfECAClusters = 0;
505 for(int i=0; i<fGeom->GetNumberOfSuperModules(); i++) { // put non empty rec.points to container
506 if(recPoints[i]) fRecPoints->AddAt(recPoints[i], fNumberOfECAClusters++);
508 AliDebug(1,Form(" Number of PC %d ", fNumberOfECAClusters));
512 // Now do real clusters
515 double e = 0.0, ehs = 0.0;
516 TIter nextdigitC(digitsC);
518 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // clean up digits
519 e = Calibrate(digit->GetAmp(), digit->GetId());
520 AliEMCALHistoUtilities::FillH1(fHists, 10, digit->GetAmp());
521 AliEMCALHistoUtilities::FillH1(fHists, 11, e);
522 if ( e < fMinECut || digit->GetTimeR() > fTimeCut )
523 digitsC->Remove(digit);
527 AliDebug(1,Form("MakeClusters: Number of digits %d -> (e %f), ehs %d\n",
528 fDigitsArr->GetEntries(),fMinECut,ehs));
532 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
533 TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
535 if(fGeom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmp(), digit->GetId()) > fECAClusteringThreshold ) ){
536 // start a new Tower RecPoint
537 if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1) ;
538 AliEMCALRecPoint *recPoint = new AliEMCALRecPoint("") ;
539 fRecPoints->AddAt(recPoint, fNumberOfECAClusters) ;
540 recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters)) ;
541 fNumberOfECAClusters++ ;
543 recPoint->SetClusterType(AliESDCaloCluster::kClusterv1);
545 recPoint->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId())) ;
546 TObjArray clusterDigits;
547 clusterDigits.AddLast(digit);
548 digitsC->Remove(digit) ;
550 AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(),
551 Calibrate(digit->GetAmp(),digit->GetId()), fECAClusteringThreshold));
553 // Grow cluster by finding neighbours
554 TIter nextClusterDigit(&clusterDigits);
555 while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) ) { // scan over digits in cluster
556 TIter nextdigitN(digitsC);
557 AliEMCALDigit *digitN = 0; // digi neighbor
558 while ( (digitN = (AliEMCALDigit *)nextdigitN()) ) { // scan over all digits to look for neighbours
559 if (AreNeighbours(digit, digitN)==1) { // call (digit,digitN) in THAT oder !!!!!
560 recPoint->AddDigit(*digitN, Calibrate(digitN->GetAmp(),digitN->GetId()) ) ;
561 clusterDigits.AddLast(digitN) ;
562 digitsC->Remove(digitN) ;
564 } // scan over digits
565 } // scan over digits already in cluster
567 AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy()));
573 AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast()));
576 void AliEMCALClusterizerv1::MakeUnfolding() const
578 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
581 //____________________________________________________________________________
582 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
584 // Shape of the shower (see EMCAL TDR)
585 // If you change this function, change also the gradient evaluation in ChiSquare()
587 Double_t r4 = r*r*r*r ;
588 Double_t r295 = TMath::Power(r, 2.95) ;
589 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
593 //____________________________________________________________________________
594 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * /*iniTower*/,
596 AliEMCALDigit ** /*maxAt*/,
597 Float_t * /*maxAtEnergy*/) const
599 // Performs the unfolding of a cluster with nMax overlapping showers
601 Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
605 //_____________________________________________________________________________
606 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
608 Double_t * /*x*/, Int_t /*iflag*/)
610 // Calculates the Chi square for the cluster unfolding minimization
611 // Number of parameters, Gradient, Chi squared, parameters, what to do
613 ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
615 //____________________________________________________________________________
616 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
618 // Print clusterizer parameters
620 TString message("\n") ;
622 if( strcmp(GetName(), "") !=0 ){
626 TString taskName(Version()) ;
628 printf("--------------- ");
629 printf(taskName.Data()) ;
631 printf("Clusterizing digits: ");
632 printf("\n ECA Local Maximum cut = %f", fECALocMaxCut);
633 printf("\n ECA Logarithmic weight = %f", fECAW0);
635 printf("\nUnfolding on\n");
637 printf("\nUnfolding off\n");
639 printf("------------------------------------------------------------------");
642 printf("AliEMCALClusterizerv1 not initialized ") ;
645 //____________________________________________________________________________
646 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
648 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
649 if(strstr(option,"deb")) {
650 printf("PrintRecPoints: Clusterization result:") ;
652 printf(" Found %d ECA Rec Points\n ",
653 fRecPoints->GetEntriesFast()) ;
656 if(strstr(option,"all")) {
657 if(strstr(option,"deb")) {
658 printf("\n-----------------------------------------------------------------------\n") ;
659 printf("Clusters in ECAL section\n") ;
660 printf("Index Ene(GeV) Multi Module GX GY GZ lX lY lZ Dispersion Lambda 1 Lambda 2 # of prim Primaries list\n") ;
668 AliEMCALHistoUtilities::FillH1(fHists, 12, double(fRecPoints->GetEntries()));
670 for (index = 0 ; index < fRecPoints->GetEntries() ; index++) {
671 AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(fRecPoints->At(index)) ;
673 //rp->GetGlobalPosition(globalpos);
675 rp->GetLocalPosition(localpos);
677 rp->GetElipsAxis(lambda);
680 primaries = rp->GetPrimaries(nprimaries);
681 if(strstr(option,"deb"))
682 printf("\n%6d %8.4f %3d %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4f %4f %2d : ",
683 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
684 globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(),
685 rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ;
687 if(rp->GetEnergy()>maxE){
688 maxE=rp->GetEnergy();
691 maxDis=rp->GetDispersion();
693 fPointE->Fill(rp->GetEnergy());
694 fPointL1->Fill(lambda[0]);
695 fPointL2->Fill(lambda[1]);
696 fPointDis->Fill(rp->GetDispersion());
697 fPointMult->Fill(rp->GetMultiplicity());
699 if(strstr(option,"deb")){
700 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
701 printf("%d ", primaries[iprimary] ) ;
709 fMaxDis->Fill(maxDis);
711 if(strstr(option,"deb"))
712 printf("\n-----------------------------------------------------------------------\n");
715 TList* AliEMCALClusterizerv1::BookHists()
717 //set up histograms for monitoring clusterizer performance
721 fPointE = new TH1F("00_pointE","point energy", 2000, 0.0, 150.);
722 fPointL1 = new TH1F("01_pointL1","point L1", 1000, 0.0, 3.);
723 fPointL2 = new TH1F("02_pointL2","point L2", 1000, 0.0, 3.);
724 fPointDis = new TH1F("03_pointDisp","point dispersion", 1000, 0.0, 10.);
725 fPointMult = new TH1F("04_pointMult","#cell in point(cluster)", 101, -0.5, 100.5);
726 fDigitAmp = new TH1F("05_digitAmp","Digit Amplitude", 2000, 0.0, 5000.);
727 fMaxE = new TH1F("06_maxE","Max point energy", 2000, 0.0, 150.);
728 fMaxL1 = new TH1F("07_maxL1","Largest (first) of eigenvalue of covariance matrix", 1000, 0.0, 3.);
729 fMaxL2 = new TH1F("08_maxL2","Smalest (second) of eigenvalue of covariace matrix", 1000, 0.0, 3.);
730 fMaxDis = new TH1F("09_maxDis","Point dispersion", 1000, 0.0, 10.); // 9
732 new TH1F("10_adcOfDigits","adc of digits(threshold control)", 1001, -0.5, 1000.5); // 10
733 new TH1F("11_energyOfDigits","energy of digits(threshold control)", 1000, 0.0, 1.); // 11
734 new TH1F("12_numberOfPoints","number of points(clusters)", 101, -0.5, 100.5); // 12
736 return AliEMCALHistoUtilities::MoveHistsToList("EmcalClusterizerv1ControlHists", kFALSE);
739 void AliEMCALClusterizerv1::SaveHists(const char *fn)
741 AliEMCALHistoUtilities::SaveListOfHists(fHists, fn, kTRUE);
744 void AliEMCALClusterizerv1::PrintRecoInfo()
746 printf(" AliEMCALClusterizerv1::PrintRecoInfo() : version %s \n", Version() );
747 TH1F *h = (TH1F*)fHists->At(12);
749 printf(" ## Multiplicity of RecPoints ## \n");
750 for(int i=1; i<=h->GetNbinsX(); i++) {
751 int nbin = int((*h)[i]);
752 int mult = int(h->GetBinCenter(i));
753 if(nbin > 0) printf(" %i : %5.5i %6.3f %% \n", mult, nbin, 100.*nbin/h->GetEntries());
758 void AliEMCALClusterizerv1::DrawLambdasHists()
762 if(fMaxL2) fMaxL2->Draw("same");
764 fMaxDis->Draw("same");