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")
41 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 // //reads gAlice from header file "..."
43 // root [1] cl->ExecuteTask()
44 // //finds RecPoints in all events stored in galice.root
45 // root [2] cl->SetDigitsBranch("digits2")
46 // //sets another title for Digitis (input) branch
47 // root [3] cl->SetRecPointsBranch("recp2")
48 // //sets another title four output branches
49 // root [4] cl->SetEmcLocalMaxCut(0.03)
50 // //set clusterization parameters
51 // root [5] cl->ExecuteTask("deb all time")
52 // //once more finds RecPoints options are
53 // // deb - print number of found rec points
54 // // deb all - print number of found RecPoints and some their characteristics
55 // // time - print benchmarking results
57 // --- ROOT system ---
66 #include "TBenchmark.h"
68 // --- Standard library ---
73 // --- AliRoot header files ---
75 #include "AliPHOSClusterizerv1.h"
76 #include "AliPHOSCpvRecPoint.h"
77 #include "AliPHOSDigit.h"
78 #include "AliPHOSDigitizer.h"
79 #include "AliPHOSEmcRecPoint.h"
81 #include "AliPHOSPpsdRecPoint.h"
82 #include "AliPHOSGetter.h"
85 ClassImp(AliPHOSClusterizerv1)
87 //____________________________________________________________________________
88 AliPHOSClusterizerv1::AliPHOSClusterizerv1() : AliPHOSClusterizer()
90 // default ctor (to be used mainly by Streamer)
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 ;
105 fHeaderFileName = "" ;
106 fDigitsBranchTitle = "" ;
107 fRecPointsInRun = 0 ;
110 //____________________________________________________________________________
111 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const char* headerFile,const char* name)
112 :AliPHOSClusterizer(headerFile, name)
114 // ctor with the indication of the file where header Tree and digits Tree are stored
117 fNumberOfCpvClusters = 0 ;
118 fNumberOfEmcClusters = 0 ;
120 fCpvClusteringThreshold = 0.0;
121 fEmcClusteringThreshold = 0.2;
122 fPpsdClusteringThreshold = 0.0000002 ;
124 fEmcLocMaxCut = 0.03 ;
125 fCpvLocMaxCut = 0.03 ;
132 fHeaderFileName = GetTitle() ;
133 fDigitsBranchTitle = GetName() ;
135 TString tempo(GetName()) ;
137 tempo.Append(Version()) ;
138 SetName(tempo.Data()) ;
139 fRecPointsInRun = 0 ;
145 //____________________________________________________________________________
146 void AliPHOSClusterizerv1::Exec(Option_t * option)
150 if( strcmp(GetName(), "")== 0 )
153 if(strstr(option,"tim"))
154 gBenchmark->Start("PHOSClusterizer");
156 if(strstr(option,"print"))
159 //check, if the branch with name of this" already exits?
160 gAlice->GetEvent(0) ;
161 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
163 TBranch * branch = 0 ;
164 Bool_t phosemcfound = kFALSE, phoscpvfound = kFALSE, clusterizerfound = kFALSE ;
166 TString taskName(GetName()) ;
167 taskName.Remove(taskName.Index(Version()) -1) ;
169 while ( (branch = (TBranch*)next()) && (!phosemcfound || !phoscpvfound || !clusterizerfound) ) {
170 if ( (strcmp(branch->GetName(), "PHOSEmcRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
171 phosemcfound = kTRUE ;
173 else if ( (strcmp(branch->GetName(), "PHOSCpvRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
174 phoscpvfound = kTRUE ;
176 else if ( (strcmp(branch->GetName(), "AliPHOSClusterizer")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
177 clusterizerfound = kTRUE ;
180 if ( phoscpvfound || phosemcfound || clusterizerfound ) {
181 cerr << "WARNING: AliPHOSClusterizer::Exec -> Emc(Cpv)RecPoints and/or Clusterizer branch with name "
182 << taskName.Data() << " already exits" << endl ;
186 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
187 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
190 for(ievent = 0; ievent < nevents; ievent++){
193 fPedestal = gime->Digitizer(taskName)->GetPedestal() ;
194 fSlope = gime->Digitizer(taskName)->GetSlope() ;
197 fNumberOfEmcClusters = 0 ;
198 fNumberOfCpvClusters = 0 ;
200 gime->Event(ievent,"D") ;
201 // if(!ReadDigits(ievent)) //reads digits for event fEvent
209 WriteRecPoints(ievent) ;
211 if(strstr(option,"deb"))
212 PrintRecPoints(option) ;
215 if(strstr(option,"tim")){
216 gBenchmark->Stop("PHOSClusterizer");
217 cout << "AliPHOSClusterizer:" << endl ;
218 cout << " took " << gBenchmark->GetCpuTime("PHOSClusterizer") << " seconds for Clusterizing "
219 << gBenchmark->GetCpuTime("PHOSClusterizer")/nevents << " seconds per event " << endl ;
225 //____________________________________________________________________________
226 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, int * maxAt, Float_t * maxAtEnergy,
227 Int_t nPar, Float_t * fitparameters) const
229 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
230 // The initial values for fitting procedure are set equal to the positions of local maxima.
231 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
233 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
234 TClonesArray * digits = gime->Digits() ;
237 gMinuit->mncler(); // Reset Minuit's list of paramters
238 gMinuit->SetPrintLevel(-1) ; // No Printout
239 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
240 // To set the address of the minimization function
242 TList * toMinuit = new TList();
243 toMinuit->AddAt(emcRP,0) ;
244 toMinuit->AddAt(digits,1) ;
246 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
248 // filling initial values for fit parameters
249 AliPHOSDigit * digit ;
253 Int_t nDigits = (Int_t) nPar / 3 ;
257 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
259 for(iDigit = 0; iDigit < nDigits; iDigit++){
260 digit = (AliPHOSDigit *) maxAt[iDigit];
265 geom->AbsToRelNumbering(digit->GetId(), relid) ;
266 geom->RelPosInModule(relid, x, z) ;
268 Float_t energy = maxAtEnergy[iDigit] ;
270 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
273 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : x = " << x << endl ;
276 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
279 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : z = " << z << endl ;
282 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
285 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : energy = " << energy << endl ;
290 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
295 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
296 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
297 gMinuit->SetMaxIterations(5);
298 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
300 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
302 if(ierflg == 4){ // Minimum not found
303 cout << "PHOS Unfolding> Fit not converged, cluster abandoned "<< endl ;
306 for(index = 0; index < nPar; index++){
309 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
310 fitparameters[index] = val ;
318 //____________________________________________________________________________
319 void AliPHOSClusterizerv1::Init()
321 // Make all memory allocations which can not be done in default constructor.
322 // Attach the Clusterizer task to the list of PHOS tasks
324 if ( strcmp(GetTitle(), "") == 0 )
325 SetTitle("galice.root") ;
327 TString taskName(GetName()) ;
328 taskName.Remove(taskName.Index(Version()) -1) ;
330 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName.Data()) ;
332 cerr << "ERROR: AliPHOSClusterizerv1::Init -> Could not obtain the Getter object !" << endl ;
337 gMinuit = new TMinuit(100) ;
339 //add Task to //YSAlice/tasks/Reconstructioner/PHOS
340 gime->PostClusterizer(this) ;
341 // create a folder on the white board //YSAlice/WhiteBoard/RecPoints/PHOS/recpointsName
342 gime->PostRecPoints(taskName.Data() ) ;
344 gime->PostDigits(taskName.Data()) ;
345 gime->PostDigitizer(taskName.Data()) ;
349 //____________________________________________________________________________
350 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
352 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
354 // = 2 are not neighbour but do not continue searching
355 // neighbours are defined as digits having at least a common vertex
356 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
357 // which is compared to a digit (d2) not yet in a cluster
359 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
364 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
367 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
369 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module and the same PPSD Module
370 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
371 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
373 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
377 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
378 rv = 2; // Difference in row numbers is too large to look further
384 if( (relid1[0] < relid2[0]) || (relid1[1] < relid2[1]) )
389 //Do NOT clusterize upper PPSD
390 if( IsInPpsd(d1) && IsInPpsd(d2) &&
392 relid1[1] < geom->GetNumberOfPadsPhi()*geom->GetNumberOfPadsPhi() ) rv = 2 ;
398 //____________________________________________________________________________
399 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
401 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
404 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
407 geom->AbsToRelNumbering(digit->GetId(), relid) ;
409 if ( relid[1] == 0 ) rv = kTRUE;
414 //____________________________________________________________________________
415 Bool_t AliPHOSClusterizerv1::IsInPpsd(AliPHOSDigit * digit) const
417 // Tells if (true) or not (false) the digit is in a PHOS-PPSD module
420 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
423 geom->AbsToRelNumbering(digit->GetId(), relid) ;
425 if ( relid[1] > 0 && relid[0] > geom->GetNCPVModules() ) rv = kTRUE;
430 //____________________________________________________________________________
431 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
433 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
436 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
439 geom->AbsToRelNumbering(digit->GetId(), relid) ;
441 if ( relid[1] > 0 && relid[0] <= geom->GetNCPVModules() ) rv = kTRUE;
445 //____________________________________________________________________________
446 Bool_t AliPHOSClusterizerv1::ReadDigits(Int_t event)
448 // reads digitis with specified title from TreeD
451 if ( gAlice->TreeD()==0) {
452 cerr << "ERROR: AliPHOSClusterizerv1::ReadDigits There is no Digit Tree" << endl;
456 //set address of the Digits and Digitizer
457 TBranch * digitsBranch = 0;
458 TBranch * digitizerBranch = 0;
459 TObjArray * lob = (TObjArray*)gAlice->TreeD()->GetListOfBranches() ;
461 TBranch * branch = 0 ;
462 Bool_t phosfound = kFALSE, digitizerfound = kFALSE ;
464 TString taskName(GetName()) ;
465 taskName.ReplaceAll(Version(), "") ;
467 while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
468 if ( (strcmp(branch->GetName(), "PHOS")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
470 digitsBranch = branch ;
473 else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
474 digitizerfound = kTRUE ;
475 digitizerBranch = branch ;
478 if ( !phosfound || !digitizerfound ) {
479 cerr << "WARNING: AliPHOSClusterizerv1::ReadDigits -> Digits and/or Digitizer branch with name " << taskName.Data()
480 << " not found" << endl ;
484 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
486 TClonesArray * digits = gime->Digits() ;
488 digitsBranch->SetAddress(&digits) ;
490 AliPHOSDigitizer * digitizer = gime->Digitizer() ;
491 digitizerBranch->SetAddress(&digitizer) ;
493 digitsBranch->GetEntry(0) ;
494 digitizerBranch->GetEntry(0) ;
496 fPedestal = digitizer->GetPedestal() ;
497 fSlope = digitizer->GetSlope() ;
502 //____________________________________________________________________________
503 void AliPHOSClusterizerv1::WriteRecPoints(Int_t event)
506 // Creates new branches with given title
507 // fills and writes into TreeR.
508 TString branchName(GetName()) ;
509 branchName.Remove(branchName.Index(Version())-1) ;
511 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
512 TObjArray * emcRecPoints = gime->EmcRecPoints(branchName) ;
513 TObjArray * cpvRecPoints = gime->CpvRecPoints(branchName) ;
514 TClonesArray * digits = gime->Digits(branchName) ;
517 //Evaluate position, dispersion and other RecPoint properties...
518 for(index = 0; index < emcRecPoints->GetEntries(); index++)
519 ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->EvalAll(fW0,digits) ;
521 emcRecPoints->Sort() ;
523 for(index = 0; index < emcRecPoints->GetEntries(); index++)
524 ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->SetIndexInList(index) ;
526 emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ;
528 //Now the same for CPV
529 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
530 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->EvalAll(fW0CPV,digits) ;
532 cpvRecPoints->Sort() ;
534 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
535 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->SetIndexInList(index) ;
537 cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
539 gAlice->GetEvent(event) ;
540 if(gAlice->TreeR()==0)
541 gAlice->MakeTree("R") ;
543 //Make branches in TreeR for RecPoints and Clusterizer
545 if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){ //generating file name
546 filename = new char[strlen(gAlice->GetBaseFile())+20] ;
547 sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ;
551 TDirectory *cwd = gDirectory;
554 Int_t bufferSize = 32000 ;
555 Int_t splitlevel = 0 ;
558 TBranch * emcBranch = gAlice->TreeR()->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
559 emcBranch->SetTitle(branchName);
561 emcBranch->SetFile(filename);
562 TIter next( emcBranch->GetListOfBranches());
564 while ((sb=(TBranch*)next())) {
565 sb->SetFile(filename);
572 TBranch * cpvBranch = gAlice->TreeR()->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
573 cpvBranch->SetTitle(branchName);
575 cpvBranch->SetFile(filename);
576 TIter next( cpvBranch->GetListOfBranches());
578 while ((sb=(TBranch*)next())) {
579 sb->SetFile(filename);
584 //And Finally clusterizer branch
585 AliPHOSClusterizerv1 * cl = (AliPHOSClusterizerv1*)gime->Clusterizer(branchName) ;
586 TBranch * clusterizerBranch = gAlice->TreeR()->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
587 &cl,bufferSize,splitlevel);
588 clusterizerBranch->SetTitle(branchName);
590 clusterizerBranch->SetFile(filename);
591 TIter next( clusterizerBranch->GetListOfBranches());
593 while ((sb=(TBranch*)next())) {
594 sb->SetFile(filename);
600 clusterizerBranch->Fill() ;
602 gAlice->TreeR()->Write(0,kOverwrite) ;
606 //____________________________________________________________________________
607 void AliPHOSClusterizerv1::MakeClusters()
609 // Steering method to construct the clusters stored in a list of Reconstructed Points
610 // A cluster is defined as a list of neighbour digits
611 TString branchName(GetName()) ;
612 branchName.Remove(branchName.Index(Version())-1) ;
614 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
616 TObjArray * emcRecPoints = gime->EmcRecPoints(branchName) ;
617 TObjArray * cpvRecPoints = gime->CpvRecPoints(branchName) ;
618 emcRecPoints->Delete() ;
619 cpvRecPoints->Delete() ;
621 TClonesArray * digits = gime->Digits(branchName) ;
622 TClonesArray * digitsC = (TClonesArray*)digits->Clone() ;
625 // Clusterization starts
627 TIter nextdigit(digitsC) ;
628 AliPHOSDigit * digit ;
629 Bool_t notremoved = kTRUE ;
631 while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digitsC
632 AliPHOSRecPoint * clu = 0 ;
634 TArrayI clusterdigitslist(1500) ;
637 if (( IsInEmc (digit) && Calibrate(digit->GetAmp()) > fEmcClusteringThreshold ) ||
638 ( IsInPpsd(digit) && Calibrate(digit->GetAmp()) > fPpsdClusteringThreshold ) ||
639 ( IsInCpv (digit) && Calibrate(digit->GetAmp()) > fCpvClusteringThreshold ) ) {
641 Int_t iDigitInCluster = 0 ;
643 if ( IsInEmc(digit) ) {
644 // start a new EMC RecPoint
645 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
646 emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
648 emcRecPoints->AddAt(new AliPHOSEmcRecPoint(), fNumberOfEmcClusters) ;
649 clu = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters) ;
650 fNumberOfEmcClusters++ ;
651 clu->AddDigit(*digit, Calibrate(digit->GetAmp())) ;
652 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
654 digitsC->Remove(digit) ;
658 // start a new PPSD/CPV cluster
659 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
660 cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
663 cpvRecPoints->AddAt(new AliPHOSPpsdRecPoint(),fNumberOfCpvClusters) ;
665 cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(), fNumberOfCpvClusters) ;
667 clu = (AliPHOSPpsdRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters) ;
668 fNumberOfCpvClusters++ ;
669 clu->AddDigit(*digit, Calibrate(digit->GetAmp()) ) ;
670 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
672 digitsC->Remove(digit) ;
675 // Here we remove remaining EMC digits, which cannot make a cluster
678 while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
680 digitsC->Remove(digit) ;
684 notremoved = kFALSE ;
691 AliPHOSDigit * digitN ;
693 while (index < iDigitInCluster){ // scan over digits already in cluster
694 digit = (AliPHOSDigit*)digits->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 digitsC->Remove(digitN) ;
707 case 2 : // too far from each other
716 } // loop over cluster
727 //____________________________________________________________________________
728 void AliPHOSClusterizerv1::MakeUnfolding()
730 // Unfolds clusters using the shape of an ElectroMagnetic shower
731 // Performs unfolding of all EMC/CPV but NOT ppsd clusters
733 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
735 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
736 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
737 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
738 TClonesArray * digits = gime->Digits() ;
740 // Unfold first EMC clusters
741 if(fNumberOfEmcClusters > 0){
743 Int_t nModulesToUnfold = geom->GetNModules() ;
745 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
747 for(index = 0 ; index < numberofNotUnfolded ; index++){
749 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) emcRecPoints->At(index) ;
750 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
753 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
754 Int_t * maxAt = new Int_t[nMultipl] ;
755 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
756 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
758 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
759 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
760 emcRecPoints->Remove(emcRecPoint);
761 emcRecPoints->Compress() ;
763 fNumberOfEmcClusters -- ;
764 numberofNotUnfolded-- ;
768 delete[] maxAtEnergy ;
771 // Unfolding of EMC clusters finished
774 // Unfold now CPV clusters
775 if(fNumberOfCpvClusters > 0){
777 Int_t nModulesToUnfold = geom->GetNCPVModules() ;
779 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
781 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
783 AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) cpvRecPoints->At(index) ;
785 if(recPoint->GetPHOSMod()> nModulesToUnfold)
788 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ;
790 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
791 Int_t * maxAt = new Int_t[nMultipl] ;
792 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
793 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,digits) ;
795 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
796 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
797 cpvRecPoints->Remove(emcRecPoint);
798 cpvRecPoints->Compress() ;
800 numberofCpvNotUnfolded-- ;
801 fNumberOfCpvClusters-- ;
805 delete[] maxAtEnergy ;
808 //Unfolding of Cpv clusters finished
812 //____________________________________________________________________________
813 Double_t AliPHOSClusterizerv1::ShowerShape(Double_t r)
815 // Shape of the shower (see PHOS TDR)
816 // If you change this function, change also the gradient evaluation in ChiSquare()
818 Double_t r4 = r*r*r*r ;
819 Double_t r295 = TMath::Power(r, 2.95) ;
820 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
824 //____________________________________________________________________________
825 void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
828 Float_t * maxAtEnergy)
830 // Performs the unfolding of a cluster with nMax overlapping showers
832 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
833 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
834 TClonesArray * digits = gime->Digits() ;
835 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
836 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
838 Int_t nPar = 3 * nMax ;
839 Float_t * fitparameters = new Float_t[nPar] ;
841 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
843 // Fit failed, return and remove cluster
844 delete[] fitparameters ;
848 // create ufolded rec points and fill them with new energy lists
849 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
850 // and later correct this number in acordance with actual energy deposition
852 Int_t nDigits = iniEmc->GetMultiplicity() ;
853 Float_t * efit = new Float_t[nDigits] ;
854 Float_t xDigit=0.,zDigit=0.,distance=0. ;
855 Float_t xpar=0.,zpar=0.,epar=0. ;
857 AliPHOSDigit * digit = 0 ;
858 Int_t * emcDigits = iniEmc->GetDigitsList() ;
862 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
863 digit = (AliPHOSDigit*) digits->At(emcDigits[iDigit] ) ;
864 geom->AbsToRelNumbering(digit->GetId(), relid) ;
865 geom->RelPosInModule(relid, xDigit, zDigit) ;
869 while(iparam < nPar ){
870 xpar = fitparameters[iparam] ;
871 zpar = fitparameters[iparam+1] ;
872 epar = fitparameters[iparam+2] ;
874 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
875 distance = TMath::Sqrt(distance) ;
876 efit[iDigit] += epar * ShowerShape(distance) ;
881 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
882 // so that energy deposited in each cell is distributed betwin new clusters proportionally
883 // to its contribution to efit
885 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
889 while(iparam < nPar ){
890 xpar = fitparameters[iparam] ;
891 zpar = fitparameters[iparam+1] ;
892 epar = fitparameters[iparam+2] ;
895 AliPHOSEmcRecPoint * emcRP = 0 ;
897 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
899 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
900 emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
902 (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint() ;
903 emcRP = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters);
904 fNumberOfEmcClusters++ ;
906 else{//create new entries in fCpvRecPoints
907 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
908 cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
910 (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint() ;
911 emcRP = (AliPHOSEmcRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters);
912 fNumberOfCpvClusters++ ;
916 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
917 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
918 geom->AbsToRelNumbering(digit->GetId(), relid) ;
919 geom->RelPosInModule(relid, xDigit, zDigit) ;
920 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
921 distance = TMath::Sqrt(distance) ;
922 ratio = epar * ShowerShape(distance) / efit[iDigit] ;
923 eDigit = emcEnergies[iDigit] * ratio ;
924 emcRP->AddDigit( *digit, eDigit ) ;
928 delete[] fitparameters ;
933 //_____________________________________________________________________________
934 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
936 // Calculates the Chi square for the cluster unfolding minimization
937 // Number of parameters, Gradient, Chi squared, parameters, what to do
940 TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
942 AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint*) toMinuit->At(0) ;
943 TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
947 // AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint *) gMinuit->GetObjectFit() ; // EmcRecPoint to fit
949 Int_t * emcDigits = emcRP->GetDigitsList() ;
951 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
953 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
955 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
960 for(iparam = 0 ; iparam < nPar ; iparam++)
961 Grad[iparam] = 0 ; // Will evaluate gradient
965 AliPHOSDigit * digit ;
968 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
970 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
976 geom->AbsToRelNumbering(digit->GetId(), relid) ;
978 geom->RelPosInModule(relid, xDigit, zDigit) ;
980 if(iflag == 2){ // calculate gradient
983 while(iParam < nPar ){
984 Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
986 distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
987 distance = TMath::Sqrt( distance ) ;
989 efit += x[iParam] * ShowerShape(distance) ;
992 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
994 while(iParam < nPar ){
995 Double_t xpar = x[iParam] ;
996 Double_t zpar = x[iParam+1] ;
997 Double_t epar = x[iParam+2] ;
998 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
999 Double_t shape = sum * ShowerShape(dr) ;
1000 Double_t r4 = dr*dr*dr*dr ;
1001 Double_t r295 = TMath::Power(dr,2.95) ;
1002 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
1003 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
1005 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
1007 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
1009 Grad[iParam] += shape ; // Derivative over energy
1016 while(iparam < nPar ){
1017 Double_t xpar = x[iparam] ;
1018 Double_t zpar = x[iparam+1] ;
1019 Double_t epar = x[iparam+2] ;
1021 Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1022 distance = TMath::Sqrt(distance) ;
1023 efit += epar * ShowerShape(distance) ;
1026 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
1027 // Here we assume, that sigma = sqrt(E)
1032 //____________________________________________________________________________
1033 void AliPHOSClusterizerv1::Print(Option_t * option)const
1035 // Print clusterizer parameters
1037 if( strcmp(GetName(), "") !=0 ){
1041 TString taskName(GetName()) ;
1042 taskName.ReplaceAll(Version(), "") ;
1044 cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl
1045 << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl
1046 << " Branch: " << fDigitsBranchTitle.Data() << endl
1048 << " EMC Clustering threshold = " << fEmcClusteringThreshold << endl
1049 << " EMC Local Maximum cut = " << fEmcLocMaxCut << endl
1050 << " EMC Logarothmic weight = " << fW0 << endl
1052 << " CPV Clustering threshold = " << fCpvClusteringThreshold << endl
1053 << " CPV Local Maximum cut = " << fCpvLocMaxCut << endl
1054 << " CPV Logarothmic weight = " << fW0CPV << endl
1056 << " PPSD Clustering threshold = " << fPpsdClusteringThreshold << endl;
1058 cout << " Unfolding on " << endl ;
1060 cout << " Unfolding off " << endl ;
1062 cout << "------------------------------------------------------------------" <<endl ;
1065 cout << " AliPHOSClusterizerv1 not initialized " << endl ;
1067 //____________________________________________________________________________
1068 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1070 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1072 TObjArray * emcRecPoints = AliPHOSGetter::GetInstance()->EmcRecPoints() ;
1073 TObjArray * cpvRecPoints = AliPHOSGetter::GetInstance()->CpvRecPoints() ;
1075 cout << "AliPHOSClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1076 cout << " Found "<< emcRecPoints->GetEntriesFast() << " EMC Rec Points and "
1077 << cpvRecPoints->GetEntriesFast() << " CPV RecPoints" << endl ;
1079 fRecPointsInRun += emcRecPoints->GetEntriesFast() ;
1080 fRecPointsInRun += cpvRecPoints->GetEntriesFast() ;
1082 if(strstr(option,"all")) {
1083 cout << "EMC clusters " << endl ;
1095 << " Primaries list " << endl;
1098 for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1099 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ;
1101 cout << setw(6) << rp->GetIndexInList() << " ";
1102 cout << setw(6) << rp->GetEnergy() << " ";
1103 cout << setw(6) << rp->GetMultiplicity()<< " ";
1104 cout << setw(6) << rp->GetPHOSMod() << " ";
1107 rp->GetLocalPosition(locpos);
1108 cout << setw(8) << locpos.X() << " ";
1109 cout << setw(8) << locpos.Y() << " ";
1110 cout << setw(8) << locpos.Z() << " ";
1113 rp->GetElipsAxis(lambda);
1114 cout << setw(10)<< lambda[0] << " ";
1115 cout << setw(10)<< lambda[1] << " ";
1120 primaries = rp->GetPrimaries(nprimaries);
1121 cout << setw(8) << primaries << " ";
1123 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1124 cout << setw(4) << primaries[iprimary] << " ";
1129 //Now plot CPV/PPSD recPoints
1130 cout << "EMC clusters " << endl ;
1139 << " Primaries list " << endl;
1141 for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1142 AliPHOSRecPoint * rp = (AliPHOSRecPoint * )cpvRecPoints->At(index) ;
1143 cout << setw(6) << rp->GetIndexInList() << " ";
1144 cout << setw(6) << rp->GetPHOSMod() << " ";
1146 if( (strcmp(rp->ClassName() , "AliPHOSPpsdRecPoint" )) == 0){
1147 AliPHOSPpsdRecPoint * ppsd = (AliPHOSPpsdRecPoint*) rp ;
1157 rp->GetLocalPosition(locpos);
1158 cout << setw(8) << locpos.X() << " ";
1159 cout << setw(8) << locpos.Y() << " ";
1160 cout << setw(8) << locpos.Z() << " ";
1164 primaries = rp->GetPrimaries(nprimaries);
1165 cout << setw(8) << primaries << " ";
1167 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1168 cout << setw(4) << primaries[iprimary] << " ";
1173 cout << "-------------------------------------------------"<<endl ;