]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - PHOS/AliPHOSClusterizerv1.cxx
Fatal if no magnetic field set for the reconstruction (Y.Belikov)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSClusterizerv1.cxx
... / ...
CommitLineData
1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/* $Id$ */
17
18/* $Log:
19 1 October 2000. Yuri Kharlov:
20 AreNeighbours()
21 PPSD upper layer is considered if number of layers>1
22
23 18 October 2000. Yuri Kharlov:
24 AliPHOSClusterizerv1()
25 CPV clusterizing parameters added
26
27 MakeClusters()
28 After first PPSD digit remove EMC digits only once
29*/
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
38// standalone mode.
39// Use Case:
40// root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root", "recpointsname", "digitsname")
41// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42// // reads gAlice from header file "galice.root", uses digits stored in the branch names "digitsname" (default = "Default")
43// // and saves recpoints in branch named "recpointsname" (default = "digitsname")
44// root [1] cl->ExecuteTask()
45// //finds RecPoints in all events stored in galice.root
46// root [2] cl->SetDigitsBranch("digits2")
47// //sets another title for Digitis (input) branch
48// root [3] cl->SetRecPointsBranch("recp2")
49// //sets another title four output branches
50// root [4] cl->SetEmcLocalMaxCut(0.03)
51// //set clusterization parameters
52// root [5] cl->ExecuteTask("deb all time")
53// //once more finds RecPoints options are
54// // deb - print number of found rec points
55// // deb all - print number of found RecPoints and some their characteristics
56// // time - print benchmarking results
57
58// --- ROOT system ---
59
60#include "TROOT.h"
61#include "TFile.h"
62#include "TFolder.h"
63#include "TMath.h"
64#include "TMinuit.h"
65#include "TTree.h"
66#include "TSystem.h"
67#include "TBenchmark.h"
68
69// --- Standard library ---
70
71#include <iostream.h>
72#include <iomanip.h>
73
74// --- AliRoot header files ---
75#include "AliPHOSCalibrationDB.h"
76#include "AliPHOSClusterizerv1.h"
77#include "AliPHOSCpvRecPoint.h"
78#include "AliPHOSDigit.h"
79#include "AliPHOSDigitizer.h"
80#include "AliPHOSEmcRecPoint.h"
81#include "AliPHOS.h"
82#include "AliPHOSGetter.h"
83#include "AliRun.h"
84
85ClassImp(AliPHOSClusterizerv1)
86
87//____________________________________________________________________________
88 AliPHOSClusterizerv1::AliPHOSClusterizerv1() : AliPHOSClusterizer()
89{
90 // default ctor (to be used mainly by Streamer)
91
92 InitParameters() ;
93 fDefaultInit = kTRUE ;
94}
95
96//____________________________________________________________________________
97AliPHOSClusterizerv1::AliPHOSClusterizerv1(const char* headerFile,const char* name, const Bool_t toSplit)
98:AliPHOSClusterizer(headerFile, name, toSplit)
99{
100 // ctor with the indication of the file where header Tree and digits Tree are stored
101
102 InitParameters() ;
103 Init() ;
104 fDefaultInit = kFALSE ;
105
106}
107
108//____________________________________________________________________________
109 AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
110{
111 // dtor
112 fSplitFile = 0 ;
113
114}
115
116//____________________________________________________________________________
117const TString AliPHOSClusterizerv1::BranchName() const
118{
119 TString branchName(GetName() ) ;
120 branchName.Remove(branchName.Index(Version())-1) ;
121 return branchName ;
122}
123
124//____________________________________________________________________________
125Float_t AliPHOSClusterizerv1::Calibrate(Int_t amp, Int_t absId) const
126{ //To be replased later by the method, reading individual parameters from the database
127 if(fCalibrationDB)
128 return fCalibrationDB->Calibrate(amp,absId) ;
129 else{ //simulation
130 if(absId <= fEmcCrystals) //calibrate as EMC
131 return fADCpedestalEmc + amp*fADCchanelEmc ;
132 else //calibrate as CPV
133 return fADCpedestalCpv+ amp*fADCchanelCpv ;
134 }
135}
136
137//____________________________________________________________________________
138void AliPHOSClusterizerv1::Exec(Option_t * option)
139{
140 // Steering method
141
142 if( strcmp(GetName(), "")== 0 )
143 Init() ;
144
145 if(strstr(option,"tim"))
146 gBenchmark->Start("PHOSClusterizer");
147
148 if(strstr(option,"print"))
149 Print("") ;
150
151 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
152 if(gime->BranchExists("RecPoints"))
153 return ;
154 Int_t nevents = gime->MaxEvent() ;
155 Int_t ievent ;
156
157 for(ievent = 0; ievent < nevents; ievent++){
158
159 gime->Event(ievent,"D") ;
160
161 if(ievent == 0)
162 GetCalibrationParameters() ;
163
164 fNumberOfEmcClusters = fNumberOfCpvClusters = 0 ;
165
166 MakeClusters() ;
167
168 if(fToUnfold)
169 MakeUnfolding() ;
170
171 WriteRecPoints(ievent) ;
172
173 if(strstr(option,"deb"))
174 PrintRecPoints(option) ;
175
176 //increment the total number of digits per run
177 fRecPointsInRun += gime->EmcRecPoints()->GetEntriesFast() ;
178 fRecPointsInRun += gime->CpvRecPoints()->GetEntriesFast() ;
179 }
180
181 if(strstr(option,"tim")){
182 gBenchmark->Stop("PHOSClusterizer");
183 cout << "AliPHOSClusterizer:" << endl ;
184 cout << " took " << gBenchmark->GetCpuTime("PHOSClusterizer") << " seconds for Clusterizing "
185 << gBenchmark->GetCpuTime("PHOSClusterizer")/nevents << " seconds per event " << endl ;
186 cout << endl ;
187 }
188
189}
190
191//____________________________________________________________________________
192Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
193 Int_t nPar, Float_t * fitparameters) const
194{
195 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
196 // The initial values for fitting procedure are set equal to the positions of local maxima.
197 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
198
199 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
200 TClonesArray * digits = gime->Digits() ;
201
202
203 gMinuit->mncler(); // Reset Minuit's list of paramters
204 gMinuit->SetPrintLevel(-1) ; // No Printout
205 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
206 // To set the address of the minimization function
207
208 TList * toMinuit = new TList();
209 toMinuit->AddAt(emcRP,0) ;
210 toMinuit->AddAt(digits,1) ;
211
212 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
213
214 // filling initial values for fit parameters
215 AliPHOSDigit * digit ;
216
217 Int_t ierflg = 0;
218 Int_t index = 0 ;
219 Int_t nDigits = (Int_t) nPar / 3 ;
220
221 Int_t iDigit ;
222
223 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
224
225 for(iDigit = 0; iDigit < nDigits; iDigit++){
226 digit = maxAt[iDigit];
227
228 Int_t relid[4] ;
229 Float_t x = 0.;
230 Float_t z = 0.;
231 geom->AbsToRelNumbering(digit->GetId(), relid) ;
232 geom->RelPosInModule(relid, x, z) ;
233
234 Float_t energy = maxAtEnergy[iDigit] ;
235
236 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
237 index++ ;
238 if(ierflg != 0){
239 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : x = " << x << endl ;
240 return kFALSE;
241 }
242 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
243 index++ ;
244 if(ierflg != 0){
245 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : z = " << z << endl ;
246 return kFALSE;
247 }
248 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
249 index++ ;
250 if(ierflg != 0){
251 cout << "PHOS Unfolding> Unable to set initial value for fit procedure : energy = " << energy << endl ;
252 return kFALSE;
253 }
254 }
255
256 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
257 // depends on it.
258 Double_t p1 = 1.0 ;
259 Double_t p2 = 0.0 ;
260
261 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
262 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
263 gMinuit->SetMaxIterations(5);
264 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
265
266 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
267
268 if(ierflg == 4){ // Minimum not found
269 cout << "PHOS Unfolding> Fit not converged, cluster abandoned "<< endl ;
270 return kFALSE ;
271 }
272 for(index = 0; index < nPar; index++){
273 Double_t err ;
274 Double_t val ;
275 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
276 fitparameters[index] = val ;
277 }
278
279 delete toMinuit ;
280 return kTRUE;
281
282}
283
284//____________________________________________________________________________
285void AliPHOSClusterizerv1::GetCalibrationParameters()
286{
287 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
288
289 const TTask * task = gime->Digitizer(BranchName()) ;
290 if(strcmp(task->IsA()->GetName(),"AliPHOSDigitizer")==0){
291 const AliPHOSDigitizer * dig = static_cast<const AliPHOSDigitizer *>(task) ;
292
293 fADCchanelEmc = dig->GetEMCchannel() ;
294 fADCpedestalEmc = dig->GetEMCpedestal();
295
296 fADCchanelCpv = dig->GetCPVchannel() ;
297 fADCpedestalCpv = dig->GetCPVpedestal() ;
298 }
299 else{
300 fCalibrationDB = gime->CalibrationDB();
301 }
302}
303
304//____________________________________________________________________________
305void AliPHOSClusterizerv1::Init()
306{
307 // Make all memory allocations which can not be done in default constructor.
308 // Attach the Clusterizer task to the list of PHOS tasks
309
310 if ( strcmp(GetTitle(), "") == 0 )
311 SetTitle("galice.root") ;
312
313 TString branchname = GetName() ;
314 branchname.Remove(branchname.Index(Version())-1) ;
315
316 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(),branchname.Data(), fToSplit ) ;
317 if ( gime == 0 ) {
318 cerr << "ERROR: AliPHOSClusterizerv1::Init -> Could not obtain the Getter object !" << endl ;
319 return ;
320 }
321
322 fSplitFile = 0 ;
323 if(fToSplit){
324 // construct the name of the file as /path/EMCAL.SDigits.root
325 //First - extract full path if necessary
326 TString fileName(GetTitle()) ;
327 Ssiz_t islash = fileName.Last('/') ;
328 if(islash<fileName.Length())
329 fileName.Remove(islash+1,fileName.Length()) ;
330 else
331 fileName="" ;
332 // Next - append the file name
333 fileName+="PHOS.RecData." ;
334 if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
335 fileName+=branchname ;
336 fileName+="." ;
337 }
338 fileName+="root" ;
339 // Finally - check if the file already opened or open the file
340 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
341 if(!fSplitFile)
342 fSplitFile = TFile::Open(fileName.Data(),"update") ;
343 }
344
345
346
347 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
348 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
349
350 if(!gMinuit)
351 gMinuit = new TMinuit(100) ;
352
353 gime->PostClusterizer(this) ;
354 gime->PostRecPoints(branchname) ;
355
356}
357
358//____________________________________________________________________________
359void AliPHOSClusterizerv1::InitParameters()
360{
361
362 fNumberOfCpvClusters = 0 ;
363 fNumberOfEmcClusters = 0 ;
364
365 fCpvClusteringThreshold = 0.0;
366 fEmcClusteringThreshold = 0.2;
367
368 fEmcLocMaxCut = 0.03 ;
369 fCpvLocMaxCut = 0.03 ;
370
371 fW0 = 4.5 ;
372 fW0CPV = 4.0 ;
373
374 fEmcTimeGate = 1.e-8 ;
375
376 fToUnfold = kTRUE ;
377
378 TString clusterizerName( GetName()) ;
379 if (clusterizerName.IsNull() )
380 clusterizerName = "Default" ;
381 clusterizerName.Append(":") ;
382 clusterizerName.Append(Version()) ;
383 SetName(clusterizerName) ;
384 fRecPointsInRun = 0 ;
385 fCalibrationDB = 0 ;
386
387
388}
389
390//____________________________________________________________________________
391Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
392{
393 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
394 // = 1 are neighbour
395 // = 2 are not neighbour but do not continue searching
396 // neighbours are defined as digits having at least a common vertex
397 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
398 // which is compared to a digit (d2) not yet in a cluster
399
400 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
401
402 Int_t rv = 0 ;
403
404 Int_t relid1[4] ;
405 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
406
407 Int_t relid2[4] ;
408 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
409
410 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module
411 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
412 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
413
414 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
415 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
416 rv = 1 ;
417 }
418 else {
419 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
420 rv = 2; // Difference in row numbers is too large to look further
421 }
422
423 }
424 else {
425
426 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
427 rv=2 ;
428
429 }
430
431 return rv ;
432}
433
434
435//____________________________________________________________________________
436Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
437{
438 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
439
440 Bool_t rv = kFALSE ;
441 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
442
443 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
444
445 if(digit->GetId() <= nEMC ) rv = kTRUE;
446
447 return rv ;
448}
449
450//____________________________________________________________________________
451Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
452{
453 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
454
455 Bool_t rv = kFALSE ;
456 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
457
458 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
459
460 if(digit->GetId() > nEMC ) rv = kTRUE;
461
462 return rv ;
463}
464
465//____________________________________________________________________________
466void AliPHOSClusterizerv1::WriteRecPoints(Int_t event)
467{
468
469 // Creates new branches with given title
470 // fills and writes into TreeR.
471
472
473 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
474 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
475 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
476 TClonesArray * digits = gime->Digits() ;
477 TTree * treeR ;
478
479 if(fToSplit){
480 if(!fSplitFile)
481 return ;
482 fSplitFile->cd() ;
483 TString name("TreeR") ;
484 name += event ;
485 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
486 }
487 else{
488 treeR = gAlice->TreeR();
489 }
490
491 if(!treeR){
492 if(fSplitFile)
493 gAlice->MakeTree("R", fSplitFile);
494 else
495 gAlice->MakeTree("R", gROOT->GetFile(GetTitle()));
496 treeR = gAlice->TreeR() ;
497 }
498
499 Int_t index ;
500 //Evaluate position, dispersion and other RecPoint properties...
501 for(index = 0; index < emcRecPoints->GetEntries(); index++)
502 ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->EvalAll(fW0,digits) ;
503
504 emcRecPoints->Sort() ;
505 for(index = 0; index < emcRecPoints->GetEntries(); index++)
506 ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->SetIndexInList(index) ;
507
508 emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ;
509
510 //Now the same for CPV
511 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
512 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->EvalAll(fW0CPV,digits) ;
513
514 cpvRecPoints->Sort() ;
515
516 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
517 ((AliPHOSRecPoint *)cpvRecPoints->At(index))->SetIndexInList(index) ;
518
519 cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
520
521 Int_t bufferSize = 32000 ;
522 Int_t splitlevel = 0 ;
523
524 //First EMC
525 TBranch * emcBranch = treeR->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
526 emcBranch->SetTitle(BranchName());
527
528 //Now CPV branch
529 TBranch * cpvBranch = treeR->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
530 cpvBranch->SetTitle(BranchName());
531
532 //And Finally clusterizer branch
533 AliPHOSClusterizerv1 * cl = this ;
534 TBranch * clusterizerBranch = treeR->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
535 &cl,bufferSize,splitlevel);
536 clusterizerBranch->SetTitle(BranchName());
537
538 emcBranch ->Fill() ;
539 cpvBranch ->Fill() ;
540 clusterizerBranch->Fill() ;
541
542 treeR->AutoSave() ;
543 if(gAlice->TreeR()!=treeR)
544 treeR->Delete();
545}
546
547//____________________________________________________________________________
548void AliPHOSClusterizerv1::MakeClusters()
549{
550 // Steering method to construct the clusters stored in a list of Reconstructed Points
551 // A cluster is defined as a list of neighbour digits
552
553 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
554
555 TObjArray * emcRecPoints = gime->EmcRecPoints(BranchName()) ;
556 TObjArray * cpvRecPoints = gime->CpvRecPoints(BranchName()) ;
557 emcRecPoints->Delete() ;
558 cpvRecPoints->Delete() ;
559
560 TClonesArray * digits = gime->Digits() ;
561 if ( !digits ) {
562 cerr << "ERROR: AliPHOSClusterizerv1::MakeClusters -> Digits with name "
563 << GetName() << " not found ! " << endl ;
564 abort() ;
565 }
566 TClonesArray * digitsC = (TClonesArray*)digits->Clone() ;
567
568
569 // Clusterization starts
570
571 TIter nextdigit(digitsC) ;
572 AliPHOSDigit * digit ;
573 Bool_t notremoved = kTRUE ;
574
575 while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digitsC
576 AliPHOSRecPoint * clu = 0 ;
577
578 TArrayI clusterdigitslist(1500) ;
579 Int_t index ;
580
581 if (( IsInEmc (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fEmcClusteringThreshold ) ||
582 ( IsInCpv (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fCpvClusteringThreshold ) ) {
583 Int_t iDigitInCluster = 0 ;
584
585 if ( IsInEmc(digit) ) {
586 // start a new EMC RecPoint
587 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
588 emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
589
590 emcRecPoints->AddAt(new AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
591 clu = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters) ;
592 fNumberOfEmcClusters++ ;
593 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId())) ;
594 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
595 iDigitInCluster++ ;
596 digitsC->Remove(digit) ;
597
598 } else {
599
600 // start a new CPV cluster
601 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
602 cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
603
604 cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
605
606 clu = (AliPHOSCpvRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters) ;
607 fNumberOfCpvClusters++ ;
608 clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId()) ) ;
609 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
610 iDigitInCluster++ ;
611 digitsC->Remove(digit) ;
612 nextdigit.Reset() ;
613
614 // Here we remove remaining EMC digits, which cannot make a cluster
615
616 if( notremoved ) {
617 while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
618 if( IsInEmc(digit) )
619 digitsC->Remove(digit) ;
620 else
621 break ;
622 }
623 notremoved = kFALSE ;
624 }
625
626 } // else
627
628 nextdigit.Reset() ;
629
630 AliPHOSDigit * digitN ;
631 index = 0 ;
632 while (index < iDigitInCluster){ // scan over digits already in cluster
633 digit = (AliPHOSDigit*)digits->At(clusterdigitslist[index]) ;
634 index++ ;
635 while ( (digitN = (AliPHOSDigit *)nextdigit()) ) { // scan over the reduced list of digits
636 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
637 switch (ineb ) {
638 case 0 : // not a neighbour
639 break ;
640 case 1 : // are neighbours
641 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->GetId() ) ) ;
642 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
643 iDigitInCluster++ ;
644 digitsC->Remove(digitN) ;
645 break ;
646 case 2 : // too far from each other
647 goto endofloop;
648 } // switch
649
650 } // while digitN
651
652 endofloop: ;
653 nextdigit.Reset() ;
654
655 } // loop over cluster
656
657 } // energy theshold
658
659
660 } // while digit
661
662 delete digitsC ;
663
664}
665
666//____________________________________________________________________________
667void AliPHOSClusterizerv1::MakeUnfolding()
668{
669 // Unfolds clusters using the shape of an ElectroMagnetic shower
670 // Performs unfolding of all EMC/CPV clusters
671
672 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
673
674
675 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
676 TObjArray * emcRecPoints = gime->EmcRecPoints(BranchName()) ;
677 TObjArray * cpvRecPoints = gime->CpvRecPoints(BranchName()) ;
678 TClonesArray * digits = gime->Digits() ;
679
680 // Unfold first EMC clusters
681 if(fNumberOfEmcClusters > 0){
682
683 Int_t nModulesToUnfold = geom->GetNModules() ;
684
685 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
686 Int_t index ;
687 for(index = 0 ; index < numberofNotUnfolded ; index++){
688
689 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) emcRecPoints->At(index) ;
690 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
691 break ;
692
693 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
694 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
695 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
696 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
697
698 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
699 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
700 emcRecPoints->Remove(emcRecPoint);
701 emcRecPoints->Compress() ;
702 index-- ;
703 fNumberOfEmcClusters -- ;
704 numberofNotUnfolded-- ;
705 }
706
707 delete[] maxAt ;
708 delete[] maxAtEnergy ;
709 }
710 }
711 // Unfolding of EMC clusters finished
712
713
714 // Unfold now CPV clusters
715 if(fNumberOfCpvClusters > 0){
716
717 Int_t nModulesToUnfold = geom->GetNModules() ;
718
719 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
720 Int_t index ;
721 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
722
723 AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) cpvRecPoints->At(index) ;
724
725 if(recPoint->GetPHOSMod()> nModulesToUnfold)
726 break ;
727
728 AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ;
729
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,fCpvLocMaxCut,digits) ;
734
735 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
736 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
737 cpvRecPoints->Remove(emcRecPoint);
738 cpvRecPoints->Compress() ;
739 index-- ;
740 numberofCpvNotUnfolded-- ;
741 fNumberOfCpvClusters-- ;
742 }
743
744 delete[] maxAt ;
745 delete[] maxAtEnergy ;
746 }
747 }
748 //Unfolding of Cpv clusters finished
749
750}
751
752//____________________________________________________________________________
753Double_t AliPHOSClusterizerv1::ShowerShape(Double_t r)
754{
755 // Shape of the shower (see PHOS TDR)
756 // If you change this function, change also the gradient evaluation in ChiSquare()
757
758 Double_t r4 = r*r*r*r ;
759 Double_t r295 = TMath::Power(r, 2.95) ;
760 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
761 return shape ;
762}
763
764//____________________________________________________________________________
765void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
766 Int_t nMax,
767 AliPHOSDigit ** maxAt,
768 Float_t * maxAtEnergy)
769{
770 // Performs the unfolding of a cluster with nMax overlapping showers
771
772 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
773
774 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
775 const TClonesArray * digits = gime->Digits() ;
776 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
777 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
778
779 Int_t nPar = 3 * nMax ;
780 Float_t * fitparameters = new Float_t[nPar] ;
781
782 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
783 if( !rv ) {
784 // Fit failed, return and remove cluster
785 delete[] fitparameters ;
786 return ;
787 }
788
789 // create ufolded rec points and fill them with new energy lists
790 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
791 // and later correct this number in acordance with actual energy deposition
792
793 Int_t nDigits = iniEmc->GetMultiplicity() ;
794 Float_t * efit = new Float_t[nDigits] ;
795 Float_t xDigit=0.,zDigit=0.,distance=0. ;
796 Float_t xpar=0.,zpar=0.,epar=0. ;
797 Int_t relid[4] ;
798 AliPHOSDigit * digit = 0 ;
799 Int_t * emcDigits = iniEmc->GetDigitsList() ;
800
801 Int_t iparam ;
802 Int_t iDigit ;
803 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
804 digit = (AliPHOSDigit*) digits->At(emcDigits[iDigit] ) ;
805 geom->AbsToRelNumbering(digit->GetId(), relid) ;
806 geom->RelPosInModule(relid, xDigit, zDigit) ;
807 efit[iDigit] = 0;
808
809 iparam = 0 ;
810 while(iparam < nPar ){
811 xpar = fitparameters[iparam] ;
812 zpar = fitparameters[iparam+1] ;
813 epar = fitparameters[iparam+2] ;
814 iparam += 3 ;
815 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
816 distance = TMath::Sqrt(distance) ;
817 efit[iDigit] += epar * ShowerShape(distance) ;
818 }
819 }
820
821
822 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
823 // so that energy deposited in each cell is distributed betwin new clusters proportionally
824 // to its contribution to efit
825
826 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
827 Float_t ratio ;
828
829 iparam = 0 ;
830 while(iparam < nPar ){
831 xpar = fitparameters[iparam] ;
832 zpar = fitparameters[iparam+1] ;
833 epar = fitparameters[iparam+2] ;
834 iparam += 3 ;
835
836 AliPHOSEmcRecPoint * emcRP = 0 ;
837
838 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
839
840 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
841 emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
842
843 (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
844 emcRP = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters);
845 fNumberOfEmcClusters++ ;
846 }
847 else{//create new entries in fCpvRecPoints
848 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
849 cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
850
851 (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
852 emcRP = (AliPHOSEmcRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters);
853 fNumberOfCpvClusters++ ;
854 }
855
856 Float_t eDigit ;
857 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
858 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
859 geom->AbsToRelNumbering(digit->GetId(), relid) ;
860 geom->RelPosInModule(relid, xDigit, zDigit) ;
861 distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
862 distance = TMath::Sqrt(distance) ;
863 ratio = epar * ShowerShape(distance) / efit[iDigit] ;
864 eDigit = emcEnergies[iDigit] * ratio ;
865 emcRP->AddDigit( *digit, eDigit ) ;
866 }
867 }
868
869 delete[] fitparameters ;
870 delete[] efit ;
871
872}
873
874//_____________________________________________________________________________
875void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
876{
877 // Calculates the Chi square for the cluster unfolding minimization
878 // Number of parameters, Gradient, Chi squared, parameters, what to do
879
880 TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
881
882 AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint*) toMinuit->At(0) ;
883 TClonesArray * digits = (TClonesArray*)toMinuit->At(1) ;
884
885
886
887 // AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint *) gMinuit->GetObjectFit() ; // EmcRecPoint to fit
888
889 Int_t * emcDigits = emcRP->GetDigitsList() ;
890
891 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
892
893 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
894
895 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
896 fret = 0. ;
897 Int_t iparam ;
898
899 if(iflag == 2)
900 for(iparam = 0 ; iparam < nPar ; iparam++)
901 Grad[iparam] = 0 ; // Will evaluate gradient
902
903 Double_t efit ;
904
905 AliPHOSDigit * digit ;
906 Int_t iDigit ;
907
908 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
909
910 digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ;
911
912 Int_t relid[4] ;
913 Float_t xDigit ;
914 Float_t zDigit ;
915
916 geom->AbsToRelNumbering(digit->GetId(), relid) ;
917
918 geom->RelPosInModule(relid, xDigit, zDigit) ;
919
920 if(iflag == 2){ // calculate gradient
921 Int_t iParam = 0 ;
922 efit = 0 ;
923 while(iParam < nPar ){
924 Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
925 iParam++ ;
926 distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ;
927 distance = TMath::Sqrt( distance ) ;
928 iParam++ ;
929 efit += x[iParam] * ShowerShape(distance) ;
930 iParam++ ;
931 }
932 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
933 iParam = 0 ;
934 while(iParam < nPar ){
935 Double_t xpar = x[iParam] ;
936 Double_t zpar = x[iParam+1] ;
937 Double_t epar = x[iParam+2] ;
938 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
939 Double_t shape = sum * ShowerShape(dr) ;
940 Double_t r4 = dr*dr*dr*dr ;
941 Double_t r295 = TMath::Power(dr,2.95) ;
942 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
943 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
944
945 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
946 iParam++ ;
947 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
948 iParam++ ;
949 Grad[iParam] += shape ; // Derivative over energy
950 iParam++ ;
951 }
952 }
953 efit = 0;
954 iparam = 0 ;
955
956 while(iparam < nPar ){
957 Double_t xpar = x[iparam] ;
958 Double_t zpar = x[iparam+1] ;
959 Double_t epar = x[iparam+2] ;
960 iparam += 3 ;
961 Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) ;
962 distance = TMath::Sqrt(distance) ;
963 efit += epar * ShowerShape(distance) ;
964 }
965
966 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
967 // Here we assume, that sigma = sqrt(E)
968 }
969
970}
971
972//____________________________________________________________________________
973void AliPHOSClusterizerv1::Print(Option_t * option)const
974{
975 // Print clusterizer parameters
976
977 if( strcmp(GetName(), "") !=0 ){
978
979 // Print parameters
980
981 TString taskName(GetName()) ;
982 taskName.ReplaceAll(Version(), "") ;
983
984 cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl
985 << "Clusterizing digits from the file: " << taskName.Data() << endl
986 << " Branch: " << GetName() << endl
987 << endl
988 << " EMC Clustering threshold = " << fEmcClusteringThreshold << endl
989 << " EMC Local Maximum cut = " << fEmcLocMaxCut << endl
990 << " EMC Logarothmic weight = " << fW0 << endl
991 << endl
992 << " CPV Clustering threshold = " << fCpvClusteringThreshold << endl
993 << " CPV Local Maximum cut = " << fCpvLocMaxCut << endl
994 << " CPV Logarothmic weight = " << fW0CPV << endl
995 << endl ;
996 if(fToUnfold)
997 cout << " Unfolding on " << endl ;
998 else
999 cout << " Unfolding off " << endl ;
1000
1001 cout << "------------------------------------------------------------------" <<endl ;
1002 }
1003 else
1004 cout << " AliPHOSClusterizerv1 not initialized " << endl ;
1005}
1006//____________________________________________________________________________
1007void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1008{
1009 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1010
1011 TObjArray * emcRecPoints = AliPHOSGetter::GetInstance()->EmcRecPoints(BranchName()) ;
1012 TObjArray * cpvRecPoints = AliPHOSGetter::GetInstance()->CpvRecPoints(BranchName()) ;
1013
1014 cout << "AliPHOSClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1015 cout << " Found "<< emcRecPoints->GetEntriesFast() << " EMC Rec Points and "
1016 << cpvRecPoints->GetEntriesFast() << " CPV RecPoints" << endl ;
1017
1018 fRecPointsInRun += emcRecPoints->GetEntriesFast() ;
1019 fRecPointsInRun += cpvRecPoints->GetEntriesFast() ;
1020
1021 if(strstr(option,"all")) {
1022 cout << "EMC clusters " << endl ;
1023 cout << " Index Ene(MeV) Multi Module X Y Z Lambda 1 Lambda 2 # of prim Primaries list " << endl;
1024
1025 Int_t index ;
1026 for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1027 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ;
1028 TVector3 locpos;
1029 rp->GetLocalPosition(locpos);
1030 Float_t lambda[2];
1031 rp->GetElipsAxis(lambda);
1032 Int_t * primaries;
1033 Int_t nprimaries;
1034 primaries = rp->GetPrimaries(nprimaries);
1035
1036 cout << setw(4) << rp->GetIndexInList() << " "
1037 << setw(7) << setprecision(3) << rp->GetEnergy() << " "
1038 << setw(3) << rp->GetMultiplicity() << " "
1039 << setw(1) << rp->GetPHOSMod() << " "
1040 << setw(6) << setprecision(2) << locpos.X() << " "
1041 << setw(6) << setprecision(2) << locpos.Y() << " "
1042 << setw(6) << setprecision(2) << locpos.Z() << " "
1043 << setw(4) << setprecision(2) << lambda[0] << " "
1044 << setw(4) << setprecision(2) << lambda[1] << " "
1045 << setw(2) << nprimaries << " " ;
1046
1047 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1048 cout << setw(4) << primaries[iprimary] << " " ;
1049 cout << endl ;
1050 }
1051
1052 //Now plot CPV recPoints
1053 cout << "EMC clusters " << endl ;
1054 cout << " Index "
1055 << " Multi "
1056 << " Module "
1057 << " X "
1058 << " Y "
1059 << " Z "
1060 << " # of prim "
1061 << " Primaries list " << endl;
1062
1063 for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1064 AliPHOSRecPoint * rp = (AliPHOSRecPoint * )cpvRecPoints->At(index) ;
1065 cout << setw(6) << rp->GetIndexInList() << " ";
1066 cout << setw(6) << rp->GetPHOSMod() << " CPV ";
1067
1068 TVector3 locpos;
1069 rp->GetLocalPosition(locpos);
1070 cout << setw(6) << locpos.X() << " ";
1071 cout << setw(6) << locpos.Y() << " ";
1072 cout << setw(6) << locpos.Z() << " ";
1073
1074 Int_t * primaries;
1075 Int_t nprimaries ;
1076 primaries = rp->GetPrimaries(nprimaries);
1077 cout << setw(6) << nprimaries << " ";
1078
1079 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1080 cout << setw(4) << primaries[iprimary] << " ";
1081 cout << endl;
1082 }
1083
1084
1085 cout << "-----------------------------------------------------------------------"<<endl ;
1086 }
1087}
1088