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