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