1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 /* History of cvs commits:
21 * Revision 1.89 2006/03/13 14:05:42 kharlov
22 * Calibration objects for EMC and CPV
24 * Revision 1.88 2006/01/11 08:54:52 hristov
25 * Additional protection in case no calibration entry was found
27 * Revision 1.87 2005/11/22 08:46:43 kharlov
28 * Updated to new CDB (Boris Polichtchouk)
30 * Revision 1.86 2005/11/14 21:52:43 hristov
33 * Revision 1.85 2005/09/27 16:08:08 hristov
34 * New version of CDB storage framework (A.Colla)
36 * Revision 1.84 2005/09/21 10:02:47 kharlov
37 * Reading calibration from CDB (Boris Polichtchouk)
39 * Revision 1.82 2005/09/02 15:43:13 kharlov
40 * Add comments in GetCalibrationParameters and Calibrate
42 * Revision 1.81 2005/09/02 14:32:07 kharlov
43 * Calibration of raw data
45 * Revision 1.80 2005/08/24 15:31:36 kharlov
46 * Setting raw digits flag
48 * Revision 1.79 2005/07/25 15:53:53 kharlov
51 * Revision 1.78 2005/05/28 14:19:04 schutz
52 * Compilation warnings fixed by T.P.
56 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
57 //////////////////////////////////////////////////////////////////////////////
58 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
59 // unfolds the clusters having several local maxima.
60 // Results are stored in TreeR#, branches PHOSEmcRP (EMC recPoints),
61 // PHOSCpvRP (CPV RecPoints) and AliPHOSClusterizer (Clusterizer with all
62 // parameters including input digits branch title, thresholds etc.)
63 // This TTask is normally called from Reconstructioner, but can as well be used in
66 // root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root", "recpointsname", "digitsname")
67 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
68 // // reads gAlice from header file "galice.root", uses digits stored in the branch names "digitsname" (default = "Default")
69 // // and saves recpoints in branch named "recpointsname" (default = "digitsname")
70 // root [1] cl->ExecuteTask()
71 // //finds RecPoints in all events stored in galice.root
72 // root [2] cl->SetDigitsBranch("digits2")
73 // //sets another title for Digitis (input) branch
74 // root [3] cl->SetRecPointsBranch("recp2")
75 // //sets another title four output branches
76 // root [4] cl->SetEmcLocalMaxCut(0.03)
77 // //set clusterization parameters
78 // root [5] cl->ExecuteTask("deb all time")
79 // //once more finds RecPoints options are
80 // // deb - print number of found rec points
81 // // deb all - print number of found RecPoints and some their characteristics
82 // // time - print benchmarking results
84 // --- ROOT system ---
89 #include "TBenchmark.h"
91 // --- Standard library ---
93 // --- AliRoot header files ---
95 #include "AliPHOSGetter.h"
96 #include "AliPHOSGeometry.h"
97 #include "AliPHOSClusterizerv1.h"
98 #include "AliPHOSEmcRecPoint.h"
99 #include "AliPHOSCpvRecPoint.h"
100 #include "AliPHOSDigit.h"
101 #include "AliPHOSDigitizer.h"
102 #include "AliPHOSCalibrationDB.h"
103 #include "AliCDBManager.h"
104 #include "AliCDBStorage.h"
105 #include "AliCDBEntry.h"
107 ClassImp(AliPHOSClusterizerv1)
109 //____________________________________________________________________________
110 AliPHOSClusterizerv1::AliPHOSClusterizerv1() : AliPHOSClusterizer()
112 // default ctor (to be used mainly by Streamer)
115 fDefaultInit = kTRUE ;
118 //____________________________________________________________________________
119 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const TString alirunFileName, const TString eventFolderName)
120 :AliPHOSClusterizer(alirunFileName, eventFolderName)
122 // ctor with the indication of the file where header Tree and digits Tree are stored
126 fDefaultInit = kFALSE ;
129 //____________________________________________________________________________
130 AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
135 //____________________________________________________________________________
136 const TString AliPHOSClusterizerv1::BranchName() const
141 //____________________________________________________________________________
142 Float_t AliPHOSClusterizerv1::Calibrate(Int_t amp, Int_t absId)
144 // Convert digitized amplitude into energy.
145 // Calibration parameters are taken from calibration data base for raw data,
146 // or from digitizer parameters for simulated data.
150 AliPHOSGetter *gime = AliPHOSGetter::Instance();
151 gime->PHOSGeometry()->AbsToRelNumbering(absId,relId) ;
152 Int_t module = relId[0];
153 Int_t column = relId[3];
154 Int_t row = relId[2];
155 if(absId <= fEmcCrystals) { //calibrate as EMC
156 fADCchanelEmc = fCalibData->GetADCchannelEmc (module,column,row);
157 fADCpedestalEmc = fCalibData->GetADCpedestalEmc(module,column,row);
158 return fADCpedestalEmc + amp*fADCchanelEmc ;
160 else { //calibrate as CPV
161 fADCchanelCpv = fCalibData->GetADCchannelCpv (module,column,row);
162 fADCpedestalCpv = fCalibData->GetADCpedestalCpv(module,column,row);
163 return fADCpedestalCpv + amp*fADCchanelCpv ;
167 if(absId <= fEmcCrystals) //calibrate as EMC
168 return fADCpedestalEmc + amp*fADCchanelEmc ;
169 else //calibrate as CPV
170 return fADCpedestalCpv+ amp*fADCchanelCpv ;
174 //____________________________________________________________________________
175 void AliPHOSClusterizerv1::Exec(Option_t *option)
177 // Steering method to perform clusterization for events
178 // in the range from fFirstEvent to fLastEvent.
179 // This range is optionally set by SetEventRange().
180 // if fLastEvent=-1 (by default), then process events until the end.
182 if(strstr(option,"tim"))
183 gBenchmark->Start("PHOSClusterizer");
185 if(strstr(option,"print")) {
190 GetCalibrationParameters() ;
192 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
194 gime->SetRawDigits(kFALSE);
196 gime->SetRawDigits(kTRUE);
198 if (fLastEvent == -1)
199 fLastEvent = gime->MaxEvent() - 1 ;
201 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // one event at the time
202 Int_t nEvents = fLastEvent - fFirstEvent + 1;
205 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
207 gime->Event(ievent ,"D"); // Read digits from simulated data
209 gime->Event(fRawReader,"W"); // Read digits from raw data
211 fNumberOfEmcClusters = fNumberOfCpvClusters = 0 ;
220 if(strstr(option,"deb"))
221 PrintRecPoints(option) ;
223 //increment the total number of recpoints per run
224 fRecPointsInRun += gime->EmcRecPoints()->GetEntriesFast() ;
225 fRecPointsInRun += gime->CpvRecPoints()->GetEntriesFast() ;
228 if(fWrite) //do not unload in "on flight" mode
231 if(strstr(option,"tim")){
232 gBenchmark->Stop("PHOSClusterizer");
233 AliInfo(Form("took %f seconds for Clusterizing %f seconds per event \n",
234 gBenchmark->GetCpuTime("PHOSClusterizer"),
235 gBenchmark->GetCpuTime("PHOSClusterizer")/nEvents )) ;
239 //____________________________________________________________________________
240 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
241 Int_t nPar, Float_t * fitparameters) const
243 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
244 // The initial values for fitting procedure are set equal to the positions of local maxima.
245 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
248 AliPHOSGetter * gime = AliPHOSGetter::Instance();
249 TClonesArray * digits = gime->Digits();
252 gMinuit->mncler(); // Reset Minuit's list of paramters
253 gMinuit->SetPrintLevel(-1) ; // No Printout
254 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
255 // To set the address of the minimization function
257 TList * toMinuit = new TList();
258 toMinuit->AddAt(emcRP,0) ;
259 toMinuit->AddAt(digits,1) ;
261 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
263 // filling initial values for fit parameters
264 AliPHOSDigit * digit ;
268 Int_t nDigits = (Int_t) nPar / 3 ;
272 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
274 for(iDigit = 0; iDigit < nDigits; iDigit++){
275 digit = maxAt[iDigit];
280 geom->AbsToRelNumbering(digit->GetId(), relid) ;
281 geom->RelPosInModule(relid, x, z) ;
283 Float_t energy = maxAtEnergy[iDigit] ;
285 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
288 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : x = %f\n", x ) ;
291 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
294 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : z =%f\n", z ) ;
297 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
300 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : energy = %f\n", energy ) ;
305 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
310 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
311 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
312 gMinuit->SetMaxIterations(5);
313 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
315 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
317 if(ierflg == 4){ // Minimum not found
318 Warning("FindFit", "PHOS Unfolding fit not converged, cluster abandoned\n" );
321 for(index = 0; index < nPar; index++){
324 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
325 fitparameters[index] = val ;
333 //____________________________________________________________________________
334 void AliPHOSClusterizerv1::GetCalibrationParameters()
336 // Set calibration parameters:
337 // if calibration database exists, they are read from database,
338 // otherwise, they are taken from digitizer.
340 // It is a user responsilibity to open CDB before reconstruction, for example:
341 // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
343 AliPHOSGetter * gime = AliPHOSGetter::Instance();
344 // fCalibData = new AliPHOSCalibData(gAlice->GetRunNumber()); //original
346 fCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
350 if ( !gime->Digitizer() )
351 gime->LoadDigitizer();
352 AliPHOSDigitizer * dig = gime->Digitizer();
353 fADCchanelEmc = dig->GetEMCchannel() ;
354 fADCpedestalEmc = dig->GetEMCpedestal();
356 fADCchanelCpv = dig->GetCPVchannel() ;
357 fADCpedestalCpv = dig->GetCPVpedestal() ;
361 //____________________________________________________________________________
362 void AliPHOSClusterizerv1::Init()
364 // Make all memory allocations which can not be done in default constructor.
365 // Attach the Clusterizer task to the list of PHOS tasks
367 AliPHOSGetter* gime = AliPHOSGetter::Instance() ;
369 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
371 AliPHOSGeometry * geom = gime->PHOSGeometry();
373 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
376 gMinuit = new TMinuit(100);
378 if ( !gime->Clusterizer() ) {
379 gime->PostClusterizer(this);
383 //____________________________________________________________________________
384 void AliPHOSClusterizerv1::InitParameters()
387 fNumberOfCpvClusters = 0 ;
388 fNumberOfEmcClusters = 0 ;
390 fCpvClusteringThreshold = 0.0;
391 fEmcClusteringThreshold = 0.2;
393 fEmcLocMaxCut = 0.03 ;
394 fCpvLocMaxCut = 0.03 ;
402 fEmcTimeGate = 1.e-8 ;
406 fRecPointsInRun = 0 ;
412 SetEventRange(0,-1) ;
415 //____________________________________________________________________________
416 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
418 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
420 // = 2 are not neighbour but do not continue searching
421 // neighbours are defined as digits having at least a common vertex
422 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
423 // which is compared to a digit (d2) not yet in a cluster
425 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
430 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
433 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
435 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module
436 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
437 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
439 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
440 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
444 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
445 rv = 2; // Difference in row numbers is too large to look further
451 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
458 //____________________________________________________________________________
459 void AliPHOSClusterizerv1::CleanDigits(TClonesArray * digits){
460 for(Int_t i=0; i<digits->GetEntriesFast(); i++){
461 AliPHOSDigit * digit = static_cast<AliPHOSDigit*>(digits->At(i)) ;
462 Float_t cut = IsInEmc(digit) ? fEmcMinE : fCpvMinE ;
463 if(Calibrate(digit->GetAmp(),digit->GetId()) < cut)
464 digits->RemoveAt(i) ;
467 for (Int_t i = 0 ; i < digits->GetEntriesFast() ; i++) {
468 AliPHOSDigit *digit = static_cast<AliPHOSDigit*>( digits->At(i) ) ;
469 digit->SetIndexInList(i) ;
473 //____________________________________________________________________________
474 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
476 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
479 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
481 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
483 if(digit->GetId() <= nEMC ) rv = kTRUE;
488 //____________________________________________________________________________
489 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
491 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
495 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
497 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
499 if(digit->GetId() > nEMC ) rv = kTRUE;
504 //____________________________________________________________________________
505 void AliPHOSClusterizerv1::Unload()
507 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
508 gime->PhosLoader()->UnloadDigits() ;
509 gime->PhosLoader()->UnloadRecPoints() ;
512 //____________________________________________________________________________
513 void AliPHOSClusterizerv1::WriteRecPoints()
516 // Creates new branches with given title
517 // fills and writes into TreeR.
519 AliPHOSGetter * gime = AliPHOSGetter::Instance();
521 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
522 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
523 TClonesArray * digits = gime->Digits() ;
527 //Evaluate position, dispersion and other RecPoint properties..
528 Int_t nEmc = emcRecPoints->GetEntriesFast();
529 for(index = 0; index < nEmc; index++){
530 AliPHOSEmcRecPoint * rp = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) );
531 rp->Purify(fEmcMinE) ;
532 if(rp->GetMultiplicity()>0.) //If this RP is not empty
533 rp->EvalAll(fW0,digits) ;
535 emcRecPoints->RemoveAt(index) ;
539 emcRecPoints->Compress() ;
540 emcRecPoints->Sort() ;
541 // emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ;
542 for(index = 0; index < emcRecPoints->GetEntries(); index++){
543 dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) )->SetIndexInList(index) ;
546 //Now the same for CPV
547 for(index = 0; index < cpvRecPoints->GetEntries(); index++){
548 AliPHOSCpvRecPoint * rp = dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(index) );
549 rp->EvalAll(fW0CPV,digits) ;
551 cpvRecPoints->Sort() ;
553 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
554 dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(index) )->SetIndexInList(index) ;
556 cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
558 if(fWrite){ //We write TreeR
559 TTree * treeR = gime->TreeR();
561 Int_t bufferSize = 32000 ;
562 Int_t splitlevel = 0 ;
565 TBranch * emcBranch = treeR->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
566 emcBranch->SetTitle(BranchName());
569 TBranch * cpvBranch = treeR->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
570 cpvBranch->SetTitle(BranchName());
575 gime->WriteRecPoints("OVERWRITE");
576 gime->WriteClusterizer("OVERWRITE");
580 //____________________________________________________________________________
581 void AliPHOSClusterizerv1::MakeClusters()
583 // Steering method to construct the clusters stored in a list of Reconstructed Points
584 // A cluster is defined as a list of neighbour digits
587 AliPHOSGetter * gime = AliPHOSGetter::Instance();
589 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
590 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
592 emcRecPoints->Delete() ;
593 cpvRecPoints->Delete() ;
595 TClonesArray * digits = gime->Digits() ;
597 //Remove digits below threshold
598 CleanDigits(digits) ;
601 TClonesArray * digitsC = static_cast<TClonesArray*>( digits->Clone() ) ;
604 // Clusterization starts
606 TIter nextdigit(digitsC) ;
607 AliPHOSDigit * digit ;
608 Bool_t notremoved = kTRUE ;
610 while ( (digit = dynamic_cast<AliPHOSDigit *>( nextdigit()) ) ) { // scan over the list of digitsC
613 AliPHOSRecPoint * clu = 0 ;
615 TArrayI clusterdigitslist(1500) ;
618 if (( IsInEmc (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fEmcClusteringThreshold ) ||
619 ( IsInCpv (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fCpvClusteringThreshold ) ) {
620 Int_t iDigitInCluster = 0 ;
622 if ( IsInEmc(digit) ) {
623 // start a new EMC RecPoint
624 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
625 emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
627 emcRecPoints->AddAt(new AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
628 clu = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(fNumberOfEmcClusters) ) ;
629 fNumberOfEmcClusters++ ;
630 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId())) ;
631 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
633 digitsC->Remove(digit) ;
637 // start a new CPV cluster
638 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
639 cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
641 cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
643 clu = dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(fNumberOfCpvClusters) ) ;
644 fNumberOfCpvClusters++ ;
645 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId()) ) ;
646 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
648 digitsC->Remove(digit) ;
651 // Here we remove remaining EMC digits, which cannot make a cluster
654 while( ( digit = dynamic_cast<AliPHOSDigit *>( nextdigit() ) ) ) {
656 digitsC->Remove(digit) ;
660 notremoved = kFALSE ;
667 AliPHOSDigit * digitN ;
669 while (index < iDigitInCluster){ // scan over digits already in cluster
670 digit = dynamic_cast<AliPHOSDigit*>( digits->At(clusterdigitslist[index]) ) ;
672 while ( (digitN = dynamic_cast<AliPHOSDigit *>( nextdigit() ) ) ) { // scan over the reduced list of digits
673 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
675 case 0 : // not a neighbour
677 case 1 : // are neighbours
678 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->GetId() ) ) ;
679 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
681 digitsC->Remove(digitN) ;
683 case 2 : // too far from each other
692 } // loop over cluster
703 //____________________________________________________________________________
704 void AliPHOSClusterizerv1::MakeUnfolding()
706 // Unfolds clusters using the shape of an ElectroMagnetic shower
707 // Performs unfolding of all EMC/CPV clusters
709 AliPHOSGetter * gime = AliPHOSGetter::Instance();
711 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
713 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
714 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
715 TClonesArray * digits = gime->Digits() ;
717 // Unfold first EMC clusters
718 if(fNumberOfEmcClusters > 0){
720 Int_t nModulesToUnfold = geom->GetNModules() ;
722 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
724 for(index = 0 ; index < numberofNotUnfolded ; index++){
726 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) ) ;
727 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
730 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
731 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
732 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
733 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
735 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
736 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
737 emcRecPoints->Remove(emcRecPoint);
738 emcRecPoints->Compress() ;
740 fNumberOfEmcClusters -- ;
741 numberofNotUnfolded-- ;
744 emcRecPoint->SetNExMax(1) ; //Only one local maximum
748 delete[] maxAtEnergy ;
751 // Unfolding of EMC clusters finished
754 // Unfold now CPV clusters
755 if(fNumberOfCpvClusters > 0){
757 Int_t nModulesToUnfold = geom->GetNModules() ;
759 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
761 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
763 AliPHOSRecPoint * recPoint = dynamic_cast<AliPHOSRecPoint *>( cpvRecPoints->At(index) ) ;
765 if(recPoint->GetPHOSMod()> nModulesToUnfold)
768 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint*>(recPoint) ;
770 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
771 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
772 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
773 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,digits) ;
775 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
776 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
777 cpvRecPoints->Remove(emcRecPoint);
778 cpvRecPoints->Compress() ;
780 numberofCpvNotUnfolded-- ;
781 fNumberOfCpvClusters-- ;
785 delete[] maxAtEnergy ;
788 //Unfolding of Cpv clusters finished
792 //____________________________________________________________________________
793 Double_t AliPHOSClusterizerv1::ShowerShape(Double_t r)
795 // Shape of the shower (see PHOS TDR)
796 // If you change this function, change also the gradient evaluation in ChiSquare()
798 Double_t r4 = r*r*r*r ;
799 Double_t r295 = TMath::Power(r, 2.95) ;
800 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
804 //____________________________________________________________________________
805 void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
807 AliPHOSDigit ** maxAt,
808 Float_t * maxAtEnergy)
810 // Performs the unfolding of a cluster with nMax overlapping showers
812 AliPHOSGetter * gime = AliPHOSGetter::Instance();
814 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
816 const TClonesArray * digits = gime->Digits() ;
817 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
818 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
820 Int_t nPar = 3 * nMax ;
821 Float_t * fitparameters = new Float_t[nPar] ;
823 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
825 // Fit failed, return and remove cluster
826 iniEmc->SetNExMax(-1) ;
827 delete[] fitparameters ;
831 // create ufolded rec points and fill them with new energy lists
832 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
833 // and later correct this number in acordance with actual energy deposition
835 Int_t nDigits = iniEmc->GetMultiplicity() ;
836 Float_t * efit = new Float_t[nDigits] ;
837 Float_t xDigit=0.,zDigit=0.,distance=0. ;
838 Float_t xpar=0.,zpar=0.,epar=0. ;
840 AliPHOSDigit * digit = 0 ;
841 Int_t * emcDigits = iniEmc->GetDigitsList() ;
845 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
846 digit = dynamic_cast<AliPHOSDigit*>( digits->At(emcDigits[iDigit] ) ) ;
847 geom->AbsToRelNumbering(digit->GetId(), relid) ;
848 geom->RelPosInModule(relid, xDigit, zDigit) ;
852 while(iparam < nPar ){
853 xpar = fitparameters[iparam] ;
854 zpar = fitparameters[iparam+1] ;
855 epar = fitparameters[iparam+2] ;
857 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
858 distance = TMath::Sqrt(distance) ;
859 efit[iDigit] += epar * ShowerShape(distance) ;
864 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
865 // so that energy deposited in each cell is distributed betwin new clusters proportionally
866 // to its contribution to efit
868 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
872 while(iparam < nPar ){
873 xpar = fitparameters[iparam] ;
874 zpar = fitparameters[iparam+1] ;
875 epar = fitparameters[iparam+2] ;
878 AliPHOSEmcRecPoint * emcRP = 0 ;
880 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
882 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
883 emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
885 (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
886 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(fNumberOfEmcClusters) ) ;
887 fNumberOfEmcClusters++ ;
888 emcRP->SetNExMax((Int_t)nPar/3) ;
890 else{//create new entries in fCpvRecPoints
891 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
892 cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
894 (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
895 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( cpvRecPoints->At(fNumberOfCpvClusters) ) ;
896 fNumberOfCpvClusters++ ;
900 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
901 digit = dynamic_cast<AliPHOSDigit*>( digits->At( emcDigits[iDigit] ) ) ;
902 geom->AbsToRelNumbering(digit->GetId(), relid) ;
903 geom->RelPosInModule(relid, xDigit, zDigit) ;
904 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
905 distance = TMath::Sqrt(distance) ;
906 ratio = epar * ShowerShape(distance) / efit[iDigit] ;
907 eDigit = emcEnergies[iDigit] * ratio ;
908 emcRP->AddDigit( *digit, eDigit ) ;
912 delete[] fitparameters ;
917 //_____________________________________________________________________________
918 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
920 // Calculates the Chi square for the cluster unfolding minimization
921 // Number of parameters, Gradient, Chi squared, parameters, what to do
923 TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;
925 AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint*>( toMinuit->At(0) ) ;
926 TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) ) ;
930 // AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( gMinuit->GetObjectFit() ) ; // EmcRecPoint to fit
932 Int_t * emcDigits = emcRP->GetDigitsList() ;
934 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
936 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
938 const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
943 for(iparam = 0 ; iparam < nPar ; iparam++)
944 Grad[iparam] = 0 ; // Will evaluate gradient
948 AliPHOSDigit * digit ;
951 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
953 digit = dynamic_cast<AliPHOSDigit*>( digits->At( emcDigits[iDigit] ) );
959 geom->AbsToRelNumbering(digit->GetId(), relid) ;
961 geom->RelPosInModule(relid, xDigit, zDigit) ;
963 if(iflag == 2){ // calculate gradient
966 while(iParam < nPar ){
967 Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
969 distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
970 distance = TMath::Sqrt( distance ) ;
972 efit += x[iParam] * ShowerShape(distance) ;
975 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
977 while(iParam < nPar ){
978 Double_t xpar = x[iParam] ;
979 Double_t zpar = x[iParam+1] ;
980 Double_t epar = x[iParam+2] ;
981 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
982 Double_t shape = sum * ShowerShape(dr) ;
983 Double_t r4 = dr*dr*dr*dr ;
984 Double_t r295 = TMath::Power(dr,2.95) ;
985 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
986 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
988 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
990 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
992 Grad[iParam] += shape ; // Derivative over energy
999 while(iparam < nPar ){
1000 Double_t xpar = x[iparam] ;
1001 Double_t zpar = x[iparam+1] ;
1002 Double_t epar = x[iparam+2] ;
1004 Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1005 distance = TMath::Sqrt(distance) ;
1006 efit += epar * ShowerShape(distance) ;
1009 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
1010 // Here we assume, that sigma = sqrt(E)
1015 //____________________________________________________________________________
1016 void AliPHOSClusterizerv1::Print(const Option_t *)const
1018 // Print clusterizer parameters
1021 TString taskName(GetName()) ;
1022 taskName.ReplaceAll(Version(), "") ;
1024 if( strcmp(GetName(), "") !=0 ) {
1026 message = "\n--------------- %s %s -----------\n" ;
1027 message += "Clusterizing digits from the file: %s\n" ;
1028 message += " Branch: %s\n" ;
1029 message += " EMC Clustering threshold = %f\n" ;
1030 message += " EMC Local Maximum cut = %f\n" ;
1031 message += " EMC Logarothmic weight = %f\n" ;
1032 message += " CPV Clustering threshold = %f\n" ;
1033 message += " CPV Local Maximum cut = %f\n" ;
1034 message += " CPV Logarothmic weight = %f\n" ;
1036 message += " Unfolding on\n" ;
1038 message += " Unfolding off\n" ;
1040 message += "------------------------------------------------------------------" ;
1043 message = " AliPHOSClusterizerv1 not initialized " ;
1045 AliInfo(Form("%s, %s %s %s %s %s %s %s %s %s %s", message.Data(),
1050 fEmcClusteringThreshold,
1053 fCpvClusteringThreshold,
1059 //____________________________________________________________________________
1060 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1062 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1064 AliPHOSGetter * gime = AliPHOSGetter::Instance();
1066 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
1067 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
1069 AliInfo(Form("\nevent %d \n Found %d EMC RecPoints and %d CPV RecPoints",
1070 gAlice->GetEvNumber(),
1071 emcRecPoints->GetEntriesFast(),
1072 cpvRecPoints->GetEntriesFast() )) ;
1074 fRecPointsInRun += emcRecPoints->GetEntriesFast() ;
1075 fRecPointsInRun += cpvRecPoints->GetEntriesFast() ;
1078 if(strstr(option,"all")) {
1079 printf("\n EMC clusters \n") ;
1080 printf("Index Ene(MeV) Multi Module X Y Z Lambdas_1 Lambda_2 # of prim Primaries list\n") ;
1082 for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1083 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ;
1085 rp->GetLocalPosition(locpos);
1087 rp->GetElipsAxis(lambda);
1090 primaries = rp->GetPrimaries(nprimaries);
1091 printf("\n%6d %8.2f %3d %2d %4.1f %4.1f %4.1f %4f %4f %2d : ",
1092 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(), rp->GetPHOSMod(),
1093 locpos.X(), locpos.Y(), locpos.Z(), lambda[0], lambda[1], nprimaries) ;
1095 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
1096 printf("%d ", primaries[iprimary] ) ;
1101 //Now plot CPV recPoints
1102 printf("\n CPV clusters \n") ;
1103 printf("Index Ene(MeV) Module X Y Z \n") ;
1104 for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1105 AliPHOSCpvRecPoint * rp = (AliPHOSCpvRecPoint * )cpvRecPoints->At(index) ;
1108 rp->GetLocalPosition(locpos);
1110 printf("\n%6d %8.2f %2d %4.1f %4.1f %4.1f \n",
1111 rp->GetIndexInList(), rp->GetEnergy(), rp->GetPHOSMod(),
1112 locpos.X(), locpos.Y(), locpos.Z()) ;