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 1 October 2000. Yuri Kharlov:
20 PPSD upper layer is considered if number of layers>1
22 18 October 2000. Yuri Kharlov:
23 AliEMCALClusterizerv1()
24 CPV clusterizing parameters added
27 After first PPSD digit remove EMC digits only once
29 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
30 //////////////////////////////////////////////////////////////////////////////
31 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
32 // unfolds the clusters having several local maxima.
33 // Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
34 // EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all
35 // parameters including input digits branch title, thresholds etc.)
36 // This TTask is normally called from Reconstructioner, but can as well be used in
39 // root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")
40 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
41 // //reads gAlice from header file "..."
42 // root [1] cl->ExecuteTask()
43 // //finds RecPoints in all events stored in galice.root
44 // root [2] cl->SetDigitsBranch("digits2")
45 // //sets another title for Digitis (input) branch
46 // root [3] cl->SetRecPointsBranch("recp2")
47 // //sets another title four output branches
48 // root [4] cl->SetTowerLocalMaxCut(0.03)
49 // //set clusterization parameters
50 // root [5] cl->ExecuteTask("deb all time")
51 // //once more finds RecPoints options are
52 // // deb - print number of found rec points
53 // // deb all - print number of found RecPoints and some their characteristics
54 // // time - print benchmarking results
56 // --- ROOT system ---
65 #include "TBenchmark.h"
67 // --- Standard library ---
72 // --- AliRoot header files ---
74 #include "AliEMCALClusterizerv1.h"
75 #include "AliEMCALDigit.h"
76 #include "AliEMCALDigitizer.h"
77 #include "AliEMCALTowerRecPoint.h"
79 #include "AliEMCALGetter.h"
82 ClassImp(AliEMCALClusterizerv1)
84 //____________________________________________________________________________
85 AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
87 // default ctor (to be used mainly by Streamer)
89 fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ;
91 fPreShoClusteringThreshold = 0.0;
92 fTowerClusteringThreshold = 0.0;
94 fTowerLocMaxCut = 0.0 ;
95 fPreShoLocMaxCut = 0.0 ;
104 fHeaderFileName = "" ;
105 fRecPointsInRun = 0 ;
108 //____________________________________________________________________________
109 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const char* headerFile,const char* name)
110 :AliEMCALClusterizer(headerFile, name)
112 // ctor with the indication of the file where header Tree and digits Tree are stored
115 fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ;
119 fPreShoClusteringThreshold = 0.0001;
120 fTowerClusteringThreshold = 0.2;
122 fTowerLocMaxCut = 0.03 ;
123 fPreShoLocMaxCut = 0.03 ;
132 fHeaderFileName = GetTitle() ;
133 fDigitsBranchTitle = GetName() ;
135 TString clusterizerName( GetName()) ;
136 clusterizerName.Append(":") ;
137 clusterizerName.Append(Version()) ;
138 SetName(clusterizerName) ;
139 fRecPointsInRun = 0 ;
144 //____________________________________________________________________________
145 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
148 //____________________________________________________________________________
149 Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp, Bool_t inpresho) const
151 if ( inpresho ) // calibrate as pre shower
152 return -fADCpedestalPreSho + amp * fADCchannelPreSho ;
154 else //calibrate as tower
155 return -fADCpedestalTower + amp * fADCchannelTower ;
157 //____________________________________________________________________________
158 void AliEMCALClusterizerv1::Exec(Option_t * option)
162 if( strcmp(GetName(), "")== 0 )
165 if(strstr(option,"tim"))
166 gBenchmark->Start("EMCALClusterizer");
168 if(strstr(option,"print"))
171 gAlice->GetEvent(0) ;
173 //check, if the branch with name of this" already exits?
174 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
176 TBranch * branch = 0 ;
177 Bool_t emcaltowerfound = kFALSE, emcalpreshofound = kFALSE, clusterizerfound = kFALSE ;
179 TString branchname = GetName() ;
180 branchname.Remove(branchname.Index(Version())-1) ;
182 while ( (branch = (TBranch*)next()) && (!emcaltowerfound || !emcalpreshofound || !clusterizerfound) ) {
183 if ( (strcmp(branch->GetName(), "EMCALTowerRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) )
184 emcaltowerfound = kTRUE ;
186 else if ( (strcmp(branch->GetName(), "EMCALPreShoRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) )
187 emcalpreshofound = kTRUE ;
189 else if ((strcmp(branch->GetName(), "AliEMCALClusterizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) )
190 clusterizerfound = kTRUE ;
193 if ( emcalpreshofound || emcaltowerfound || clusterizerfound ) {
194 cerr << "WARNING: AliEMCALClusterizer::Exec -> Tower(PreSho)RecPoints and/or Clusterizer branch with name "
195 << branchname.Data() << " already exits" << endl ;
199 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
200 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
203 for(ievent = 0; ievent < nevents; ievent++){
206 GetCalibrationParameters() ;
208 fNumberOfTowerClusters = fNumberOfPreShoClusters = 0 ;
210 gime->Event(ievent,"D") ;
217 WriteRecPoints(ievent) ;
219 if(strstr(option,"deb"))
220 PrintRecPoints(option) ;
222 //increment the total number of digits per run
223 fRecPointsInRun += gime->TowerRecPoints()->GetEntriesFast() ;
224 fRecPointsInRun += gime->PreShowerRecPoints()->GetEntriesFast() ;
227 if(strstr(option,"tim")){
228 gBenchmark->Stop("EMCALClusterizer");
229 cout << "AliEMCALClusterizer:" << endl ;
230 cout << " took " << gBenchmark->GetCpuTime("EMCALClusterizer") << " seconds for Clusterizing "
231 << gBenchmark->GetCpuTime("EMCALClusterizer")/nevents << " seconds per event " << endl ;
237 //____________________________________________________________________________
238 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
239 Int_t nPar, Float_t * fitparameters) const
241 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
242 // The initial values for fitting procedure are set equal to the positions of local maxima.
243 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
245 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
246 TClonesArray * digits = gime->Digits() ;
249 gMinuit->mncler(); // Reset Minuit's list of paramters
250 gMinuit->SetPrintLevel(-1) ; // No Printout
251 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
252 // To set the address of the minimization function
254 TList * toMinuit = new TList();
255 toMinuit->AddAt(emcRP,0) ;
256 toMinuit->AddAt(digits,1) ;
258 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
260 // filling initial values for fit parameters
261 AliEMCALDigit * digit ;
265 Int_t nDigits = (Int_t) nPar / 3 ;
269 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
271 for(iDigit = 0; iDigit < nDigits; iDigit++){
272 digit = maxAt[iDigit];
277 geom->AbsToRelNumbering(digit->GetId(), relid) ;
278 geom->PosInAlice(relid, x, z) ;
280 Float_t energy = maxAtEnergy[iDigit] ;
282 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
285 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : x = " << x << endl ;
288 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
291 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : z = " << z << endl ;
294 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
297 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : energy = " << energy << endl ;
302 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
307 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
308 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
309 gMinuit->SetMaxIterations(5);
310 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
312 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
314 if(ierflg == 4){ // Minimum not found
315 cout << "EMCAL Unfolding> Fit not converged, cluster abandoned "<< endl ;
318 for(index = 0; index < nPar; index++){
321 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
322 fitparameters[index] = val ;
330 //____________________________________________________________________________
331 void AliEMCALClusterizerv1::GetCalibrationParameters()
333 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
334 TString branchname = GetName() ;
335 branchname.Remove(branchname.Index(Version())-1) ;
337 AliEMCALDigitizer * dig = gime->Digitizer(branchname) ;
339 fADCchannelTower = dig->GetTowerchannel() ;
340 fADCpedestalTower = dig->GetTowerpedestal();
342 fADCchannelPreSho = dig->GetPreShochannel() ;
343 fADCpedestalPreSho = dig->GetPreShopedestal() ;
346 //____________________________________________________________________________
347 void AliEMCALClusterizerv1::Init()
349 // Make all memory allocations which can not be done in default constructor.
350 // Attach the Clusterizer task to the list of EMCAL tasks
352 if ( strcmp(GetTitle(), "") == 0 )
353 SetTitle("galice.root") ;
355 TString branchname = GetName() ;
356 branchname.Remove(branchname.Index(Version())-1) ;
358 AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), branchname, "update") ;
360 cerr << "ERROR: AliEMCALClusterizerv1::Init -> Could not obtain the Getter object !" << endl ;
364 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
365 fNTowers = geom->GetNZ() * geom->GetNPhi() ;
368 gMinuit = new TMinuit(100) ;
370 gime->PostClusterizer(this) ;
371 // create a folder on the white board
372 gime->PostRecPoints(branchname ) ;
374 gime->PostDigits(branchname) ;
375 gime->PostDigitizer(branchname) ;
379 //____________________________________________________________________________
380 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
382 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
384 // = 2 are not neighbour but do not continue searching
385 // neighbours are defined as digits having at least a common vertex
386 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
387 // which is compared to a digit (d2) not yet in a cluster
389 AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
394 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
397 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
399 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same EMCAL Arm
400 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
401 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
403 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
404 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fTimeGate))
408 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
409 rv = 2; // Difference in row numbers is too large to look further
415 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
424 //____________________________________________________________________________
425 Bool_t AliEMCALClusterizerv1::IsInTower(AliEMCALDigit * digit) const
427 // Tells if (true) or not (false) the digit is in a EMCAL-Tower
430 if (!digit->IsInPreShower())
435 //____________________________________________________________________________
436 Bool_t AliEMCALClusterizerv1::IsInPreShower(AliEMCALDigit * digit) const
438 // Tells if (true) or not (false) the digit is in a EMCAL-PreShower
441 if (digit->IsInPreShower())
446 //____________________________________________________________________________
447 void AliEMCALClusterizerv1::WriteRecPoints(Int_t event)
450 // Creates new branches with given title
451 // fills and writes into TreeR.
453 TString branchName(GetName() ) ;
454 branchName.Remove(branchName.Index(Version())-1) ;
456 AliEMCALGetter *gime = AliEMCALGetter::GetInstance() ;
457 TObjArray * towerRecPoints = gime->TowerRecPoints(branchName) ;
458 TObjArray * preshoRecPoints = gime->PreShowerRecPoints(branchName) ;
459 TClonesArray * digits = gime->Digits(branchName) ;
462 //Evaluate position, dispersion and other RecPoint properties...
463 for(index = 0; index < towerRecPoints->GetEntries(); index++)
464 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->EvalAll(fW0,digits) ;
466 towerRecPoints->Sort() ;
468 for(index = 0; index < towerRecPoints->GetEntries(); index++)
469 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->SetIndexInList(index) ;
471 towerRecPoints->Expand(towerRecPoints->GetEntriesFast()) ;
473 //Now the same for CPV
474 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
475 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->EvalAll(fW0CPV,digits) ;
477 preshoRecPoints->Sort() ;
479 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
480 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->SetIndexInList(index) ;
482 preshoRecPoints->Expand(preshoRecPoints->GetEntriesFast()) ;
484 //Make branches in TreeR for RecPoints and Clusterizer
486 if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){ //generating file name
487 filename = new char[strlen(gAlice->GetBaseFile())+20] ;
488 sprintf(filename,"%s/EMCAL.Reco.root",gAlice->GetBaseFile()) ;
492 TDirectory *cwd = gDirectory;
495 Int_t bufferSize = 32000 ;
496 Int_t splitlevel = 0 ;
499 TBranch * emcBranch = gAlice->TreeR()->Branch("EMCALTowerRP","TObjArray",&towerRecPoints,bufferSize,splitlevel);
500 emcBranch->SetTitle(branchName);
502 emcBranch->SetFile(filename);
503 TIter next( emcBranch->GetListOfBranches());
505 while ((sb=(TBranch*)next())) {
506 sb->SetFile(filename);
513 TBranch * cpvBranch = gAlice->TreeR()->Branch("EMCALPreShoRP","TObjArray",&preshoRecPoints,bufferSize,splitlevel);
514 cpvBranch->SetTitle(branchName);
516 cpvBranch->SetFile(filename);
517 TIter next( cpvBranch->GetListOfBranches());
519 while ((sb=(TBranch*)next())) {
520 sb->SetFile(filename);
525 //And Finally clusterizer branch
526 AliEMCALClusterizerv1 * cl = (AliEMCALClusterizerv1*)gime->Clusterizer(branchName) ;
527 TBranch * clusterizerBranch = gAlice->TreeR()->Branch("AliEMCALClusterizer","AliEMCALClusterizerv1",
528 &cl,bufferSize,splitlevel);
529 clusterizerBranch->SetTitle(branchName);
531 clusterizerBranch->SetFile(filename);
532 TIter next( clusterizerBranch->GetListOfBranches());
534 while ((sb=(TBranch*)next())) {
535 sb->SetFile(filename);
541 clusterizerBranch->Fill() ;
543 gAlice->TreeR()->Write(0,kOverwrite) ;
547 //____________________________________________________________________________
548 void AliEMCALClusterizerv1::MakeClusters()
550 // Steering method to construct the clusters stored in a list of Reconstructed Points
551 // A cluster is defined as a list of neighbour digits
553 TString branchName(GetName()) ;
554 branchName.Remove(branchName.Index(Version())-1) ;
556 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
558 TObjArray * towerRecPoints = gime->TowerRecPoints(branchName) ;
559 TObjArray * preshoRecPoints = gime->PreShowerRecPoints(branchName) ;
560 towerRecPoints->Delete() ;
561 preshoRecPoints->Delete() ;
563 TClonesArray * digits = gime->Digits(branchName) ;
564 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
567 // Clusterization starts
569 TIter nextdigit(digitsC) ;
570 AliEMCALDigit * digit ;
571 Bool_t notremoved = kTRUE ;
573 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
574 AliEMCALRecPoint * clu = 0 ;
576 TArrayI clusterdigitslist(1500) ;
579 if (( IsInTower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fTowerClusteringThreshold ) ||
580 ( IsInPreShower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fPreShoClusteringThreshold ) ) {
582 Int_t iDigitInCluster = 0 ;
584 if ( IsInTower(digit) ) {
585 // start a new Tower RecPoint
586 if(fNumberOfTowerClusters >= towerRecPoints->GetSize())
587 towerRecPoints->Expand(2*fNumberOfTowerClusters+1) ;
589 towerRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfTowerClusters) ;
590 clu = dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(fNumberOfTowerClusters)) ;
591 fNumberOfTowerClusters++ ;
592 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower())) ;
593 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
595 digitsC->Remove(digit) ;
599 // start a new Pre Shower cluster
600 if(fNumberOfPreShoClusters >= preshoRecPoints->GetSize())
601 preshoRecPoints->Expand(2*fNumberOfPreShoClusters+1);
603 preshoRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfPreShoClusters) ;
605 clu = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(fNumberOfPreShoClusters)) ;
606 fNumberOfPreShoClusters++ ;
607 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower() ) );
608 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
610 digitsC->Remove(digit) ;
613 // Here we remove remaining Tower digits, which cannot make a cluster
616 while( ( digit = dynamic_cast<AliEMCALDigit *>(nextdigit()) ) ) {
617 if( IsInTower(digit) )
618 digitsC->Remove(digit) ;
622 notremoved = kFALSE ;
629 AliEMCALDigit * digitN ;
631 while (index < iDigitInCluster){ // scan over digits already in cluster
632 digit = (AliEMCALDigit*)digits->At(clusterdigitslist[index]) ;
634 while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits
635 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
637 case 0 : // not a neighbour
639 case 1 : // are neighbours
640 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->IsInPreShower() ) ) ;
641 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
643 digitsC->Remove(digitN) ;
645 case 2 : // too far from each other
654 } // loop over cluster
665 //____________________________________________________________________________
666 void AliEMCALClusterizerv1::MakeUnfolding()
668 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
670 // // Unfolds clusters using the shape of an ElectroMagnetic shower
671 // // Performs unfolding of all EMC/CPV clusters
673 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
675 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
676 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
677 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
678 // TClonesArray * digits = gime->Digits() ;
680 // // Unfold first EMC clusters
681 // if(fNumberOfTowerClusters > 0){
683 // Int_t nModulesToUnfold = geom->GetNModules() ;
685 // Int_t numberofNotUnfolded = fNumberOfTowerClusters ;
687 // for(index = 0 ; index < numberofNotUnfolded ; index++){
689 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint *) emcRecPoints->At(index) ;
690 // if(emcRecPoint->GetEMCALMod()> nModulesToUnfold)
693 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
694 // Int_t * maxAt = new Int_t[nMultipl] ;
695 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
696 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fTowerLocMaxCut,digits) ;
698 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
699 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
700 // emcRecPoints->Remove(emcRecPoint);
701 // emcRecPoints->Compress() ;
703 // fNumberOfTowerClusters -- ;
704 // numberofNotUnfolded-- ;
708 // delete[] maxAtEnergy ;
711 // // Unfolding of EMC clusters finished
714 // // Unfold now CPV clusters
715 // if(fNumberOfPreShoClusters > 0){
717 // Int_t nModulesToUnfold = geom->GetNModules() ;
719 // Int_t numberofPreShoNotUnfolded = fNumberOfPreShoClusters ;
721 // for(index = 0 ; index < numberofPreShoNotUnfolded ; index++){
723 // AliEMCALRecPoint * recPoint = (AliEMCALRecPoint *) cpvRecPoints->At(index) ;
725 // if(recPoint->GetEMCALMod()> nModulesToUnfold)
728 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint*) recPoint ;
730 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
731 // Int_t * maxAt = new Int_t[nMultipl] ;
732 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
733 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fPreShoLocMaxCut,digits) ;
735 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
736 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
737 // cpvRecPoints->Remove(emcRecPoint);
738 // cpvRecPoints->Compress() ;
740 // numberofPreShoNotUnfolded-- ;
741 // fNumberOfPreShoClusters-- ;
745 // delete[] maxAtEnergy ;
748 // //Unfolding of PreSho clusters finished
752 //____________________________________________________________________________
753 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
755 // Shape of the shower (see EMCAL TDR)
756 // If you change this function, change also the gradient evaluation in ChiSquare()
758 Double_t r4 = r*r*r*r ;
759 Double_t r295 = TMath::Power(r, 2.95) ;
760 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
764 //____________________________________________________________________________
765 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * iniTower,
767 AliEMCALDigit ** maxAt,
768 Float_t * maxAtEnergy)
770 // Performs the unfolding of a cluster with nMax overlapping showers
772 Fatal("AliEMCALClusterizerv1::UnfoldCluster", "--> Unfolding not implemented") ;
774 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
775 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
776 // const TClonesArray * digits = gime->Digits() ;
777 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
778 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
780 // Int_t nPar = 3 * nMax ;
781 // Float_t * fitparameters = new Float_t[nPar] ;
783 // Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
785 // // Fit failed, return and remove cluster
786 // delete[] fitparameters ;
790 // // create ufolded rec points and fill them with new energy lists
791 // // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
792 // // and later correct this number in acordance with actual energy deposition
794 // Int_t nDigits = iniTower->GetMultiplicity() ;
795 // Float_t * efit = new Float_t[nDigits] ;
796 // Float_t xDigit=0.,zDigit=0.,distance=0. ;
797 // Float_t xpar=0.,zpar=0.,epar=0. ;
799 // AliEMCALDigit * digit = 0 ;
800 // Int_t * emcDigits = iniTower->GetDigitsList() ;
804 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
805 // digit = (AliEMCALDigit*) digits->At(emcDigits[iDigit] ) ;
806 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
807 // geom->RelPosInModule(relid, xDigit, zDigit) ;
811 // while(iparam < nPar ){
812 // xpar = fitparameters[iparam] ;
813 // zpar = fitparameters[iparam+1] ;
814 // epar = fitparameters[iparam+2] ;
816 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
817 // distance = TMath::Sqrt(distance) ;
818 // efit[iDigit] += epar * ShowerShape(distance) ;
823 // // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
824 // // so that energy deposited in each cell is distributed betwin new clusters proportionally
825 // // to its contribution to efit
827 // Float_t * emcEnergies = iniTower->GetEnergiesList() ;
831 // while(iparam < nPar ){
832 // xpar = fitparameters[iparam] ;
833 // zpar = fitparameters[iparam+1] ;
834 // epar = fitparameters[iparam+2] ;
837 // AliEMCALTowerRecPoint * emcRP = 0 ;
839 // if(iniTower->IsTower()){ //create new entries in fTowerRecPoints...
841 // if(fNumberOfTowerClusters >= emcRecPoints->GetSize())
842 // emcRecPoints->Expand(2*fNumberOfTowerClusters) ;
844 // (*emcRecPoints)[fNumberOfTowerClusters] = new AliEMCALTowerRecPoint("") ;
845 // emcRP = (AliEMCALTowerRecPoint *) emcRecPoints->At(fNumberOfTowerClusters);
846 // fNumberOfTowerClusters++ ;
848 // else{//create new entries in fPreShoRecPoints
849 // if(fNumberOfPreShoClusters >= cpvRecPoints->GetSize())
850 // cpvRecPoints->Expand(2*fNumberOfPreShoClusters) ;
852 // (*cpvRecPoints)[fNumberOfPreShoClusters] = new AliEMCALPreShoRecPoint("") ;
853 // emcRP = (AliEMCALTowerRecPoint *) cpvRecPoints->At(fNumberOfPreShoClusters);
854 // fNumberOfPreShoClusters++ ;
858 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
859 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
860 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
861 // geom->RelPosInModule(relid, xDigit, zDigit) ;
862 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
863 // distance = TMath::Sqrt(distance) ;
864 // ratio = epar * ShowerShape(distance) / efit[iDigit] ;
865 // eDigit = emcEnergies[iDigit] * ratio ;
866 // emcRP->AddDigit( *digit, eDigit ) ;
870 // delete[] fitparameters ;
875 //_____________________________________________________________________________
876 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
878 // Calculates the Chi square for the cluster unfolding minimization
879 // Number of parameters, Gradient, Chi squared, parameters, what to do
882 // Fatal("AliEMCALClusterizerv1::UnfoldingChiSquare","-->Unfolding not implemented") ;
884 // TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
886 // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint*) toMinuit->At(0) ;
887 // TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
891 // // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint *) gMinuit->GetObjectFit() ; // TowerRecPoint to fit
893 // Int_t * emcDigits = emcRP->GetDigitsList() ;
895 // Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
897 // Float_t * emcEnergies = emcRP->GetEnergiesList() ;
899 // const AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
904 // for(iparam = 0 ; iparam < nPar ; iparam++)
905 // Grad[iparam] = 0 ; // Will evaluate gradient
909 // AliEMCALDigit * digit ;
912 // for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
914 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
920 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
922 // geom->RelPosInModule(relid, xDigit, zDigit) ;
924 // if(iflag == 2){ // calculate gradient
925 // Int_t iParam = 0 ;
927 // while(iParam < nPar ){
928 // Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
930 // distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
931 // distance = TMath::Sqrt( distance ) ;
933 // efit += x[iParam] * ShowerShape(distance) ;
936 // Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
938 // while(iParam < nPar ){
939 // Double_t xpar = x[iParam] ;
940 // Double_t zpar = x[iParam+1] ;
941 // Double_t epar = x[iParam+2] ;
942 // Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
943 // Double_t shape = sum * ShowerShape(dr) ;
944 // Double_t r4 = dr*dr*dr*dr ;
945 // Double_t r295 = TMath::Power(dr,2.95) ;
946 // Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
947 // 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
949 // Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
951 // Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
953 // Grad[iParam] += shape ; // Derivative over energy
960 // while(iparam < nPar ){
961 // Double_t xpar = x[iparam] ;
962 // Double_t zpar = x[iparam+1] ;
963 // Double_t epar = x[iparam+2] ;
965 // Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
966 // distance = TMath::Sqrt(distance) ;
967 // efit += epar * ShowerShape(distance) ;
970 // fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
971 // // Here we assume, that sigma = sqrt(E)
976 //____________________________________________________________________________
977 void AliEMCALClusterizerv1::Print(Option_t * option)const
979 // Print clusterizer parameters
981 if( strcmp(GetName(), "") !=0 ){
985 TString taskName(GetName()) ;
986 taskName.ReplaceAll(Version(), "") ;
988 cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl
989 << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl
990 << " Branch: " << fDigitsBranchTitle.Data() << endl
992 << " EMC Clustering threshold = " << fTowerClusteringThreshold << endl
993 << " EMC Local Maximum cut = " << fTowerLocMaxCut << endl
994 << " EMC Logarothmic weight = " << fW0 << endl
996 << " CPV Clustering threshold = " << fPreShoClusteringThreshold << endl
997 << " CPV Local Maximum cut = " << fPreShoLocMaxCut << endl
998 << " CPV Logarothmic weight = " << fW0CPV << endl
1001 cout << " Unfolding on " << endl ;
1003 cout << " Unfolding off " << endl ;
1005 cout << "------------------------------------------------------------------" <<endl ;
1008 cout << " AliEMCALClusterizerv1 not initialized " << endl ;
1010 //____________________________________________________________________________
1011 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
1013 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
1015 TObjArray * towerRecPoints = AliEMCALGetter::GetInstance()->TowerRecPoints() ;
1016 TObjArray * preshoRecPoints = AliEMCALGetter::GetInstance()->PreShowerRecPoints() ;
1018 cout << "AliEMCALClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1019 cout << " Found "<< towerRecPoints->GetEntriesFast() << " TOWER Rec Points and "
1020 << preshoRecPoints->GetEntriesFast() << " PRE SHOWER RecPoints" << endl ;
1022 fRecPointsInRun += towerRecPoints->GetEntriesFast() ;
1023 fRecPointsInRun += preshoRecPoints->GetEntriesFast() ;
1025 if(strstr(option,"all")) {
1027 cout << "Tower clusters " << endl ;
1028 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
1031 for (index = 0 ; index < towerRecPoints->GetEntries() ; index++) {
1032 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint * >(towerRecPoints->At(index)) ;
1034 rp->GetGlobalPosition(globalpos);
1036 rp->GetElipsAxis(lambda);
1039 primaries = rp->GetPrimaries(nprimaries);
1041 cout << setw(4) << rp->GetIndexInList() << " "
1042 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
1043 << setw(3) << rp->GetMultiplicity() << " "
1044 << setw(1) << rp->GetEMCALArm() << " "
1045 << setw(5) << setprecision(4) << globalpos.X() << " "
1046 << setw(5) << setprecision(4) << globalpos.Y() << " "
1047 << setw(5) << setprecision(4) << globalpos.Z() << " "
1048 << setw(4) << setprecision(2) << lambda[0] << " "
1049 << setw(4) << setprecision(2) << lambda[1] << " "
1050 << setw(2) << nprimaries << " " ;
1052 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1053 cout << setw(4) << primaries[iprimary] << " " ;
1057 //Now plot Pre shower recPoints
1059 cout << "-----------------------------------------------------------------------"<<endl ;
1061 cout << "PreShower clusters " << endl ;
1062 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
1064 for (index = 0 ; index < preshoRecPoints->GetEntries() ; index++) {
1065 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(index)) ;
1067 rp->GetGlobalPosition(globalpos);
1069 rp->GetElipsAxis(lambda);
1072 primaries = rp->GetPrimaries(nprimaries);
1074 cout << setw(4) << rp->GetIndexInList() << " "
1075 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
1076 << setw(3) << rp->GetMultiplicity() << " "
1077 << setw(1) << rp->GetEMCALArm() << " "
1078 << setw(5) << setprecision(4) << globalpos.X() << " "
1079 << setw(5) << setprecision(4) << globalpos.Y() << " "
1080 << setw(5) << setprecision(4) << globalpos.Z() << " "
1081 << setw(4) << setprecision(2) << lambda[0] << " "
1082 << setw(4) << setprecision(2) << lambda[1] << " "
1083 << setw(2) << nprimaries << " " ;
1085 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1086 cout << setw(4) << primaries[iprimary] << " " ;
1090 cout << "-----------------------------------------------------------------------"<<endl ;