2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
19 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
20 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
21 // of new IO (à la PHOS)
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 ---
57 #include "TBenchmark.h"
59 // --- Standard library ---
62 // --- AliRoot header files ---
63 #include "AliEMCALGetter.h"
64 #include "AliEMCALClusterizerv1.h"
65 #include "AliEMCALRecPoint.h"
66 #include "AliEMCALDigit.h"
67 #include "AliEMCALDigitizer.h"
69 #include "AliEMCALGeometry.h"
71 ClassImp(AliEMCALClusterizerv1)
73 //____________________________________________________________________________
74 AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
76 // default ctor (to be used mainly by Streamer)
79 fDefaultInit = kTRUE ;
82 //____________________________________________________________________________
83 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const TString alirunFileName, const TString eventFolderName)
84 :AliEMCALClusterizer(alirunFileName, eventFolderName)
86 // ctor with the indication of the file where header Tree and digits Tree are stored
90 fDefaultInit = kFALSE ;
94 //____________________________________________________________________________
95 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
101 //____________________________________________________________________________
102 const TString AliEMCALClusterizerv1::BranchName() const
108 //____________________________________________________________________________
109 Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp) const
111 //To be replased later by the method, reading individual parameters from the database
112 return -fADCpedestalECA + amp * fADCchannelECA ;
115 //____________________________________________________________________________
116 void AliEMCALClusterizerv1::Exec(Option_t * option)
118 // Steering method to perform clusterization for events
119 // in the range from fFirstEvent to fLastEvent.
120 // This range is optionally set by SetEventRange().
121 // if fLastEvent=-1 (by default), then process events until the end.
123 if(strstr(option,"tim"))
124 gBenchmark->Start("EMCALClusterizer");
126 if(strstr(option,"print"))
129 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
131 if (fLastEvent == -1)
132 fLastEvent = gime->MaxEvent() - 1;
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(AliEMCALRecPoint * 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();
280 gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());
282 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
284 fNTowers = geom->GetNZ() * geom->GetNPhi() ;
286 gMinuit = new TMinuit(100) ;
288 if ( !gime->Clusterizer() )
289 gime->PostClusterizer(this);
292 //____________________________________________________________________________
293 void AliEMCALClusterizerv1::InitParameters()
295 // Initializes the parameters for the Clusterizer
296 fNumberOfECAClusters = 0;
297 fECAClusteringThreshold = 1.0; // must be adjusted according to the noise level set by digitizer
298 fECALocMaxCut = 0.03 ;
302 fRecPointsInRun = 0 ;
306 //____________________________________________________________________________
307 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
309 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
311 // = 2 are not neighbour but do not continue searching
312 // neighbours are defined as digits having at least a common vertex
313 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
314 // which is compared to a digit (d2) not yet in a cluster
316 AliEMCALGeometry * geom = AliEMCALGetter::Instance()->EMCALGeometry() ;
321 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
324 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
327 Int_t rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;
328 Int_t coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;
330 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
334 if((relid2[0] > relid1[0]) && (relid2[1] > relid1[1]+1))
335 rv = 2; // Difference in row numbers is too large to look further
339 printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d ",
340 rv, d1->GetId(), relid1[0], relid1[1],
341 d2->GetId(), relid2[0], relid2[1]) ;
346 //____________________________________________________________________________
347 void AliEMCALClusterizerv1::Unload()
349 // Unloads the Digits and RecPoints
350 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
351 gime->EmcalLoader()->UnloadDigits() ;
352 gime->EmcalLoader()->UnloadRecPoints() ;
355 //____________________________________________________________________________
356 void AliEMCALClusterizerv1::WriteRecPoints()
359 // Creates new branches with given title
360 // fills and writes into TreeR.
362 AliEMCALGetter *gime = AliEMCALGetter::Instance() ;
364 TObjArray * aECARecPoints = gime->ECARecPoints() ;
366 TClonesArray * digits = gime->Digits() ;
367 TTree * treeR = gime->TreeR(); ;
371 //Evaluate position, dispersion and other RecPoint properties for EC section
372 for(index = 0; index < aECARecPoints->GetEntries(); index++)
373 (dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(index)))->EvalAll(fECAW0,digits) ;
375 aECARecPoints->Sort() ;
377 for(index = 0; index < aECARecPoints->GetEntries(); index++)
378 (dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(index)))->SetIndexInList(index) ;
380 aECARecPoints->Expand(aECARecPoints->GetEntriesFast()) ;
382 Int_t bufferSize = 32000 ;
383 Int_t splitlevel = 0 ;
386 TBranch * branchECA = treeR->Branch("EMCALECARP","TObjArray",&aECARecPoints,bufferSize,splitlevel);
387 branchECA->SetTitle(BranchName());
391 gime->WriteRecPoints("OVERWRITE");
392 gime->WriteClusterizer("OVERWRITE");
395 //____________________________________________________________________________
396 void AliEMCALClusterizerv1::MakeClusters()
398 // Steering method to construct the clusters stored in a list of Reconstructed Points
399 // A cluster is defined as a list of neighbour digits
401 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
403 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
405 TObjArray * aECARecPoints = gime->ECARecPoints() ;
407 aECARecPoints->Delete() ;
409 TClonesArray * digits = gime->Digits() ;
410 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
412 // Clusterization starts
413 TIter nextdigit(digitsC) ;
414 AliEMCALDigit * digit;
416 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
417 AliEMCALRecPoint * clu = 0 ;
419 TArrayI clusterECAdigitslist(50000);
422 printf("MakeClusters: id = %d, ene = %f , thre = %f", digit->GetId(),Calibrate(digit->GetAmp()),
423 fECAClusteringThreshold) ;
426 if ( geom->IsInECA(digit->GetId()) && (Calibrate(digit->GetAmp()) > fECAClusteringThreshold ) ){
427 Int_t iDigitInECACluster = 0;
428 // start a new Tower RecPoint
429 if(fNumberOfECAClusters >= aECARecPoints->GetSize())
430 aECARecPoints->Expand(2*fNumberOfECAClusters+1) ;
431 AliEMCALRecPoint * rp = new AliEMCALRecPoint("") ;
432 aECARecPoints->AddAt(rp, fNumberOfECAClusters) ;
433 clu = dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(fNumberOfECAClusters)) ;
434 fNumberOfECAClusters++ ;
435 clu->AddDigit(*digit, Calibrate(digit->GetAmp())) ;
436 clusterECAdigitslist[iDigitInECACluster] = digit->GetIndexInList() ;
437 iDigitInECACluster++ ;
438 digitsC->Remove(digit) ;
440 printf("MakeClusters: OK id = %d, ene = %f , thre = %f ", digit->GetId(),Calibrate(digit->GetAmp()), fECAClusteringThreshold) ;
443 AliEMCALDigit * digitN ;
446 // Find the neighbours
447 while (index < iDigitInECACluster){ // scan over digits already in cluster
448 digit = (AliEMCALDigit*)digits->At(clusterECAdigitslist[index]) ;
450 while ( (digitN = (AliEMCALDigit *)nextdigit())) { // scan over the reduced list of digits
451 // check that the digit is above the min E Cut
452 if(Calibrate(digitN->GetAmp()) < fMinECut ) digitsC->Remove(digitN);
453 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
455 case 0 : // not a neighbour
457 case 1 : // are neighbours
458 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp()) ) ;
459 clusterECAdigitslist[iDigitInECACluster] = digitN->GetIndexInList() ;
460 iDigitInECACluster++ ;
461 digitsC->Remove(digitN) ;
463 case 2 : // too far from each other
470 } // loop over ECA cluster
471 } // energy threshold
472 else if(Calibrate(digit->GetAmp()) < fMinECut ){
473 digitsC->Remove(digit);
479 //____________________________________________________________________________
480 void AliEMCALClusterizerv1::MakeUnfolding() const
482 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
486 //____________________________________________________________________________
487 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
489 // Shape of the shower (see EMCAL TDR)
490 // If you change this function, change also the gradient evaluation in ChiSquare()
492 Double_t r4 = r*r*r*r ;
493 Double_t r295 = TMath::Power(r, 2.95) ;
494 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
498 //____________________________________________________________________________
499 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * /*iniTower*/,
501 AliEMCALDigit ** /*maxAt*/,
502 Float_t * /*maxAtEnergy*/) const
504 // Performs the unfolding of a cluster with nMax overlapping showers
506 Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
510 //_____________________________________________________________________________
511 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
513 Double_t * /*x*/, Int_t /*iflag*/)
515 // Calculates the Chi square for the cluster unfolding minimization
516 // Number of parameters, Gradient, Chi squared, parameters, what to do
518 ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
520 //____________________________________________________________________________
521 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
523 // Print clusterizer parameters
525 TString message("\n") ;
527 if( strcmp(GetName(), "") !=0 ){
531 TString taskName(GetName()) ;
532 taskName.ReplaceAll(Version(), "") ;
534 printf("--------------- ");
535 printf(taskName.Data()) ;
538 printf("-----------\n");
539 printf("Clusterizing digits from the file: ");
540 printf(taskName.Data());
541 printf("\n Branch: ");
543 printf("\n ECA Local Maximum cut = %f", fECALocMaxCut);
544 printf("\n ECA Logarithmic weight = %f", fECAW0);
546 printf("\nUnfolding on\n");
548 printf("\nUnfolding off\n");
550 printf("------------------------------------------------------------------");
553 printf("AliEMCALClusterizerv1 not initialized ") ;
556 //____________________________________________________________________________
557 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
559 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
561 TObjArray * aECARecPoints = AliEMCALGetter::Instance()->ECARecPoints() ;
562 printf("PrintRecPoints: Clusterization result:") ;
564 printf("event # %d\n", gAlice->GetEvNumber() ) ;
565 printf(" Found %d ECA Rec Points\n ",
566 aECARecPoints->GetEntriesFast()) ;
568 fRecPointsInRun += aECARecPoints->GetEntriesFast() ;
570 if(strstr(option,"all")) {
572 printf("\n-----------------------------------------------------------------------\n") ;
573 printf("Clusters in ECAL section\n") ;
574 printf("Index Ene(GeV) Multi Module phi r theta X Y Z Dispersion Lambda 1 Lambda 2 # of prim Primaries list\n") ;
576 for (index = 0 ; index < aECARecPoints->GetEntries() ; index++) {
577 AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(aECARecPoints->At(index)) ;
579 rp->GetGlobalPosition(globalpos);
581 rp->GetLocalPosition(localpos);
583 rp->GetElipsAxis(lambda);
586 primaries = rp->GetPrimaries(nprimaries);
587 printf("\n%6d %8.4f %3d %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4f %4f %2d : ",
588 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
589 globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(),
590 rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ;
591 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
592 printf("%d ", primaries[iprimary] ) ;
595 printf("\n-----------------------------------------------------------------------\n");