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