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 "AliPHOSClusterizerv1.h"
73 #include "AliPHOSCpvRecPoint.h"
74 #include "AliPHOSDigit.h"
75 #include "AliPHOSDigitizer.h"
76 #include "AliPHOSEmcRecPoint.h"
78 #include "AliPHOSGetter.h"
80 #include "AliPHOSCalibrManager.h"
81 #include "AliPHOSCalibrationData.h"
83 ClassImp(AliPHOSClusterizerv1)
85 //____________________________________________________________________________
86 AliPHOSClusterizerv1::AliPHOSClusterizerv1() : AliPHOSClusterizer()
88 // default ctor (to be used mainly by Streamer)
91 fDefaultInit = kTRUE ;
94 //____________________________________________________________________________
95 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const char* headerFile,const char* name, const Bool_t toSplit)
96 :AliPHOSClusterizer(headerFile, name, toSplit)
98 // ctor with the indication of the file where header Tree and digits Tree are stored
102 fDefaultInit = kFALSE ;
106 //____________________________________________________________________________
107 AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
118 //____________________________________________________________________________
119 const TString AliPHOSClusterizerv1::BranchName() const
121 // Returns branch name of the AliPHOSClusterizer
123 TString branchName(GetName() ) ;
124 branchName.Remove(branchName.Index(Version())-1) ;
128 //____________________________________________________________________________
129 Float_t AliPHOSClusterizerv1::Calibrate(Int_t amp, Int_t absId) const
131 // Returns energy value in the cell absId with the digit amplitude amp
133 if(fPedestals || fGains ){ //use calibration data
134 if(!fPedestals || !fGains ){
135 Error("Calibrate","Either Pedestals of Gains not set!") ;
138 Float_t en=(amp - fPedestals->Data(absId))*fGains->Data(absId) ;
145 if(absId <= fEmcCrystals) //calibrate as EMC
146 return fADCpedestalEmc + amp*fADCchanelEmc ;
147 else //calibrate as CPV
148 return fADCpedestalCpv+ amp*fADCchanelCpv ;
152 //____________________________________________________________________________
153 void AliPHOSClusterizerv1::Exec(Option_t * option)
156 if( strcmp(GetName(), "")== 0 )
159 if(strstr(option,"tim"))
160 gBenchmark->Start("PHOSClusterizer");
162 if(strstr(option,"print"))
165 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
166 if(gime->BranchExists("RecPoints"))
169 //test, if this is real data:
170 TBranch * eh = gAlice->TreeE()->GetBranch("AliPHOSBeamTestEvent") ;
171 if(eh){ //Read CalibrationData
172 AliPHOSCalibrManager * cmngr=AliPHOSCalibrManager::GetInstance() ;
173 if(!cmngr){ //not defined yet
174 cmngr=AliPHOSCalibrManager::GetInstance("PHOSBTCalibration.root") ;
178 fPedestals = new AliPHOSCalibrationData("Pedestals",fCalibrVersion) ;
179 cmngr->ReadFromRoot(*fPedestals,fCalibrRun) ;
182 fGains = new AliPHOSCalibrationData("Gains",fCalibrVersion) ;
183 cmngr->ReadFromRoot(*fGains,fCalibrRun) ;
186 Int_t nevents = gime->MaxEvent() ;
189 for(ievent = 0; ievent < nevents; ievent++){
190 gime->Event(ievent,"D") ;
191 Int_t pattern = gime->EventPattern() ;
193 printf("pattern %d \n",pattern) ;
195 Error("Exec","Patterns for reconstruction of events are not set, use SetPatterns() ") ;
198 Bool_t skip = kTRUE ;
199 Int_t npat = fPatterns->GetSize() ;
200 for(Int_t i = 0; i<npat;i++)
201 if(pattern==fPatterns->At(i)){
210 GetCalibrationParameters() ;
212 fNumberOfEmcClusters = fNumberOfCpvClusters = 0 ;
219 WriteRecPoints(ievent) ;
221 if(strstr(option,"deb"))
222 PrintRecPoints(option) ;
224 //increment the total number of digits per run
225 fRecPointsInRun += gime->EmcRecPoints()->GetEntriesFast() ;
226 fRecPointsInRun += gime->CpvRecPoints()->GetEntriesFast() ;
229 if(strstr(option,"tim")){
230 gBenchmark->Stop("PHOSClusterizer");
231 Info("Exec", " took %f seconds for Clusterizing %f seconds per event \n",
232 gBenchmark->GetCpuTime("PHOSClusterizer"),
233 gBenchmark->GetCpuTime("PHOSClusterizer")/nevents ) ;
237 //____________________________________________________________________________
238 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
239 Int_t nPar, Float_t * fitparameters) const
241 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
242 // The initial values for fitting procedure are set equal to the positions of local maxima.
243 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
245 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
246 TClonesArray * digits = gime->Digits() ;
249 gMinuit->mncler(); // Reset Minuit's list of paramters
250 gMinuit->SetPrintLevel(-1) ; // No Printout
251 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
252 // To set the address of the minimization function
254 TList * toMinuit = new TList();
255 toMinuit->AddAt(emcRP,0) ;
256 toMinuit->AddAt(digits,1) ;
258 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
260 // filling initial values for fit parameters
261 AliPHOSDigit * digit ;
265 Int_t nDigits = (Int_t) nPar / 3 ;
269 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
271 for(iDigit = 0; iDigit < nDigits; iDigit++){
272 digit = maxAt[iDigit];
277 geom->AbsToRelNumbering(digit->GetId(), relid) ;
278 geom->RelPosInModule(relid, x, z) ;
280 Float_t energy = maxAtEnergy[iDigit] ;
282 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
285 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : x = %f\n", x ) ;
288 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
291 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : z =%f\n", z ) ;
294 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
297 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : energy = %f\n", energy ) ;
302 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
307 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
308 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
309 gMinuit->SetMaxIterations(5);
310 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
312 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
314 if(ierflg == 4){ // Minimum not found
315 Warning("FindFit", "PHOS Unfolding fit not converged, cluster abandoned\n" );
318 for(index = 0; index < nPar; index++){
321 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
322 fitparameters[index] = val ;
330 //____________________________________________________________________________
331 void AliPHOSClusterizerv1::GetCalibrationParameters()
333 // Get calibration parameters from AliPHOSDigitizer
335 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
337 const TTask * task = gime->Digitizer(BranchName()) ;
338 if(strcmp(task->IsA()->GetName(),"AliPHOSDigitizer")==0){
339 const AliPHOSDigitizer * dig = static_cast<const AliPHOSDigitizer *>(task) ;
341 fADCchanelEmc = dig->GetEMCchannel() ;
342 fADCpedestalEmc = dig->GetEMCpedestal();
344 fADCchanelCpv = dig->GetCPVchannel() ;
345 fADCpedestalCpv = dig->GetCPVpedestal() ;
349 // fCalibrationDB = gime->CalibrationDB();
353 //____________________________________________________________________________
354 void AliPHOSClusterizerv1::Init()
356 // Make all memory allocations which can not be done in default constructor.
357 // Attach the Clusterizer task to the list of PHOS tasks
359 if ( strcmp(GetTitle(), "") == 0 )
360 SetTitle("galice.root") ;
362 TString branchname = GetName() ;
363 branchname.Remove(branchname.Index(Version())-1) ;
365 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(),branchname.Data(), fToSplit ) ;
367 Error("Init", "Could not obtain the Getter object !") ;
373 // construct the name of the file as /path/EMCAL.SDigits.root
374 //First - extract full path if necessary
375 TString fileName(GetTitle()) ;
376 Ssiz_t islash = fileName.Last('/') ;
377 if(islash<fileName.Length())
378 fileName.Remove(islash+1,fileName.Length()) ;
381 // Next - append the file name
382 fileName+="PHOS.RecData." ;
383 if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
384 fileName+=branchname ;
388 // Finally - check if the file already opened or open the file
389 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
391 fSplitFile = TFile::Open(fileName.Data(),"update") ;
396 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
397 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
400 gMinuit = new TMinuit(100) ;
402 gime->PostClusterizer(this) ;
403 gime->PostRecPoints(branchname) ;
407 //____________________________________________________________________________
408 void AliPHOSClusterizerv1::InitParameters()
410 // Set initial parameters for clusterization
412 fNumberOfCpvClusters = 0 ;
413 fNumberOfEmcClusters = 0 ;
415 fCpvClusteringThreshold = 0.0;
416 fEmcClusteringThreshold = 0.2;
418 fEmcLocMaxCut = 0.03 ;
419 fCpvLocMaxCut = 0.03 ;
424 fEmcTimeGate = 1.e-8 ;
428 fPurifyThreshold = 0.012 ;
432 fCalibrVersion= "v1" ;
436 TString clusterizerName( GetName()) ;
437 if (clusterizerName.IsNull() )
438 clusterizerName = "Default" ;
439 clusterizerName.Append(":") ;
440 clusterizerName.Append(Version()) ;
441 SetName(clusterizerName) ;
442 fRecPointsInRun = 0 ;
447 //____________________________________________________________________________
448 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
450 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
452 // = 2 are not neighbour but do not continue searching
453 // neighbours are defined as digits having at least a common vertex
454 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
455 // which is compared to a digit (d2) not yet in a cluster
457 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
462 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
465 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
467 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module
468 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
469 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
471 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
472 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
476 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
477 rv = 2; // Difference in row numbers is too large to look further
483 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
492 //____________________________________________________________________________
493 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
495 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
498 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
500 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
502 if(digit->GetId() <= nEMC ) rv = kTRUE;
507 //____________________________________________________________________________
508 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
510 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
513 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
515 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
517 if(digit->GetId() > nEMC ) rv = kTRUE;
522 //____________________________________________________________________________
523 void AliPHOSClusterizerv1::WriteRecPoints(Int_t event)
526 // Creates new branches with given title
527 // fills and writes into TreeR.
530 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
531 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
532 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
533 TClonesArray * digits = gime->Digits() ;
540 TString name("TreeR") ;
542 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
545 treeR = gAlice->TreeR();
550 gAlice->MakeTree("R", fSplitFile);
552 gAlice->MakeTree("R", gROOT->GetFile(GetTitle()));
553 treeR = gAlice->TreeR() ;
557 //Evaluate position, dispersion and other RecPoint properties...
558 for(index = 0; index < emcRecPoints->GetEntries(); index++){
559 AliPHOSEmcRecPoint * emcrp = static_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(index)) ;
560 emcrp->Purify(fPurifyThreshold) ;
561 if(emcrp->GetMultiplicity())
562 emcrp->EvalAll(fW0,digits) ;
564 emcRecPoints->Remove(emcrp) ;
567 emcRecPoints->Compress() ;
568 emcRecPoints->Sort() ;
569 for(index = 0; index < emcRecPoints->GetEntries(); index++)
570 ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->SetIndexInList(index) ;
572 emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ;
574 //Now the same for CPV
575 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
576 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->EvalAll(fW0CPV,digits) ;
578 cpvRecPoints->Sort() ;
580 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
581 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->SetIndexInList(index) ;
583 cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
585 Int_t bufferSize = 32000 ;
586 Int_t splitlevel = 0 ;
589 TBranch * emcBranch = treeR->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
590 emcBranch->SetTitle(BranchName());
593 TBranch * cpvBranch = treeR->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
594 cpvBranch->SetTitle(BranchName());
596 //And Finally clusterizer branch
597 AliPHOSClusterizerv1 * cl = this ;
598 TBranch * clusterizerBranch = treeR->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
599 &cl,bufferSize,splitlevel);
600 clusterizerBranch->SetTitle(BranchName());
604 clusterizerBranch->Fill() ;
607 if(gAlice->TreeR()!=treeR)
611 //____________________________________________________________________________
612 void AliPHOSClusterizerv1::MakeClusters()
614 // Steering method to construct the clusters stored in a list of Reconstructed Points
615 // A cluster is defined as a list of neighbour digits
617 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
619 TObjArray * emcRecPoints = gime->EmcRecPoints(BranchName()) ;
620 TObjArray * cpvRecPoints = gime->CpvRecPoints(BranchName()) ;
621 emcRecPoints->Delete() ;
622 cpvRecPoints->Delete() ;
624 TClonesArray * digits = gime->Digits() ;
626 Fatal("MakeClusters", "Digits with name %s not found !", BranchName().Data() ) ;
628 TClonesArray * digitsC = (TClonesArray*)digits->Clone() ;
631 // Clusterization starts
633 TIter nextdigit(digitsC) ;
634 AliPHOSDigit * digit ;
635 Bool_t notremoved = kTRUE ;
637 while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digitsC
638 AliPHOSRecPoint * clu = 0 ;
640 TArrayI clusterdigitslist(1500) ;
642 // cout << "digit " << digit << " " << Calibrate(digit->GetAmp(),digit->GetId()) << endl ;
644 if (( IsInEmc (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fEmcClusteringThreshold ) ||
645 ( IsInCpv (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fCpvClusteringThreshold ) ) {
646 Int_t iDigitInCluster = 0 ;
648 if ( IsInEmc(digit) ) {
649 // start a new EMC RecPoint
650 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
651 emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
653 emcRecPoints->AddAt(new AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
654 clu = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters) ;
655 fNumberOfEmcClusters++ ;
656 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId())) ;
657 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
659 digitsC->Remove(digit) ;
663 // start a new CPV cluster
664 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
665 cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
667 cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
669 clu = (AliPHOSCpvRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters) ;
670 fNumberOfCpvClusters++ ;
671 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId()) ) ;
672 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
674 digitsC->Remove(digit) ;
677 // Here we remove remaining EMC digits, which cannot make a cluster
680 while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
682 digitsC->Remove(digit) ;
686 notremoved = kFALSE ;
693 AliPHOSDigit * digitN ;
695 while (index < iDigitInCluster){ // scan over digits already in cluster
696 digit = (AliPHOSDigit*)digits->At(clusterdigitslist[index]) ;
698 while ( (digitN = (AliPHOSDigit *)nextdigit()) ) { // scan over the reduced list of digits
699 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
701 case 0 : // not a neighbour
703 case 1 : // are neighbours
704 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->GetId() ) ) ;
705 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
707 digitsC->Remove(digitN) ;
709 case 2 : // too far from each other
718 } // loop over cluster
729 //____________________________________________________________________________
730 void AliPHOSClusterizerv1::MakeUnfolding()
732 // Unfolds clusters using the shape of an ElectroMagnetic shower
733 // Performs unfolding of all EMC/CPV clusters
735 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
738 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
739 TObjArray * emcRecPoints = gime->EmcRecPoints(BranchName()) ;
740 TObjArray * cpvRecPoints = gime->CpvRecPoints(BranchName()) ;
741 TClonesArray * digits = gime->Digits() ;
743 // Unfold first EMC clusters
744 if(fNumberOfEmcClusters > 0){
746 Int_t nModulesToUnfold = geom->GetNModules() ;
748 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
750 for(index = 0 ; index < numberofNotUnfolded ; index++){
752 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) emcRecPoints->At(index) ;
753 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
756 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
757 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
758 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
759 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
761 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
762 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
763 if(emcRecPoint->GetNExMax()!=-1){ //Do not remove clusters which failed to unfold
764 emcRecPoints->Remove(emcRecPoint);
765 emcRecPoints->Compress() ;
767 fNumberOfEmcClusters -- ;
768 numberofNotUnfolded-- ;
772 emcRecPoint->SetNExMax(1) ; //Only one local maximum
775 delete[] maxAtEnergy ;
778 // Unfolding of EMC clusters finished
781 // Unfold now CPV clusters
782 if(fNumberOfCpvClusters > 0){
784 Int_t nModulesToUnfold = geom->GetNModules() ;
786 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
788 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
790 AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) cpvRecPoints->At(index) ;
792 if(recPoint->GetPHOSMod()> nModulesToUnfold)
795 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ;
797 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
798 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
799 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
800 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,digits) ;
802 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
803 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
804 cpvRecPoints->Remove(emcRecPoint);
805 cpvRecPoints->Compress() ;
807 numberofCpvNotUnfolded-- ;
808 fNumberOfCpvClusters-- ;
812 delete[] maxAtEnergy ;
815 //Unfolding of Cpv clusters finished
819 //____________________________________________________________________________
820 Double_t AliPHOSClusterizerv1::ShowerShape(Double_t r)
822 // Shape of the shower (see PHOS TDR)
823 // If you change this function, change also the gradient evaluation in ChiSquare()
825 Double_t r4 = r*r*r*r ;
826 Double_t r295 = TMath::Power(r, 2.95) ;
827 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
831 //____________________________________________________________________________
832 void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
834 AliPHOSDigit ** maxAt,
835 Float_t * maxAtEnergy)
837 // Performs the unfolding of a cluster with nMax overlapping showers
839 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
841 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
842 const TClonesArray * digits = gime->Digits() ;
843 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
844 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
846 Int_t nPar = 3 * nMax ;
847 Float_t * fitparameters = new Float_t[nPar] ;
849 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
851 // Fit failed, return
852 iniEmc->SetNExMax(-1) ;
853 delete[] fitparameters ;
857 // create ufolded rec points and fill them with new energy lists
858 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
859 // and later correct this number in acordance with actual energy deposition
861 Int_t nDigits = iniEmc->GetMultiplicity() ;
862 Float_t * efit = new Float_t[nDigits] ;
863 Float_t xDigit=0.,zDigit=0.,distance=0. ;
864 Float_t xpar=0.,zpar=0.,epar=0. ;
866 AliPHOSDigit * digit = 0 ;
867 Int_t * emcDigits = iniEmc->GetDigitsList() ;
871 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
872 digit = (AliPHOSDigit*) digits->At(emcDigits[iDigit] ) ;
873 geom->AbsToRelNumbering(digit->GetId(), relid) ;
874 geom->RelPosInModule(relid, xDigit, zDigit) ;
878 while(iparam < nPar ){
879 xpar = fitparameters[iparam] ;
880 zpar = fitparameters[iparam+1] ;
881 epar = fitparameters[iparam+2] ;
883 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
884 distance = TMath::Sqrt(distance) ;
885 efit[iDigit] += epar * ShowerShape(distance) ;
890 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
891 // so that energy deposited in each cell is distributed betwin new clusters proportionally
892 // to its contribution to efit
894 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
898 while(iparam < nPar ){
899 xpar = fitparameters[iparam] ;
900 zpar = fitparameters[iparam+1] ;
901 epar = fitparameters[iparam+2] ;
904 AliPHOSEmcRecPoint * emcRP = 0 ;
906 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
908 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
909 emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
911 (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
912 emcRP = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters);
913 fNumberOfEmcClusters++ ;
914 emcRP->SetNExMax((Int_t)nPar/3) ;
916 else{//create new entries in fCpvRecPoints
917 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
918 cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
920 (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
921 emcRP = (AliPHOSEmcRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters);
922 fNumberOfCpvClusters++ ;
926 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
927 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
928 geom->AbsToRelNumbering(digit->GetId(), relid) ;
929 geom->RelPosInModule(relid, xDigit, zDigit) ;
930 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
931 distance = TMath::Sqrt(distance) ;
932 ratio = epar * ShowerShape(distance) / efit[iDigit] ;
933 eDigit = emcEnergies[iDigit] * ratio ;
934 emcRP->AddDigit( *digit, eDigit ) ;
938 delete[] fitparameters ;
943 //_____________________________________________________________________________
944 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
946 // Calculates the Chi square for the cluster unfolding minimization
947 // Number of parameters, Gradient, Chi squared, parameters, what to do
949 TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
951 AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint*) toMinuit->At(0) ;
952 TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
956 // AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint *) gMinuit->GetObjectFit() ; // EmcRecPoint to fit
958 Int_t * emcDigits = emcRP->GetDigitsList() ;
960 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
962 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
964 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
969 for(iparam = 0 ; iparam < nPar ; iparam++)
970 Grad[iparam] = 0 ; // Will evaluate gradient
974 AliPHOSDigit * digit ;
977 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
979 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
985 geom->AbsToRelNumbering(digit->GetId(), relid) ;
987 geom->RelPosInModule(relid, xDigit, zDigit) ;
989 if(iflag == 2){ // calculate gradient
992 while(iParam < nPar ){
993 Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
995 distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
996 distance = TMath::Sqrt( distance ) ;
998 efit += x[iParam] * ShowerShape(distance) ;
1002 if(emcEnergies[iDigit] > 0)
1003 sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ;
1004 // Here we assume, that sigma = sqrt(E)
1006 while(iParam < nPar ){
1007 Double_t xpar = x[iParam] ;
1008 Double_t zpar = x[iParam+1] ;
1009 Double_t epar = x[iParam+2] ;
1010 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
1011 Double_t shape = sum * ShowerShape(dr) ;
1012 Double_t r4 = dr*dr*dr*dr ;
1013 Double_t r295 = TMath::Power(dr,2.95) ;
1014 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
1015 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
1017 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
1019 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
1021 Grad[iParam] += shape ; // Derivative over energy
1028 while(iparam < nPar ){
1029 Double_t xpar = x[iparam] ;
1030 Double_t zpar = x[iparam+1] ;
1031 Double_t epar = x[iparam+2] ;
1033 Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1034 distance = TMath::Sqrt(distance) ;
1035 efit += epar * ShowerShape(distance) ;
1038 if(emcEnergies[iDigit] > 0) //for the case if rounding error
1039 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/
1040 emcEnergies[iDigit] ;
1041 // Here we assume, that sigma = sqrt(E)
1046 //____________________________________________________________________________
1047 void AliPHOSClusterizerv1::Print(Option_t * option)const
1049 // Print clusterizer parameters
1052 TString taskName(GetName()) ;
1053 taskName.ReplaceAll(Version(), "") ;
1055 if( strcmp(GetName(), "") !=0 ) {
1057 message = "\n--------------- %s %s -----------\n" ;
1058 message += "Clusterizing digits from the file: %s\n" ;
1059 message += " Branch: %s\n" ;
1060 message += " EMC Clustering threshold = %f\n" ;
1061 message += " EMC Local Maximum cut = %f\n" ;
1062 message += " EMC Logarothmic weight = %f\n" ;
1063 message += " CPV Clustering threshold = %f\n" ;
1064 message += " CPV Local Maximum cut = %f\n" ;
1065 message += " CPV Logarothmic weight = %f\n" ;
1067 message += " Unfolding on\n" ;
1069 message += " Unfolding off\n" ;
1071 message += "------------------------------------------------------------------" ;
1074 message = " AliPHOSClusterizerv1 not initialized " ;
1076 Info("Print", message.Data(),
1081 fEmcClusteringThreshold,
1084 fCpvClusteringThreshold,
1090 //____________________________________________________________________________
1091 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1093 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1095 TObjArray * emcRecPoints = AliPHOSGetter::GetInstance()->EmcRecPoints(BranchName()) ;
1096 TObjArray * cpvRecPoints = AliPHOSGetter::GetInstance()->CpvRecPoints(BranchName()) ;
1098 if(strstr(option,"deb")){
1099 printf("event %d, found %d EmcRecPoints and %d CPV RecPoints \n",
1100 gAlice->GetEvNumber(),
1101 emcRecPoints->GetEntriesFast(),
1102 cpvRecPoints->GetEntriesFast()) ;
1106 if(strstr(option,"all")) {
1107 printf("\nEMC clusters\n" );
1108 printf("Index Ene(MeV) Multi Module X Y Z Lambda 1 Lambda 2 # of prim Primaries list\n" );
1110 for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1111 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ;
1113 rp->GetLocalPosition(locpos);
1115 rp->GetElipsAxis(lambda);
1118 primaries = rp->GetPrimaries(nprimaries);
1119 printf("\n%6d %8.2f %3d %2d %4.1f %4.1f %4.1f %4f %4f %2d : ",
1120 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(), rp->GetPHOSMod(),
1121 locpos.X(), locpos.Y(), locpos.Z(), lambda[0], lambda[1], nprimaries) ;
1123 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
1124 printf("%d ", primaries[iprimary] ) ;
1127 //Now plot CPV recPoints
1128 printf("\n CPV RecPoints \n") ;
1129 printf("Index Ene(MeV) Module X Y Z # of prim Primaries list\n") ;
1130 for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1131 AliPHOSRecPoint * rp = (AliPHOSRecPoint * )cpvRecPoints->At(index) ;
1134 rp->GetLocalPosition(locpos);
1138 primaries = rp->GetPrimaries(nprimaries);
1139 printf("\n%6d %8.2f %2d %4.1f %4.1f %4.1f %2d : ",
1140 rp->GetIndexInList(), rp->GetEnergy(), rp->GetPHOSMod(),
1141 locpos.X(), locpos.Y(), locpos.Z(), nprimaries) ;
1142 primaries = rp->GetPrimaries(nprimaries);
1143 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
1144 printf("%d ", primaries[iprimary] ) ;