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