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(GetTitle()) ;
130 if (fLastEvent == -1)
131 fLastEvent = gime->MaxEvent() - 1 ;
133 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent());
134 Int_t nEvents = fLastEvent - fFirstEvent + 1;
137 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
138 gime->Event(ievent,"D") ;
140 GetCalibrationParameters() ;
142 fNumberOfECAClusters = 0;
151 if(strstr(option,"deb"))
152 PrintRecPoints(option) ;
154 //increment the total number of recpoints per run
155 fRecPointsInRun += gime->ECARecPoints()->GetEntriesFast() ;
160 if(strstr(option,"tim")){
161 gBenchmark->Stop("EMCALClusterizer");
162 printf("Exec took %f seconds for Clusterizing %f seconds per event",
163 gBenchmark->GetCpuTime("EMCALClusterizer"), gBenchmark->GetCpuTime("EMCALClusterizer")/nEvents ) ;
167 //____________________________________________________________________________
168 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
169 Int_t nPar, Float_t * fitparameters) const
171 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
172 // The initial values for fitting procedure are set equal to the positions of local maxima.
173 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
175 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
176 TClonesArray * digits = gime->Digits() ;
178 gMinuit->mncler(); // Reset Minuit's list of paramters
179 gMinuit->SetPrintLevel(-1) ; // No Printout
180 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
181 // To set the address of the minimization function
182 TList * toMinuit = new TList();
183 toMinuit->AddAt(emcRP,0) ;
184 toMinuit->AddAt(digits,1) ;
186 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
188 // filling initial values for fit parameters
189 AliEMCALDigit * digit ;
193 Int_t nDigits = (Int_t) nPar / 3 ;
197 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
199 for(iDigit = 0; iDigit < nDigits; iDigit++){
200 digit = maxAt[iDigit];
205 geom->AbsToRelNumbering(digit->GetId(), relid) ;
206 geom->PosInAlice(relid, x, z) ;
208 Float_t energy = maxAtEnergy[iDigit] ;
210 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
213 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : x = %f", x ) ;
216 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
219 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : z = %f", z) ;
222 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
225 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : energy = %f", energy) ;
230 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
235 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
236 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
237 gMinuit->SetMaxIterations(5);
238 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
240 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
242 if(ierflg == 4){ // Minimum not found
243 Error("FindFit", "EMCAL Unfolding Fit not converged, cluster abandoned " ) ;
246 for(index = 0; index < nPar; index++){
249 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
250 fitparameters[index] = val ;
258 //____________________________________________________________________________
259 void AliEMCALClusterizerv1::GetCalibrationParameters()
261 // Gets the parameters for the calibration from the digitizer
262 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
264 if ( !gime->Digitizer() )
265 gime->LoadDigitizer();
266 AliEMCALDigitizer * dig = gime->Digitizer();
268 fADCchannelECA = dig->GetECAchannel() ;
269 fADCpedestalECA = dig->GetECApedestal();
272 //____________________________________________________________________________
273 void AliEMCALClusterizerv1::Init()
275 // Make all memory allocations which can not be done in default constructor.
276 // Attach the Clusterizer task to the list of EMCAL tasks
278 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());
280 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
282 fNTowers = geom->GetNZ() * geom->GetNPhi() ;
284 gMinuit = new TMinuit(100) ;
286 if ( !gime->Clusterizer() )
287 gime->PostClusterizer(this);
290 //____________________________________________________________________________
291 void AliEMCALClusterizerv1::InitParameters()
293 // Initializes the parameters for the Clusterizer
294 fNumberOfECAClusters = 0;
295 fECAClusteringThreshold = 0.0045; // must be adjusted according to the noise leve set by digitizer
296 fECALocMaxCut = 0.03 ;
300 fRecPointsInRun = 0 ;
303 //____________________________________________________________________________
304 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
306 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
308 // = 2 are not neighbour but do not continue searching
309 // neighbours are defined as digits having at least a common vertex
310 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
311 // which is compared to a digit (d2) not yet in a cluster
313 AliEMCALGeometry * geom = AliEMCALGetter::Instance()->EMCALGeometry() ;
318 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
321 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
323 if ( (relid1[0] == relid2[0])){ // inside the same EMCAL Arm
324 Int_t rowdiff = TMath::Abs( relid1[1] - relid2[1] ) ;
325 Int_t coldiff = TMath::Abs( relid1[2] - relid2[2] ) ;
327 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
328 if(TMath::Abs(d1->GetTime() - d2->GetTime() ) < fTimeGate)
332 if((relid2[1] > relid1[1]) && (relid2[2] > relid1[2]+1))
333 rv = 2; // Difference in row numbers is too large to look further
339 if(relid1[0] < relid2[0])
344 printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d,%d \n id2=%d, relid2=%d,%d,%d ",
345 rv, d1->GetId(), relid1[0], relid1[1], relid1[2],
346 d2->GetId(), relid2[0], relid2[1], relid2[2]) ;
351 //____________________________________________________________________________
352 void AliEMCALClusterizerv1::Unload()
354 // Unloads the Digits and RecPoints
355 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
356 gime->EmcalLoader()->UnloadDigits() ;
357 gime->EmcalLoader()->UnloadRecPoints() ;
360 //____________________________________________________________________________
361 void AliEMCALClusterizerv1::WriteRecPoints()
364 // Creates new branches with given title
365 // fills and writes into TreeR.
367 AliEMCALGetter *gime = AliEMCALGetter::Instance() ;
369 TObjArray * aECARecPoints = gime->ECARecPoints() ;
371 TClonesArray * digits = gime->Digits() ;
372 TTree * treeR = gime->TreeR(); ;
376 //Evaluate position, dispersion and other RecPoint properties for EC section
377 for(index = 0; index < aECARecPoints->GetEntries(); index++)
378 (dynamic_cast<AliEMCALTowerRecPoint *>(aECARecPoints->At(index)))->EvalAll(fECAW0,digits) ;
380 aECARecPoints->Sort() ;
382 for(index = 0; index < aECARecPoints->GetEntries(); index++)
383 (dynamic_cast<AliEMCALTowerRecPoint *>(aECARecPoints->At(index)))->SetIndexInList(index) ;
385 aECARecPoints->Expand(aECARecPoints->GetEntriesFast()) ;
387 Int_t bufferSize = 32000 ;
388 Int_t splitlevel = 0 ;
391 TBranch * branchECA = treeR->Branch("EMCALECARP","TObjArray",&aECARecPoints,bufferSize,splitlevel);
392 branchECA->SetTitle(BranchName());
396 gime->WriteRecPoints("OVERWRITE");
397 gime->WriteClusterizer("OVERWRITE");
400 //____________________________________________________________________________
401 void AliEMCALClusterizerv1::MakeClusters()
403 // Steering method to construct the clusters stored in a list of Reconstructed Points
404 // A cluster is defined as a list of neighbour digits
406 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
408 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
410 TObjArray * aECARecPoints = gime->ECARecPoints() ;
412 aECARecPoints->Delete() ;
414 TClonesArray * digits = gime->Digits() ;
417 AliEMCALDigit * digit ;
420 // count the number of digits in ECA section
421 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) { // scan over the list of digits
422 if (geom->IsInECA(digit->GetId()))
425 Error("MakeClusters", "id = %d is a wrong ID!", digit->GetId()) ;
429 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
430 // Clusterization starts
432 TIter nextdigit(digitsC) ;
434 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
435 AliEMCALRecPoint * clu = 0 ;
437 TArrayI clusterECAdigitslist(50);
439 Bool_t inECA = kFALSE;
440 if( geom->IsInECA(digit->GetId()) ) {
445 printf("MakeClusters: id = %d, ene = %f , thre = %f",
446 digit->GetId(),Calibrate(digit->GetAmp()), fECAClusteringThreshold) ;
448 if (inECA && (Calibrate(digit->GetAmp()) > fECAClusteringThreshold ) ){
450 Int_t iDigitInECACluster = 0;
452 if( geom->IsInECA(digit->GetId()) ) {
453 // start a new Tower RecPoint
454 if(fNumberOfECAClusters >= aECARecPoints->GetSize())
455 aECARecPoints->Expand(2*fNumberOfECAClusters+1) ;
456 AliEMCALTowerRecPoint * rp = new AliEMCALTowerRecPoint("") ;
458 aECARecPoints->AddAt(rp, fNumberOfECAClusters) ;
459 clu = dynamic_cast<AliEMCALTowerRecPoint *>(aECARecPoints->At(fNumberOfECAClusters)) ;
460 fNumberOfECAClusters++ ;
461 clu->AddDigit(*digit, Calibrate(digit->GetAmp())) ;
462 clusterECAdigitslist[iDigitInECACluster] = digit->GetIndexInList() ;
463 iDigitInECACluster++ ;
464 digitsC->Remove(digit) ;
466 printf("MakeClusters: OK id = %d, ene = %f , thre = %f ", digit->GetId(),Calibrate(digit->GetAmp()), fECAClusteringThreshold) ;
471 AliEMCALDigit * digitN ;
476 while (index < iDigitInECACluster){ // scan over digits already in cluster
477 digit = (AliEMCALDigit*)digits->At(clusterECAdigitslist[index]) ;
479 while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits
480 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
482 case 0 : // not a neighbour
484 case 1 : // are neighbours
485 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp()) ) ;
486 clusterECAdigitslist[iDigitInECACluster] = digitN->GetIndexInList() ;
487 iDigitInECACluster++ ;
488 digitsC->Remove(digitN) ;
490 case 2 : // too far from each other
498 } // loop over ECA cluster
504 //____________________________________________________________________________
505 void AliEMCALClusterizerv1::MakeUnfolding() const
507 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
511 //____________________________________________________________________________
512 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
514 // Shape of the shower (see EMCAL TDR)
515 // If you change this function, change also the gradient evaluation in ChiSquare()
517 Double_t r4 = r*r*r*r ;
518 Double_t r295 = TMath::Power(r, 2.95) ;
519 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
523 //____________________________________________________________________________
524 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * /*iniTower*/,
526 AliEMCALDigit ** /*maxAt*/,
527 Float_t * /*maxAtEnergy*/) const
529 // Performs the unfolding of a cluster with nMax overlapping showers
531 Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
535 //_____________________________________________________________________________
536 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
538 Double_t * /*x*/, Int_t /*iflag*/)
540 // Calculates the Chi square for the cluster unfolding minimization
541 // Number of parameters, Gradient, Chi squared, parameters, what to do
543 ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
545 //____________________________________________________________________________
546 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
548 // Print clusterizer parameters
550 TString message("\n") ;
552 if( strcmp(GetName(), "") !=0 ){
556 TString taskName(GetName()) ;
557 taskName.ReplaceAll(Version(), "") ;
559 printf("--------------- ");
560 printf(taskName.Data()) ;
563 printf("-----------\n");
564 printf("Clusterizing digits from the file: ");
565 printf(taskName.Data());
566 printf("\n Branch: ");
568 printf("\n ECA Local Maximum cut = %f", fECALocMaxCut);
569 printf("\n ECA Logarithmic weight = %f", fECAW0);
571 printf("\nUnfolding on\n");
573 printf("\nUnfolding off\n");
575 printf("------------------------------------------------------------------");
578 printf("AliEMCALClusterizerv1 not initialized ") ;
581 //____________________________________________________________________________
582 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
584 // 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");