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 //////////////////////////////////////////////////////////////////////////////
22 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
23 // unfolds the clusters having several local maxima.
24 // Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
25 // EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all
26 // parameters including input digits branch title, thresholds etc.)
27 // This TTask is normally called from Reconstructioner, but can as well be used in
30 // root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")
31 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
32 // //reads gAlice from header file "..."
33 // root [1] cl->ExecuteTask()
34 // //finds RecPoints in all events stored in galice.root
35 // root [2] cl->SetDigitsBranch("digits2")
36 // //sets another title for Digitis (input) branch
37 // root [3] cl->SetRecPointsBranch("recp2")
38 // //sets another title four output branches
39 // root [4] cl->SetTowerLocalMaxCut(0.03)
40 // //set clusterization parameters
41 // root [5] cl->ExecuteTask("deb all time")
42 // //once more finds RecPoints options are
43 // // deb - print number of found rec points
44 // // deb all - print number of found RecPoints and some their characteristics
45 // // time - print benchmarking results
47 // --- ROOT system ---
56 #include "TBenchmark.h"
58 // --- Standard library ---
61 // --- AliRoot header files ---
62 #include "AliEMCALGetter.h"
63 #include "AliEMCALClusterizerv1.h"
64 #include "AliEMCALTowerRecPoint.h"
65 #include "AliEMCALDigit.h"
66 #include "AliEMCALDigitizer.h"
68 #include "AliEMCALGeometry.h"
70 ClassImp(AliEMCALClusterizerv1)
72 //____________________________________________________________________________
73 AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
75 // default ctor (to be used mainly by Streamer)
78 fDefaultInit = kTRUE ;
81 //____________________________________________________________________________
82 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const TString alirunFileName, const TString eventFolderName)
83 :AliEMCALClusterizer(alirunFileName, eventFolderName)
85 // ctor with the indication of the file where header Tree and digits Tree are stored
89 fDefaultInit = kFALSE ;
93 //____________________________________________________________________________
94 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
100 //____________________________________________________________________________
101 const TString AliEMCALClusterizerv1::BranchName() const
107 //____________________________________________________________________________
108 Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp) const
110 //To be replased later by the method, reading individual parameters from the database
111 return -fADCpedestalECA + amp * fADCchannelECA ;
114 //____________________________________________________________________________
115 void AliEMCALClusterizerv1::Exec(Option_t * option)
117 // Steering method to perform clusterization for events
118 // in the range from fFirstEvent to fLastEvent.
119 // This range is optionally set by SetEventRange().
120 // if fLastEvent=-1 (by default), then process events until the end.
122 if(strstr(option,"tim"))
123 gBenchmark->Start("EMCALClusterizer");
125 if(strstr(option,"print"))
128 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
130 if (fLastEvent == -1)
131 fLastEvent = gime->MaxEvent() - 1 ;
133 fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent());
134 Int_t nEvents = fLastEvent - fFirstEvent + 1;
138 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
139 gime->Event(ievent,"D") ;
141 GetCalibrationParameters() ;
143 fNumberOfECAClusters = 0;
152 if(strstr(option,"deb"))
153 PrintRecPoints(option) ;
155 //increment the total number of recpoints per run
156 fRecPointsInRun += gime->ECARecPoints()->GetEntriesFast() ;
161 if(strstr(option,"tim")){
162 gBenchmark->Stop("EMCALClusterizer");
163 printf("Exec took %f seconds for Clusterizing %f seconds per event",
164 gBenchmark->GetCpuTime("EMCALClusterizer"), gBenchmark->GetCpuTime("EMCALClusterizer")/nEvents ) ;
168 //____________________________________________________________________________
169 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
170 Int_t nPar, Float_t * fitparameters) const
172 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
173 // The initial values for fitting procedure are set equal to the positions of local maxima.
174 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
176 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
177 TClonesArray * digits = gime->Digits() ;
179 gMinuit->mncler(); // Reset Minuit's list of paramters
180 gMinuit->SetPrintLevel(-1) ; // No Printout
181 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
182 // To set the address of the minimization function
183 TList * toMinuit = new TList();
184 toMinuit->AddAt(emcRP,0) ;
185 toMinuit->AddAt(digits,1) ;
187 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
189 // filling initial values for fit parameters
190 AliEMCALDigit * digit ;
194 Int_t nDigits = (Int_t) nPar / 3 ;
198 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
200 for(iDigit = 0; iDigit < nDigits; iDigit++){
201 digit = maxAt[iDigit];
206 geom->AbsToRelNumbering(digit->GetId(), relid) ;
207 geom->PosInAlice(relid, x, z) ;
209 Float_t energy = maxAtEnergy[iDigit] ;
211 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
214 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : x = %f", x ) ;
217 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
220 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : z = %f", z) ;
223 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
226 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : energy = %f", energy) ;
231 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
236 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
237 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
238 gMinuit->SetMaxIterations(5);
239 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
241 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
243 if(ierflg == 4){ // Minimum not found
244 Error("FindFit", "EMCAL Unfolding Fit not converged, cluster abandoned " ) ;
247 for(index = 0; index < nPar; index++){
250 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
251 fitparameters[index] = val ;
259 //____________________________________________________________________________
260 void AliEMCALClusterizerv1::GetCalibrationParameters()
262 // Gets the parameters for the calibration from the digitizer
263 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
265 if ( !gime->Digitizer() )
266 gime->LoadDigitizer();
267 AliEMCALDigitizer * dig = gime->Digitizer();
269 fADCchannelECA = dig->GetECAchannel() ;
270 fADCpedestalECA = dig->GetECApedestal();
273 //____________________________________________________________________________
274 void AliEMCALClusterizerv1::Init()
276 // Make all memory allocations which can not be done in default constructor.
277 // Attach the Clusterizer task to the list of EMCAL tasks
279 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());
281 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
283 fNTowers = geom->GetNZ() * geom->GetNPhi() ;
285 gMinuit = new TMinuit(100) ;
287 if ( !gime->Clusterizer() )
288 gime->PostClusterizer(this);
291 //____________________________________________________________________________
292 void AliEMCALClusterizerv1::InitParameters()
294 // Initializes the parameters for the Clusterizer
295 fNumberOfECAClusters = 0;
296 fECAClusteringThreshold = 0.0045; // must be adjusted according to the noise leve set by digitizer
297 fECALocMaxCut = 0.03 ;
301 fRecPointsInRun = 0 ;
304 //____________________________________________________________________________
305 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
307 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
309 // = 2 are not neighbour but do not continue searching
310 // neighbours are defined as digits having at least a common vertex
311 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
312 // which is compared to a digit (d2) not yet in a cluster
314 AliEMCALGeometry * geom = AliEMCALGetter::Instance()->EMCALGeometry() ;
319 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
322 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
324 if ( (relid1[0] == relid2[0])){ // inside the same EMCAL Arm
325 Int_t rowdiff = TMath::Abs( relid1[1] - relid2[1] ) ;
326 Int_t coldiff = TMath::Abs( relid1[2] - relid2[2] ) ;
328 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
329 if(TMath::Abs(d1->GetTime() - d2->GetTime() ) < fTimeGate)
333 if((relid2[1] > relid1[1]) && (relid2[2] > relid1[2]+1))
334 rv = 2; // Difference in row numbers is too large to look further
340 if(relid1[0] < relid2[0])
345 printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d,%d \n id2=%d, relid2=%d,%d,%d ",
346 rv, d1->GetId(), relid1[0], relid1[1], relid1[2],
347 d2->GetId(), relid2[0], relid2[1], relid2[2]) ;
352 //____________________________________________________________________________
353 void AliEMCALClusterizerv1::Unload()
355 // Unloads the Digits and RecPoints
356 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
357 gime->EmcalLoader()->UnloadDigits() ;
358 gime->EmcalLoader()->UnloadRecPoints() ;
361 //____________________________________________________________________________
362 void AliEMCALClusterizerv1::WriteRecPoints()
365 // Creates new branches with given title
366 // fills and writes into TreeR.
368 AliEMCALGetter *gime = AliEMCALGetter::Instance() ;
370 TObjArray * aECARecPoints = gime->ECARecPoints() ;
372 TClonesArray * digits = gime->Digits() ;
373 TTree * treeR = gime->TreeR(); ;
377 //Evaluate position, dispersion and other RecPoint properties for EC section
378 for(index = 0; index < aECARecPoints->GetEntries(); index++)
379 (dynamic_cast<AliEMCALTowerRecPoint *>(aECARecPoints->At(index)))->EvalAll(fECAW0,digits) ;
381 aECARecPoints->Sort() ;
383 for(index = 0; index < aECARecPoints->GetEntries(); index++)
384 (dynamic_cast<AliEMCALTowerRecPoint *>(aECARecPoints->At(index)))->SetIndexInList(index) ;
386 aECARecPoints->Expand(aECARecPoints->GetEntriesFast()) ;
388 Int_t bufferSize = 32000 ;
389 Int_t splitlevel = 0 ;
392 TBranch * branchECA = treeR->Branch("EMCALECARP","TObjArray",&aECARecPoints,bufferSize,splitlevel);
393 branchECA->SetTitle(BranchName());
397 gime->WriteRecPoints("OVERWRITE");
398 gime->WriteClusterizer("OVERWRITE");
401 //____________________________________________________________________________
402 void AliEMCALClusterizerv1::MakeClusters()
404 // Steering method to construct the clusters stored in a list of Reconstructed Points
405 // A cluster is defined as a list of neighbour digits
407 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
409 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
411 TObjArray * aECARecPoints = gime->ECARecPoints() ;
413 aECARecPoints->Delete() ;
415 TClonesArray * digits = gime->Digits() ;
418 AliEMCALDigit * digit ;
421 // count the number of digits in ECA section
422 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) { // scan over the list of digits
423 if (geom->IsInECA(digit->GetId()))
426 Error("MakeClusters", "id = %d is a wrong ID!", digit->GetId()) ;
430 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
431 // Clusterization starts
433 TIter nextdigit(digitsC) ;
435 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
436 AliEMCALRecPoint * clu = 0 ;
438 TArrayI clusterECAdigitslist(50);
440 Bool_t inECA = kFALSE;
441 if( geom->IsInECA(digit->GetId()) ) {
446 printf("MakeClusters: id = %d, ene = %f , thre = %f",
447 digit->GetId(),Calibrate(digit->GetAmp()), fECAClusteringThreshold) ;
449 if (inECA && (Calibrate(digit->GetAmp()) > fECAClusteringThreshold ) ){
451 Int_t iDigitInECACluster = 0;
453 if( geom->IsInECA(digit->GetId()) ) {
454 // start a new Tower RecPoint
455 if(fNumberOfECAClusters >= aECARecPoints->GetSize())
456 aECARecPoints->Expand(2*fNumberOfECAClusters+1) ;
457 AliEMCALTowerRecPoint * rp = new AliEMCALTowerRecPoint("") ;
459 aECARecPoints->AddAt(rp, fNumberOfECAClusters) ;
460 clu = dynamic_cast<AliEMCALTowerRecPoint *>(aECARecPoints->At(fNumberOfECAClusters)) ;
461 fNumberOfECAClusters++ ;
462 clu->AddDigit(*digit, Calibrate(digit->GetAmp())) ;
463 clusterECAdigitslist[iDigitInECACluster] = digit->GetIndexInList() ;
464 iDigitInECACluster++ ;
465 digitsC->Remove(digit) ;
467 printf("MakeClusters: OK id = %d, ene = %f , thre = %f ", digit->GetId(),Calibrate(digit->GetAmp()), fECAClusteringThreshold) ;
472 AliEMCALDigit * digitN ;
477 while (index < iDigitInECACluster){ // scan over digits already in cluster
478 digit = (AliEMCALDigit*)digits->At(clusterECAdigitslist[index]) ;
480 while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits
481 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
483 case 0 : // not a neighbour
485 case 1 : // are neighbours
486 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp()) ) ;
487 clusterECAdigitslist[iDigitInECACluster] = digitN->GetIndexInList() ;
488 iDigitInECACluster++ ;
489 digitsC->Remove(digitN) ;
491 case 2 : // too far from each other
499 } // loop over ECA cluster
505 //____________________________________________________________________________
506 void AliEMCALClusterizerv1::MakeUnfolding() const
508 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
512 //____________________________________________________________________________
513 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
515 // Shape of the shower (see EMCAL TDR)
516 // If you change this function, change also the gradient evaluation in ChiSquare()
518 Double_t r4 = r*r*r*r ;
519 Double_t r295 = TMath::Power(r, 2.95) ;
520 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
524 //____________________________________________________________________________
525 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * /*iniTower*/,
527 AliEMCALDigit ** /*maxAt*/,
528 Float_t * /*maxAtEnergy*/) const
530 // Performs the unfolding of a cluster with nMax overlapping showers
532 Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
536 //_____________________________________________________________________________
537 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
539 Double_t * /*x*/, Int_t /*iflag*/)
541 // Calculates the Chi square for the cluster unfolding minimization
542 // Number of parameters, Gradient, Chi squared, parameters, what to do
544 ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
546 //____________________________________________________________________________
547 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
549 // Print clusterizer parameters
551 TString message("\n") ;
553 if( strcmp(GetName(), "") !=0 ){
557 TString taskName(GetName()) ;
558 taskName.ReplaceAll(Version(), "") ;
560 printf("--------------- ");
561 printf(taskName.Data()) ;
564 printf("-----------\n");
565 printf("Clusterizing digits from the file: ");
566 printf(taskName.Data());
567 printf("\n Branch: ");
569 printf("\n ECA Local Maximum cut = %f", fECALocMaxCut);
570 printf("\n ECA Logarothmic weight = %f", fECAW0);
572 printf("\nUnfolding on\n");
574 printf("\nUnfolding off\n");
576 printf("------------------------------------------------------------------");
579 printf("AliEMCALClusterizerv1 not initialized ") ;
582 //____________________________________________________________________________
583 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
585 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
586 TObjArray * aECARecPoints = AliEMCALGetter::Instance()->ECARecPoints() ;
587 printf("PrintRecPoints: Clusterization result:") ;
589 printf("event # %d\n", gAlice->GetEvNumber() ) ;
590 printf(" Found %d ECA Rec Points\n ",
591 aECARecPoints->GetEntriesFast()) ;
593 fRecPointsInRun += aECARecPoints->GetEntriesFast() ;
595 if(strstr(option,"all")) {
597 printf("\n-----------------------------------------------------------------------\n") ;
598 printf("Clusters in ECAL section\n") ;
599 printf("Index Ene(GeV) Multi Module phi r theta X Y Z Dispersion Lambda 1 Lambda 2 # of prim Primaries list\n") ;
601 for (index = 0 ; index < aECARecPoints->GetEntries() ; index++) {
602 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint * >(aECARecPoints->At(index)) ;
604 rp->GetGlobalPosition(globalpos);
606 rp->GetLocalPosition(localpos);
608 rp->GetElipsAxis(lambda);
611 primaries = rp->GetPrimaries(nprimaries);
612 printf("\n%6d %8.4f %3d %2d %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4f %4f %2d : ",
613 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(), rp->GetEMCALArm(),
614 globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(),
615 rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ;
616 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
617 printf("%d ", primaries[iprimary] ) ;
620 printf("\n-----------------------------------------------------------------------\n");