15-feb-2006 NvE Class IcePandel introduced. In the /macros directory an example job
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALClusterizerv1.cxx
CommitLineData
4aa81449 1
483b0559 2/**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
7 * *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
173558f2 16
483b0559 17/* $Id$ */
803d1ab0 18
483b0559 19//*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
05a92d59 20// August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
21// of new IO (à la PHOS)
483b0559 22//////////////////////////////////////////////////////////////////////////////
23// Clusterization class. Performs clusterization (collects neighbouring active cells) and
24// unfolds the clusters having several local maxima.
25// Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
26// EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all
27// parameters including input digits branch title, thresholds etc.)
28// This TTask is normally called from Reconstructioner, but can as well be used in
29// standalone mode.
30// Use Case:
31// root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")
32// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
33// //reads gAlice from header file "..."
34// root [1] cl->ExecuteTask()
35// //finds RecPoints in all events stored in galice.root
36// root [2] cl->SetDigitsBranch("digits2")
37// //sets another title for Digitis (input) branch
38// root [3] cl->SetRecPointsBranch("recp2")
39// //sets another title four output branches
40// root [4] cl->SetTowerLocalMaxCut(0.03)
41// //set clusterization parameters
42// root [5] cl->ExecuteTask("deb all time")
43// //once more finds RecPoints options are
44// // deb - print number of found rec points
45// // deb all - print number of found RecPoints and some their characteristics
46// // time - print benchmarking results
47
48// --- ROOT system ---
49
50#include "TROOT.h"
51#include "TFile.h"
52#include "TFolder.h"
53#include "TMath.h"
54#include "TMinuit.h"
55#include "TTree.h"
56#include "TSystem.h"
57#include "TBenchmark.h"
58
59// --- Standard library ---
60
173558f2 61
483b0559 62// --- AliRoot header files ---
88cb7938 63#include "AliEMCALGetter.h"
483b0559 64#include "AliEMCALClusterizerv1.h"
70a93198 65#include "AliEMCALRecPoint.h"
483b0559 66#include "AliEMCALDigit.h"
67#include "AliEMCALDigitizer.h"
483b0559 68#include "AliEMCAL.h"
05a92d59 69#include "AliEMCALGeometry.h"
483b0559 70
71ClassImp(AliEMCALClusterizerv1)
1963b290 72
73Int_t addOn[20][60][60];
74
483b0559 75//____________________________________________________________________________
76 AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
77{
78 // default ctor (to be used mainly by Streamer)
79
839828a6 80 InitParameters() ;
1963b290 81 fDefaultInit = kTRUE ;
cf24bdc4 82 for(Int_t is=0;is<20;is++){
83 for(Int_t i=0;i<60;i++){
84 for(Int_t j=0;j<60;j++){
85 addOn[is][i][j]=0;
86 }
87 }
88 }
89//PH cout<<"file to read 1"<<endl;
1963b290 90 ReadFile();
cf24bdc4 91//PH cout<<"file read 1"<<endl;
483b0559 92}
93
94//____________________________________________________________________________
88cb7938 95AliEMCALClusterizerv1::AliEMCALClusterizerv1(const TString alirunFileName, const TString eventFolderName)
96:AliEMCALClusterizer(alirunFileName, eventFolderName)
483b0559 97{
98 // ctor with the indication of the file where header Tree and digits Tree are stored
99
839828a6 100 InitParameters() ;
483b0559 101 Init() ;
05a92d59 102 fDefaultInit = kFALSE ;
1963b290 103 for(Int_t is=0;is<20;is++){
cf24bdc4 104 for(Int_t i=0;i<60;i++){
105 for(Int_t j=0;j<60;j++){
106 addOn[is][i][j]=0;
107 }
108 }
109 }
110//PH cout<<"file to read 2"<<endl;
1963b290 111 ReadFile();
cf24bdc4 112//PH cout<<"file read 2"<<endl;
483b0559 113
114}
05a92d59 115
483b0559 116//____________________________________________________________________________
117 AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
118{
ef305168 119 // dtor
839828a6 120
ef305168 121}
122
123//____________________________________________________________________________
124const TString AliEMCALClusterizerv1::BranchName() const
9859bfc0 125{
88cb7938 126 return GetName();
127
483b0559 128}
839828a6 129
483b0559 130//____________________________________________________________________________
fdebddeb 131Float_t AliEMCALClusterizerv1::Calibrate(Int_t amp) const
9859bfc0 132{
fdebddeb 133 //To be replased later by the method, reading individual parameters from the database
134 return -fADCpedestalECA + amp * fADCchannelECA ;
483b0559 135}
05a92d59 136
483b0559 137//____________________________________________________________________________
138void AliEMCALClusterizerv1::Exec(Option_t * option)
139{
12c4dd95 140 // Steering method to perform clusterization for events
141 // in the range from fFirstEvent to fLastEvent.
142 // This range is optionally set by SetEventRange().
143 // if fLastEvent=-1 (by default), then process events until the end.
483b0559 144
483b0559 145 if(strstr(option,"tim"))
146 gBenchmark->Start("EMCALClusterizer");
147
148 if(strstr(option,"print"))
149 Print("") ;
150
3c9d4047 151 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
88cb7938 152
12c4dd95 153 if (fLastEvent == -1)
cf24bdc4 154 fLastEvent = gime->MaxEvent() - 1;
12c4dd95 155 Int_t nEvents = fLastEvent - fFirstEvent + 1;
483b0559 156
12c4dd95 157 Int_t ievent ;
12c4dd95 158 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
05a92d59 159 gime->Event(ievent,"D") ;
88cb7938 160 GetCalibrationParameters() ;
483b0559 161
fdebddeb 162 fNumberOfECAClusters = 0;
05a92d59 163
483b0559 164 MakeClusters() ;
ed611565 165
483b0559 166 if(fToUnfold)
167 MakeUnfolding() ;
168
9e5d2067 169 WriteRecPoints() ;
483b0559 170
171 if(strstr(option,"deb"))
172 PrintRecPoints(option) ;
173
fdebddeb 174 //increment the total number of recpoints per run
88cb7938 175 fRecPointsInRun += gime->ECARecPoints()->GetEntriesFast() ;
88cb7938 176 }
483b0559 177
88cb7938 178 Unload();
179
483b0559 180 if(strstr(option,"tim")){
181 gBenchmark->Stop("EMCALClusterizer");
fdebddeb 182 printf("Exec took %f seconds for Clusterizing %f seconds per event",
4fbfb21c 183 gBenchmark->GetCpuTime("EMCALClusterizer"), gBenchmark->GetCpuTime("EMCALClusterizer")/nEvents ) ;
483b0559 184 }
1963b290 185
186 SaveHists();
187
483b0559 188}
189
190//____________________________________________________________________________
70a93198 191Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
483b0559 192 Int_t nPar, Float_t * fitparameters) const
193{
194 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
195 // The initial values for fitting procedure are set equal to the positions of local maxima.
196 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
197
88cb7938 198 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
483b0559 199 TClonesArray * digits = gime->Digits() ;
200
483b0559 201 gMinuit->mncler(); // Reset Minuit's list of paramters
202 gMinuit->SetPrintLevel(-1) ; // No Printout
203 gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
204 // To set the address of the minimization function
483b0559 205 TList * toMinuit = new TList();
206 toMinuit->AddAt(emcRP,0) ;
207 toMinuit->AddAt(digits,1) ;
208
209 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
210
211 // filling initial values for fit parameters
212 AliEMCALDigit * digit ;
213
214 Int_t ierflg = 0;
215 Int_t index = 0 ;
216 Int_t nDigits = (Int_t) nPar / 3 ;
217
218 Int_t iDigit ;
219
220 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
221
222 for(iDigit = 0; iDigit < nDigits; iDigit++){
a0636361 223 digit = maxAt[iDigit];
483b0559 224
70a93198 225 Int_t relid[2] ;
483b0559 226 Float_t x = 0.;
227 Float_t z = 0.;
228 geom->AbsToRelNumbering(digit->GetId(), relid) ;
229 geom->PosInAlice(relid, x, z) ;
230
231 Float_t energy = maxAtEnergy[iDigit] ;
232
233 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
234 index++ ;
235 if(ierflg != 0){
9859bfc0 236 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : x = %f", x ) ;
483b0559 237 return kFALSE;
238 }
239 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
240 index++ ;
241 if(ierflg != 0){
9859bfc0 242 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : z = %f", z) ;
483b0559 243 return kFALSE;
244 }
245 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
246 index++ ;
247 if(ierflg != 0){
9859bfc0 248 Error("FindFit", "EMCAL Unfolding Unable to set initial value for fit procedure : energy = %f", energy) ;
483b0559 249 return kFALSE;
250 }
251 }
252
253 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
254 // depends on it.
255 Double_t p1 = 1.0 ;
256 Double_t p2 = 0.0 ;
257
258 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
259 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
260 gMinuit->SetMaxIterations(5);
261 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
262
263 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
264
265 if(ierflg == 4){ // Minimum not found
9859bfc0 266 Error("FindFit", "EMCAL Unfolding Fit not converged, cluster abandoned " ) ;
483b0559 267 return kFALSE ;
268 }
269 for(index = 0; index < nPar; index++){
270 Double_t err ;
271 Double_t val ;
272 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
273 fitparameters[index] = val ;
274 }
275
276 delete toMinuit ;
277 return kTRUE;
278
279}
280
281//____________________________________________________________________________
282void AliEMCALClusterizerv1::GetCalibrationParameters()
283{
fdebddeb 284 // Gets the parameters for the calibration from the digitizer
88cb7938 285 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
483b0559 286
88cb7938 287 if ( !gime->Digitizer() )
288 gime->LoadDigitizer();
fdebddeb 289 AliEMCALDigitizer * dig = gime->Digitizer();
ed611565 290
88cb7938 291 fADCchannelECA = dig->GetECAchannel() ;
292 fADCpedestalECA = dig->GetECApedestal();
cf24bdc4 293//PH cout<<"ChannelECA, peds "<<fADCchannelECA<<" "<<fADCpedestalECA<<endl;
483b0559 294}
05a92d59 295
483b0559 296//____________________________________________________________________________
297void AliEMCALClusterizerv1::Init()
298{
299 // Make all memory allocations which can not be done in default constructor.
300 // Attach the Clusterizer task to the list of EMCAL tasks
301
7c193632 302 AliEMCALGetter * gime = AliEMCALGetter::Instance();
303 if(!gime)
304 gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());
483b0559 305
88cb7938 306 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
cf24bdc4 307//PH cout<<"gime,geom "<<gime<<" "<<geom<<endl;
05a92d59 308
1963b290 309//Sub fNTowers = geom->GetNZ() * geom->GetNPhi() ;
310 fNTowers =400;
483b0559 311 if(!gMinuit)
312 gMinuit = new TMinuit(100) ;
1963b290 313 //Sub if ( !gime->Clusterizer() )
314 //Sub gime->PostClusterizer(this);
315 BookHists();
cf24bdc4 316//PH cout<<"hists booked "<<endl;
483b0559 317}
318
319//____________________________________________________________________________
839828a6 320void AliEMCALClusterizerv1::InitParameters()
fdebddeb 321{
322 // Initializes the parameters for the Clusterizer
323 fNumberOfECAClusters = 0;
1963b290 324
325 fECAClusteringThreshold = 0.5; // value obtained from Alexei
326 fECALocMaxCut = 0.03;
327
88cb7938 328 fECAW0 = 4.5 ;
839828a6 329 fTimeGate = 1.e-8 ;
839828a6 330 fToUnfold = kFALSE ;
fdebddeb 331 fRecPointsInRun = 0 ;
0b7f0b26 332 fMinECut = 0.3;
839828a6 333}
334
335//____________________________________________________________________________
483b0559 336Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
337{
338 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
339 // = 1 are neighbour
340 // = 2 are not neighbour but do not continue searching
341 // neighbours are defined as digits having at least a common vertex
342 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
343 // which is compared to a digit (d2) not yet in a cluster
344
88cb7938 345 AliEMCALGeometry * geom = AliEMCALGetter::Instance()->EMCALGeometry() ;
483b0559 346
347 Int_t rv = 0 ;
348
70a93198 349 Int_t relid1[2] ;
1963b290 350 //Sub geom->AbsToRelNumbering(d1->GetId(), relid1) ;
351 Int_t nSupMod=0;
352 Int_t nTower=0;
353 Int_t nIphi=0;
354 Int_t nIeta=0;
355 Int_t iphi=0;
356 Int_t ieta=0;
357 geom->GetCellIndex(d1->GetId(), nSupMod,nTower,nIphi,nIeta);
d87bd045 358 geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
1963b290 359 relid1[0]=ieta;
360 relid1[1]=iphi;
361
362
363 Int_t nSupMod1=0;
364 Int_t nTower1=0;
365 Int_t nIphi1=0;
366 Int_t nIeta1=0;
367 Int_t iphi1=0;
368 Int_t ieta1=0;
369 geom->GetCellIndex(d2->GetId(), nSupMod1,nTower1,nIphi1,nIeta1);
d87bd045 370 geom->GetCellPhiEtaIndexInSModule(nSupMod1, nTower1,nIphi1,nIeta1, iphi1,ieta1);
1963b290 371 Int_t relid2[2] ;
372 relid2[0]=ieta1;
373 relid2[1]=iphi1;
374
375 //Sub geom->AbsToRelNumbering(d2->GetId(), relid2) ;
ed611565 376
483b0559 377
70a93198 378 Int_t rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;
379 Int_t coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;
380
381 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
70a93198 382 rv = 1 ;
383 }
483b0559 384 else {
70a93198 385 if((relid2[0] > relid1[0]) && (relid2[1] > relid1[1]+1))
386 rv = 2; // Difference in row numbers is too large to look further
483b0559 387 }
70a93198 388
ed611565 389 if (gDebug == 2 )
1963b290 390if(rv==1)printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n",
70a93198 391 rv, d1->GetId(), relid1[0], relid1[1],
392 d2->GetId(), relid2[0], relid2[1]) ;
ed611565 393
483b0559 394 return rv ;
395}
396
88cb7938 397//____________________________________________________________________________
398void AliEMCALClusterizerv1::Unload()
399{
fdebddeb 400 // Unloads the Digits and RecPoints
88cb7938 401 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
402 gime->EmcalLoader()->UnloadDigits() ;
403 gime->EmcalLoader()->UnloadRecPoints() ;
404}
483b0559 405
483b0559 406//____________________________________________________________________________
9e5d2067 407void AliEMCALClusterizerv1::WriteRecPoints()
483b0559 408{
409
410 // Creates new branches with given title
411 // fills and writes into TreeR.
412
88cb7938 413 AliEMCALGetter *gime = AliEMCALGetter::Instance() ;
ed611565 414
88cb7938 415 TObjArray * aECARecPoints = gime->ECARecPoints() ;
ed611565 416
05a92d59 417 TClonesArray * digits = gime->Digits() ;
88cb7938 418 TTree * treeR = gime->TreeR(); ;
05a92d59 419
05a92d59 420 Int_t index ;
ed611565 421
ed611565 422 //Evaluate position, dispersion and other RecPoint properties for EC section
88cb7938 423 for(index = 0; index < aECARecPoints->GetEntries(); index++)
70a93198 424 (dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(index)))->EvalAll(fECAW0,digits) ;
ed611565 425
88cb7938 426 aECARecPoints->Sort() ;
483b0559 427
88cb7938 428 for(index = 0; index < aECARecPoints->GetEntries(); index++)
70a93198 429 (dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(index)))->SetIndexInList(index) ;
483b0559 430
88cb7938 431 aECARecPoints->Expand(aECARecPoints->GetEntriesFast()) ;
ed611565 432
483b0559 433 Int_t bufferSize = 32000 ;
9859bfc0 434 Int_t splitlevel = 0 ;
ed611565 435
436 //EC section branch
88cb7938 437 TBranch * branchECA = treeR->Branch("EMCALECARP","TObjArray",&aECARecPoints,bufferSize,splitlevel);
438 branchECA->SetTitle(BranchName());
ed611565 439
88cb7938 440 branchECA->Fill() ;
483b0559 441
cf24bdc4 442 gime->WriteRecPoints("OVERWRITE");
443 gime->WriteClusterizer("OVERWRITE");
483b0559 444}
445
446//____________________________________________________________________________
447void AliEMCALClusterizerv1::MakeClusters()
448{
449 // Steering method to construct the clusters stored in a list of Reconstructed Points
450 // A cluster is defined as a list of neighbour digits
05a92d59 451
88cb7938 452 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
ed611565 453
454 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
455
88cb7938 456 TObjArray * aECARecPoints = gime->ECARecPoints() ;
ed611565 457
88cb7938 458 aECARecPoints->Delete() ;
ed611565 459
05a92d59 460 TClonesArray * digits = gime->Digits() ;
483b0559 461 TClonesArray * digitsC = dynamic_cast<TClonesArray*>(digits->Clone()) ;
4aa81449 462
463 // Clusterization starts
464 TIter nextdigit(digitsC) ;
465 AliEMCALDigit * digit;
1963b290 466
467//just for hist
468/*
469 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
470 digitAmp->Fill(digit->GetAmp());
471 }
472 */
473///////////
474
483b0559 475 while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
476 AliEMCALRecPoint * clu = 0 ;
477
b048087c 478 TArrayI clusterECAdigitslist(50000);
4aa81449 479
1963b290 480//Sub if (gDebug == 2) {
481 // printf("MakeClusters: id = %d, ene = %f , thre = %f \n", digit->GetId(),Calibrate(digit->GetAmp()),
482// fECAClusteringThreshold) ;
483// }
484////////////////////////temp solution, embedding///////////////////
485 int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
486 int iphi=0, ieta=0;
487 geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
d87bd045 488 geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
1963b290 489
490 ///////////////////////////////////
491
492//cout<<ieta<<" "<<iphi<<endl;
493
494// if ( geom->IsInECA(digit->GetId()) && (Calibrate(digit->GetAmp()) > fECAClusteringThreshold ) ){
495 if (geom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) > fECAClusteringThreshold ) ){
cf24bdc4 496 //if(addOn[nSupMod-1][ieta-1][iphi-1]>0)cout<<"11 digit, add "<<ieta<<" "<<iphi<<" "<<addOn[nSupMod-1][ieta-1][iphi-1]<<" "<<digit->GetAmp()<<endl;
1963b290 497// cout<<"crossed the threshold "<<endl;
4aa81449 498 Int_t iDigitInECACluster = 0;
499 // start a new Tower RecPoint
500 if(fNumberOfECAClusters >= aECARecPoints->GetSize())
88cb7938 501 aECARecPoints->Expand(2*fNumberOfECAClusters+1) ;
4aa81449 502 AliEMCALRecPoint * rp = new AliEMCALRecPoint("") ;
503 aECARecPoints->AddAt(rp, fNumberOfECAClusters) ;
504 clu = dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(fNumberOfECAClusters)) ;
505 fNumberOfECAClusters++ ;
1963b290 506 clu->AddDigit(*digit, Calibrate(digit->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1])) ;
4aa81449 507 clusterECAdigitslist[iDigitInECACluster] = digit->GetIndexInList() ;
508 iDigitInECACluster++ ;
1963b290 509// cout<<" 1st setp:cluno, digNo "<<fNumberOfECAClusters<<" "<<iDigitInECACluster<<endl;
4aa81449 510 digitsC->Remove(digit) ;
511 if (gDebug == 2 )
1963b290 512 printf("MakeClusters: OK id = %d, ene = %f , thre = %f \n", digit->GetId(),Calibrate(digit->GetAmp()), fECAClusteringThreshold) ;
483b0559 513 nextdigit.Reset() ;
514
515 AliEMCALDigit * digitN ;
ed611565 516 Int_t index = 0 ;
517
4aa81449 518 // Find the neighbours
88cb7938 519 while (index < iDigitInECACluster){ // scan over digits already in cluster
520 digit = (AliEMCALDigit*)digits->At(clusterECAdigitslist[index]) ;
483b0559 521 index++ ;
c8be63ae 522 while ( (digitN = (AliEMCALDigit *)nextdigit())) { // scan over the reduced list of digits
1963b290 523// cout<<"we have new digit "<<endl;
c8be63ae 524 // check that the digit is above the min E Cut
1963b290 525 ////////////////////
526 geom->GetCellIndex(digitN->GetId(), nSupMod,nTower,nIphi,nIeta);
d87bd045 527 geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
1963b290 528 ////////////////
529 if(Calibrate(digitN->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) < fMinECut ) digitsC->Remove(digitN);
530// cout<<" new digit above ECut "<<endl;
fdebddeb 531 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
1963b290 532// cout<<" new digit neighbour?? "<<ineb<<endl;
483b0559 533 switch (ineb ) {
534 case 0 : // not a neighbour
535 break ;
536 case 1 : // are neighbours
cf24bdc4 537 //if(addOn[nSupMod-1][ieta-1][iphi-1]>0)cout<<"22 digit, add "<<nSupMod<<" "<<ieta<<" "<<iphi<<" "<<addOn[nSupMod-1][ieta-1][iphi-1]<<" "<<digit->GetAmp()<<endl;
1963b290 538 clu->AddDigit(*digitN, Calibrate( digitN->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) ) ;
88cb7938 539 clusterECAdigitslist[iDigitInECACluster] = digitN->GetIndexInList() ;
540 iDigitInECACluster++ ;
1963b290 541// cout<<"when neighbour: cluno, digNo "<<digit->GetId()<<" "<<fNumberOfECAClusters<<" "<<iDigitInECACluster<<endl;
483b0559 542 digitsC->Remove(digitN) ;
1963b290 543// break ;
544// case 2 : // too far from each other
545//Subh goto endofloop1;
546// cout<<"earlier go to end of loop"<<endl;
483b0559 547 } // switch
1963b290 548 //cout<<"in nextDigit loop "<<fNumberOfECAClusters<<" "<<iDigitInECACluster<<endl;
483b0559 549 } // while digitN
550
1963b290 551//Sub endofloop1: ;
483b0559 552 nextdigit.Reset() ;
88cb7938 553 } // loop over ECA cluster
c8be63ae 554 } // energy threshold
1963b290 555 else if(Calibrate(digit->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) < fMinECut ){
cf24bdc4 556 //if(addOn[nSupMod-1][ieta-1][iphi-1]>0)cout<<"33 digit, add "<<ieta<<" "<<iphi<<" "<<addOn[nSupMod-1][ieta-1][iphi-1]<<" "<<digit->GetAmp()<<endl;
4aa81449 557 digitsC->Remove(digit);
c8be63ae 558 }
1963b290 559 //cout<<"after endofloop: cluno, digNo "<<fNumberOfECAClusters<<endl;
9859bfc0 560 } // while digit
483b0559 561 delete digitsC ;
cf24bdc4 562cout<<"total no of clusters "<<fNumberOfECAClusters<<" from "<<digits->GetEntriesFast()<<" digits"<<endl;
483b0559 563}
564
565//____________________________________________________________________________
fdebddeb 566void AliEMCALClusterizerv1::MakeUnfolding() const
483b0559 567{
568 Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
9859bfc0 569
483b0559 570}
571
572//____________________________________________________________________________
573Double_t AliEMCALClusterizerv1::ShowerShape(Double_t r)
574{
575 // Shape of the shower (see EMCAL TDR)
576 // If you change this function, change also the gradient evaluation in ChiSquare()
577
578 Double_t r4 = r*r*r*r ;
579 Double_t r295 = TMath::Power(r, 2.95) ;
580 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
581 return shape ;
582}
583
584//____________________________________________________________________________
70a93198 585void AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * /*iniTower*/,
9e5d2067 586 Int_t /*nMax*/,
587 AliEMCALDigit ** /*maxAt*/,
fdebddeb 588 Float_t * /*maxAtEnergy*/) const
483b0559 589{
590 // Performs the unfolding of a cluster with nMax overlapping showers
483b0559 591
9859bfc0 592 Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
173558f2 593
483b0559 594}
595
596//_____________________________________________________________________________
9e5d2067 597void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
598 Double_t & /*fret*/,
599 Double_t * /*x*/, Int_t /*iflag*/)
483b0559 600{
601 // Calculates the Chi square for the cluster unfolding minimization
602 // Number of parameters, Gradient, Chi squared, parameters, what to do
483b0559 603
9859bfc0 604 ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
483b0559 605}
483b0559 606//____________________________________________________________________________
9e5d2067 607void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
483b0559 608{
609 // Print clusterizer parameters
610
9859bfc0 611 TString message("\n") ;
612
483b0559 613 if( strcmp(GetName(), "") !=0 ){
614
615 // Print parameters
616
fdebddeb 617 TString taskName(GetName()) ;
483b0559 618 taskName.ReplaceAll(Version(), "") ;
9859bfc0 619
fdebddeb 620 printf("--------------- ");
621 printf(taskName.Data()) ;
622 printf(" ");
623 printf(GetTitle()) ;
624 printf("-----------\n");
625 printf("Clusterizing digits from the file: ");
626 printf(taskName.Data());
627 printf("\n Branch: ");
628 printf(GetName());
629 printf("\n ECA Local Maximum cut = %f", fECALocMaxCut);
b481a360 630 printf("\n ECA Logarithmic weight = %f", fECAW0);
483b0559 631 if(fToUnfold)
fdebddeb 632 printf("\nUnfolding on\n");
483b0559 633 else
fdebddeb 634 printf("\nUnfolding off\n");
483b0559 635
fdebddeb 636 printf("------------------------------------------------------------------");
483b0559 637 }
638 else
fdebddeb 639 printf("AliEMCALClusterizerv1 not initialized ") ;
483b0559 640}
173558f2 641
483b0559 642//____________________________________________________________________________
643void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
644{
b481a360 645 // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
646
88cb7938 647 TObjArray * aECARecPoints = AliEMCALGetter::Instance()->ECARecPoints() ;
fdebddeb 648 printf("PrintRecPoints: Clusterization result:") ;
ed611565 649
650 printf("event # %d\n", gAlice->GetEvNumber() ) ;
fdebddeb 651 printf(" Found %d ECA Rec Points\n ",
652 aECARecPoints->GetEntriesFast()) ;
483b0559 653
88cb7938 654 fRecPointsInRun += aECARecPoints->GetEntriesFast() ;
11f9c5ff 655
483b0559 656 if(strstr(option,"all")) {
fdebddeb 657 Int_t index =0;
ed611565 658 printf("\n-----------------------------------------------------------------------\n") ;
659 printf("Clusters in ECAL section\n") ;
660 printf("Index Ene(GeV) Multi Module phi r theta X Y Z Dispersion Lambda 1 Lambda 2 # of prim Primaries list\n") ;
1963b290 661 Float_t maxE=0;
662 Float_t maxL1=0;
663 Float_t maxL2=0;
664 Float_t maxDis=0;
665
88cb7938 666 for (index = 0 ; index < aECARecPoints->GetEntries() ; index++) {
70a93198 667 AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(aECARecPoints->At(index)) ;
483b0559 668 TVector3 globalpos;
1963b290 669 //rp->GetGlobalPosition(globalpos);
ed611565 670 TVector3 localpos;
671 rp->GetLocalPosition(localpos);
483b0559 672 Float_t lambda[2];
673 rp->GetElipsAxis(lambda);
674 Int_t * primaries;
675 Int_t nprimaries;
676 primaries = rp->GetPrimaries(nprimaries);
70a93198 677 printf("\n%6d %8.4f %3d %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4.1f %4f %4f %2d : ",
678 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
ed611565 679 globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(),
680 rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ;
1963b290 681 /////////////
682 if(rp->GetEnergy()>maxE){
683 maxE=rp->GetEnergy();
684 maxL1=lambda[0];
685 maxL2=lambda[1];
686 maxDis=rp->GetDispersion();
687 }
688 pointE->Fill(rp->GetEnergy());
689 pointL1->Fill(lambda[0]);
690 pointL2->Fill(lambda[1]);
691 pointDis->Fill(rp->GetDispersion());
692 pointMult->Fill(rp->GetMultiplicity());
693 /////////////
9859bfc0 694 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
ed611565 695 printf("%d ", primaries[iprimary] ) ;
9859bfc0 696 }
483b0559 697 }
1963b290 698
699 MaxE->Fill(maxE);
700 MaxL1->Fill(maxL1);
701 MaxL2->Fill(maxL2);
702 MaxDis->Fill(maxDis);
703
704
ed611565 705 printf("\n-----------------------------------------------------------------------\n");
483b0559 706 }
707}
1963b290 708 void AliEMCALClusterizerv1::BookHists()
709{
710 pointE = new TH1F("pointE","point energy", 2000, 0.0, 150.);
711 pointL1 = new TH1F("pointL1","point L1", 1000, 0.0, 3.);
712 pointL2 = new TH1F("pointL2","point L2", 1000, 0.0, 3.);
713 pointDis = new TH1F("pointDis","point Dis", 1000, 0.0, 3.);
714 pointMult = new TH1F("pointMult","point Mult", 100, 0.0, 100.);
715 digitAmp = new TH1F("digitAmp","Digit Amplitude", 2000, 0.0, 5000.);
716 MaxE = new TH1F("maxE","Max point energy", 2000, 0.0, 150.);
717 MaxL1 = new TH1F("maxL1","Max point L1", 1000, 0.0, 3.);
718 MaxL2 = new TH1F("maxL2","Max point L2", 1000, 0.0, 3.);
719 MaxDis = new TH1F("maxDis","Max point Dis", 1000, 0.0, 3.);
720}
721void AliEMCALClusterizerv1::SaveHists()
722{
723 recofile=new TFile("reco.root","RECREATE");
724 pointE->Write();
725 pointL1->Write();
726 pointL2->Write();
727 pointDis->Write();
728 pointMult->Write();
729 digitAmp->Write();
730 MaxE->Write();
731 MaxL1->Write();
732 MaxL2->Write();
733 MaxDis->Write();
734recofile->Close();
735}
736
737void AliEMCALClusterizerv1::ReadFile()
738{
739 return; // 3-jan-05
740 FILE *fp = fopen("hijing1.dat","r");
741 for(Int_t line=0;line<9113;line++){
742 Int_t eg,l1,l2,sm;
743 Int_t ncols0;
744 ncols0 = fscanf(fp,"%d %d %d %d",&sm,&l1,&l2,&eg);
745 // cout<<eg<<" "<<l1<<" "<<l2<<endl;
746 addOn[sm-1][l1-1][l2-1]=eg;
747 //addOn[sm-1][l1-1][l2-1]=0;
748 }
749}
750