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)
90 fDefaultInit = kTRUE ;
93 //____________________________________________________________________________
94 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const char* headerFile,const char* name)
95 :AliEMCALClusterizer(headerFile, name)
97 // ctor with the indication of the file where header Tree and digits Tree are stored
100 fDefaultInit = kFALSE ;
104 //____________________________________________________________________________
105 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
108 // fDefaultInit = kTRUE if Clusterizer created by default ctor (to get just the parameters)
111 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
113 // remove the task from the folder list
114 gime->RemoveTask("C",GetName()) ;
116 // remove the RecPoints from the folder list
117 TString name(GetName()) ;
118 name.Remove(name.Index(":")) ;
119 gime->RemoveObjects("D", name) ; // Digits
120 gime->RemoveObjects("RT", name) ; // TowerRecPoints
121 gime->RemoveObjects("RP", name) ; // PreShoRecPoints
129 //____________________________________________________________________________
130 const TString AliEMCALClusterizerv1::BranchName() const
132 TString branchName(GetName() ) ;
133 branchName.Remove(branchName.Index(Version())-1) ;
137 //____________________________________________________________________________
138 Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp, Bool_t inpresho) const
140 if ( inpresho ) // calibrate as pre shower
141 return -fADCpedestalPreSho + amp * fADCchannelPreSho ;
143 else //calibrate as tower
144 return -fADCpedestalTower + amp * fADCchannelTower ;
146 //____________________________________________________________________________
147 void AliEMCALClusterizerv1::Exec(Option_t * option)
151 if( strcmp(GetName(), "")== 0 )
154 if(strstr(option,"tim"))
155 gBenchmark->Start("EMCALClusterizer");
157 if(strstr(option,"print"))
160 gAlice->GetEvent(0) ;
162 //check, if the branch with name of this" already exits?
163 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
165 TBranch * branch = 0 ;
166 Bool_t emcaltowerfound = kFALSE, emcalpreshofound = kFALSE, clusterizerfound = kFALSE ;
168 TString branchname = GetName() ;
169 branchname.Remove(branchname.Index(Version())-1) ;
171 while ( (branch = (TBranch*)next()) && (!emcaltowerfound || !emcalpreshofound || !clusterizerfound) ) {
172 if ( (strcmp(branch->GetName(), "EMCALTowerRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) )
173 emcaltowerfound = kTRUE ;
175 else if ( (strcmp(branch->GetName(), "EMCALPreShoRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) )
176 emcalpreshofound = kTRUE ;
178 else if ((strcmp(branch->GetName(), "AliEMCALClusterizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) )
179 clusterizerfound = kTRUE ;
182 if ( emcalpreshofound || emcaltowerfound || clusterizerfound ) {
183 cerr << "WARNING: AliEMCALClusterizer::Exec -> Tower(PreSho)RecPoints and/or Clusterizer branch with name "
184 << branchname.Data() << " already exits" << endl ;
188 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
189 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
192 for(ievent = 0; ievent < nevents; ievent++){
195 GetCalibrationParameters() ;
197 fNumberOfTowerClusters = fNumberOfPreShoClusters = 0 ;
199 gime->Event(ievent,"D") ;
206 WriteRecPoints(ievent) ;
208 if(strstr(option,"deb"))
209 PrintRecPoints(option) ;
211 //increment the total number of digits per run
212 fRecPointsInRun += gime->TowerRecPoints()->GetEntriesFast() ;
213 fRecPointsInRun += gime->PreShowerRecPoints()->GetEntriesFast() ;
216 if(strstr(option,"tim")){
217 gBenchmark->Stop("EMCALClusterizer");
218 cout << "AliEMCALClusterizer:" << endl ;
219 cout << " took " << gBenchmark->GetCpuTime("EMCALClusterizer") << " seconds for Clusterizing "
220 << gBenchmark->GetCpuTime("EMCALClusterizer")/nevents << " seconds per event " << endl ;
226 //____________________________________________________________________________
227 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
228 Int_t nPar, Float_t * fitparameters) const
230 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
231 // The initial values for fitting procedure are set equal to the positions of local maxima.
232 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
234 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
235 TClonesArray * digits = gime->Digits() ;
238 gMinuit->mncler(); // Reset Minuit's list of paramters
239 gMinuit->SetPrintLevel(-1) ; // No Printout
240 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
241 // To set the address of the minimization function
243 TList * toMinuit = new TList();
244 toMinuit->AddAt(emcRP,0) ;
245 toMinuit->AddAt(digits,1) ;
247 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
249 // filling initial values for fit parameters
250 AliEMCALDigit * digit ;
254 Int_t nDigits = (Int_t) nPar / 3 ;
258 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
260 for(iDigit = 0; iDigit < nDigits; iDigit++){
261 digit = maxAt[iDigit];
266 geom->AbsToRelNumbering(digit->GetId(), relid) ;
267 geom->PosInAlice(relid, x, z) ;
269 Float_t energy = maxAtEnergy[iDigit] ;
271 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
274 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : x = " << x << endl ;
277 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
280 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : z = " << z << endl ;
283 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
286 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : energy = " << energy << endl ;
291 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
296 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
297 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
298 gMinuit->SetMaxIterations(5);
299 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
301 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
303 if(ierflg == 4){ // Minimum not found
304 cout << "EMCAL Unfolding> Fit not converged, cluster abandoned "<< endl ;
307 for(index = 0; index < nPar; index++){
310 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
311 fitparameters[index] = val ;
319 //____________________________________________________________________________
320 void AliEMCALClusterizerv1::GetCalibrationParameters()
322 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
323 TString branchname = GetName() ;
324 branchname.Remove(branchname.Index(Version())-1) ;
326 AliEMCALDigitizer * dig = gime->Digitizer(branchname) ;
328 fADCchannelTower = dig->GetTowerchannel() ;
329 fADCpedestalTower = dig->GetTowerpedestal();
331 fADCchannelPreSho = dig->GetPreShochannel() ;
332 fADCpedestalPreSho = dig->GetPreShopedestal() ;
335 //____________________________________________________________________________
336 void AliEMCALClusterizerv1::Init()
338 // Make all memory allocations which can not be done in default constructor.
339 // Attach the Clusterizer task to the list of EMCAL tasks
341 if ( strcmp(GetTitle(), "") == 0 )
342 SetTitle("galice.root") ;
344 TString branchname = GetName() ;
345 branchname.Remove(branchname.Index(Version())-1) ;
347 AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), branchname, "update") ;
349 cerr << "ERROR: AliEMCALClusterizerv1::Init -> Could not obtain the Getter object !" << endl ;
353 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
354 fNTowers = geom->GetNZ() * geom->GetNPhi() ;
357 gMinuit = new TMinuit(100) ;
359 gime->PostClusterizer(this) ;
360 // create a folder on the white board
361 gime->PostRecPoints(branchname ) ;
363 gime->PostDigits(branchname) ;
364 gime->PostDigitizer(branchname) ;
368 //____________________________________________________________________________
369 void AliEMCALClusterizerv1::InitParameters()
371 fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ;
375 fPreShoClusteringThreshold = 0.0001;
376 fTowerClusteringThreshold = 0.2;
378 fTowerLocMaxCut = 0.03 ;
379 fPreShoLocMaxCut = 0.03 ;
388 fHeaderFileName = GetTitle() ;
389 fDigitsBranchTitle = GetName() ;
391 TString clusterizerName( GetName()) ;
392 if (clusterizerName.IsNull() )
393 clusterizerName = "Default" ;
394 clusterizerName.Append(":") ;
395 clusterizerName.Append(Version()) ;
396 SetName(clusterizerName) ;
397 fRecPointsInRun = 0 ;
400 //____________________________________________________________________________
401 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
403 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
405 // = 2 are not neighbour but do not continue searching
406 // neighbours are defined as digits having at least a common vertex
407 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
408 // which is compared to a digit (d2) not yet in a cluster
410 AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
415 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
418 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
420 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same EMCAL Arm
421 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
422 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
424 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
425 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fTimeGate))
429 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
430 rv = 2; // Difference in row numbers is too large to look further
436 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
445 //____________________________________________________________________________
446 Bool_t AliEMCALClusterizerv1::IsInTower(AliEMCALDigit * digit) const
448 // Tells if (true) or not (false) the digit is in a EMCAL-Tower
451 if (!digit->IsInPreShower())
456 //____________________________________________________________________________
457 Bool_t AliEMCALClusterizerv1::IsInPreShower(AliEMCALDigit * digit) const
459 // Tells if (true) or not (false) the digit is in a EMCAL-PreShower
462 if (digit->IsInPreShower())
467 //____________________________________________________________________________
468 void AliEMCALClusterizerv1::WriteRecPoints(Int_t event)
471 // Creates new branches with given title
472 // fills and writes into TreeR.
474 AliEMCALGetter *gime = AliEMCALGetter::GetInstance() ;
475 TObjArray * towerRecPoints = gime->TowerRecPoints(BranchName()) ;
476 TObjArray * preshoRecPoints = gime->PreShowerRecPoints(BranchName()) ;
477 TClonesArray * digits = gime->Digits(BranchName()) ;
480 if (!gAlice->TreeR() )
481 gAlice->MakeTree("R", fSplitFile);
482 treeR = gAlice->TreeR() ;
485 //Evaluate position, dispersion and other RecPoint properties...
486 for(index = 0; index < towerRecPoints->GetEntries(); index++)
487 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->EvalAll(fW0,digits) ;
489 towerRecPoints->Sort() ;
491 for(index = 0; index < towerRecPoints->GetEntries(); index++)
492 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->SetIndexInList(index) ;
494 towerRecPoints->Expand(towerRecPoints->GetEntriesFast()) ;
496 //Now the same for pre shower
497 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
498 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->EvalAll(fW0CPV,digits) ;
500 preshoRecPoints->Sort() ;
502 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
503 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->SetIndexInList(index) ;
505 preshoRecPoints->Expand(preshoRecPoints->GetEntriesFast()) ;
507 //Make branches in TreeR for RecPoints and Clusterizer
510 Int_t bufferSize = 32000 ;
511 Int_t splitlevel = 0 ;
514 TBranch * emcBranch = treeR->Branch("EMCALTowerRP","TObjArray",&towerRecPoints,bufferSize,splitlevel);
515 emcBranch->SetTitle(BranchName());
518 //Now Pre Shower branch
519 TBranch * cpvBranch = treeR->Branch("EMCALPreShoRP","TObjArray",&preshoRecPoints,bufferSize,splitlevel);
520 cpvBranch->SetTitle(BranchName());
523 //And Finally clusterizer branch
524 AliEMCALClusterizerv1 * cl = (AliEMCALClusterizerv1*)gime->Clusterizer(BranchName()) ;
525 TBranch * clusterizerBranch = treeR->Branch("AliEMCALClusterizer","AliEMCALClusterizerv1",
526 &cl,bufferSize,splitlevel);
527 clusterizerBranch->SetTitle(BranchName());
531 clusterizerBranch->Fill() ;
533 treeR->AutoSave() ; //Write(0,kOverwrite) ;
537 //____________________________________________________________________________
538 void AliEMCALClusterizerv1::MakeClusters()
540 // Steering method to construct the clusters stored in a list of Reconstructed Points
541 // A cluster is defined as a list of neighbour digits
543 TString branchName(GetName()) ;
544 branchName.Remove(branchName.Index(Version())-1) ;
546 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
548 TObjArray * towerRecPoints = gime->TowerRecPoints(branchName) ;
549 TObjArray * preshoRecPoints = gime->PreShowerRecPoints(branchName) ;
550 towerRecPoints->Delete() ;
551 preshoRecPoints->Delete() ;
553 TClonesArray * digits = gime->Digits(branchName) ;
554 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
557 // Clusterization starts
559 TIter nextdigit(digitsC) ;
560 AliEMCALDigit * digit ;
561 Bool_t notremoved = kTRUE ;
563 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
564 AliEMCALRecPoint * clu = 0 ;
566 TArrayI clusterdigitslist(1500) ;
569 if (( IsInTower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fTowerClusteringThreshold ) ||
570 ( IsInPreShower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fPreShoClusteringThreshold ) ) {
572 Int_t iDigitInCluster = 0 ;
574 if ( IsInTower(digit) ) {
575 // start a new Tower RecPoint
576 if(fNumberOfTowerClusters >= towerRecPoints->GetSize())
577 towerRecPoints->Expand(2*fNumberOfTowerClusters+1) ;
579 towerRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfTowerClusters) ;
580 clu = dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(fNumberOfTowerClusters)) ;
581 fNumberOfTowerClusters++ ;
582 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower())) ;
583 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
585 digitsC->Remove(digit) ;
589 // start a new Pre Shower cluster
590 if(fNumberOfPreShoClusters >= preshoRecPoints->GetSize())
591 preshoRecPoints->Expand(2*fNumberOfPreShoClusters+1);
593 preshoRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfPreShoClusters) ;
595 clu = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(fNumberOfPreShoClusters)) ;
596 fNumberOfPreShoClusters++ ;
597 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower() ) );
598 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
600 digitsC->Remove(digit) ;
603 // Here we remove remaining Tower digits, which cannot make a cluster
606 while( ( digit = dynamic_cast<AliEMCALDigit *>(nextdigit()) ) ) {
607 if( IsInTower(digit) )
608 digitsC->Remove(digit) ;
612 notremoved = kFALSE ;
619 AliEMCALDigit * digitN ;
621 while (index < iDigitInCluster){ // scan over digits already in cluster
622 digit = (AliEMCALDigit*)digits->At(clusterdigitslist[index]) ;
624 while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits
625 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
627 case 0 : // not a neighbour
629 case 1 : // are neighbours
630 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->IsInPreShower() ) ) ;
631 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
633 digitsC->Remove(digitN) ;
635 case 2 : // too far from each other
644 } // loop over cluster
655 //____________________________________________________________________________
656 void AliEMCALClusterizerv1::MakeUnfolding()
658 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
660 // // Unfolds clusters using the shape of an ElectroMagnetic shower
661 // // Performs unfolding of all EMC/CPV clusters
663 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
665 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
666 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
667 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
668 // TClonesArray * digits = gime->Digits() ;
670 // // Unfold first EMC clusters
671 // if(fNumberOfTowerClusters > 0){
673 // Int_t nModulesToUnfold = geom->GetNModules() ;
675 // Int_t numberofNotUnfolded = fNumberOfTowerClusters ;
677 // for(index = 0 ; index < numberofNotUnfolded ; index++){
679 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint *) emcRecPoints->At(index) ;
680 // if(emcRecPoint->GetEMCALMod()> nModulesToUnfold)
683 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
684 // Int_t * maxAt = new Int_t[nMultipl] ;
685 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
686 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fTowerLocMaxCut,digits) ;
688 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
689 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
690 // emcRecPoints->Remove(emcRecPoint);
691 // emcRecPoints->Compress() ;
693 // fNumberOfTowerClusters -- ;
694 // numberofNotUnfolded-- ;
698 // delete[] maxAtEnergy ;
701 // // Unfolding of EMC clusters finished
704 // // Unfold now CPV clusters
705 // if(fNumberOfPreShoClusters > 0){
707 // Int_t nModulesToUnfold = geom->GetNModules() ;
709 // Int_t numberofPreShoNotUnfolded = fNumberOfPreShoClusters ;
711 // for(index = 0 ; index < numberofPreShoNotUnfolded ; index++){
713 // AliEMCALRecPoint * recPoint = (AliEMCALRecPoint *) cpvRecPoints->At(index) ;
715 // if(recPoint->GetEMCALMod()> nModulesToUnfold)
718 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint*) recPoint ;
720 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
721 // Int_t * maxAt = new Int_t[nMultipl] ;
722 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
723 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fPreShoLocMaxCut,digits) ;
725 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
726 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
727 // cpvRecPoints->Remove(emcRecPoint);
728 // cpvRecPoints->Compress() ;
730 // numberofPreShoNotUnfolded-- ;
731 // fNumberOfPreShoClusters-- ;
735 // delete[] maxAtEnergy ;
738 // //Unfolding of PreSho clusters finished
742 //____________________________________________________________________________
743 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
745 // Shape of the shower (see EMCAL TDR)
746 // If you change this function, change also the gradient evaluation in ChiSquare()
748 Double_t r4 = r*r*r*r ;
749 Double_t r295 = TMath::Power(r, 2.95) ;
750 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
754 //____________________________________________________________________________
755 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * iniTower,
757 AliEMCALDigit ** maxAt,
758 Float_t * maxAtEnergy)
760 // Performs the unfolding of a cluster with nMax overlapping showers
762 Fatal("AliEMCALClusterizerv1::UnfoldCluster", "--> Unfolding not implemented") ;
764 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
765 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
766 // const TClonesArray * digits = gime->Digits() ;
767 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
768 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
770 // Int_t nPar = 3 * nMax ;
771 // Float_t * fitparameters = new Float_t[nPar] ;
773 // Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
775 // // Fit failed, return and remove cluster
776 // delete[] fitparameters ;
780 // // create ufolded rec points and fill them with new energy lists
781 // // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
782 // // and later correct this number in acordance with actual energy deposition
784 // Int_t nDigits = iniTower->GetMultiplicity() ;
785 // Float_t * efit = new Float_t[nDigits] ;
786 // Float_t xDigit=0.,zDigit=0.,distance=0. ;
787 // Float_t xpar=0.,zpar=0.,epar=0. ;
789 // AliEMCALDigit * digit = 0 ;
790 // Int_t * emcDigits = iniTower->GetDigitsList() ;
794 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
795 // digit = (AliEMCALDigit*) digits->At(emcDigits[iDigit] ) ;
796 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
797 // geom->RelPosInModule(relid, xDigit, zDigit) ;
801 // while(iparam < nPar ){
802 // xpar = fitparameters[iparam] ;
803 // zpar = fitparameters[iparam+1] ;
804 // epar = fitparameters[iparam+2] ;
806 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
807 // distance = TMath::Sqrt(distance) ;
808 // efit[iDigit] += epar * ShowerShape(distance) ;
813 // // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
814 // // so that energy deposited in each cell is distributed betwin new clusters proportionally
815 // // to its contribution to efit
817 // Float_t * emcEnergies = iniTower->GetEnergiesList() ;
821 // while(iparam < nPar ){
822 // xpar = fitparameters[iparam] ;
823 // zpar = fitparameters[iparam+1] ;
824 // epar = fitparameters[iparam+2] ;
827 // AliEMCALTowerRecPoint * emcRP = 0 ;
829 // if(iniTower->IsTower()){ //create new entries in fTowerRecPoints...
831 // if(fNumberOfTowerClusters >= emcRecPoints->GetSize())
832 // emcRecPoints->Expand(2*fNumberOfTowerClusters) ;
834 // (*emcRecPoints)[fNumberOfTowerClusters] = new AliEMCALTowerRecPoint("") ;
835 // emcRP = (AliEMCALTowerRecPoint *) emcRecPoints->At(fNumberOfTowerClusters);
836 // fNumberOfTowerClusters++ ;
838 // else{//create new entries in fPreShoRecPoints
839 // if(fNumberOfPreShoClusters >= cpvRecPoints->GetSize())
840 // cpvRecPoints->Expand(2*fNumberOfPreShoClusters) ;
842 // (*cpvRecPoints)[fNumberOfPreShoClusters] = new AliEMCALPreShoRecPoint("") ;
843 // emcRP = (AliEMCALTowerRecPoint *) cpvRecPoints->At(fNumberOfPreShoClusters);
844 // fNumberOfPreShoClusters++ ;
848 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
849 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
850 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
851 // geom->RelPosInModule(relid, xDigit, zDigit) ;
852 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
853 // distance = TMath::Sqrt(distance) ;
854 // ratio = epar * ShowerShape(distance) / efit[iDigit] ;
855 // eDigit = emcEnergies[iDigit] * ratio ;
856 // emcRP->AddDigit( *digit, eDigit ) ;
860 // delete[] fitparameters ;
865 //_____________________________________________________________________________
866 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
868 // Calculates the Chi square for the cluster unfolding minimization
869 // Number of parameters, Gradient, Chi squared, parameters, what to do
872 // Fatal("AliEMCALClusterizerv1::UnfoldingChiSquare","-->Unfolding not implemented") ;
874 // TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
876 // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint*) toMinuit->At(0) ;
877 // TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
881 // // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint *) gMinuit->GetObjectFit() ; // TowerRecPoint to fit
883 // Int_t * emcDigits = emcRP->GetDigitsList() ;
885 // Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
887 // Float_t * emcEnergies = emcRP->GetEnergiesList() ;
889 // const AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
894 // for(iparam = 0 ; iparam < nPar ; iparam++)
895 // Grad[iparam] = 0 ; // Will evaluate gradient
899 // AliEMCALDigit * digit ;
902 // for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
904 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
910 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
912 // geom->RelPosInModule(relid, xDigit, zDigit) ;
914 // if(iflag == 2){ // calculate gradient
915 // Int_t iParam = 0 ;
917 // while(iParam < nPar ){
918 // Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
920 // distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
921 // distance = TMath::Sqrt( distance ) ;
923 // efit += x[iParam] * ShowerShape(distance) ;
926 // Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
928 // while(iParam < nPar ){
929 // Double_t xpar = x[iParam] ;
930 // Double_t zpar = x[iParam+1] ;
931 // Double_t epar = x[iParam+2] ;
932 // Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
933 // Double_t shape = sum * ShowerShape(dr) ;
934 // Double_t r4 = dr*dr*dr*dr ;
935 // Double_t r295 = TMath::Power(dr,2.95) ;
936 // Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
937 // 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
939 // Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
941 // Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
943 // Grad[iParam] += shape ; // Derivative over energy
950 // while(iparam < nPar ){
951 // Double_t xpar = x[iparam] ;
952 // Double_t zpar = x[iparam+1] ;
953 // Double_t epar = x[iparam+2] ;
955 // Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
956 // distance = TMath::Sqrt(distance) ;
957 // efit += epar * ShowerShape(distance) ;
960 // fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
961 // // Here we assume, that sigma = sqrt(E)
966 //____________________________________________________________________________
967 void AliEMCALClusterizerv1::Print(Option_t * option)const
969 // Print clusterizer parameters
971 if( strcmp(GetName(), "") !=0 ){
975 TString taskName(GetName()) ;
976 taskName.ReplaceAll(Version(), "") ;
978 cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl
979 << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl
980 << " Branch: " << fDigitsBranchTitle.Data() << endl
982 << " EMC Clustering threshold = " << fTowerClusteringThreshold << endl
983 << " EMC Local Maximum cut = " << fTowerLocMaxCut << endl
984 << " EMC Logarothmic weight = " << fW0 << endl
986 << " CPV Clustering threshold = " << fPreShoClusteringThreshold << endl
987 << " CPV Local Maximum cut = " << fPreShoLocMaxCut << endl
988 << " CPV Logarothmic weight = " << fW0CPV << endl
991 cout << " Unfolding on " << endl ;
993 cout << " Unfolding off " << endl ;
995 cout << "------------------------------------------------------------------" <<endl ;
998 cout << " AliEMCALClusterizerv1 not initialized " << endl ;
1000 //____________________________________________________________________________
1001 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
1003 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
1005 TObjArray * towerRecPoints = AliEMCALGetter::GetInstance()->TowerRecPoints() ;
1006 TObjArray * preshoRecPoints = AliEMCALGetter::GetInstance()->PreShowerRecPoints() ;
1008 cout << "AliEMCALClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1009 cout << " Found "<< towerRecPoints->GetEntriesFast() << " TOWER Rec Points and "
1010 << preshoRecPoints->GetEntriesFast() << " PRE SHOWER RecPoints" << endl ;
1012 fRecPointsInRun += towerRecPoints->GetEntriesFast() ;
1013 fRecPointsInRun += preshoRecPoints->GetEntriesFast() ;
1015 if(strstr(option,"all")) {
1017 cout << "Tower clusters " << endl ;
1018 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
1021 for (index = 0 ; index < towerRecPoints->GetEntries() ; index++) {
1022 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint * >(towerRecPoints->At(index)) ;
1024 rp->GetGlobalPosition(globalpos);
1026 rp->GetElipsAxis(lambda);
1029 primaries = rp->GetPrimaries(nprimaries);
1031 cout << setw(4) << rp->GetIndexInList() << " "
1032 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
1033 << setw(3) << rp->GetMultiplicity() << " "
1034 << setw(1) << rp->GetEMCALArm() << " "
1035 << setw(5) << setprecision(4) << globalpos.X() << " "
1036 << setw(5) << setprecision(4) << globalpos.Y() << " "
1037 << setw(5) << setprecision(4) << globalpos.Z() << " "
1038 << setw(4) << setprecision(2) << lambda[0] << " "
1039 << setw(4) << setprecision(2) << lambda[1] << " "
1040 << setw(2) << nprimaries << " " ;
1042 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1043 cout << setw(4) << primaries[iprimary] << " " ;
1047 //Now plot Pre shower recPoints
1049 cout << "-----------------------------------------------------------------------"<<endl ;
1051 cout << "PreShower clusters " << endl ;
1052 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
1054 for (index = 0 ; index < preshoRecPoints->GetEntries() ; index++) {
1055 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(index)) ;
1057 rp->GetGlobalPosition(globalpos);
1059 rp->GetElipsAxis(lambda);
1062 primaries = rp->GetPrimaries(nprimaries);
1064 cout << setw(4) << rp->GetIndexInList() << " "
1065 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
1066 << setw(3) << rp->GetMultiplicity() << " "
1067 << setw(1) << rp->GetEMCALArm() << " "
1068 << setw(5) << setprecision(4) << globalpos.X() << " "
1069 << setw(5) << setprecision(4) << globalpos.Y() << " "
1070 << setw(5) << setprecision(4) << globalpos.Z() << " "
1071 << setw(4) << setprecision(2) << lambda[0] << " "
1072 << setw(4) << setprecision(2) << lambda[1] << " "
1073 << setw(2) << nprimaries << " " ;
1075 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1076 cout << setw(4) << primaries[iprimary] << " " ;
1080 cout << "-----------------------------------------------------------------------"<<endl ;