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.92 2006/04/29 20:26:46 hristov
22 * Separate EMC and CPV calibration (Yu.Kharlov)
24 * Revision 1.91 2006/04/22 10:30:17 hristov
25 * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
27 * Revision 1.90 2006/04/11 15:22:59 hristov
28 * run number in query set to -1: forces AliCDBManager to use its run number (A.Colla)
30 * Revision 1.89 2006/03/13 14:05:42 kharlov
31 * Calibration objects for EMC and CPV
33 * Revision 1.88 2006/01/11 08:54:52 hristov
34 * Additional protection in case no calibration entry was found
36 * Revision 1.87 2005/11/22 08:46:43 kharlov
37 * Updated to new CDB (Boris Polichtchouk)
39 * Revision 1.86 2005/11/14 21:52:43 hristov
42 * Revision 1.85 2005/09/27 16:08:08 hristov
43 * New version of CDB storage framework (A.Colla)
45 * Revision 1.84 2005/09/21 10:02:47 kharlov
46 * Reading calibration from CDB (Boris Polichtchouk)
48 * Revision 1.82 2005/09/02 15:43:13 kharlov
49 * Add comments in GetCalibrationParameters and Calibrate
51 * Revision 1.81 2005/09/02 14:32:07 kharlov
52 * Calibration of raw data
54 * Revision 1.80 2005/08/24 15:31:36 kharlov
55 * Setting raw digits flag
57 * Revision 1.79 2005/07/25 15:53:53 kharlov
60 * Revision 1.78 2005/05/28 14:19:04 schutz
61 * Compilation warnings fixed by T.P.
65 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
66 //////////////////////////////////////////////////////////////////////////////
67 // Clusterization class. Performs clusterization (collects neighbouring active cells) and
68 // unfolds the clusters having several local maxima.
69 // Results are stored in TreeR#, branches PHOSEmcRP (EMC recPoints),
70 // PHOSCpvRP (CPV RecPoints) and AliPHOSClusterizer (Clusterizer with all
71 // parameters including input digits branch title, thresholds etc.)
72 // This TTask is normally called from Reconstructioner, but can as well be used in
75 // root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root", "recpointsname", "digitsname")
76 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
77 // // reads gAlice from header file "galice.root", uses digits stored in the branch names "digitsname" (default = "Default")
78 // // and saves recpoints in branch named "recpointsname" (default = "digitsname")
79 // root [1] cl->ExecuteTask()
80 // //finds RecPoints in all events stored in galice.root
81 // root [2] cl->SetDigitsBranch("digits2")
82 // //sets another title for Digitis (input) branch
83 // root [3] cl->SetRecPointsBranch("recp2")
84 // //sets another title four output branches
85 // root [4] cl->SetEmcLocalMaxCut(0.03)
86 // //set clusterization parameters
87 // root [5] cl->ExecuteTask("deb all time")
88 // //once more finds RecPoints options are
89 // // deb - print number of found rec points
90 // // deb all - print number of found RecPoints and some their characteristics
91 // // time - print benchmarking results
93 // --- ROOT system ---
98 #include "TBenchmark.h"
100 // --- Standard library ---
102 // --- AliRoot header files ---
104 #include "AliPHOSGetter.h"
105 #include "AliPHOSGeometry.h"
106 #include "AliPHOSClusterizerv1.h"
107 #include "AliPHOSEmcRecPoint.h"
108 #include "AliPHOSCpvRecPoint.h"
109 #include "AliPHOSDigit.h"
110 #include "AliPHOSDigitizer.h"
111 #include "AliPHOSCalibrationDB.h"
112 #include "AliCDBManager.h"
113 #include "AliCDBStorage.h"
114 #include "AliCDBEntry.h"
116 ClassImp(AliPHOSClusterizerv1)
118 //____________________________________________________________________________
119 AliPHOSClusterizerv1::AliPHOSClusterizerv1() : AliPHOSClusterizer()
121 // default ctor (to be used mainly by Streamer)
124 fDefaultInit = kTRUE ;
127 //____________________________________________________________________________
128 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const TString alirunFileName, const TString eventFolderName)
129 :AliPHOSClusterizer(alirunFileName, eventFolderName)
131 // ctor with the indication of the file where header Tree and digits Tree are stored
135 fDefaultInit = kFALSE ;
138 //____________________________________________________________________________
139 AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
144 //____________________________________________________________________________
145 const TString AliPHOSClusterizerv1::BranchName() const
150 //____________________________________________________________________________
151 Float_t AliPHOSClusterizerv1::CalibrateEMC(Float_t amp, Int_t absId)
153 // Convert EMC measured amplitude into real energy.
154 // Calibration parameters are taken from calibration data base for raw data,
155 // or from digitizer parameters for simulated data.
159 AliPHOSGetter *gime = AliPHOSGetter::Instance();
160 gime->PHOSGeometry()->AbsToRelNumbering(absId,relId) ;
161 Int_t module = relId[0];
162 Int_t column = relId[3];
163 Int_t row = relId[2];
164 if(absId <= fEmcCrystals) { // this is EMC
165 fADCchanelEmc = fCalibData->GetADCchannelEmc (module,column,row);
166 return amp*fADCchanelEmc ;
170 if(absId <= fEmcCrystals) // this is EMC
171 return fADCpedestalEmc + amp*fADCchanelEmc ;
176 //____________________________________________________________________________
177 Float_t AliPHOSClusterizerv1::CalibrateCPV(Int_t amp, Int_t absId)
179 // Convert digitized CPV amplitude into charge.
180 // Calibration parameters are taken from calibration data base for raw data,
181 // or from digitizer parameters for simulated data.
185 AliPHOSGetter *gime = AliPHOSGetter::Instance();
186 gime->PHOSGeometry()->AbsToRelNumbering(absId,relId) ;
187 Int_t module = relId[0];
188 Int_t column = relId[3];
189 Int_t row = relId[2];
190 if(absId > fEmcCrystals) { // this is CPV
191 fADCchanelCpv = fCalibData->GetADCchannelCpv (module,column,row);
192 fADCpedestalCpv = fCalibData->GetADCpedestalCpv(module,column,row);
193 return fADCpedestalCpv + amp*fADCchanelCpv ;
197 if(absId > fEmcCrystals) // this is CPV
198 return fADCpedestalCpv+ amp*fADCchanelCpv ;
203 //____________________________________________________________________________
204 void AliPHOSClusterizerv1::Exec(Option_t *option)
206 // Steering method to perform clusterization for events
207 // in the range from fFirstEvent to fLastEvent.
208 // This range is optionally set by SetEventRange().
209 // if fLastEvent=-1 (by default), then process events until the end.
211 if(strstr(option,"tim"))
212 gBenchmark->Start("PHOSClusterizer");
214 if(strstr(option,"print")) {
219 GetCalibrationParameters() ;
221 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
223 gime->SetRawDigits(kFALSE);
225 gime->SetRawDigits(kTRUE);
227 if (fLastEvent == -1)
228 fLastEvent = gime->MaxEvent() - 1 ;
230 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // one event at the time
231 Int_t nEvents = fLastEvent - fFirstEvent + 1;
234 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
236 gime->Event(ievent ,"D"); // Read digits from simulated data
238 gime->Event(fRawReader,"W",fIsOldRCUFormat); // Read digits from raw data
240 fNumberOfEmcClusters = fNumberOfCpvClusters = 0 ;
249 if(strstr(option,"deb"))
250 PrintRecPoints(option) ;
252 //increment the total number of recpoints per run
253 fRecPointsInRun += gime->EmcRecPoints()->GetEntriesFast() ;
254 fRecPointsInRun += gime->CpvRecPoints()->GetEntriesFast() ;
255 if (fRawReader != 0) {
256 AliRunLoader * rl = AliRunLoader::GetRunLoader(gime->PhosLoader()->GetTitle());
257 Int_t iEvent = ievent;
258 rl->SetEventNumber(++iEvent);
262 if(fWrite) //do not unload in "on flight" mode
265 if(strstr(option,"tim")){
266 gBenchmark->Stop("PHOSClusterizer");
267 AliInfo(Form("took %f seconds for Clusterizing %f seconds per event \n",
268 gBenchmark->GetCpuTime("PHOSClusterizer"),
269 gBenchmark->GetCpuTime("PHOSClusterizer")/nEvents )) ;
273 //____________________________________________________________________________
274 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
275 Int_t nPar, Float_t * fitparameters) const
277 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
278 // The initial values for fitting procedure are set equal to the positions of local maxima.
279 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
282 AliPHOSGetter * gime = AliPHOSGetter::Instance();
283 TClonesArray * digits = gime->Digits();
286 gMinuit->mncler(); // Reset Minuit's list of paramters
287 gMinuit->SetPrintLevel(-1) ; // No Printout
288 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
289 // To set the address of the minimization function
291 TList * toMinuit = new TList();
292 toMinuit->AddAt(emcRP,0) ;
293 toMinuit->AddAt(digits,1) ;
295 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
297 // filling initial values for fit parameters
298 AliPHOSDigit * digit ;
302 Int_t nDigits = (Int_t) nPar / 3 ;
306 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
308 for(iDigit = 0; iDigit < nDigits; iDigit++){
309 digit = maxAt[iDigit];
314 geom->AbsToRelNumbering(digit->GetId(), relid) ;
315 geom->RelPosInModule(relid, x, z) ;
317 Float_t energy = maxAtEnergy[iDigit] ;
319 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
322 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : x = %f\n", x ) ;
325 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
328 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : z =%f\n", z ) ;
331 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
334 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : energy = %f\n", energy ) ;
339 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
344 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
345 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
346 gMinuit->SetMaxIterations(5);
347 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
349 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
351 if(ierflg == 4){ // Minimum not found
352 Warning("FindFit", "PHOS Unfolding fit not converged, cluster abandoned\n" );
355 for(index = 0; index < nPar; index++){
358 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
359 fitparameters[index] = val ;
367 //____________________________________________________________________________
368 void AliPHOSClusterizerv1::GetCalibrationParameters()
370 // Set calibration parameters:
371 // if calibration database exists, they are read from database,
372 // otherwise, they are taken from digitizer.
374 // It is a user responsilibity to open CDB before reconstruction, for example:
375 // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
377 AliPHOSGetter * gime = AliPHOSGetter::Instance();
378 // fCalibData = new AliPHOSCalibData(gAlice->GetRunNumber()); //original
380 fCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
384 if ( !gime->Digitizer() )
385 gime->LoadDigitizer();
386 AliPHOSDigitizer * dig = gime->Digitizer();
387 fADCchanelEmc = dig->GetEMCchannel() ;
388 fADCpedestalEmc = dig->GetEMCpedestal();
390 fADCchanelCpv = dig->GetCPVchannel() ;
391 fADCpedestalCpv = dig->GetCPVpedestal() ;
395 //____________________________________________________________________________
396 void AliPHOSClusterizerv1::Init()
398 // Make all memory allocations which can not be done in default constructor.
399 // Attach the Clusterizer task to the list of PHOS tasks
401 AliPHOSGetter* gime = AliPHOSGetter::Instance() ;
403 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
405 AliPHOSGeometry * geom = gime->PHOSGeometry();
407 AliError("Could not find PHOS geometry! Loading the default one !");
408 geom = AliPHOSGeometry::GetInstance("IHEP","");
411 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
414 gMinuit = new TMinuit(100);
416 if ( !gime->Clusterizer() ) {
417 gime->PostClusterizer(this);
421 //____________________________________________________________________________
422 void AliPHOSClusterizerv1::InitParameters()
425 fNumberOfCpvClusters = 0 ;
426 fNumberOfEmcClusters = 0 ;
428 fCpvClusteringThreshold = 0.0;
429 fEmcClusteringThreshold = 0.2;
431 fEmcLocMaxCut = 0.03 ;
432 fCpvLocMaxCut = 0.03 ;
440 fEmcTimeGate = 1.e-8 ;
444 fRecPointsInRun = 0 ;
450 SetEventRange(0,-1) ;
452 fIsOldRCUFormat = kFALSE;
455 //____________________________________________________________________________
456 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
458 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
460 // = 2 are not neighbour but do not continue searching
461 // neighbours are defined as digits having at least a common vertex
462 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
463 // which is compared to a digit (d2) not yet in a cluster
465 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
470 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
473 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
475 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module
476 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
477 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
479 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
480 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
484 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
485 rv = 2; // Difference in row numbers is too large to look further
491 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
498 //____________________________________________________________________________
499 void AliPHOSClusterizerv1::CleanDigits(TClonesArray * digits)
501 // Remove digits with amplitudes below threshold
503 for(Int_t i=0; i<digits->GetEntriesFast(); i++){
504 AliPHOSDigit * digit = static_cast<AliPHOSDigit*>(digits->At(i)) ;
505 if ( (IsInEmc(digit) && CalibrateEMC(digit->GetEnergy(),digit->GetId()) < fEmcMinE) ||
506 (IsInCpv(digit) && CalibrateCPV(digit->GetAmp() ,digit->GetId()) < fCpvMinE) )
507 digits->RemoveAt(i) ;
510 for (Int_t i = 0 ; i < digits->GetEntriesFast() ; i++) {
511 AliPHOSDigit *digit = static_cast<AliPHOSDigit*>( digits->At(i) ) ;
512 digit->SetIndexInList(i) ;
516 //____________________________________________________________________________
517 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
519 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
522 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
524 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
526 if(digit->GetId() <= nEMC ) rv = kTRUE;
531 //____________________________________________________________________________
532 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
534 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
538 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
540 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
542 if(digit->GetId() > nEMC ) rv = kTRUE;
547 //____________________________________________________________________________
548 void AliPHOSClusterizerv1::Unload()
550 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
551 gime->PhosLoader()->UnloadDigits() ;
552 gime->PhosLoader()->UnloadRecPoints() ;
555 //____________________________________________________________________________
556 void AliPHOSClusterizerv1::WriteRecPoints()
559 // Creates new branches with given title
560 // fills and writes into TreeR.
562 AliPHOSGetter * gime = AliPHOSGetter::Instance();
564 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
565 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
566 TClonesArray * digits = gime->Digits() ;
570 //Evaluate position, dispersion and other RecPoint properties..
571 Int_t nEmc = emcRecPoints->GetEntriesFast();
572 for(index = 0; index < nEmc; index++){
573 AliPHOSEmcRecPoint * rp = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) );
574 rp->Purify(fEmcMinE) ;
575 if(rp->GetMultiplicity()>0.) //If this RP is not empty
576 rp->EvalAll(fW0,digits) ;
578 emcRecPoints->RemoveAt(index) ;
582 emcRecPoints->Compress() ;
583 emcRecPoints->Sort() ;
584 // emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ;
585 for(index = 0; index < emcRecPoints->GetEntries(); index++){
586 dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) )->SetIndexInList(index) ;
589 //Now the same for CPV
590 for(index = 0; index < cpvRecPoints->GetEntries(); index++){
591 AliPHOSCpvRecPoint * rp = dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(index) );
592 rp->EvalAll(fW0CPV,digits) ;
594 cpvRecPoints->Sort() ;
596 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
597 dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(index) )->SetIndexInList(index) ;
599 cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
601 if(fWrite){ //We write TreeR
602 TTree * treeR = gime->TreeR();
604 Int_t bufferSize = 32000 ;
605 Int_t splitlevel = 0 ;
608 TBranch * emcBranch = treeR->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
609 emcBranch->SetTitle(BranchName());
612 TBranch * cpvBranch = treeR->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
613 cpvBranch->SetTitle(BranchName());
618 gime->WriteRecPoints("OVERWRITE");
619 gime->WriteClusterizer("OVERWRITE");
623 //____________________________________________________________________________
624 void AliPHOSClusterizerv1::MakeClusters()
626 // Steering method to construct the clusters stored in a list of Reconstructed Points
627 // A cluster is defined as a list of neighbour digits
630 AliPHOSGetter * gime = AliPHOSGetter::Instance();
632 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
633 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
635 emcRecPoints->Delete() ;
636 cpvRecPoints->Delete() ;
638 TClonesArray * digits = gime->Digits() ;
640 //Remove digits below threshold
641 CleanDigits(digits) ;
643 TClonesArray * digitsC = static_cast<TClonesArray*>( digits->Clone() ) ;
645 // Clusterization starts
647 TIter nextdigit(digitsC) ;
648 AliPHOSDigit * digit ;
649 Bool_t notremoved = kTRUE ;
651 while ( (digit = dynamic_cast<AliPHOSDigit *>( nextdigit()) ) ) { // scan over the list of digitsC
654 AliPHOSRecPoint * clu = 0 ;
656 TArrayI clusterdigitslist(1500) ;
659 if (( IsInEmc (digit) &&
660 CalibrateEMC(digit->GetEnergy(),digit->GetId()) > fEmcClusteringThreshold ) ||
662 CalibrateCPV(digit->GetAmp() ,digit->GetId()) > fCpvClusteringThreshold ) ) {
663 Int_t iDigitInCluster = 0 ;
665 if ( IsInEmc(digit) ) {
666 // start a new EMC RecPoint
667 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
668 emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
670 emcRecPoints->AddAt(new AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
671 clu = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(fNumberOfEmcClusters) ) ;
672 fNumberOfEmcClusters++ ;
673 clu->AddDigit(*digit, CalibrateEMC(digit->GetEnergy(),digit->GetId())) ;
674 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
676 digitsC->Remove(digit) ;
680 // start a new CPV cluster
681 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
682 cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
684 cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
686 clu = dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(fNumberOfCpvClusters) ) ;
687 fNumberOfCpvClusters++ ;
688 clu->AddDigit(*digit, CalibrateCPV(digit->GetAmp(),digit->GetId()) ) ;
689 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
691 digitsC->Remove(digit) ;
694 // Here we remove remaining EMC digits, which cannot make a cluster
697 while( ( digit = dynamic_cast<AliPHOSDigit *>( nextdigit() ) ) ) {
699 digitsC->Remove(digit) ;
703 notremoved = kFALSE ;
710 AliPHOSDigit * digitN ;
712 while (index < iDigitInCluster){ // scan over digits already in cluster
713 digit = dynamic_cast<AliPHOSDigit*>( digits->At(clusterdigitslist[index]) ) ;
715 while ( (digitN = dynamic_cast<AliPHOSDigit *>( nextdigit() ) ) ) { // scan over the reduced list of digits
716 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
718 case 0 : // not a neighbour
720 case 1 : // are neighbours
721 if (IsInEmc (digitN))
722 clu->AddDigit(*digitN, CalibrateEMC( digitN->GetEnergy(), digitN->GetId() ) );
724 clu->AddDigit(*digitN, CalibrateCPV( digitN->GetAmp() , digitN->GetId() ) );
726 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
728 digitsC->Remove(digitN) ;
730 case 2 : // too far from each other
739 } // loop over cluster
750 //____________________________________________________________________________
751 void AliPHOSClusterizerv1::MakeUnfolding()
753 // Unfolds clusters using the shape of an ElectroMagnetic shower
754 // Performs unfolding of all EMC/CPV clusters
756 AliPHOSGetter * gime = AliPHOSGetter::Instance();
758 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
760 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
761 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
762 TClonesArray * digits = gime->Digits() ;
764 // Unfold first EMC clusters
765 if(fNumberOfEmcClusters > 0){
767 Int_t nModulesToUnfold = geom->GetNModules() ;
769 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
771 for(index = 0 ; index < numberofNotUnfolded ; index++){
773 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) ) ;
774 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
777 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
778 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
779 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
780 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
782 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
783 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
784 emcRecPoints->Remove(emcRecPoint);
785 emcRecPoints->Compress() ;
787 fNumberOfEmcClusters -- ;
788 numberofNotUnfolded-- ;
791 emcRecPoint->SetNExMax(1) ; //Only one local maximum
795 delete[] maxAtEnergy ;
798 // Unfolding of EMC clusters finished
801 // Unfold now CPV clusters
802 if(fNumberOfCpvClusters > 0){
804 Int_t nModulesToUnfold = geom->GetNModules() ;
806 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
808 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
810 AliPHOSRecPoint * recPoint = dynamic_cast<AliPHOSRecPoint *>( cpvRecPoints->At(index) ) ;
812 if(recPoint->GetPHOSMod()> nModulesToUnfold)
815 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint*>(recPoint) ;
817 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
818 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
819 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
820 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,digits) ;
822 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
823 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
824 cpvRecPoints->Remove(emcRecPoint);
825 cpvRecPoints->Compress() ;
827 numberofCpvNotUnfolded-- ;
828 fNumberOfCpvClusters-- ;
832 delete[] maxAtEnergy ;
835 //Unfolding of Cpv clusters finished
839 //____________________________________________________________________________
840 Double_t AliPHOSClusterizerv1::ShowerShape(Double_t r)
842 // Shape of the shower (see PHOS TDR)
843 // If you change this function, change also the gradient evaluation in ChiSquare()
845 Double_t r4 = r*r*r*r ;
846 Double_t r295 = TMath::Power(r, 2.95) ;
847 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
851 //____________________________________________________________________________
852 void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
854 AliPHOSDigit ** maxAt,
855 Float_t * maxAtEnergy)
857 // Performs the unfolding of a cluster with nMax overlapping showers
859 AliPHOSGetter * gime = AliPHOSGetter::Instance();
861 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
863 const TClonesArray * digits = gime->Digits() ;
864 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
865 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
867 Int_t nPar = 3 * nMax ;
868 Float_t * fitparameters = new Float_t[nPar] ;
870 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
872 // Fit failed, return and remove cluster
873 iniEmc->SetNExMax(-1) ;
874 delete[] fitparameters ;
878 // create ufolded rec points and fill them with new energy lists
879 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
880 // and later correct this number in acordance with actual energy deposition
882 Int_t nDigits = iniEmc->GetMultiplicity() ;
883 Float_t * efit = new Float_t[nDigits] ;
884 Float_t xDigit=0.,zDigit=0.,distance=0. ;
885 Float_t xpar=0.,zpar=0.,epar=0. ;
887 AliPHOSDigit * digit = 0 ;
888 Int_t * emcDigits = iniEmc->GetDigitsList() ;
892 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
893 digit = dynamic_cast<AliPHOSDigit*>( digits->At(emcDigits[iDigit] ) ) ;
894 geom->AbsToRelNumbering(digit->GetId(), relid) ;
895 geom->RelPosInModule(relid, xDigit, zDigit) ;
899 while(iparam < nPar ){
900 xpar = fitparameters[iparam] ;
901 zpar = fitparameters[iparam+1] ;
902 epar = fitparameters[iparam+2] ;
904 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
905 distance = TMath::Sqrt(distance) ;
906 efit[iDigit] += epar * ShowerShape(distance) ;
911 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
912 // so that energy deposited in each cell is distributed betwin new clusters proportionally
913 // to its contribution to efit
915 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
919 while(iparam < nPar ){
920 xpar = fitparameters[iparam] ;
921 zpar = fitparameters[iparam+1] ;
922 epar = fitparameters[iparam+2] ;
925 AliPHOSEmcRecPoint * emcRP = 0 ;
927 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
929 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
930 emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
932 (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
933 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(fNumberOfEmcClusters) ) ;
934 fNumberOfEmcClusters++ ;
935 emcRP->SetNExMax((Int_t)nPar/3) ;
937 else{//create new entries in fCpvRecPoints
938 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
939 cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
941 (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
942 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( cpvRecPoints->At(fNumberOfCpvClusters) ) ;
943 fNumberOfCpvClusters++ ;
947 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
948 digit = dynamic_cast<AliPHOSDigit*>( digits->At( emcDigits[iDigit] ) ) ;
949 geom->AbsToRelNumbering(digit->GetId(), relid) ;
950 geom->RelPosInModule(relid, xDigit, zDigit) ;
951 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
952 distance = TMath::Sqrt(distance) ;
953 ratio = epar * ShowerShape(distance) / efit[iDigit] ;
954 eDigit = emcEnergies[iDigit] * ratio ;
955 emcRP->AddDigit( *digit, eDigit ) ;
959 delete[] fitparameters ;
964 //_____________________________________________________________________________
965 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
967 // Calculates the Chi square for the cluster unfolding minimization
968 // Number of parameters, Gradient, Chi squared, parameters, what to do
970 TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;
972 AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint*>( toMinuit->At(0) ) ;
973 TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) ) ;
977 // AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( gMinuit->GetObjectFit() ) ; // EmcRecPoint to fit
979 Int_t * emcDigits = emcRP->GetDigitsList() ;
981 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
983 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
985 const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
990 for(iparam = 0 ; iparam < nPar ; iparam++)
991 Grad[iparam] = 0 ; // Will evaluate gradient
995 AliPHOSDigit * digit ;
998 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
1000 digit = dynamic_cast<AliPHOSDigit*>( digits->At( emcDigits[iDigit] ) );
1006 geom->AbsToRelNumbering(digit->GetId(), relid) ;
1008 geom->RelPosInModule(relid, xDigit, zDigit) ;
1010 if(iflag == 2){ // calculate gradient
1013 while(iParam < nPar ){
1014 Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
1016 distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
1017 distance = TMath::Sqrt( distance ) ;
1019 efit += x[iParam] * ShowerShape(distance) ;
1022 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
1024 while(iParam < nPar ){
1025 Double_t xpar = x[iParam] ;
1026 Double_t zpar = x[iParam+1] ;
1027 Double_t epar = x[iParam+2] ;
1028 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
1029 Double_t shape = sum * ShowerShape(dr) ;
1030 Double_t r4 = dr*dr*dr*dr ;
1031 Double_t r295 = TMath::Power(dr,2.95) ;
1032 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
1033 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
1035 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
1037 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
1039 Grad[iParam] += shape ; // Derivative over energy
1046 while(iparam < nPar ){
1047 Double_t xpar = x[iparam] ;
1048 Double_t zpar = x[iparam+1] ;
1049 Double_t epar = x[iparam+2] ;
1051 Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
1052 distance = TMath::Sqrt(distance) ;
1053 efit += epar * ShowerShape(distance) ;
1056 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
1057 // Here we assume, that sigma = sqrt(E)
1062 //____________________________________________________________________________
1063 void AliPHOSClusterizerv1::Print(const Option_t *)const
1065 // Print clusterizer parameters
1068 TString taskName(GetName()) ;
1069 taskName.ReplaceAll(Version(), "") ;
1071 if( strcmp(GetName(), "") !=0 ) {
1073 message = "\n--------------- %s %s -----------\n" ;
1074 message += "Clusterizing digits from the file: %s\n" ;
1075 message += " Branch: %s\n" ;
1076 message += " EMC Clustering threshold = %f\n" ;
1077 message += " EMC Local Maximum cut = %f\n" ;
1078 message += " EMC Logarothmic weight = %f\n" ;
1079 message += " CPV Clustering threshold = %f\n" ;
1080 message += " CPV Local Maximum cut = %f\n" ;
1081 message += " CPV Logarothmic weight = %f\n" ;
1083 message += " Unfolding on\n" ;
1085 message += " Unfolding off\n" ;
1087 message += "------------------------------------------------------------------" ;
1090 message = " AliPHOSClusterizerv1 not initialized " ;
1092 AliInfo(Form("%s, %s %s %s %s %s %s %s %s %s %s", message.Data(),
1097 fEmcClusteringThreshold,
1100 fCpvClusteringThreshold,
1106 //____________________________________________________________________________
1107 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1109 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1111 AliPHOSGetter * gime = AliPHOSGetter::Instance();
1113 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
1114 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
1116 AliInfo(Form("\nevent %d \n Found %d EMC RecPoints and %d CPV RecPoints",
1117 gAlice->GetEvNumber(),
1118 emcRecPoints->GetEntriesFast(),
1119 cpvRecPoints->GetEntriesFast() )) ;
1121 fRecPointsInRun += emcRecPoints->GetEntriesFast() ;
1122 fRecPointsInRun += cpvRecPoints->GetEntriesFast() ;
1125 if(strstr(option,"all")) {
1126 printf("\n EMC clusters \n") ;
1127 printf("Index Ene(MeV) Multi Module X Y Z Lambdas_1 Lambda_2 # of prim Primaries list\n") ;
1129 for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1130 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ;
1132 rp->GetLocalPosition(locpos);
1134 rp->GetElipsAxis(lambda);
1137 primaries = rp->GetPrimaries(nprimaries);
1138 printf("\n%6d %8.2f %3d %2d %4.1f %4.1f %4.1f %4f %4f %2d : ",
1139 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(), rp->GetPHOSMod(),
1140 locpos.X(), locpos.Y(), locpos.Z(), lambda[0], lambda[1], nprimaries) ;
1142 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
1143 printf("%d ", primaries[iprimary] ) ;
1148 //Now plot CPV recPoints
1149 printf("\n CPV clusters \n") ;
1150 printf("Index Ene(MeV) Module X Y Z \n") ;
1151 for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1152 AliPHOSCpvRecPoint * rp = (AliPHOSCpvRecPoint * )cpvRecPoints->At(index) ;
1155 rp->GetLocalPosition(locpos);
1157 printf("\n%6d %8.2f %2d %4.1f %4.1f %4.1f \n",
1158 rp->GetIndexInList(), rp->GetEnergy(), rp->GetPHOSMod(),
1159 locpos.X(), locpos.Y(), locpos.Z()) ;