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 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
31 // of new IO (à la PHOS)
32 //////////////////////////////////////////////////////////////////////////////
33 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
34 // unfolds the clusters having several local maxima.
35 // Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
36 // EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all
37 // parameters including input digits branch title, thresholds etc.)
38 // This TTask is normally called from Reconstructioner, but can as well be used in
41 // root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")
42 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
43 // //reads gAlice from header file "..."
44 // root [1] cl->ExecuteTask()
45 // //finds RecPoints in all events stored in galice.root
46 // root [2] cl->SetDigitsBranch("digits2")
47 // //sets another title for Digitis (input) branch
48 // root [3] cl->SetRecPointsBranch("recp2")
49 // //sets another title four output branches
50 // root [4] cl->SetTowerLocalMaxCut(0.03)
51 // //set clusterization parameters
52 // root [5] cl->ExecuteTask("deb all time")
53 // //once more finds RecPoints options are
54 // // deb - print number of found rec points
55 // // deb all - print number of found RecPoints and some their characteristics
56 // // time - print benchmarking results
58 // --- ROOT system ---
67 #include "TBenchmark.h"
69 // --- Standard library ---
74 // --- AliRoot header files ---
76 #include "AliEMCALClusterizerv1.h"
77 #include "AliEMCALDigit.h"
78 #include "AliEMCALDigitizer.h"
79 #include "AliEMCALTowerRecPoint.h"
81 #include "AliEMCALGetter.h"
82 #include "AliEMCALGeometry.h"
85 ClassImp(AliEMCALClusterizerv1)
87 //____________________________________________________________________________
88 AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
90 // default ctor (to be used mainly by Streamer)
93 fDefaultInit = kTRUE ;
96 //____________________________________________________________________________
97 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const char* headerFile, const char* name, const Bool_t toSplit)
98 :AliEMCALClusterizer(headerFile, name, toSplit)
100 // ctor with the indication of the file where header Tree and digits Tree are stored
104 fDefaultInit = kFALSE ;
108 //____________________________________________________________________________
109 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
116 //____________________________________________________________________________
117 const TString AliEMCALClusterizerv1::BranchName() const
119 TString branchName(GetName() ) ;
120 branchName.Remove(branchName.Index(Version())-1) ;
124 //____________________________________________________________________________
125 Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp, Bool_t inpresho) const
126 {//To be replased later by the method, reading individual parameters from the database
128 if ( inpresho ) // calibrate as pre shower
129 return -fADCpedestalPreSho + amp * fADCchannelPreSho ;
131 else //calibrate as tower
132 return -fADCpedestalTower + amp * fADCchannelTower ;
135 //____________________________________________________________________________
136 void AliEMCALClusterizerv1::Exec(Option_t * option)
140 if( strcmp(GetName(), "")== 0 )
143 if(strstr(option,"tim"))
144 gBenchmark->Start("EMCALClusterizer");
146 if(strstr(option,"print"))
149 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
150 if(gime->BranchExists("RecPoints"))
152 Int_t nevents = gime->MaxEvent() ;
155 for(ievent = 0; ievent < nevents; ievent++){
157 gime->Event(ievent,"D") ;
160 GetCalibrationParameters() ;
162 fNumberOfTowerClusters = fNumberOfPreShoClusters = 0 ;
169 WriteRecPoints(ievent) ;
171 if(strstr(option,"deb"))
172 PrintRecPoints(option) ;
174 //increment the total number of digits per run
175 fRecPointsInRun += gime->TowerRecPoints()->GetEntriesFast() ;
176 fRecPointsInRun += gime->PreShowerRecPoints()->GetEntriesFast() ;
179 if(strstr(option,"tim")){
180 gBenchmark->Stop("EMCALClusterizer");
181 cout << "AliEMCALClusterizer:" << endl ;
182 cout << " took " << gBenchmark->GetCpuTime("EMCALClusterizer") << " seconds for Clusterizing "
183 << gBenchmark->GetCpuTime("EMCALClusterizer")/nevents << " seconds per event " << endl ;
189 //____________________________________________________________________________
190 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
191 Int_t nPar, Float_t * fitparameters) const
193 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
194 // The initial values for fitting procedure are set equal to the positions of local maxima.
195 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
197 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
198 TClonesArray * digits = gime->Digits() ;
201 gMinuit->mncler(); // Reset Minuit's list of paramters
202 gMinuit->SetPrintLevel(-1) ; // No Printout
203 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
204 // To set the address of the minimization function
206 TList * toMinuit = new TList();
207 toMinuit->AddAt(emcRP,0) ;
208 toMinuit->AddAt(digits,1) ;
210 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
212 // filling initial values for fit parameters
213 AliEMCALDigit * digit ;
217 Int_t nDigits = (Int_t) nPar / 3 ;
221 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
223 for(iDigit = 0; iDigit < nDigits; iDigit++){
224 digit = maxAt[iDigit];
229 geom->AbsToRelNumbering(digit->GetId(), relid) ;
230 geom->PosInAlice(relid, x, z) ;
232 Float_t energy = maxAtEnergy[iDigit] ;
234 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
237 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : x = " << x << endl ;
240 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
243 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : z = " << z << endl ;
246 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
249 cout << "EMCAL Unfolding> Unable to set initial value for fit procedure : energy = " << energy << endl ;
254 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
259 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
260 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
261 gMinuit->SetMaxIterations(5);
262 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
264 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
266 if(ierflg == 4){ // Minimum not found
267 cout << "EMCAL Unfolding> Fit not converged, cluster abandoned "<< endl ;
270 for(index = 0; index < nPar; index++){
273 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
274 fitparameters[index] = val ;
282 //____________________________________________________________________________
283 void AliEMCALClusterizerv1::GetCalibrationParameters()
285 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
286 const AliEMCALDigitizer * dig = gime->Digitizer(BranchName()) ;
288 fADCchannelTower = dig->GetTowerchannel() ;
289 fADCpedestalTower = dig->GetTowerpedestal();
291 fADCchannelPreSho = dig->GetPreShochannel() ;
292 fADCpedestalPreSho = dig->GetPreShopedestal() ;
296 //____________________________________________________________________________
297 void AliEMCALClusterizerv1::Init()
299 // Make all memory allocations which can not be done in default constructor.
300 // Attach the Clusterizer task to the list of EMCAL tasks
302 if ( strcmp(GetTitle(), "") == 0 )
303 SetTitle("galice.root") ;
305 TString branchname = GetName() ;
306 branchname.Remove(branchname.Index(Version())-1) ;
308 AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), branchname.Data(), fToSplit ) ;
310 cerr << "ERROR: AliEMCALClusterizerv1::Init -> Could not obtain the Getter object !" << endl ;
316 // construct the name of the file as /path/EMCAL.SDigits.root
317 //First - extract full path if necessary
318 TString fileName(GetTitle()) ;
319 Ssiz_t islash = fileName.Last('/') ;
320 if(islash<fileName.Length())
321 fileName.Remove(islash+1,fileName.Length()) ;
324 // Next - append the file name
325 fileName+="EMCAL.RecData." ;
326 if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
327 fileName+=branchname ;
331 // Finally - check if the file already opened or open the file
332 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
334 fSplitFile = TFile::Open(fileName.Data(),"update") ;
338 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
339 fNTowers = geom->GetNZ() * geom->GetNPhi() ;
342 gMinuit = new TMinuit(100) ;
344 gime->PostClusterizer(this) ;
345 gime->PostRecPoints(branchname ) ;
349 //____________________________________________________________________________
350 void AliEMCALClusterizerv1::InitParameters()
352 fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ;
356 fPreShoClusteringThreshold = 0.0001;
357 fTowerClusteringThreshold = 0.2;
359 fTowerLocMaxCut = 0.03 ;
360 fPreShoLocMaxCut = 0.03 ;
369 TString clusterizerName( GetName()) ;
370 if (clusterizerName.IsNull() )
371 clusterizerName = "Default" ;
372 clusterizerName.Append(":") ;
373 clusterizerName.Append(Version()) ;
374 SetName(clusterizerName) ;
375 fRecPointsInRun = 0 ;
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 AliEMCALGetter *gime = AliEMCALGetter::GetInstance() ;
454 TObjArray * towerRecPoints = gime->TowerRecPoints() ;
455 TObjArray * preshoRecPoints = gime->PreShowerRecPoints() ;
456 TClonesArray * digits = gime->Digits() ;
463 TString name("TreeR") ;
465 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
468 treeR = gAlice->TreeR();
472 gAlice->MakeTree("R", fSplitFile);
473 treeR = gAlice->TreeR() ;
477 //Evaluate position, dispersion and other RecPoint properties...
478 for(index = 0; index < towerRecPoints->GetEntries(); index++)
479 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->EvalAll(fW0,digits) ;
481 towerRecPoints->Sort() ;
483 for(index = 0; index < towerRecPoints->GetEntries(); index++)
484 (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->SetIndexInList(index) ;
486 towerRecPoints->Expand(towerRecPoints->GetEntriesFast()) ;
488 //Now the same for pre shower
489 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
490 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->EvalAll(fW0CPV,digits) ;
492 preshoRecPoints->Sort() ;
494 for(index = 0; index < preshoRecPoints->GetEntries(); index++)
495 (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->SetIndexInList(index) ;
497 preshoRecPoints->Expand(preshoRecPoints->GetEntriesFast()) ;
499 Int_t bufferSize = 32000 ;
500 Int_t splitlevel = 0 ;
503 TBranch * towerBranch = treeR->Branch("EMCALTowerRP","TObjArray",&towerRecPoints,bufferSize,splitlevel);
504 towerBranch->SetTitle(BranchName());
506 //Now Pre Shower branch
507 TBranch * preshoBranch = treeR->Branch("EMCALPreShoRP","TObjArray",&preshoRecPoints,bufferSize,splitlevel);
508 preshoBranch->SetTitle(BranchName());
510 //And Finally clusterizer branch
511 AliEMCALClusterizerv1 * cl = (AliEMCALClusterizerv1*)gime->Clusterizer(BranchName()) ;
512 TBranch * clusterizerBranch = treeR->Branch("AliEMCALClusterizer","AliEMCALClusterizerv1",
513 &cl,bufferSize,splitlevel);
514 clusterizerBranch->SetTitle(BranchName());
516 towerBranch ->Fill() ;
517 preshoBranch ->Fill() ;
518 clusterizerBranch->Fill() ;
520 treeR->AutoSave() ; //Write(0,kOverwrite) ;
521 if(gAlice->TreeR()!=treeR)
525 //____________________________________________________________________________
526 void AliEMCALClusterizerv1::MakeClusters()
528 // Steering method to construct the clusters stored in a list of Reconstructed Points
529 // A cluster is defined as a list of neighbour digits
531 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
533 TObjArray * towerRecPoints = gime->TowerRecPoints(BranchName()) ;
534 TObjArray * preshoRecPoints = gime->PreShowerRecPoints(BranchName()) ;
535 towerRecPoints->Delete() ;
536 preshoRecPoints->Delete() ;
538 TClonesArray * digits = gime->Digits() ;
540 cerr << "ERROR: AliEMCALClusterizerv1::MakeClusters -> Digits with name "
541 << GetName() << " not found ! " << endl ;
544 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
547 // Clusterization starts
549 TIter nextdigit(digitsC) ;
550 AliEMCALDigit * digit ;
551 Bool_t notremoved = kTRUE ;
553 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
554 AliEMCALRecPoint * clu = 0 ;
556 TArrayI clusterdigitslist(1500) ;
559 if (( IsInTower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fTowerClusteringThreshold ) ||
560 ( IsInPreShower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fPreShoClusteringThreshold ) ) {
562 Int_t iDigitInCluster = 0 ;
564 if ( IsInTower(digit) ) {
565 // start a new Tower RecPoint
566 if(fNumberOfTowerClusters >= towerRecPoints->GetSize())
567 towerRecPoints->Expand(2*fNumberOfTowerClusters+1) ;
569 towerRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfTowerClusters) ;
570 clu = dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(fNumberOfTowerClusters)) ;
571 fNumberOfTowerClusters++ ;
572 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower())) ;
573 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
575 digitsC->Remove(digit) ;
579 // start a new Pre Shower cluster
580 if(fNumberOfPreShoClusters >= preshoRecPoints->GetSize())
581 preshoRecPoints->Expand(2*fNumberOfPreShoClusters+1);
583 preshoRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfPreShoClusters) ;
585 clu = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(fNumberOfPreShoClusters)) ;
586 fNumberOfPreShoClusters++ ;
587 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower() ) );
588 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
590 digitsC->Remove(digit) ;
593 // Here we remove remaining Tower digits, which cannot make a cluster
596 while( ( digit = dynamic_cast<AliEMCALDigit *>(nextdigit()) ) ) {
597 if( IsInTower(digit) )
598 digitsC->Remove(digit) ;
602 notremoved = kFALSE ;
609 AliEMCALDigit * digitN ;
611 while (index < iDigitInCluster){ // scan over digits already in cluster
612 digit = (AliEMCALDigit*)digits->At(clusterdigitslist[index]) ;
614 while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits
615 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
617 case 0 : // not a neighbour
619 case 1 : // are neighbours
620 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->IsInPreShower() ) ) ;
621 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
623 digitsC->Remove(digitN) ;
625 case 2 : // too far from each other
634 } // loop over cluster
645 //____________________________________________________________________________
646 void AliEMCALClusterizerv1::MakeUnfolding()
648 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
650 // // Unfolds clusters using the shape of an ElectroMagnetic shower
651 // // Performs unfolding of all EMC/CPV clusters
653 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
655 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
656 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
657 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
658 // TClonesArray * digits = gime->Digits() ;
660 // // Unfold first EMC clusters
661 // if(fNumberOfTowerClusters > 0){
663 // Int_t nModulesToUnfold = geom->GetNModules() ;
665 // Int_t numberofNotUnfolded = fNumberOfTowerClusters ;
667 // for(index = 0 ; index < numberofNotUnfolded ; index++){
669 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint *) emcRecPoints->At(index) ;
670 // if(emcRecPoint->GetEMCALMod()> nModulesToUnfold)
673 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
674 // Int_t * maxAt = new Int_t[nMultipl] ;
675 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
676 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fTowerLocMaxCut,digits) ;
678 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
679 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
680 // emcRecPoints->Remove(emcRecPoint);
681 // emcRecPoints->Compress() ;
683 // fNumberOfTowerClusters -- ;
684 // numberofNotUnfolded-- ;
688 // delete[] maxAtEnergy ;
691 // // Unfolding of EMC clusters finished
694 // // Unfold now CPV clusters
695 // if(fNumberOfPreShoClusters > 0){
697 // Int_t nModulesToUnfold = geom->GetNModules() ;
699 // Int_t numberofPreShoNotUnfolded = fNumberOfPreShoClusters ;
701 // for(index = 0 ; index < numberofPreShoNotUnfolded ; index++){
703 // AliEMCALRecPoint * recPoint = (AliEMCALRecPoint *) cpvRecPoints->At(index) ;
705 // if(recPoint->GetEMCALMod()> nModulesToUnfold)
708 // AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint*) recPoint ;
710 // Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
711 // Int_t * maxAt = new Int_t[nMultipl] ;
712 // Float_t * maxAtEnergy = new Float_t[nMultipl] ;
713 // Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fPreShoLocMaxCut,digits) ;
715 // if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
716 // UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
717 // cpvRecPoints->Remove(emcRecPoint);
718 // cpvRecPoints->Compress() ;
720 // numberofPreShoNotUnfolded-- ;
721 // fNumberOfPreShoClusters-- ;
725 // delete[] maxAtEnergy ;
728 // //Unfolding of PreSho clusters finished
732 //____________________________________________________________________________
733 Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
735 // Shape of the shower (see EMCAL TDR)
736 // If you change this function, change also the gradient evaluation in ChiSquare()
738 Double_t r4 = r*r*r*r ;
739 Double_t r295 = TMath::Power(r, 2.95) ;
740 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
744 //____________________________________________________________________________
745 void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * iniTower,
747 AliEMCALDigit ** maxAt,
748 Float_t * maxAtEnergy)
750 // Performs the unfolding of a cluster with nMax overlapping showers
752 Fatal("AliEMCALClusterizerv1::UnfoldCluster", "--> Unfolding not implemented") ;
754 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
755 // const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
756 // const TClonesArray * digits = gime->Digits() ;
757 // TObjArray * emcRecPoints = gime->TowerRecPoints() ;
758 // TObjArray * cpvRecPoints = gime->PreShoRecPoints() ;
760 // Int_t nPar = 3 * nMax ;
761 // Float_t * fitparameters = new Float_t[nPar] ;
763 // Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
765 // // Fit failed, return and remove cluster
766 // delete[] fitparameters ;
770 // // create ufolded rec points and fill them with new energy lists
771 // // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
772 // // and later correct this number in acordance with actual energy deposition
774 // Int_t nDigits = iniTower->GetMultiplicity() ;
775 // Float_t * efit = new Float_t[nDigits] ;
776 // Float_t xDigit=0.,zDigit=0.,distance=0. ;
777 // Float_t xpar=0.,zpar=0.,epar=0. ;
779 // AliEMCALDigit * digit = 0 ;
780 // Int_t * emcDigits = iniTower->GetDigitsList() ;
784 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
785 // digit = (AliEMCALDigit*) digits->At(emcDigits[iDigit] ) ;
786 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
787 // geom->RelPosInModule(relid, xDigit, zDigit) ;
791 // while(iparam < nPar ){
792 // xpar = fitparameters[iparam] ;
793 // zpar = fitparameters[iparam+1] ;
794 // epar = fitparameters[iparam+2] ;
796 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
797 // distance = TMath::Sqrt(distance) ;
798 // efit[iDigit] += epar * ShowerShape(distance) ;
803 // // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
804 // // so that energy deposited in each cell is distributed betwin new clusters proportionally
805 // // to its contribution to efit
807 // Float_t * emcEnergies = iniTower->GetEnergiesList() ;
811 // while(iparam < nPar ){
812 // xpar = fitparameters[iparam] ;
813 // zpar = fitparameters[iparam+1] ;
814 // epar = fitparameters[iparam+2] ;
817 // AliEMCALTowerRecPoint * emcRP = 0 ;
819 // if(iniTower->IsTower()){ //create new entries in fTowerRecPoints...
821 // if(fNumberOfTowerClusters >= emcRecPoints->GetSize())
822 // emcRecPoints->Expand(2*fNumberOfTowerClusters) ;
824 // (*emcRecPoints)[fNumberOfTowerClusters] = new AliEMCALTowerRecPoint("") ;
825 // emcRP = (AliEMCALTowerRecPoint *) emcRecPoints->At(fNumberOfTowerClusters);
826 // fNumberOfTowerClusters++ ;
828 // else{//create new entries in fPreShoRecPoints
829 // if(fNumberOfPreShoClusters >= cpvRecPoints->GetSize())
830 // cpvRecPoints->Expand(2*fNumberOfPreShoClusters) ;
832 // (*cpvRecPoints)[fNumberOfPreShoClusters] = new AliEMCALPreShoRecPoint("") ;
833 // emcRP = (AliEMCALTowerRecPoint *) cpvRecPoints->At(fNumberOfPreShoClusters);
834 // fNumberOfPreShoClusters++ ;
838 // for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
839 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
840 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
841 // geom->RelPosInModule(relid, xDigit, zDigit) ;
842 // distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
843 // distance = TMath::Sqrt(distance) ;
844 // ratio = epar * ShowerShape(distance) / efit[iDigit] ;
845 // eDigit = emcEnergies[iDigit] * ratio ;
846 // emcRP->AddDigit( *digit, eDigit ) ;
850 // delete[] fitparameters ;
855 //_____________________________________________________________________________
856 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
858 // Calculates the Chi square for the cluster unfolding minimization
859 // Number of parameters, Gradient, Chi squared, parameters, what to do
862 // Fatal("AliEMCALClusterizerv1::UnfoldingChiSquare","-->Unfolding not implemented") ;
864 // TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
866 // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint*) toMinuit->At(0) ;
867 // TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
871 // // AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint *) gMinuit->GetObjectFit() ; // TowerRecPoint to fit
873 // Int_t * emcDigits = emcRP->GetDigitsList() ;
875 // Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
877 // Float_t * emcEnergies = emcRP->GetEnergiesList() ;
879 // const AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
884 // for(iparam = 0 ; iparam < nPar ; iparam++)
885 // Grad[iparam] = 0 ; // Will evaluate gradient
889 // AliEMCALDigit * digit ;
892 // for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
894 // digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ;
900 // geom->AbsToRelNumbering(digit->GetId(), relid) ;
902 // geom->RelPosInModule(relid, xDigit, zDigit) ;
904 // if(iflag == 2){ // calculate gradient
905 // Int_t iParam = 0 ;
907 // while(iParam < nPar ){
908 // Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
910 // distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
911 // distance = TMath::Sqrt( distance ) ;
913 // efit += x[iParam] * ShowerShape(distance) ;
916 // Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
918 // while(iParam < nPar ){
919 // Double_t xpar = x[iParam] ;
920 // Double_t zpar = x[iParam+1] ;
921 // Double_t epar = x[iParam+2] ;
922 // Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
923 // Double_t shape = sum * ShowerShape(dr) ;
924 // Double_t r4 = dr*dr*dr*dr ;
925 // Double_t r295 = TMath::Power(dr,2.95) ;
926 // Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
927 // 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
929 // Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
931 // Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
933 // Grad[iParam] += shape ; // Derivative over energy
940 // while(iparam < nPar ){
941 // Double_t xpar = x[iparam] ;
942 // Double_t zpar = x[iparam+1] ;
943 // Double_t epar = x[iparam+2] ;
945 // Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
946 // distance = TMath::Sqrt(distance) ;
947 // efit += epar * ShowerShape(distance) ;
950 // fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
951 // // Here we assume, that sigma = sqrt(E)
956 //____________________________________________________________________________
957 void AliEMCALClusterizerv1::Print(Option_t * option)const
959 // Print clusterizer parameters
961 if( strcmp(GetName(), "") !=0 ){
965 TString taskName(GetName()) ;
966 taskName.ReplaceAll(Version(), "") ;
968 cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl
969 << "Clusterizing digits from the file: " << taskName.Data() << endl
970 << " Branch: " << GetName() << endl
972 << " EMC Clustering threshold = " << fTowerClusteringThreshold << endl
973 << " EMC Local Maximum cut = " << fTowerLocMaxCut << endl
974 << " EMC Logarothmic weight = " << fW0 << endl
976 << " CPV Clustering threshold = " << fPreShoClusteringThreshold << endl
977 << " CPV Local Maximum cut = " << fPreShoLocMaxCut << endl
978 << " CPV Logarothmic weight = " << fW0CPV << endl
981 cout << " Unfolding on " << endl ;
983 cout << " Unfolding off " << endl ;
985 cout << "------------------------------------------------------------------" <<endl ;
988 cout << " AliEMCALClusterizerv1 not initialized " << endl ;
990 //____________________________________________________________________________
991 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
993 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
995 TObjArray * towerRecPoints = AliEMCALGetter::GetInstance()->TowerRecPoints() ;
996 TObjArray * preshoRecPoints = AliEMCALGetter::GetInstance()->PreShowerRecPoints() ;
998 cout << "AliEMCALClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
999 cout << " Found "<< towerRecPoints->GetEntriesFast() << " TOWER Rec Points and "
1000 << preshoRecPoints->GetEntriesFast() << " PRE SHOWER RecPoints" << endl ;
1002 fRecPointsInRun += towerRecPoints->GetEntriesFast() ;
1003 fRecPointsInRun += preshoRecPoints->GetEntriesFast() ;
1005 if(strstr(option,"all")) {
1007 cout << "Tower clusters " << endl ;
1008 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
1011 for (index = 0 ; index < towerRecPoints->GetEntries() ; index++) {
1012 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint * >(towerRecPoints->At(index)) ;
1014 rp->GetGlobalPosition(globalpos);
1016 rp->GetElipsAxis(lambda);
1019 primaries = rp->GetPrimaries(nprimaries);
1021 cout << setw(4) << rp->GetIndexInList() << " "
1022 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
1023 << setw(3) << rp->GetMultiplicity() << " "
1024 << setw(1) << rp->GetEMCALArm() << " "
1025 << setw(5) << setprecision(4) << globalpos.X() << " "
1026 << setw(5) << setprecision(4) << globalpos.Y() << " "
1027 << setw(5) << setprecision(4) << globalpos.Z() << " "
1028 << setw(4) << setprecision(2) << lambda[0] << " "
1029 << setw(4) << setprecision(2) << lambda[1] << " "
1030 << setw(2) << nprimaries << " " ;
1032 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1033 cout << setw(4) << primaries[iprimary] << " " ;
1037 //Now plot Pre shower recPoints
1039 cout << "-----------------------------------------------------------------------"<<endl ;
1041 cout << "PreShower clusters " << endl ;
1042 cout << " Index Ene(MeV) Multi Module phi r theta Lambda 1 Lambda 2 # of prim Primaries list " << endl;
1044 for (index = 0 ; index < preshoRecPoints->GetEntries() ; index++) {
1045 AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(index)) ;
1047 rp->GetGlobalPosition(globalpos);
1049 rp->GetElipsAxis(lambda);
1052 primaries = rp->GetPrimaries(nprimaries);
1054 cout << setw(4) << rp->GetIndexInList() << " "
1055 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
1056 << setw(3) << rp->GetMultiplicity() << " "
1057 << setw(1) << rp->GetEMCALArm() << " "
1058 << setw(5) << setprecision(4) << globalpos.X() << " "
1059 << setw(5) << setprecision(4) << globalpos.Y() << " "
1060 << setw(5) << setprecision(4) << globalpos.Z() << " "
1061 << setw(4) << setprecision(2) << lambda[0] << " "
1062 << setw(4) << setprecision(2) << lambda[1] << " "
1063 << setw(2) << nprimaries << " " ;
1065 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1066 cout << setw(4) << primaries[iprimary] << " " ;
1070 cout << "-----------------------------------------------------------------------"<<endl ;