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 **************************************************************************/
19 1 October 2000. Yuri Kharlov:
21 PPSD upper layer is considered if number of layers>1
23 18 October 2000. Yuri Kharlov:
24 AliPHOSClusterizerv1()
25 CPV clusterizing parameters added
28 After first PPSD digit remove EMC digits only once
30 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
31 //////////////////////////////////////////////////////////////////////////////
32 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
33 // unfolds the clusters having several local maxima.
34 // Results are stored in TreeR#, branches PHOSEmcRP (EMC recPoints),
35 // PHOSCpvRP (CPV RecPoints) and AliPHOSClusterizer (Clusterizer with all
36 // parameters including input digits branch title, thresholds etc.)
37 // This TTask is normally called from Reconstructioner, but can as well be used in
40 // root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root", "recpointsname", "digitsname")
41 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 // // reads gAlice from header file "galice.root", uses digits stored in the branch names "digitsname" (default = "Default")
43 // // and saves recpoints in branch named "recpointsname" (default = "digitsname")
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->SetEmcLocalMaxCut(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 ---
71 // --- AliRoot header files ---
72 #include "AliPHOSCalibrationDB.h"
73 #include "AliPHOSClusterizerv1.h"
74 #include "AliPHOSCpvRecPoint.h"
75 #include "AliPHOSDigit.h"
76 #include "AliPHOSDigitizer.h"
77 #include "AliPHOSEmcRecPoint.h"
79 #include "AliPHOSGetter.h"
82 ClassImp(AliPHOSClusterizerv1)
84 //____________________________________________________________________________
85 AliPHOSClusterizerv1::AliPHOSClusterizerv1() : AliPHOSClusterizer()
87 // default ctor (to be used mainly by Streamer)
90 fDefaultInit = kTRUE ;
93 //____________________________________________________________________________
94 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const char* headerFile,const char* name, const Bool_t toSplit)
95 :AliPHOSClusterizer(headerFile, name, toSplit)
97 // ctor with the indication of the file where header Tree and digits Tree are stored
101 fDefaultInit = kFALSE ;
105 //____________________________________________________________________________
106 AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
113 //____________________________________________________________________________
114 const TString AliPHOSClusterizerv1::BranchName() const
116 TString branchName(GetName() ) ;
117 branchName.Remove(branchName.Index(Version())-1) ;
121 //____________________________________________________________________________
122 Float_t AliPHOSClusterizerv1::Calibrate(Int_t amp, Int_t absId) const
123 { //To be replased later by the method, reading individual parameters from the database
125 return fCalibrationDB->Calibrate(amp,absId) ;
127 if(absId <= fEmcCrystals) //calibrate as EMC
128 return fADCpedestalEmc + amp*fADCchanelEmc ;
129 else //calibrate as CPV
130 return fADCpedestalCpv+ amp*fADCchanelCpv ;
134 //____________________________________________________________________________
135 void AliPHOSClusterizerv1::Exec(Option_t * option)
139 if( strcmp(GetName(), "")== 0 )
142 if(strstr(option,"tim"))
143 gBenchmark->Start("PHOSClusterizer");
145 if(strstr(option,"print"))
148 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
149 if(gime->BranchExists("RecPoints"))
151 Int_t nevents = gime->MaxEvent() ;
154 for(ievent = 0; ievent < nevents; ievent++){
156 gime->Event(ievent,"D") ;
159 GetCalibrationParameters() ;
161 fNumberOfEmcClusters = fNumberOfCpvClusters = 0 ;
168 WriteRecPoints(ievent) ;
170 if(strstr(option,"deb"))
171 PrintRecPoints(option) ;
173 //increment the total number of digits per run
174 fRecPointsInRun += gime->EmcRecPoints()->GetEntriesFast() ;
175 fRecPointsInRun += gime->CpvRecPoints()->GetEntriesFast() ;
178 if(strstr(option,"tim")){
179 gBenchmark->Stop("PHOSClusterizer");
180 Info("Exec", " took %f seconds for Clusterizing %f seconds per event \n",
181 gBenchmark->GetCpuTime("PHOSClusterizer"),
182 gBenchmark->GetCpuTime("PHOSClusterizer")/nevents ) ;
186 //____________________________________________________________________________
187 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
188 Int_t nPar, Float_t * fitparameters) const
190 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
191 // The initial values for fitting procedure are set equal to the positions of local maxima.
192 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
194 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
195 TClonesArray * digits = gime->Digits() ;
198 gMinuit->mncler(); // Reset Minuit's list of paramters
199 gMinuit->SetPrintLevel(-1) ; // No Printout
200 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
201 // To set the address of the minimization function
203 TList * toMinuit = new TList();
204 toMinuit->AddAt(emcRP,0) ;
205 toMinuit->AddAt(digits,1) ;
207 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
209 // filling initial values for fit parameters
210 AliPHOSDigit * digit ;
214 Int_t nDigits = (Int_t) nPar / 3 ;
218 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
220 for(iDigit = 0; iDigit < nDigits; iDigit++){
221 digit = maxAt[iDigit];
226 geom->AbsToRelNumbering(digit->GetId(), relid) ;
227 geom->RelPosInModule(relid, x, z) ;
229 Float_t energy = maxAtEnergy[iDigit] ;
231 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
234 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : x = %f\n", x ) ;
237 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
240 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : z =%f\n", z ) ;
243 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
246 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : energy = %f\n", energy ) ;
251 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
256 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
257 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
258 gMinuit->SetMaxIterations(5);
259 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
261 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
263 if(ierflg == 4){ // Minimum not found
264 Warning("FindFit", "PHOS Unfolding fit not converged, cluster abandoned\n" );
267 for(index = 0; index < nPar; index++){
270 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
271 fitparameters[index] = val ;
279 //____________________________________________________________________________
280 void AliPHOSClusterizerv1::GetCalibrationParameters()
282 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
284 const TTask * task = gime->Digitizer(BranchName()) ;
285 if(strcmp(task->IsA()->GetName(),"AliPHOSDigitizer")==0){
286 const AliPHOSDigitizer * dig = static_cast<const AliPHOSDigitizer *>(task) ;
288 fADCchanelEmc = dig->GetEMCchannel() ;
289 fADCpedestalEmc = dig->GetEMCpedestal();
291 fADCchanelCpv = dig->GetCPVchannel() ;
292 fADCpedestalCpv = dig->GetCPVpedestal() ;
295 fCalibrationDB = gime->CalibrationDB();
299 //____________________________________________________________________________
300 void AliPHOSClusterizerv1::Init()
302 // Make all memory allocations which can not be done in default constructor.
303 // Attach the Clusterizer task to the list of PHOS tasks
305 if ( strcmp(GetTitle(), "") == 0 )
306 SetTitle("galice.root") ;
308 TString branchname = GetName() ;
309 branchname.Remove(branchname.Index(Version())-1) ;
311 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(),branchname.Data(), fToSplit ) ;
313 Error("Init", "Could not obtain the Getter object !") ;
319 // construct the name of the file as /path/EMCAL.SDigits.root
320 //First - extract full path if necessary
321 TString fileName(GetTitle()) ;
322 Ssiz_t islash = fileName.Last('/') ;
323 if(islash<fileName.Length())
324 fileName.Remove(islash+1,fileName.Length()) ;
327 // Next - append the file name
328 fileName+="PHOS.RecData." ;
329 if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
330 fileName+=branchname ;
334 // Finally - check if the file already opened or open the file
335 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
337 fSplitFile = TFile::Open(fileName.Data(),"update") ;
342 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
343 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
346 gMinuit = new TMinuit(100) ;
348 gime->PostClusterizer(this) ;
349 gime->PostRecPoints(branchname) ;
353 //____________________________________________________________________________
354 void AliPHOSClusterizerv1::InitParameters()
357 fNumberOfCpvClusters = 0 ;
358 fNumberOfEmcClusters = 0 ;
360 fCpvClusteringThreshold = 0.0;
361 fEmcClusteringThreshold = 0.2;
363 fEmcLocMaxCut = 0.03 ;
364 fCpvLocMaxCut = 0.03 ;
369 fEmcTimeGate = 1.e-8 ;
373 TString clusterizerName( GetName()) ;
374 if (clusterizerName.IsNull() )
375 clusterizerName = "Default" ;
376 clusterizerName.Append(":") ;
377 clusterizerName.Append(Version()) ;
378 SetName(clusterizerName) ;
379 fRecPointsInRun = 0 ;
385 //____________________________________________________________________________
386 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
388 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
390 // = 2 are not neighbour but do not continue searching
391 // neighbours are defined as digits having at least a common vertex
392 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
393 // which is compared to a digit (d2) not yet in a cluster
395 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
400 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
403 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
405 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module
406 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
407 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
409 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
410 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
414 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
415 rv = 2; // Difference in row numbers is too large to look further
421 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
430 //____________________________________________________________________________
431 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
433 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
436 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
438 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
440 if(digit->GetId() <= nEMC ) rv = kTRUE;
445 //____________________________________________________________________________
446 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
448 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
451 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
453 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
455 if(digit->GetId() > nEMC ) rv = kTRUE;
460 //____________________________________________________________________________
461 void AliPHOSClusterizerv1::WriteRecPoints(Int_t event)
464 // Creates new branches with given title
465 // fills and writes into TreeR.
468 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
469 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
470 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
471 TClonesArray * digits = gime->Digits() ;
478 TString name("TreeR") ;
480 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
483 treeR = gAlice->TreeR();
488 gAlice->MakeTree("R", fSplitFile);
490 gAlice->MakeTree("R", gROOT->GetFile(GetTitle()));
491 treeR = gAlice->TreeR() ;
495 //Evaluate position, dispersion and other RecPoint properties...
496 for(index = 0; index < emcRecPoints->GetEntries(); index++)
497 ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->EvalAll(fW0,digits) ;
499 emcRecPoints->Sort() ;
500 for(index = 0; index < emcRecPoints->GetEntries(); index++)
501 ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->SetIndexInList(index) ;
503 emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ;
505 //Now the same for CPV
506 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
507 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->EvalAll(fW0CPV,digits) ;
509 cpvRecPoints->Sort() ;
511 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
512 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->SetIndexInList(index) ;
514 cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
516 Int_t bufferSize = 32000 ;
517 Int_t splitlevel = 0 ;
520 TBranch * emcBranch = treeR->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
521 emcBranch->SetTitle(BranchName());
524 TBranch * cpvBranch = treeR->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
525 cpvBranch->SetTitle(BranchName());
527 //And Finally clusterizer branch
528 AliPHOSClusterizerv1 * cl = this ;
529 TBranch * clusterizerBranch = treeR->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
530 &cl,bufferSize,splitlevel);
531 clusterizerBranch->SetTitle(BranchName());
535 clusterizerBranch->Fill() ;
538 if(gAlice->TreeR()!=treeR)
542 //____________________________________________________________________________
543 void AliPHOSClusterizerv1::MakeClusters()
545 // Steering method to construct the clusters stored in a list of Reconstructed Points
546 // A cluster is defined as a list of neighbour digits
548 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
550 TObjArray * emcRecPoints = gime->EmcRecPoints(BranchName()) ;
551 TObjArray * cpvRecPoints = gime->CpvRecPoints(BranchName()) ;
552 emcRecPoints->Delete() ;
553 cpvRecPoints->Delete() ;
555 TClonesArray * digits = gime->Digits() ;
557 Fatal("MakeClusters", "Digits with name %s not found !", BranchName().Data() ) ;
559 TClonesArray * digitsC = (TClonesArray*)digits->Clone() ;
562 // Clusterization starts
564 TIter nextdigit(digitsC) ;
565 AliPHOSDigit * digit ;
566 Bool_t notremoved = kTRUE ;
568 while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digitsC
569 AliPHOSRecPoint * clu = 0 ;
571 TArrayI clusterdigitslist(1500) ;
574 if (( IsInEmc (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fEmcClusteringThreshold ) ||
575 ( IsInCpv (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fCpvClusteringThreshold ) ) {
576 Int_t iDigitInCluster = 0 ;
578 if ( IsInEmc(digit) ) {
579 // start a new EMC RecPoint
580 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
581 emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
583 emcRecPoints->AddAt(new AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
584 clu = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters) ;
585 fNumberOfEmcClusters++ ;
586 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId())) ;
587 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
589 digitsC->Remove(digit) ;
593 // start a new CPV cluster
594 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
595 cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
597 cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
599 clu = (AliPHOSCpvRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters) ;
600 fNumberOfCpvClusters++ ;
601 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId()) ) ;
602 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
604 digitsC->Remove(digit) ;
607 // Here we remove remaining EMC digits, which cannot make a cluster
610 while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
612 digitsC->Remove(digit) ;
616 notremoved = kFALSE ;
623 AliPHOSDigit * digitN ;
625 while (index < iDigitInCluster){ // scan over digits already in cluster
626 digit = (AliPHOSDigit*)digits->At(clusterdigitslist[index]) ;
628 while ( (digitN = (AliPHOSDigit *)nextdigit()) ) { // scan over the reduced list of digits
629 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
631 case 0 : // not a neighbour
633 case 1 : // are neighbours
634 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->GetId() ) ) ;
635 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
637 digitsC->Remove(digitN) ;
639 case 2 : // too far from each other
648 } // loop over cluster
659 //____________________________________________________________________________
660 void AliPHOSClusterizerv1::MakeUnfolding()
662 // Unfolds clusters using the shape of an ElectroMagnetic shower
663 // Performs unfolding of all EMC/CPV clusters
665 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
668 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
669 TObjArray * emcRecPoints = gime->EmcRecPoints(BranchName()) ;
670 TObjArray * cpvRecPoints = gime->CpvRecPoints(BranchName()) ;
671 TClonesArray * digits = gime->Digits() ;
673 // Unfold first EMC clusters
674 if(fNumberOfEmcClusters > 0){
676 Int_t nModulesToUnfold = geom->GetNModules() ;
678 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
680 for(index = 0 ; index < numberofNotUnfolded ; index++){
682 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) emcRecPoints->At(index) ;
683 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
686 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
687 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
688 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
689 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
691 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
692 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
693 emcRecPoints->Remove(emcRecPoint);
694 emcRecPoints->Compress() ;
696 fNumberOfEmcClusters -- ;
697 numberofNotUnfolded-- ;
701 delete[] maxAtEnergy ;
704 // Unfolding of EMC clusters finished
707 // Unfold now CPV clusters
708 if(fNumberOfCpvClusters > 0){
710 Int_t nModulesToUnfold = geom->GetNModules() ;
712 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
714 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
716 AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) cpvRecPoints->At(index) ;
718 if(recPoint->GetPHOSMod()> nModulesToUnfold)
721 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ;
723 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
724 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
725 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
726 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,digits) ;
728 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
729 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
730 cpvRecPoints->Remove(emcRecPoint);
731 cpvRecPoints->Compress() ;
733 numberofCpvNotUnfolded-- ;
734 fNumberOfCpvClusters-- ;
738 delete[] maxAtEnergy ;
741 //Unfolding of Cpv clusters finished
745 //____________________________________________________________________________
746 Double_t AliPHOSClusterizerv1::ShowerShape(Double_t r)
748 // Shape of the shower (see PHOS TDR)
749 // If you change this function, change also the gradient evaluation in ChiSquare()
751 Double_t r4 = r*r*r*r ;
752 Double_t r295 = TMath::Power(r, 2.95) ;
753 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
757 //____________________________________________________________________________
758 void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
760 AliPHOSDigit ** maxAt,
761 Float_t * maxAtEnergy)
763 // Performs the unfolding of a cluster with nMax overlapping showers
765 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
767 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
768 const TClonesArray * digits = gime->Digits() ;
769 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
770 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
772 Int_t nPar = 3 * nMax ;
773 Float_t * fitparameters = new Float_t[nPar] ;
775 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
777 // Fit failed, return and remove cluster
778 delete[] fitparameters ;
782 // create ufolded rec points and fill them with new energy lists
783 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
784 // and later correct this number in acordance with actual energy deposition
786 Int_t nDigits = iniEmc->GetMultiplicity() ;
787 Float_t * efit = new Float_t[nDigits] ;
788 Float_t xDigit=0.,zDigit=0.,distance=0. ;
789 Float_t xpar=0.,zpar=0.,epar=0. ;
791 AliPHOSDigit * digit = 0 ;
792 Int_t * emcDigits = iniEmc->GetDigitsList() ;
796 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
797 digit = (AliPHOSDigit*) digits->At(emcDigits[iDigit] ) ;
798 geom->AbsToRelNumbering(digit->GetId(), relid) ;
799 geom->RelPosInModule(relid, xDigit, zDigit) ;
803 while(iparam < nPar ){
804 xpar = fitparameters[iparam] ;
805 zpar = fitparameters[iparam+1] ;
806 epar = fitparameters[iparam+2] ;
808 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
809 distance = TMath::Sqrt(distance) ;
810 efit[iDigit] += epar * ShowerShape(distance) ;
815 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
816 // so that energy deposited in each cell is distributed betwin new clusters proportionally
817 // to its contribution to efit
819 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
823 while(iparam < nPar ){
824 xpar = fitparameters[iparam] ;
825 zpar = fitparameters[iparam+1] ;
826 epar = fitparameters[iparam+2] ;
829 AliPHOSEmcRecPoint * emcRP = 0 ;
831 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
833 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
834 emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
836 (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
837 emcRP = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters);
838 fNumberOfEmcClusters++ ;
840 else{//create new entries in fCpvRecPoints
841 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
842 cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
844 (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
845 emcRP = (AliPHOSEmcRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters);
846 fNumberOfCpvClusters++ ;
850 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
851 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
852 geom->AbsToRelNumbering(digit->GetId(), relid) ;
853 geom->RelPosInModule(relid, xDigit, zDigit) ;
854 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
855 distance = TMath::Sqrt(distance) ;
856 ratio = epar * ShowerShape(distance) / efit[iDigit] ;
857 eDigit = emcEnergies[iDigit] * ratio ;
858 emcRP->AddDigit( *digit, eDigit ) ;
862 delete[] fitparameters ;
867 //_____________________________________________________________________________
868 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
870 // Calculates the Chi square for the cluster unfolding minimization
871 // Number of parameters, Gradient, Chi squared, parameters, what to do
873 TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
875 AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint*) toMinuit->At(0) ;
876 TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
880 // AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint *) gMinuit->GetObjectFit() ; // EmcRecPoint to fit
882 Int_t * emcDigits = emcRP->GetDigitsList() ;
884 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
886 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
888 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
893 for(iparam = 0 ; iparam < nPar ; iparam++)
894 Grad[iparam] = 0 ; // Will evaluate gradient
898 AliPHOSDigit * digit ;
901 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
903 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
909 geom->AbsToRelNumbering(digit->GetId(), relid) ;
911 geom->RelPosInModule(relid, xDigit, zDigit) ;
913 if(iflag == 2){ // calculate gradient
916 while(iParam < nPar ){
917 Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
919 distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
920 distance = TMath::Sqrt( distance ) ;
922 efit += x[iParam] * ShowerShape(distance) ;
925 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
927 while(iParam < nPar ){
928 Double_t xpar = x[iParam] ;
929 Double_t zpar = x[iParam+1] ;
930 Double_t epar = x[iParam+2] ;
931 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
932 Double_t shape = sum * ShowerShape(dr) ;
933 Double_t r4 = dr*dr*dr*dr ;
934 Double_t r295 = TMath::Power(dr,2.95) ;
935 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
936 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
938 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
940 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
942 Grad[iParam] += shape ; // Derivative over energy
949 while(iparam < nPar ){
950 Double_t xpar = x[iparam] ;
951 Double_t zpar = x[iparam+1] ;
952 Double_t epar = x[iparam+2] ;
954 Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
955 distance = TMath::Sqrt(distance) ;
956 efit += epar * ShowerShape(distance) ;
959 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
960 // Here we assume, that sigma = sqrt(E)
965 //____________________________________________________________________________
966 void AliPHOSClusterizerv1::Print(Option_t * option)const
968 // Print clusterizer parameters
971 TString taskName(GetName()) ;
972 taskName.ReplaceAll(Version(), "") ;
974 if( strcmp(GetName(), "") !=0 ) {
976 message = "\n--------------- %s %s -----------\n" ;
977 message += "Clusterizing digits from the file: %s\n" ;
978 message += " Branch: %s\n" ;
979 message += " EMC Clustering threshold = %f\n" ;
980 message += " EMC Local Maximum cut = %f\n" ;
981 message += " EMC Logarothmic weight = %f\n" ;
982 message += " CPV Clustering threshold = %f\n" ;
983 message += " CPV Local Maximum cut = %f\n" ;
984 message += " CPV Logarothmic weight = %f\n" ;
986 message += " Unfolding on\n" ;
988 message += " Unfolding off\n" ;
990 message += "------------------------------------------------------------------" ;
993 message = " AliPHOSClusterizerv1 not initialized " ;
995 Info("Print", message.Data(),
1000 fEmcClusteringThreshold,
1003 fCpvClusteringThreshold,
1009 //____________________________________________________________________________
1010 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1012 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1014 TObjArray * emcRecPoints = AliPHOSGetter::GetInstance()->EmcRecPoints(BranchName()) ;
1015 TObjArray * cpvRecPoints = AliPHOSGetter::GetInstance()->CpvRecPoints(BranchName()) ;
1018 message = "\nevent " ;
1019 message += gAlice->GetEvNumber() ;
1020 message += "\n Found " ;
1021 message += emcRecPoints->GetEntriesFast() ;
1022 message += "EMC RecPoints and " ;
1023 message += cpvRecPoints->GetEntriesFast() ;
1024 message += " CPV RecPoints \n" ;
1026 fRecPointsInRun += emcRecPoints->GetEntriesFast() ;
1027 fRecPointsInRun += cpvRecPoints->GetEntriesFast() ;
1029 char * tempo = new char[8192];
1031 if(strstr(option,"all")) {
1032 message += "\nEMC clusters\n" ;
1033 message += "Index Ene(MeV) Multi Module X Y Z Lambda 1 Lambda 2 # of prim Primaries list\n" ;
1035 for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1036 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ;
1038 rp->GetLocalPosition(locpos);
1040 rp->GetElipsAxis(lambda);
1043 primaries = rp->GetPrimaries(nprimaries);
1044 sprintf(tempo, "\n%6d %8.2f %3d %2d %4.1f %4.1f %4.1f %4f %4f %2d : ",
1045 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(), rp->GetPHOSMod(),
1046 locpos.X(), locpos.Y(), locpos.Z(), lambda[0], lambda[1], nprimaries) ;
1049 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
1050 sprintf(tempo, "%d ", primaries[iprimary] ) ;
1054 //Now plot CPV recPoints
1055 message += "Index Ene(MeV) Module X Y Z # of prim Primaries list\n" ;
1056 for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1057 AliPHOSRecPoint * rp = (AliPHOSRecPoint * )cpvRecPoints->At(index) ;
1060 rp->GetLocalPosition(locpos);
1064 primaries = rp->GetPrimaries(nprimaries);
1065 sprintf(tempo, "\n%6d %8.2f %2d %4.1f %4.1f %4.1f %2d : ",
1066 rp->GetIndexInList(), rp->GetEnergy(), rp->GetPHOSMod(),
1067 locpos.X(), locpos.Y(), locpos.Z(), nprimaries) ;
1069 primaries = rp->GetPrimaries(nprimaries);
1070 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
1071 sprintf(tempo, "%d ", primaries[iprimary] ) ;
1078 Info("Print", message.Data() ) ;