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 // unfolding of the clusters with 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 normally called from Reconstructioner, but as well can be used it in
39 // root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root")
40 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
41 // //reads gAlice from header file "..."
42 // root [1] cl->ExecuteTask()
43 // //finds RecPoints in all events stored in galice.root
44 // root [2] cl->SetDigitsBranch("digits2")
45 // //sets another title for Digitis (input) branch
46 // root [3] cl->SetRecPointsBranch("recp2")
47 // //sets another title four output branches
48 // root [4] cl->SetEmcLocalMaxCut(0.03)
49 // //set clusterization parameters
50 // root [5] cl->ExecuteTask("deb all time")
51 // //once more finds RecPoints options are
52 // // deb - print number of found rec points
53 // // deb all - print number of found RecPoints and some their characteristics
54 // // time - print benchmarking results
56 // --- ROOT system ---
65 #include "TBenchmark.h"
67 // --- Standard library ---
72 // --- AliRoot header files ---
74 #include "AliPHOSClusterizerv1.h"
75 #include "AliPHOSCpvRecPoint.h"
76 #include "AliPHOSDigit.h"
77 #include "AliPHOSDigitizer.h"
78 #include "AliPHOSEmcRecPoint.h"
80 #include "AliPHOSPpsdRecPoint.h"
83 ClassImp(AliPHOSClusterizerv1)
85 //____________________________________________________________________________
86 AliPHOSClusterizerv1::AliPHOSClusterizerv1():AliPHOSClusterizer()
88 // default ctor (to be used mainly by Streamer)
89 SetName("AliPHOSClusterizer");
90 SetTitle("Version 1") ;
92 fNumberOfCpvClusters = 0 ;
93 fNumberOfEmcClusters = 0 ;
95 fCpvClusteringThreshold = 0.0;
96 fEmcClusteringThreshold = 0.2;
97 fPpsdClusteringThreshold = 0.0000002 ;
99 fEmcLocMaxCut = 0.03 ;
100 fCpvLocMaxCut = 0.03 ;
112 fIsInitialized = kFALSE ;
115 //____________________________________________________________________________
116 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const char* headerFile,const char* digitsFile):AliPHOSClusterizer()
118 SetName("AliPHOSClusterizer");
119 SetTitle("Version 1") ;
121 fNumberOfCpvClusters = 0 ;
122 fNumberOfEmcClusters = 0 ;
124 fCpvClusteringThreshold = 0.0;
125 fEmcClusteringThreshold = 0.2;
126 fPpsdClusteringThreshold = 0.0000002 ;
128 fEmcLocMaxCut = 0.03 ;
129 fCpvLocMaxCut = 0.03 ;
136 fHeaderFileName = headerFile ;
137 fDigitsBranchTitle = digitsFile ;
139 TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
142 if(fHeaderFileName.Contains("rfio")) // if we read file using HPSS
143 file = TFile::Open(fHeaderFileName.Data(),"update") ;
145 file = new TFile(fHeaderFileName.Data(),"update") ;
146 gAlice = (AliRun *) file->Get("gAlice") ;
149 AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;
150 fGeom = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
152 fDigits = new TClonesArray("AliPHOSDigit",10) ;
153 fDigitizer = new AliPHOSDigitizer() ;
154 fEmcRecPoints = new TObjArray(200) ;
155 fCpvRecPoints = new TObjArray(200) ;
157 if(!gMinuit) gMinuit = new TMinuit(100) ;
159 // add Task to //root/Tasks folder
160 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
161 roottasks->Add(this) ;
164 fIsInitialized = kTRUE ;
167 //____________________________________________________________________________
168 void AliPHOSClusterizerv1::Exec(Option_t * option){
171 if(!fIsInitialized) Init() ;
173 if(strstr(option,"tim"))
174 gBenchmark->Start("PHOSClusterizer");
176 Int_t nEvents = (Int_t) gAlice->TreeE()->GetEntries() ;
178 for(fEvent = 0;fEvent< nEvents; fEvent++){
179 if(!ReadDigits()) //reads digits for event fEvent
183 if(fToUnfold) MakeUnfolding() ;
185 if(strstr(option,"deb"))
186 PrintRecPoints(option) ;
189 if(strstr(option,"tim")){
190 gBenchmark->Stop("PHOSClusterizer");
191 cout << "AliPHOSClusterizer:" << endl ;
192 cout << " took " << gBenchmark->GetCpuTime("PHOSClusterizer") << " seconds for Clusterizing "
193 << gBenchmark->GetCpuTime("PHOSClusterizer")/nEvents << " seconds per event " << endl ;
200 //____________________________________________________________________________
201 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, int * maxAt, Float_t * maxAtEnergy,
202 Int_t nPar, Float_t * fitparameters)
204 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
205 // the initial values for fitting procedure are set in the positions of local maxima.
206 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
208 gMinuit->mncler(); // Reset Minuit's list of paramters
209 gMinuit->SetPrintLevel(-1) ; // No Printout
210 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
211 // To set the address of the minimization function
213 TList * toMinuit = new TList();
214 toMinuit->AddAt(emcRP,0) ;
215 toMinuit->AddAt(fDigits,1) ;
217 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
219 // filling initial values for fit parameters
220 AliPHOSDigit * digit ;
224 Int_t nDigits = (Int_t) nPar / 3 ;
229 for(iDigit = 0; iDigit < nDigits; iDigit++){
230 digit = (AliPHOSDigit *) maxAt[iDigit];
235 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
236 fGeom->RelPosInModule(relid, x, z) ;
238 Float_t energy = maxAtEnergy[iDigit] ;
240 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
243 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : x = " << x << endl ;
246 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
249 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : z = " << z << endl ;
252 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
255 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : energy = " << energy << endl ;
260 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
265 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
266 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
267 gMinuit->SetMaxIterations(5);
268 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
270 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
272 if(ierflg == 4){ // Minimum not found
273 cout << "PHOS Unfolding> Fit not converged, cluster abandoned "<< endl ;
276 for(index = 0; index < nPar; index++){
279 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
280 fitparameters[index] = val ;
288 //____________________________________________________________________________
289 void AliPHOSClusterizerv1::Init(){
290 //Make all memory allocations which can not be done in default constructor.
292 if(fHeaderFileName.IsNull())
293 fHeaderFileName = "galice.root" ;
296 TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
299 file = new TFile(fHeaderFileName.Data(),"update") ;
300 gAlice = (AliRun *) file->Get("gAlice") ;
303 AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;
304 fGeom = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
306 fDigits = new TClonesArray("AliPHOSDigit",10) ;
307 fDigitizer = new AliPHOSDigitizer() ;
308 fEmcRecPoints = new TObjArray(200) ;
309 fCpvRecPoints = new TObjArray(200) ;
311 if(!gMinuit) gMinuit = new TMinuit(100) ;
313 // add Task to //root/Tasks folder
314 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
315 roottasks->Add(this) ;
317 fIsInitialized = kTRUE ;
320 //____________________________________________________________________________
321 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
323 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
325 // = 2 are not neighbour but do not continue searching
326 // neighbours are defined as digits having at least common vertex
327 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
328 // which is compared to a digit (d2) not yet in a cluster
333 fGeom->AbsToRelNumbering(d1->GetId(), relid1) ;
336 fGeom->AbsToRelNumbering(d2->GetId(), relid2) ;
338 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module and the same PPSD Module
339 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
340 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
342 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
346 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
347 rv = 2; // Difference in row numbers is too large to look further
353 if( (relid1[0] < relid2[0]) || (relid1[1] < relid2[1]) )
358 //Do NOT clusterize upper PPSD
359 if( IsInPpsd(d1) && IsInPpsd(d2) &&
361 relid1[1] < fGeom->GetNumberOfPadsPhi()*fGeom->GetNumberOfPadsPhi() ) rv = 2 ;
367 //____________________________________________________________________________
368 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
370 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
375 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
377 if ( relid[1] == 0 ) rv = kTRUE;
382 //____________________________________________________________________________
383 Bool_t AliPHOSClusterizerv1::IsInPpsd(AliPHOSDigit * digit) const
385 // Tells if (true) or not (false) the digit is in a PHOS-PPSD module
390 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
392 if ( relid[1] > 0 && relid[0] > fGeom->GetNCPVModules() ) rv = kTRUE;
397 //____________________________________________________________________________
398 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
400 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
405 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
407 if ( relid[1] > 0 && relid[0] <= fGeom->GetNCPVModules() ) rv = kTRUE;
411 //____________________________________________________________________________
412 Bool_t AliPHOSClusterizerv1::ReadDigits(){
413 //reads digitis with specified title from TreeD
415 fNumberOfEmcClusters = 0 ;
416 fNumberOfCpvClusters = 0 ;
418 // Get Digits Tree header from file
419 gAlice->GetEvent(fEvent) ;
420 gAlice->SetEvent(fEvent) ;
422 TTree * treeD = gAlice->TreeD() ;
426 sprintf(treeName,"TreeD%d",fEvent);
427 cout << "Error in AliPHOSClusterizerv1 : no "<<treeName << endl ;
428 cout << " Do nothing " << endl ;
432 TBranch * digitsBranch = 0;
433 TBranch * digitizerBranch = 0;
435 TObjArray * branches = treeD->GetListOfBranches() ;
437 Bool_t phosNotFound = kTRUE ;
438 Bool_t digitizerNotFound = kTRUE ;
440 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
443 digitsBranch=(TBranch *) branches->At(ibranch) ;
444 if( fDigitsBranchTitle.CompareTo(digitsBranch->GetTitle())==0 )
445 if( strcmp(digitsBranch->GetName(),"PHOS") == 0)
446 phosNotFound = kFALSE ;
449 if(digitizerNotFound){
450 digitizerBranch = (TBranch *) branches->At(ibranch) ;
451 if( fDigitsBranchTitle.CompareTo(digitizerBranch->GetTitle()) == 0)
452 if( strcmp(digitizerBranch->GetName(),"AliPHOSDigitizer") == 0)
453 digitizerNotFound = kFALSE ;
458 if(digitizerNotFound || phosNotFound){
459 cout << "ERROR in AliPHOSClusterizerv1: " << endl ;
460 cout << " Can't find Branch with digits or Digitizer "<< endl ; ;
461 cout << " Do nothing" <<endl ;
465 digitsBranch->SetAddress(&fDigits) ;
466 digitizerBranch->SetAddress(&fDigitizer) ;
468 digitsBranch->GetEntry(0) ;
469 digitizerBranch->GetEntry(0) ;
471 fPedestal = fDigitizer->GetPedestal() ;
472 fSlope = fDigitizer->GetSlope() ;
476 //____________________________________________________________________________
477 void AliPHOSClusterizerv1::WriteRecPoints(){
478 // checks, if PHOSEmcRP etc. branches with given title already exist,
479 // exits without writing, otherwise create new branches with given title
480 // fills and wrights TreeR.
483 //Evaluate poisition, dispersion and other RecPoint properties...
484 for(index = 0; index < fEmcRecPoints->GetEntries(); index++)
485 ((AliPHOSEmcRecPoint *)fEmcRecPoints->At(index))->EvalAll(fW0,fDigits) ;
487 fEmcRecPoints->Sort() ;
489 for(index = 0; index < fEmcRecPoints->GetEntries(); index++)
490 ((AliPHOSEmcRecPoint *)fEmcRecPoints->At(index))->SetIndexInList(index) ;
492 fEmcRecPoints->Expand(fEmcRecPoints->GetEntriesFast()) ;
494 //Now the same for CPV
495 for(index = 0; index < fCpvRecPoints->GetEntries(); index++)
496 ((AliPHOSRecPoint *)fCpvRecPoints->At(index))->EvalAll(fW0CPV,fDigits) ;
498 fCpvRecPoints->Sort() ;
500 for(index = 0; index < fCpvRecPoints->GetEntries(); index++)
501 ((AliPHOSRecPoint *)fCpvRecPoints->At(index))->SetIndexInList(index) ;
503 fCpvRecPoints->Expand(fCpvRecPoints->GetEntriesFast()) ;
505 if(gAlice->TreeR()==0)
506 gAlice->MakeTree("R") ;
509 //Check, if branches already exist
510 TBranch * emcBranch = 0;
511 TBranch * cpvBranch = 0;
512 TBranch * clusterizerBranch = 0;
514 TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
516 Bool_t emcNotFound = kTRUE ;
517 Bool_t cpvNotFound = kTRUE ;
518 Bool_t clusterizerNotFound = kTRUE ;
520 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
523 emcBranch=(TBranch *) branches->At(ibranch) ;
524 if( fRecPointsBranchTitle.CompareTo(emcBranch->GetTitle())==0 ){
525 if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0) {
526 emcNotFound = kFALSE ;
531 cpvBranch=(TBranch *) branches->At(ibranch) ;
532 if( fRecPointsBranchTitle.CompareTo(cpvBranch->GetTitle())==0 ){
533 if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0) {
534 cpvNotFound = kFALSE ;
539 if(clusterizerNotFound){
540 clusterizerBranch = (TBranch *) branches->At(ibranch) ;
541 if( fRecPointsBranchTitle.CompareTo(clusterizerBranch->GetTitle()) == 0){
542 if( strcmp(clusterizerBranch->GetName(),"AliPHOSClusterizer") == 0) {
543 clusterizerNotFound = kFALSE ;
550 if(!(clusterizerNotFound && emcNotFound && cpvNotFound)){
551 cout << "AliPHOSClusterizer error" << endl;
552 cout << " Branches PHOSEmcRP, PHOSCpvRP and AliPHOSClusterizer " << endl ;
553 cout << " with title '" << fRecPointsBranchTitle.Data() <<"' already exist" << endl ;
554 cout << " can not overwrite " << endl ;
558 //Make branches in TreeR for RecPoints and Clusterizer
560 if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){ //generating file name
561 filename = new char[strlen(gAlice->GetBaseFile())+20] ;
562 sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ;
566 TDirectory *cwd = gDirectory;
569 Int_t bufferSize = 32000 ;
570 Int_t splitlevel = 0 ;
571 emcBranch = gAlice->TreeR()->Branch("PHOSEmcRP","TObjArray",&fEmcRecPoints,bufferSize,splitlevel);
572 emcBranch->SetTitle(fRecPointsBranchTitle.Data());
574 emcBranch->SetFile(filename);
575 TIter next( emcBranch->GetListOfBranches());
577 while ((sb=(TBranch*)next())) {
578 sb->SetFile(filename);
584 cpvBranch = gAlice->TreeR()->Branch("PHOSCpvRP","TObjArray",&fCpvRecPoints,bufferSize,splitlevel);
585 cpvBranch->SetTitle(fRecPointsBranchTitle.Data());
587 cpvBranch->SetFile(filename);
588 TIter next( cpvBranch->GetListOfBranches());
590 while ((sb=(TBranch*)next())) {
591 sb->SetFile(filename);
596 //And Finally clusterizer branch
597 AliPHOSClusterizerv1 * cl = this ;
598 clusterizerBranch = gAlice->TreeR()->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
599 &cl,bufferSize,splitlevel);
600 clusterizerBranch->SetTitle(fRecPointsBranchTitle.Data());
602 clusterizerBranch->SetFile(filename);
603 TIter next( clusterizerBranch->GetListOfBranches());
605 while ((sb=(TBranch*)next())) {
606 sb->SetFile(filename);
614 clusterizerBranch->Fill() ;
615 gAlice->TreeR()->Write(0,kOverwrite) ;
619 //____________________________________________________________________________
620 void AliPHOSClusterizerv1::MakeClusters()
622 // Steering method to construct the clusters stored in a list of Reconstructed Points
623 // A cluster is defined as a list of neighbour digits
624 fEmcRecPoints->Clear() ;
625 fCpvRecPoints->Clear() ;
628 // Clusterization starts
629 TClonesArray * digits = (TClonesArray*)fDigits->Clone() ;
631 TIter nextdigit(digits) ;
632 AliPHOSDigit * digit ;
633 Bool_t notremoved = kTRUE ;
635 while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digits
636 AliPHOSRecPoint * clu = 0 ;
638 TArrayI clusterdigitslist(1000) ;
641 if (( IsInEmc (digit) && Calibrate(digit->GetAmp()) > fEmcClusteringThreshold ) ||
642 ( IsInPpsd(digit) && Calibrate(digit->GetAmp()) > fPpsdClusteringThreshold ) ||
643 ( IsInCpv (digit) && Calibrate(digit->GetAmp()) > fCpvClusteringThreshold ) ) {
645 Int_t iDigitInCluster = 0 ;
647 if ( IsInEmc(digit) ) {
648 // start a new EMC RecPoint
649 if(fNumberOfEmcClusters >= fEmcRecPoints->GetSize()) fEmcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
650 fEmcRecPoints->AddAt(new AliPHOSEmcRecPoint(), fNumberOfEmcClusters) ;
651 clu = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(fNumberOfEmcClusters) ;
652 fNumberOfEmcClusters++ ;
653 clu->AddDigit(*digit, Calibrate(digit->GetAmp())) ;
654 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
656 digits->Remove(digit) ;
660 // start a new PPSD/CPV cluster
661 if(fNumberOfCpvClusters >= fCpvRecPoints->GetSize()) fCpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
663 fCpvRecPoints->AddAt(new AliPHOSPpsdRecPoint(),fNumberOfCpvClusters) ;
665 fCpvRecPoints->AddAt(new AliPHOSCpvRecPoint(), fNumberOfCpvClusters) ;
666 clu = (AliPHOSPpsdRecPoint *) fCpvRecPoints->At(fNumberOfCpvClusters) ;
667 fNumberOfCpvClusters++ ;
669 clu->AddDigit(*digit, Calibrate(digit->GetAmp()) ) ;
670 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
672 digits->Remove(digit) ;
675 // Here we remove resting EMC digits, which cannot make cluster
678 while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
680 digits->Remove(digit) ;
684 notremoved = kFALSE ;
691 AliPHOSDigit * digitN ;
693 while (index < iDigitInCluster){ // scan over digits already in cluster
694 digit = (AliPHOSDigit*)fDigits->At(clusterdigitslist[index]) ;
696 while ( (digitN = (AliPHOSDigit *)nextdigit()) ) { // scan over the reduced list of digits
697 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
699 case 0 : // not a neighbour
701 case 1 : // are neighbours
702 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp() ) ) ;
703 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
705 digits->Remove(digitN) ;
707 case 2 : // too far from each other
716 } // loop over cluster
727 //____________________________________________________________________________
728 void AliPHOSClusterizerv1::MakeUnfolding(){
729 //Unfolds clusters using the shape of ElectroMagnetic shower
730 // Performs unfolding of all EMC/CPV but NOT ppsd clusters
732 //Unfold first EMC clusters
733 if(fNumberOfEmcClusters > 0){
735 Int_t nModulesToUnfold = fGeom->GetNModules() ;
737 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
739 for(index = 0 ; index < numberofNotUnfolded ; index++){
741 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(index) ;
742 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
745 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
746 Int_t * maxAt = new Int_t[nMultipl] ;
747 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
748 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,fDigits) ;
750 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
751 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
752 fEmcRecPoints->Remove(emcRecPoint);
753 fEmcRecPoints->Compress() ;
755 fNumberOfEmcClusters -- ;
756 numberofNotUnfolded-- ;
760 delete[] maxAtEnergy ;
763 //Unfolding of EMC clusters finished
766 //Unfold now CPV clusters
767 if(fNumberOfCpvClusters > 0){
769 Int_t nModulesToUnfold = fGeom->GetNCPVModules() ;
771 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
773 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
775 AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) fCpvRecPoints->At(index) ;
777 if(recPoint->GetPHOSMod()> nModulesToUnfold)
780 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ;
782 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
783 Int_t * maxAt = new Int_t[nMultipl] ;
784 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
785 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,fDigits) ;
787 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
788 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
789 fCpvRecPoints->Remove(emcRecPoint);
790 fCpvRecPoints->Compress() ;
792 numberofCpvNotUnfolded-- ;
793 fNumberOfCpvClusters-- ;
797 delete[] maxAtEnergy ;
800 //Unfolding of Cpv clusters finished
804 //____________________________________________________________________________
805 Double_t AliPHOSClusterizerv1::ShowerShape(Double_t r)
807 // Shape of the shower (see PHOS TDR)
808 // If you change this function, change also the gradien evaluation in ChiSquare()
810 Double_t r4 = r*r*r*r ;
811 Double_t r295 = TMath::Power(r, 2.95) ;
812 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
816 //____________________________________________________________________________
817 void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
820 Float_t * maxAtEnergy)
822 // Performs the unfolding of a cluster with nMax overlapping showers
824 Int_t nPar = 3 * nMax ;
825 Float_t * fitparameters = new Float_t[nPar] ;
827 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
829 // Fit failed, return and remove cluster
830 delete[] fitparameters ;
834 // create ufolded rec points and fill them with new energy lists
835 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
836 // and later correct this number in acordance with actual energy deposition
838 Int_t nDigits = iniEmc->GetMultiplicity() ;
839 Float_t * efit = new Float_t[nDigits] ;
840 Float_t xDigit,zDigit,distance ;
841 Float_t xpar,zpar,epar ;
843 AliPHOSDigit * digit ;
844 Int_t * emcDigits = iniEmc->GetDigitsList() ;
848 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
849 digit = (AliPHOSDigit*) fDigits->At(emcDigits[iDigit] ) ;
850 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
851 fGeom->RelPosInModule(relid, xDigit, zDigit) ;
855 while(iparam < nPar ){
856 xpar = fitparameters[iparam] ;
857 zpar = fitparameters[iparam+1] ;
858 epar = fitparameters[iparam+2] ;
860 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
861 distance = TMath::Sqrt(distance) ;
862 efit[iDigit] += epar * ShowerShape(distance) ;
867 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
868 // so that energy deposited in each cell is distributed betwin new clusters proportionally
869 // to its contribution to efit
871 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
875 while(iparam < nPar ){
876 xpar = fitparameters[iparam] ;
877 zpar = fitparameters[iparam+1] ;
878 epar = fitparameters[iparam+2] ;
881 AliPHOSEmcRecPoint * emcRP ;
883 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
885 if(fNumberOfEmcClusters >= fEmcRecPoints->GetSize())
886 fEmcRecPoints->Expand(2*fNumberOfEmcClusters) ;
888 (*fEmcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint() ;
889 emcRP = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(fNumberOfEmcClusters);
890 fNumberOfEmcClusters++ ;
892 else{//create new entries in fCpvRecPoints
893 if(fNumberOfCpvClusters >= fCpvRecPoints->GetSize())
894 fCpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
896 (*fCpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint() ;
897 emcRP = (AliPHOSEmcRecPoint *) fCpvRecPoints->At(fNumberOfCpvClusters);
898 fNumberOfCpvClusters++ ;
902 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
903 digit = (AliPHOSDigit*) fDigits->At( emcDigits[iDigit] ) ;
904 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
905 fGeom->RelPosInModule(relid, xDigit, zDigit) ;
906 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
907 distance = TMath::Sqrt(distance) ;
908 ratio = epar * ShowerShape(distance) / efit[iDigit] ;
909 eDigit = emcEnergies[iDigit] * ratio ;
910 emcRP->AddDigit( *digit, eDigit ) ;
914 delete[] fitparameters ;
919 //_____________________________________________________________________________
920 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
922 // Calculates th Chi square for the cluster unfolding minimization
923 // Number of parameters, Gradient, Chi squared, parameters, what to do
926 TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
928 AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint*) toMinuit->At(0) ;
929 TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
933 // AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint *) gMinuit->GetObjectFit() ; // EmcRecPoint to fit
935 Int_t * emcDigits = emcRP->GetDigitsList() ;
937 Int_t nOfDigits = emcRP->GetDigitsMultiplicity() ;
939 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
941 AliPHOSGeometry * geom = AliPHOSGeometry::GetInstance() ;
947 for(iparam = 0 ; iparam < nPar ; iparam++)
948 Grad[iparam] = 0 ; // Will evaluate gradient
952 AliPHOSDigit * digit ;
955 for( iDigit = 0 ; iDigit < nOfDigits ; iDigit++) {
957 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
963 geom->AbsToRelNumbering(digit->GetId(), relid) ;
965 geom->RelPosInModule(relid, xDigit, zDigit) ;
967 if(iflag == 2){ // calculate gradient
970 while(iParam < nPar ){
971 Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
973 distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
974 distance = TMath::Sqrt( distance ) ;
976 efit += x[iParam] * ShowerShape(distance) ;
979 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
981 while(iParam < nPar ){
982 Double_t xpar = x[iParam] ;
983 Double_t zpar = x[iParam+1] ;
984 Double_t epar = x[iParam+2] ;
985 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
986 Double_t shape = sum * ShowerShape(dr) ;
987 Double_t r4 = dr*dr*dr*dr ;
988 Double_t r295 = TMath::Power(dr,2.95) ;
989 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
990 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
992 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
994 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
996 Grad[iParam] += shape ; // Derivative over energy
1003 while(iparam < nPar ){
1004 Double_t xpar = x[iparam] ;
1005 Double_t zpar = x[iparam+1] ;
1006 Double_t epar = x[iparam+2] ;
1008 Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1009 distance = TMath::Sqrt(distance) ;
1010 efit += epar * ShowerShape(distance) ;
1013 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
1014 // Here we assume, that sigma = sqrt(E)
1019 //____________________________________________________________________________
1020 void AliPHOSClusterizerv1::Print(Option_t * option)const
1026 cout << "---------------"<< GetName() << " " << GetTitle()<< "-----------" << endl
1027 << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl
1028 << " Branch: " << fDigitsBranchTitle.Data() << endl
1030 << " EMC Clustering threshold = " << fEmcClusteringThreshold << endl
1031 << " EMC Local Maximum cut = " << fEmcLocMaxCut << endl
1032 << " EMC Logarothmic weight = " << fW0 << endl
1034 << " CPV Clustering threshold = " << fCpvClusteringThreshold << endl
1035 << " CPV Local Maximum cut = " << fCpvLocMaxCut << endl
1036 << " CPV Logarothmic weight = " << fW0CPV << endl
1038 << " PPSD Clustering threshold = " << fPpsdClusteringThreshold << endl;
1040 cout << " Unfolding on " << endl ;
1042 cout << " Unfolding off " << endl ;
1044 cout << "------------------------------------------------------------------" <<endl ;
1047 cout << " AliPHOSClusterizerv1 not initialized " << endl ;
1049 //____________________________________________________________________________
1050 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option){
1051 //Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1053 cout << "AliPHOSClusterizerv1: " << endl ;
1054 cout << " Found "<< fEmcRecPoints->GetEntriesFast() << " EMC Rec Points and "
1055 << fCpvRecPoints->GetEntriesFast() << " CPV RecPoints" << endl ;
1057 if(strstr(option,"all")) {
1058 cout << "EMC clusters " << endl ;
1070 << " Primaries list " << endl;
1073 for (index = 0 ; index < fEmcRecPoints->GetEntries() ; index++) {
1074 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )fEmcRecPoints->At(index) ;
1076 cout << setw(6) << rp->GetIndexInList() << " ";
1077 cout << setw(6) << rp->GetEnergy() << " ";
1078 cout << setw(6) << rp->GetMultiplicity()<< " ";
1079 cout << setw(6) << rp->GetPHOSMod() << " ";
1082 rp->GetLocalPosition(locpos);
1083 cout << setw(8) << locpos.X() << " ";
1084 cout << setw(8) << locpos.Y() << " ";
1085 cout << setw(8) << locpos.Z() << " ";
1088 rp->GetElipsAxis(lambda);
1089 cout << setw(10)<< lambda[0] << " ";
1090 cout << setw(10)<< lambda[1] << " ";
1095 primaries = rp->GetPrimaries(nprimaries);
1096 cout << setw(8) << primaries << " ";
1098 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1099 cout << setw(4) << primaries[iprimary] << " ";
1104 //Now plot PCV/PPSD recPoints
1105 cout << "EMC clusters " << endl ;
1114 << " Primaries list " << endl;
1116 for (index = 0 ; index < fEmcRecPoints->GetEntries() ; index++) {
1117 AliPHOSRecPoint * rp = (AliPHOSRecPoint * )fCpvRecPoints->At(index) ;
1118 cout << setw(6) << rp->GetIndexInList() << " ";
1119 cout << setw(6) << rp->GetPHOSMod() << " ";
1121 if( (strcmp(rp->ClassName() , "AliPHOSPpsdRecPoint" )) == 0){
1122 AliPHOSPpsdRecPoint * ppsd = (AliPHOSPpsdRecPoint*) rp ;
1132 rp->GetLocalPosition(locpos);
1133 cout << setw(8) << locpos.X() << " ";
1134 cout << setw(8) << locpos.Y() << " ";
1135 cout << setw(8) << locpos.Z() << " ";
1139 primaries = rp->GetPrimaries(nprimaries);
1140 cout << setw(8) << primaries << " ";
1142 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1143 cout << setw(4) << primaries[iprimary] << " ";
1148 cout << "-------------------------------------------------"<<endl ;