1 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
7 * Author: The ALICE Off-line Project. *
9 * Contributors are mentioned in the code where appropriate. *
13 * Permission to use, copy, modify and distribute this software and its *
15 * documentation strictly for non-commercial purposes is hereby granted *
17 * without fee, provided that the above copyright notice appears in all *
19 * copies and that both the copyright notice and this permission notice *
21 * appear in the supporting documentation. The authors make no claims *
23 * about the suitability of this software for any purpose. It is *
25 * provided "as is" without express or implied warranty. *
27 **************************************************************************/
35 1 October 2000. Yuri Kharlov:
39 PPSD upper layer is considered if number of layers>1
43 18 October 2000. Yuri Kharlov:
45 AliEMCALClusterizerv1()
47 CPV clusterizing parameters added
53 After first PPSD digit remove EMC digits only once
57 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
59 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
61 // of new IO (à la PHOS)
63 //////////////////////////////////////////////////////////////////////////////
65 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
67 // unfolds the clusters having several local maxima.
69 // Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
71 // EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all
73 // parameters including input digits branch title, thresholds etc.)
75 // This TTask is normally called from Reconstructioner, but can as well be used in
81 // root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")
83 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
85 // //reads gAlice from header file "..."
87 // root [1] cl->ExecuteTask()
89 // //finds RecPoints in all events stored in galice.root
91 // root [2] cl->SetDigitsBranch("digits2")
93 // //sets another title for Digitis (input) branch
95 // root [3] cl->SetRecPointsBranch("recp2")
97 // //sets another title four output branches
99 // root [4] cl->SetTowerLocalMaxCut(0.03)
101 // //set clusterization parameters
103 // root [5] cl->ExecuteTask("deb all time")
105 // //once more finds RecPoints options are
107 // // deb - print number of found rec points
109 // // deb all - print number of found RecPoints and some their characteristics
111 // // time - print benchmarking results
115 // --- ROOT system ---
133 #include "TBenchmark.h"
137 // --- Standard library ---
141 #include <iostream.h>
147 // --- AliRoot header files ---
151 #include "AliEMCALClusterizerv1.h"
153 #include "AliEMCALDigit.h"
155 #include "AliEMCALDigitizer.h"
157 #include "AliEMCALTowerRecPoint.h"
159 #include "AliEMCAL.h"
161 #include "AliEMCALGetter.h"
163 #include "AliEMCALGeometry.h"
169 ClassImp(AliEMCALClusterizerv1)
173 //____________________________________________________________________________
175 AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
179 // default ctor (to be used mainly by Streamer)
185 fDefaultInit = kTRUE ;
191 //____________________________________________________________________________
193 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const char* headerFile, const char* name, const Bool_t toSplit)
195 :AliEMCALClusterizer(headerFile, name, toSplit)
199 // ctor with the indication of the file where header Tree and digits Tree are stored
207 fDefaultInit = kFALSE ;
215 //____________________________________________________________________________
217 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
231 //____________________________________________________________________________
233 const TString AliEMCALClusterizerv1::BranchName() const
237 TString branchName(GetName() ) ;
239 branchName.Remove(branchName.Index(Version())-1) ;
247 //____________________________________________________________________________
249 Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp, Bool_t inpresho) const
251 {//To be replased later by the method, reading individual parameters from the database
255 if ( inpresho ) // calibrate as pre shower
257 return -fADCpedestalPreSho + amp * fADCchannelPreSho ;
261 else //calibrate as tower
263 return -fADCpedestalTower + amp * fADCchannelTower ;
269 //____________________________________________________________________________
271 void AliEMCALClusterizerv1::Exec(Option_t * option)
279 if( strcmp(GetName(), "")== 0 )
285 if(strstr(option,"tim"))
287 gBenchmark->Start("EMCALClusterizer");
291 if(strstr(option,"print"))
297 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
299 if(gime->BranchExists("RecPoints"))
303 Int_t nevents = gime->MaxEvent() ;
309 for(ievent = 0; ievent < nevents; ievent++){
313 gime->Event(ievent,"D") ;
319 GetCalibrationParameters() ;
323 fNumberOfTowerClusters = fNumberOfPreShoClusters = 0 ;
337 WriteRecPoints(ievent) ;
341 if(strstr(option,"deb"))
343 PrintRecPoints(option) ;
347 //increment the total number of digits per run
349 fRecPointsInRun += gime->TowerRecPoints()->GetEntriesFast() ;
351 fRecPointsInRun += gime->PreShowerRecPoints()->GetEntriesFast() ;
357 if(strstr(option,"tim")){
359 gBenchmark->Stop("EMCALClusterizer");
361 cout << "AliEMCALClusterizer:" << endl ;
363 cout << " took " << gBenchmark->GetCpuTime("EMCALClusterizer") << " seconds for Clusterizing "
365 << gBenchmark->GetCpuTime("EMCALClusterizer")/nevents << " seconds per event " << endl ;
377 //____________________________________________________________________________
379 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
381 Int_t nPar, Float_t * fitparameters) const
385 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
387 // The initial values for fitting procedure are set equal to the positions of local maxima.
389 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
393 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
395 TClonesArray * digits = gime->Digits() ;
401 gMinuit->mncler(); // Reset Minuit's list of paramters
403 gMinuit->SetPrintLevel(-1) ; // No Printout
405 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
407 // To set the address of the minimization function
411 TList * toMinuit = new TList();
413 toMinuit->AddAt(emcRP,0) ;
415 toMinuit->AddAt(digits,1) ;
419 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
423 // filling initial values for fit parameters
425 AliEMCALDigit * digit ;
433 Int_t nDigits = (Int_t) nPar / 3 ;
441 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
445 for(iDigit = 0; iDigit < nDigits; iDigit++){
447 digit = maxAt[iDigit];
457 geom->AbsToRelNumbering(digit->GetId(), relid) ;
459 geom->PosInAlice(relid, x, z) ;
463 Float_t energy = maxAtEnergy[iDigit] ;
467 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
473 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : x = " << x << endl ;
479 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
485 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : z = " << z << endl ;
491 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
497 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : energy = " << energy << endl ;
507 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
517 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
519 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
521 gMinuit->SetMaxIterations(5);
523 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
527 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
531 if(ierflg == 4){ // Minimum not found
533 cout << "EMCAL Unfolding> Fit not converged, cluster abandoned "<< endl ;
539 for(index = 0; index < nPar; index++){
545 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
547 fitparameters[index] = val ;
563 //____________________________________________________________________________
565 void AliEMCALClusterizerv1::GetCalibrationParameters()
569 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
571 const AliEMCALDigitizer * dig = gime->Digitizer(BranchName()) ;
575 fADCchannelTower = dig->GetTowerchannel() ;
577 fADCpedestalTower = dig->GetTowerpedestal();
581 fADCchannelPreSho = dig->GetPreShochannel() ;
583 fADCpedestalPreSho = dig->GetPreShopedestal() ;
591 //____________________________________________________________________________
593 void AliEMCALClusterizerv1::Init()
597 // Make all memory allocations which can not be done in default constructor.
599 // Attach the Clusterizer task to the list of EMCAL tasks
603 if ( strcmp(GetTitle(), "") == 0 )
605 SetTitle("galice.root") ;
609 TString branchname = GetName() ;
611 branchname.Remove(branchname.Index(Version())-1) ;
615 AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), branchname.Data(), fToSplit ) ;
619 cerr << "ERROR: AliEMCALClusterizerv1::Init -> Could not obtain the Getter object !" << endl ;
631 // construct the name of the file as /path/EMCAL.SDigits.root
633 //First - extract full path if necessary
635 TString fileName(GetTitle()) ;
637 Ssiz_t islash = fileName.Last('/') ;
639 if(islash<fileName.Length())
641 fileName.Remove(islash+1,fileName.Length()) ;
647 // Next - append the file name
649 fileName+="EMCAL.RecData." ;
651 if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
653 fileName+=branchname ;
661 // Finally - check if the file already opened or open the file
663 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
667 fSplitFile = TFile::Open(fileName.Data(),"update") ;
675 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
677 fNTowers = geom->GetNZ() * geom->GetNPhi() ;
683 gMinuit = new TMinuit(100) ;
687 gime->PostClusterizer(this) ;
689 gime->PostRecPoints(branchname ) ;
697 //____________________________________________________________________________
699 void AliEMCALClusterizerv1::InitParameters()
703 fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ;
711 fPreShoClusteringThreshold = 0.0001;
713 fTowerClusteringThreshold = 0.2;
717 fTowerLocMaxCut = 0.03 ;
719 fPreShoLocMaxCut = 0.03 ;
737 TString clusterizerName( GetName()) ;
739 if (clusterizerName.IsNull() )
741 clusterizerName = "Default" ;
743 clusterizerName.Append(":") ;
745 clusterizerName.Append(Version()) ;
747 SetName(clusterizerName) ;
749 fRecPointsInRun = 0 ;
757 //____________________________________________________________________________
759 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
763 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
767 // = 2 are not neighbour but do not continue searching
769 // neighbours are defined as digits having at least a common vertex
771 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
773 // which is compared to a digit (d2) not yet in a cluster
777 AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
787 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
793 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
797 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same EMCAL Arm
799 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
801 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
805 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
807 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fTimeGate))
815 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
817 rv = 2; // Difference in row numbers is too large to look further
829 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
847 //____________________________________________________________________________
849 Bool_t AliEMCALClusterizerv1::IsInTower(AliEMCALDigit * digit) const
853 // Tells if (true) or not (false) the digit is in a EMCAL-Tower
859 if (!digit->IsInPreShower())
869 //____________________________________________________________________________
871 Bool_t AliEMCALClusterizerv1::IsInPreShower(AliEMCALDigit * digit) const
875 // Tells if (true) or not (false) the digit is in a EMCAL-PreShower
881 if (digit->IsInPreShower())
891 //____________________________________________________________________________
893 void AliEMCALClusterizerv1::WriteRecPoints(Int_t event)
899 // Creates new branches with given title
901 // fills and writes into TreeR.
905 AliEMCALGetter *gime = AliEMCALGetter::GetInstance() ;
907 TObjArray * towerRecPoints = gime->TowerRecPoints() ;
909 TObjArray * preshoRecPoints = gime->PreShowerRecPoints() ;
911 TClonesArray * digits = gime->Digits() ;
925 TString name("TreeR") ;
929 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
935 treeR = gAlice->TreeR();
943 gAlice->MakeTree("R", fSplitFile);
945 treeR = gAlice->TreeR() ;
953 //Evaluate position, dispersion and other RecPoint properties...
955 for(index = 0; index < towerRecPoints->GetEntries(); index++)
957 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->EvalAll(fW0,digits) ;
961 towerRecPoints->Sort() ;
965 for(index = 0; index < towerRecPoints->GetEntries(); index++)
967 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->SetIndexInList(index) ;
971 towerRecPoints->Expand(towerRecPoints->GetEntriesFast()) ;
975 //Now the same for pre shower
977 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
979 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->EvalAll(fW0CPV,digits) ;
983 preshoRecPoints->Sort() ;
987 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
989 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->SetIndexInList(index) ;
993 preshoRecPoints->Expand(preshoRecPoints->GetEntriesFast()) ;
997 Int_t bufferSize = 32000 ;
999 Int_t splitlevel = 0 ;
1003 //First Tower branch
1005 TBranch * towerBranch = treeR->Branch("EMCALTowerRP","TObjArray",&towerRecPoints,bufferSize,splitlevel);
1007 towerBranch->SetTitle(BranchName());
1011 //Now Pre Shower branch
1013 TBranch * preshoBranch = treeR->Branch("EMCALPreShoRP","TObjArray",&preshoRecPoints,bufferSize,splitlevel);
1015 preshoBranch->SetTitle(BranchName());
1019 //And Finally clusterizer branch
1021 AliEMCALClusterizerv1 * cl = (AliEMCALClusterizerv1*)gime->Clusterizer(BranchName()) ;
1023 TBranch * clusterizerBranch = treeR->Branch("AliEMCALClusterizer","AliEMCALClusterizerv1",
1025 &cl,bufferSize,splitlevel);
1027 clusterizerBranch->SetTitle(BranchName());
1031 towerBranch ->Fill() ;
1033 preshoBranch ->Fill() ;
1035 clusterizerBranch->Fill() ;
1039 treeR->AutoSave() ; //Write(0,kOverwrite) ;
1041 if(gAlice->TreeR()!=treeR)
1049 //____________________________________________________________________________
1051 void AliEMCALClusterizerv1::MakeClusters()
1055 // Steering method to construct the clusters stored in a list of Reconstructed Points
1057 // A cluster is defined as a list of neighbour digits
1061 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
1065 TObjArray * towerRecPoints = gime->TowerRecPoints(BranchName()) ;
1067 TObjArray * preshoRecPoints = gime->PreShowerRecPoints(BranchName()) ;
1069 towerRecPoints->Delete() ;
1071 preshoRecPoints->Delete() ;
1075 TClonesArray * digits = gime->Digits() ;
1079 cerr << "ERROR: AliEMCALClusterizerv1::MakeClusters -> Digits with name "
1081 << GetName() << " not found ! " << endl ;
1087 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
1093 // Clusterization starts
1097 TIter nextdigit(digitsC) ;
1099 AliEMCALDigit * digit ;
1101 Bool_t notremoved = kTRUE ;
1105 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
1107 AliEMCALRecPoint * clu = 0 ;
1111 TArrayI clusterdigitslist(1500) ;
1117 if (( IsInTower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fTowerClusteringThreshold ) ||
1119 ( IsInPreShower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fPreShoClusteringThreshold ) ) {
1123 Int_t iDigitInCluster = 0 ;
1127 if ( IsInTower(digit) ) {
1129 // start a new Tower RecPoint
1131 if(fNumberOfTowerClusters >= towerRecPoints->GetSize())
1133 towerRecPoints->Expand(2*fNumberOfTowerClusters+1) ;
1137 towerRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfTowerClusters) ;
1139 clu = dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(fNumberOfTowerClusters)) ;
1141 fNumberOfTowerClusters++ ;
1143 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower())) ;
1145 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
1149 digitsC->Remove(digit) ;
1157 // start a new Pre Shower cluster
1159 if(fNumberOfPreShoClusters >= preshoRecPoints->GetSize())
1161 preshoRecPoints->Expand(2*fNumberOfPreShoClusters+1);
1165 preshoRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfPreShoClusters) ;
1169 clu = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(fNumberOfPreShoClusters)) ;
1171 fNumberOfPreShoClusters++ ;
1173 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower() ) );
1175 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
1179 digitsC->Remove(digit) ;
1185 // Here we remove remaining Tower digits, which cannot make a cluster
1191 while( ( digit = dynamic_cast<AliEMCALDigit *>(nextdigit()) ) ) {
1193 if( IsInTower(digit) )
1195 digitsC->Remove(digit) ;
1203 notremoved = kFALSE ;
1217 AliEMCALDigit * digitN ;
1221 while (index < iDigitInCluster){ // scan over digits already in cluster
1223 digit = (AliEMCALDigit*)digits->At(clusterdigitslist[index]) ;
1227 while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits
1229 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
1233 case 0 : // not a neighbour
1237 case 1 : // are neighbours
1239 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->IsInPreShower() ) ) ;
1241 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
1245 digitsC->Remove(digitN) ;
1249 case 2 : // too far from each other
1267 } // loop over cluster
1271 } // energy theshold
1289 //____________________________________________________________________________
1291 void AliEMCALClusterizerv1::MakeUnfolding()
1295 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
1299 // // Unfolds clusters using the shape of an ElectroMagnetic shower
1301 // // Performs unfolding of all EMC/CPV clusters
1305 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
1309 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
1311 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
1313 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
1315 // TClonesArray * digits = gime->Digits() ;
1319 // // Unfold first EMC clusters
1321 // if(fNumberOfTowerClusters > 0){
1325 // Int_t nModulesToUnfold = geom->GetNModules() ;
1329 // Int_t numberofNotUnfolded = fNumberOfTowerClusters ;
1333 // for(index = 0 ; index < numberofNotUnfolded ; index++){
1337 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint *) emcRecPoints->At(index) ;
1339 // if(emcRecPoint->GetEMCALMod()> nModulesToUnfold)
1345 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
1347 // Int_t * maxAt = new Int_t[nMultipl] ;
1349 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
1351 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fTowerLocMaxCut,digits) ;
1355 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
1357 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
1359 // emcRecPoints->Remove(emcRecPoint);
1361 // emcRecPoints->Compress() ;
1365 // fNumberOfTowerClusters -- ;
1367 // numberofNotUnfolded-- ;
1375 // delete[] maxAtEnergy ;
1381 // // Unfolding of EMC clusters finished
1387 // // Unfold now CPV clusters
1389 // if(fNumberOfPreShoClusters > 0){
1393 // Int_t nModulesToUnfold = geom->GetNModules() ;
1397 // Int_t numberofPreShoNotUnfolded = fNumberOfPreShoClusters ;
1401 // for(index = 0 ; index < numberofPreShoNotUnfolded ; index++){
1405 // AliEMCALRecPoint * recPoint = (AliEMCALRecPoint *) cpvRecPoints->At(index) ;
1409 // if(recPoint->GetEMCALMod()> nModulesToUnfold)
1415 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint*) recPoint ;
1419 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
1421 // Int_t * maxAt = new Int_t[nMultipl] ;
1423 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
1425 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fPreShoLocMaxCut,digits) ;
1429 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
1431 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
1433 // cpvRecPoints->Remove(emcRecPoint);
1435 // cpvRecPoints->Compress() ;
1439 // numberofPreShoNotUnfolded-- ;
1441 // fNumberOfPreShoClusters-- ;
1449 // delete[] maxAtEnergy ;
1455 // //Unfolding of PreSho clusters finished
1463 //____________________________________________________________________________
1465 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
1469 // Shape of the shower (see EMCAL TDR)
1471 // If you change this function, change also the gradient evaluation in ChiSquare()
1475 Double_t r4 = r*r*r*r ;
1477 Double_t r295 = TMath::Power(r, 2.95) ;
1479 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
1487 //____________________________________________________________________________
1489 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * iniTower,
1493 AliEMCALDigit ** maxAt,
1495 Float_t * maxAtEnergy)
1499 // Performs the unfolding of a cluster with nMax overlapping showers
1503 Fatal("AliEMCALClusterizerv1::UnfoldCluster", "--> Unfolding not implemented") ;
1507 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
1509 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
1511 // const TClonesArray * digits = gime->Digits() ;
1513 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
1515 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
1519 // Int_t nPar = 3 * nMax ;
1521 // Float_t * fitparameters = new Float_t[nPar] ;
1525 // Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
1529 // // Fit failed, return and remove cluster
1531 // delete[] fitparameters ;
1539 // // create ufolded rec points and fill them with new energy lists
1541 // // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
1543 // // and later correct this number in acordance with actual energy deposition
1547 // Int_t nDigits = iniTower->GetMultiplicity() ;
1549 // Float_t * efit = new Float_t[nDigits] ;
1551 // Float_t xDigit=0.,zDigit=0.,distance=0. ;
1553 // Float_t xpar=0.,zpar=0.,epar=0. ;
1557 // AliEMCALDigit * digit = 0 ;
1559 // Int_t * emcDigits = iniTower->GetDigitsList() ;
1567 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
1569 // digit = (AliEMCALDigit*) digits->At(emcDigits[iDigit] ) ;
1571 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
1573 // geom->RelPosInModule(relid, xDigit, zDigit) ;
1575 // efit[iDigit] = 0;
1581 // while(iparam < nPar ){
1583 // xpar = fitparameters[iparam] ;
1585 // zpar = fitparameters[iparam+1] ;
1587 // epar = fitparameters[iparam+2] ;
1591 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1593 // distance = TMath::Sqrt(distance) ;
1595 // efit[iDigit] += epar * ShowerShape(distance) ;
1605 // // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
1607 // // so that energy deposited in each cell is distributed betwin new clusters proportionally
1609 // // to its contribution to efit
1613 // Float_t * emcEnergies = iniTower->GetEnergiesList() ;
1621 // while(iparam < nPar ){
1623 // xpar = fitparameters[iparam] ;
1625 // zpar = fitparameters[iparam+1] ;
1627 // epar = fitparameters[iparam+2] ;
1633 // AliEMCALTowerRecPoint * emcRP = 0 ;
1637 // if(iniTower->IsTower()){ //create new entries in fTowerRecPoints...
1641 // if(fNumberOfTowerClusters >= emcRecPoints->GetSize())
1643 // emcRecPoints->Expand(2*fNumberOfTowerClusters) ;
1647 // (*emcRecPoints)[fNumberOfTowerClusters] = new AliEMCALTowerRecPoint("") ;
1649 // emcRP = (AliEMCALTowerRecPoint *) emcRecPoints->At(fNumberOfTowerClusters);
1651 // fNumberOfTowerClusters++ ;
1655 // else{//create new entries in fPreShoRecPoints
1657 // if(fNumberOfPreShoClusters >= cpvRecPoints->GetSize())
1659 // cpvRecPoints->Expand(2*fNumberOfPreShoClusters) ;
1663 // (*cpvRecPoints)[fNumberOfPreShoClusters] = new AliEMCALPreShoRecPoint("") ;
1665 // emcRP = (AliEMCALTowerRecPoint *) cpvRecPoints->At(fNumberOfPreShoClusters);
1667 // fNumberOfPreShoClusters++ ;
1675 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
1677 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
1679 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
1681 // geom->RelPosInModule(relid, xDigit, zDigit) ;
1683 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1685 // distance = TMath::Sqrt(distance) ;
1687 // ratio = epar * ShowerShape(distance) / efit[iDigit] ;
1689 // eDigit = emcEnergies[iDigit] * ratio ;
1691 // emcRP->AddDigit( *digit, eDigit ) ;
1699 // delete[] fitparameters ;
1709 //_____________________________________________________________________________
1711 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
1715 // Calculates the Chi square for the cluster unfolding minimization
1717 // Number of parameters, Gradient, Chi squared, parameters, what to do
1723 // Fatal("AliEMCALClusterizerv1::UnfoldingChiSquare","-->Unfolding not implemented") ;
1727 // TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
1731 // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint*) toMinuit->At(0) ;
1733 // TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
1741 // // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint *) gMinuit->GetObjectFit() ; // TowerRecPoint to fit
1745 // Int_t * emcDigits = emcRP->GetDigitsList() ;
1749 // Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
1753 // Float_t * emcEnergies = emcRP->GetEnergiesList() ;
1757 // const AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
1767 // for(iparam = 0 ; iparam < nPar ; iparam++)
1769 // Grad[iparam] = 0 ; // Will evaluate gradient
1777 // AliEMCALDigit * digit ;
1783 // for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
1787 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
1799 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
1803 // geom->RelPosInModule(relid, xDigit, zDigit) ;
1807 // if(iflag == 2){ // calculate gradient
1809 // Int_t iParam = 0 ;
1813 // while(iParam < nPar ){
1815 // Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
1819 // distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
1821 // distance = TMath::Sqrt( distance ) ;
1825 // efit += x[iParam] * ShowerShape(distance) ;
1831 // Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
1835 // while(iParam < nPar ){
1837 // Double_t xpar = x[iParam] ;
1839 // Double_t zpar = x[iParam+1] ;
1841 // Double_t epar = x[iParam+2] ;
1843 // Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
1845 // Double_t shape = sum * ShowerShape(dr) ;
1847 // Double_t r4 = dr*dr*dr*dr ;
1849 // Double_t r295 = TMath::Power(dr,2.95) ;
1851 // Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
1853 // 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
1857 // Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
1861 // Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
1865 // Grad[iParam] += shape ; // Derivative over energy
1879 // while(iparam < nPar ){
1881 // Double_t xpar = x[iparam] ;
1883 // Double_t zpar = x[iparam+1] ;
1885 // Double_t epar = x[iparam+2] ;
1889 // Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1891 // distance = TMath::Sqrt(distance) ;
1893 // efit += epar * ShowerShape(distance) ;
1899 // fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
1901 // // Here we assume, that sigma = sqrt(E)
1911 //____________________________________________________________________________
1913 void AliEMCALClusterizerv1::Print(Option_t * option)const
1917 // Print clusterizer parameters
1921 if( strcmp(GetName(), "") !=0 ){
1929 TString taskName(GetName()) ;
1931 taskName.ReplaceAll(Version(), "") ;
1935 cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl
1937 << "Clusterizing digits from the file: " << taskName.Data() << endl
1939 << " Branch: " << GetName() << endl
1943 << " EMC Clustering threshold = " << fTowerClusteringThreshold << endl
1945 << " EMC Local Maximum cut = " << fTowerLocMaxCut << endl
1947 << " EMC Logarothmic weight = " << fW0 << endl
1951 << " CPV Clustering threshold = " << fPreShoClusteringThreshold << endl
1953 << " CPV Local Maximum cut = " << fPreShoLocMaxCut << endl
1955 << " CPV Logarothmic weight = " << fW0CPV << endl
1961 cout << " Unfolding on " << endl ;
1965 cout << " Unfolding off " << endl ;
1969 cout << "------------------------------------------------------------------" <<endl ;
1975 cout << " AliEMCALClusterizerv1 not initialized " << endl ;
1979 //____________________________________________________________________________
1981 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
1985 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
1989 TObjArray * towerRecPoints = AliEMCALGetter::GetInstance()->TowerRecPoints() ;
1991 TObjArray * preshoRecPoints = AliEMCALGetter::GetInstance()->PreShowerRecPoints() ;
1995 cout << "AliEMCALClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1997 cout << " Found "<< towerRecPoints->GetEntriesFast() << " TOWER Rec Points and "
1999 << preshoRecPoints->GetEntriesFast() << " PRE SHOWER RecPoints" << endl ;
2003 fRecPointsInRun += towerRecPoints->GetEntriesFast() ;
2005 fRecPointsInRun += preshoRecPoints->GetEntriesFast() ;
2009 if(strstr(option,"all")) {
2013 cout << "Tower clusters " << endl ;
2015 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
2021 for (index = 0 ; index < towerRecPoints->GetEntries() ; index++) {
2023 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint * >(towerRecPoints->At(index)) ;
2027 rp->GetGlobalPosition(globalpos);
2031 rp->GetElipsAxis(lambda);
2037 primaries = rp->GetPrimaries(nprimaries);
2041 cout << setw(4) << rp->GetIndexInList() << " "
2043 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
2045 << setw(3) << rp->GetMultiplicity() << " "
2047 << setw(1) << rp->GetEMCALArm() << " "
2049 << setw(5) << setprecision(4) << globalpos.X() << " "
2051 << setw(5) << setprecision(4) << globalpos.Y() << " "
2053 << setw(5) << setprecision(4) << globalpos.Z() << " "
2055 << setw(4) << setprecision(2) << lambda[0] << " "
2057 << setw(4) << setprecision(2) << lambda[1] << " "
2059 << setw(2) << nprimaries << " " ;
2063 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
2065 cout << setw(4) << primaries[iprimary] << " " ;
2073 //Now plot Pre shower recPoints
2077 cout << "-----------------------------------------------------------------------"<<endl ;
2081 cout << "PreShower clusters " << endl ;
2083 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
2087 for (index = 0 ; index < preshoRecPoints->GetEntries() ; index++) {
2089 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(index)) ;
2093 rp->GetGlobalPosition(globalpos);
2097 rp->GetElipsAxis(lambda);
2103 primaries = rp->GetPrimaries(nprimaries);
2107 cout << setw(4) << rp->GetIndexInList() << " "
2109 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
2111 << setw(3) << rp->GetMultiplicity() << " "
2113 << setw(1) << rp->GetEMCALArm() << " "
2115 << setw(5) << setprecision(4) << globalpos.X() << " "
2117 << setw(5) << setprecision(4) << globalpos.Y() << " "
2119 << setw(5) << setprecision(4) << globalpos.Z() << " "
2121 << setw(4) << setprecision(2) << lambda[0] << " "
2123 << setw(4) << setprecision(2) << lambda[1] << " "
2125 << setw(2) << nprimaries << " " ;
2129 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
2131 cout << setw(4) << primaries[iprimary] << " " ;
2139 cout << "-----------------------------------------------------------------------"<<endl ;