]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSClusterizerv1.cxx
eff C++ warnings corrected.
[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
702ab87e 18/* History of cvs commits:
19 *
20 * $Log$
f36dfa9a 21 * Revision 1.104 2007/04/27 16:55:53 kharlov
22 * Calibration stops if PHOS CDB objects do not exist
23 *
76d78b7a 24 * Revision 1.103 2007/04/11 11:55:45 policheh
25 * SetDistancesToBadChannels() added.
26 *
af0b0570 27 * Revision 1.102 2007/03/28 19:18:15 kharlov
28 * RecPoints recalculation in TSM removed
29 *
8c1fb709 30 * Revision 1.101 2007/03/06 06:51:27 kharlov
31 * Calculation of cluster properties dep. on vertex posponed to TrackSegmentMaker
32 *
a978bcea 33 * Revision 1.100 2007/01/10 11:07:26 kharlov
34 * Raw digits writing to file (B.Polichtchouk)
35 *
cbddd97f 36 * Revision 1.99 2006/11/07 16:49:51 kharlov
37 * Corrections for next event switch in case of raw data (B.Polichtchouk)
38 *
e1c29d81 39 * Revision 1.98 2006/10/27 17:14:27 kharlov
40 * Introduce AliDebug and AliLog (B.Polichtchouk)
41 *
e1afea68 42 * Revision 1.97 2006/08/29 11:41:19 kharlov
43 * Missing implementation of ctors and = operator are added
44 *
f1133801 45 * Revision 1.96 2006/08/25 16:56:30 kharlov
46 * Compliance with Effective C++
47 *
0378398c 48 * Revision 1.95 2006/08/11 12:36:26 cvetan
49 * Update of the PHOS code needed in order to read and reconstruct the beam test raw data (i.e. without an existing galice.root)
50 *
52c5f046 51 * Revision 1.94 2006/08/07 12:27:49 hristov
52 * Removing obsolete code which affected the event numbering scheme
53 *
4d4750bd 54 * Revision 1.93 2006/08/01 12:20:17 cvetan
55 * 1. Adding a possibility to read and reconstruct an old rcu formatted raw data. This is controlled by an option of AliReconstruction and AliPHOSReconstructor. 2. In case of raw data processing (without galice.root) create the default AliPHOSGeometry object. Most likely this should be moved to the CDB
56 *
f5eaa851 57 * Revision 1.92 2006/04/29 20:26:46 hristov
58 * Separate EMC and CPV calibration (Yu.Kharlov)
59 *
e95226ae 60 * Revision 1.91 2006/04/22 10:30:17 hristov
61 * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
62 *
27a73a5d 63 * Revision 1.90 2006/04/11 15:22:59 hristov
64 * run number in query set to -1: forces AliCDBManager to use its run number (A.Colla)
65 *
28871337 66 * Revision 1.89 2006/03/13 14:05:42 kharlov
67 * Calibration objects for EMC and CPV
68 *
fc6706cb 69 * Revision 1.88 2006/01/11 08:54:52 hristov
70 * Additional protection in case no calibration entry was found
71 *
e521c5a0 72 * Revision 1.87 2005/11/22 08:46:43 kharlov
73 * Updated to new CDB (Boris Polichtchouk)
74 *
ef629623 75 * Revision 1.86 2005/11/14 21:52:43 hristov
76 * Coding conventions
77 *
fdf65bb5 78 * Revision 1.85 2005/09/27 16:08:08 hristov
79 * New version of CDB storage framework (A.Colla)
80 *
9e1ceb13 81 * Revision 1.84 2005/09/21 10:02:47 kharlov
82 * Reading calibration from CDB (Boris Polichtchouk)
83 *
8d6d4fb5 84 * Revision 1.82 2005/09/02 15:43:13 kharlov
85 * Add comments in GetCalibrationParameters and Calibrate
86 *
5c149aeb 87 * Revision 1.81 2005/09/02 14:32:07 kharlov
88 * Calibration of raw data
89 *
44ae287e 90 * Revision 1.80 2005/08/24 15:31:36 kharlov
91 * Setting raw digits flag
92 *
9735c615 93 * Revision 1.79 2005/07/25 15:53:53 kharlov
94 * Read raw data
95 *
3cf4f75f 96 * Revision 1.78 2005/05/28 14:19:04 schutz
97 * Compilation warnings fixed by T.P.
98 *
702ab87e 99 */
100
9a1398dd 101//*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
d15a28e7 102//////////////////////////////////////////////////////////////////////////////
9a1398dd 103// Clusterization class. Performs clusterization (collects neighbouring active cells) and
a4e98857 104// unfolds the clusters having several local maxima.
105// Results are stored in TreeR#, branches PHOSEmcRP (EMC recPoints),
9a1398dd 106// PHOSCpvRP (CPV RecPoints) and AliPHOSClusterizer (Clusterizer with all
f035f6ce 107// parameters including input digits branch title, thresholds etc.)
a4e98857 108// This TTask is normally called from Reconstructioner, but can as well be used in
109// standalone mode.
110// Use Case:
fc12304f 111// root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root", "recpointsname", "digitsname")
a4e98857 112// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
fc12304f 113// // reads gAlice from header file "galice.root", uses digits stored in the branch names "digitsname" (default = "Default")
114// // and saves recpoints in branch named "recpointsname" (default = "digitsname")
a4e98857 115// root [1] cl->ExecuteTask()
9a1398dd 116// //finds RecPoints in all events stored in galice.root
a4e98857 117// root [2] cl->SetDigitsBranch("digits2")
f035f6ce 118// //sets another title for Digitis (input) branch
a4e98857 119// root [3] cl->SetRecPointsBranch("recp2")
f035f6ce 120// //sets another title four output branches
a4e98857 121// root [4] cl->SetEmcLocalMaxCut(0.03)
9a1398dd 122// //set clusterization parameters
a4e98857 123// root [5] cl->ExecuteTask("deb all time")
9a1398dd 124// //once more finds RecPoints options are
125// // deb - print number of found rec points
126// // deb all - print number of found RecPoints and some their characteristics
127// // time - print benchmarking results
ed4205d8 128
d15a28e7 129// --- ROOT system ---
130
131#include "TMath.h"
9a1398dd 132#include "TMinuit.h"
133#include "TTree.h"
9a1398dd 134#include "TBenchmark.h"
d15a28e7 135
136// --- Standard library ---
137
d15a28e7 138// --- AliRoot header files ---
351dd634 139#include "AliLog.h"
e1c29d81 140#include "AliRunLoader.h"
a978bcea 141#include "AliGenerator.h"
88cb7938 142#include "AliPHOSGetter.h"
e957fea8 143#include "AliPHOSGeometry.h"
d15a28e7 144#include "AliPHOSClusterizerv1.h"
88cb7938 145#include "AliPHOSEmcRecPoint.h"
9a1398dd 146#include "AliPHOSCpvRecPoint.h"
d15a28e7 147#include "AliPHOSDigit.h"
9a1398dd 148#include "AliPHOSDigitizer.h"
ba898748 149#include "AliPHOSCalibrationDB.h"
9e1ceb13 150#include "AliCDBManager.h"
eaa5012b 151#include "AliCDBStorage.h"
fdf65bb5 152#include "AliCDBEntry.h"
d15a28e7 153
154ClassImp(AliPHOSClusterizerv1)
f035f6ce 155
d15a28e7 156//____________________________________________________________________________
0378398c 157AliPHOSClusterizerv1::AliPHOSClusterizerv1() :
158 AliPHOSClusterizer(),
159 fDefaultInit(0), fEmcCrystals(0), fToUnfold(0),
160 fWrite(0), fNumberOfEmcClusters(0), fNumberOfCpvClusters(0),
161 fCalibData(0), fADCchanelEmc(0), fADCpedestalEmc(0),
162 fADCchanelCpv(0), fADCpedestalCpv(0), fEmcClusteringThreshold(0),
163 fCpvClusteringThreshold(0), fEmcMinE(0), fCpvMinE(0),
164 fEmcLocMaxCut(0), fW0(0), fCpvLocMaxCut(0),
165 fW0CPV(0), fRecPointsInRun(0), fEmcTimeGate(0),
166 fIsOldRCUFormat(0)
d15a28e7 167{
f035f6ce 168 // default ctor (to be used mainly by Streamer)
f035f6ce 169
8d0f3f77 170 InitParameters() ;
92f521a9 171 fDefaultInit = kTRUE ;
9a1398dd 172}
7b7c1533 173
9a1398dd 174//____________________________________________________________________________
0378398c 175AliPHOSClusterizerv1::AliPHOSClusterizerv1(const TString alirunFileName, const TString eventFolderName) :
176 AliPHOSClusterizer(alirunFileName, eventFolderName),
177 fDefaultInit(0), fEmcCrystals(0), fToUnfold(0),
178 fWrite(0), fNumberOfEmcClusters(0), fNumberOfCpvClusters(0),
179 fCalibData(0), fADCchanelEmc(0), fADCpedestalEmc(0),
180 fADCchanelCpv(0), fADCpedestalCpv(0), fEmcClusteringThreshold(0),
181 fCpvClusteringThreshold(0), fEmcMinE(0), fCpvMinE(0),
182 fEmcLocMaxCut(0), fW0(0), fCpvLocMaxCut(0),
183 fW0CPV(0), fRecPointsInRun(0), fEmcTimeGate(0),
184 fIsOldRCUFormat(0)
9a1398dd 185{
a4e98857 186 // ctor with the indication of the file where header Tree and digits Tree are stored
9a1398dd 187
8d0f3f77 188 InitParameters() ;
2bd5457f 189 Init() ;
92f521a9 190 fDefaultInit = kFALSE ;
9a1398dd 191}
fc12304f 192
f1133801 193//____________________________________________________________________________
194AliPHOSClusterizerv1::AliPHOSClusterizerv1(const AliPHOSClusterizerv1 & obj) :
195 AliPHOSClusterizer(obj),
196 fDefaultInit(0), fEmcCrystals(0), fToUnfold(0),
197 fWrite(0), fNumberOfEmcClusters(0), fNumberOfCpvClusters(0),
198 fCalibData(0), fADCchanelEmc(0), fADCpedestalEmc(0),
199 fADCchanelCpv(0), fADCpedestalCpv(0), fEmcClusteringThreshold(0),
200 fCpvClusteringThreshold(0), fEmcMinE(0), fCpvMinE(0),
201 fEmcLocMaxCut(0), fW0(0), fCpvLocMaxCut(0),
202 fW0CPV(0), fRecPointsInRun(0), fEmcTimeGate(0),
203 fIsOldRCUFormat(0)
204{
205 // Copy constructor
206}
9688c1dd 207//____________________________________________________________________________
208 AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
209{
0bc3b8ed 210 // dtor
211
9688c1dd 212}
fc12304f 213//____________________________________________________________________________
214const TString AliPHOSClusterizerv1::BranchName() const
215{
88cb7938 216 return GetName();
fc12304f 217}
218
9a1398dd 219//____________________________________________________________________________
e95226ae 220Float_t AliPHOSClusterizerv1::CalibrateEMC(Float_t amp, Int_t absId)
88cb7938 221{
e95226ae 222 // Convert EMC measured amplitude into real energy.
5c149aeb 223 // Calibration parameters are taken from calibration data base for raw data,
224 // or from digitizer parameters for simulated data.
225
44ae287e 226 if(fCalibData){
227 Int_t relId[4];
228 AliPHOSGetter *gime = AliPHOSGetter::Instance();
229 gime->PHOSGeometry()->AbsToRelNumbering(absId,relId) ;
230 Int_t module = relId[0];
8d6d4fb5 231 Int_t column = relId[3];
eaa5012b 232 Int_t row = relId[2];
e95226ae 233 if(absId <= fEmcCrystals) { // this is EMC
44ae287e 234 fADCchanelEmc = fCalibData->GetADCchannelEmc (module,column,row);
e95226ae 235 return amp*fADCchanelEmc ;
44ae287e 236 }
e95226ae 237 }
238 else{ //simulation
239 if(absId <= fEmcCrystals) // this is EMC
240 return fADCpedestalEmc + amp*fADCchanelEmc ;
241 }
242 return 0;
243}
244
245//____________________________________________________________________________
246Float_t AliPHOSClusterizerv1::CalibrateCPV(Int_t amp, Int_t absId)
247{
248 // Convert digitized CPV amplitude into charge.
249 // Calibration parameters are taken from calibration data base for raw data,
250 // or from digitizer parameters for simulated data.
251
252 if(fCalibData){
253 Int_t relId[4];
254 AliPHOSGetter *gime = AliPHOSGetter::Instance();
255 gime->PHOSGeometry()->AbsToRelNumbering(absId,relId) ;
256 Int_t module = relId[0];
257 Int_t column = relId[3];
258 Int_t row = relId[2];
259 if(absId > fEmcCrystals) { // this is CPV
fc6706cb 260 fADCchanelCpv = fCalibData->GetADCchannelCpv (module,column,row);
261 fADCpedestalCpv = fCalibData->GetADCpedestalCpv(module,column,row);
262 return fADCpedestalCpv + amp*fADCchanelCpv ;
263 }
ba898748 264 }
265 else{ //simulation
e95226ae 266 if(absId > fEmcCrystals) // this is CPV
ba898748 267 return fADCpedestalCpv+ amp*fADCchanelCpv ;
268 }
e95226ae 269 return 0;
3758d9fc 270}
548f0134 271
3758d9fc 272//____________________________________________________________________________
eabde521 273void AliPHOSClusterizerv1::Exec(Option_t *option)
a4e98857 274{
eabde521 275 // Steering method to perform clusterization for events
276 // in the range from fFirstEvent to fLastEvent.
277 // This range is optionally set by SetEventRange().
278 // if fLastEvent=-1 (by default), then process events until the end.
9a1398dd 279
280 if(strstr(option,"tim"))
7d493c2c 281 gBenchmark->Start("PHOSClusterizer");
9a1398dd 282
eabde521 283 if(strstr(option,"print")) {
88cb7938 284 Print() ;
eabde521 285 return ;
286 }
bed9e3fb 287
eaa5012b 288 GetCalibrationParameters() ;
289
ba898748 290 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
9735c615 291 if (fRawReader == 0)
292 gime->SetRawDigits(kFALSE);
293 else
294 gime->SetRawDigits(kTRUE);
88cb7938 295
9cae45af 296 if (fLastEvent == -1)
297 fLastEvent = gime->MaxEvent() - 1 ;
298 else
fb43ada4 299 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // one event at the time
eabde521 300 Int_t nEvents = fLastEvent - fFirstEvent + 1;
ba898748 301
9cae45af 302 Int_t ievent ;
9cae45af 303 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
3cf4f75f 304 if (fRawReader == 0)
305 gime->Event(ievent ,"D"); // Read digits from simulated data
f5eaa851 306 else {
e1c29d81 307 AliRunLoader * rl = AliRunLoader::GetRunLoader(gime->PhosLoader()->GetTitle());
308 rl->GetEvent(ievent);
f5eaa851 309 gime->Event(fRawReader,"W",fIsOldRCUFormat); // Read digits from raw data
310 }
548f0134 311 fNumberOfEmcClusters = fNumberOfCpvClusters = 0 ;
9a1398dd 312
88cb7938 313 MakeClusters() ;
e1afea68 314
315 AliDebug(2,Form(" ---- Printing clusters (%d) of event %d.\n",
316 gime->EmcRecPoints()->GetEntries(),ievent));
317 if(AliLog::GetGlobalDebugLevel()>1)
318 gime->EmcRecPoints()->Print();
88cb7938 319
fc12304f 320 if(fToUnfold)
321 MakeUnfolding() ;
af0b0570 322
88cb7938 323 WriteRecPoints();
7b7c1533 324
94de8339 325 if(strstr(option,"deb"))
326 PrintRecPoints(option) ;
327
88cb7938 328 //increment the total number of recpoints per run
fbf811ec 329 fRecPointsInRun += gime->EmcRecPoints()->GetEntriesFast() ;
330 fRecPointsInRun += gime->CpvRecPoints()->GetEntriesFast() ;
27a73a5d 331 }
88cb7938 332
ba898748 333 if(fWrite) //do not unload in "on flight" mode
334 Unload();
9a1398dd 335
336 if(strstr(option,"tim")){
337 gBenchmark->Stop("PHOSClusterizer");
351dd634 338 AliInfo(Form("took %f seconds for Clusterizing %f seconds per event \n",
21cd0c07 339 gBenchmark->GetCpuTime("PHOSClusterizer"),
351dd634 340 gBenchmark->GetCpuTime("PHOSClusterizer")/nEvents )) ;
88cb7938 341 }
9a1398dd 342}
343
344//____________________________________________________________________________
a0636361 345Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
d1de15f5 346 Int_t nPar, Float_t * fitparameters) const
9a1398dd 347{
348 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
a4e98857 349 // The initial values for fitting procedure are set equal to the positions of local maxima.
f035f6ce 350 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
9a1398dd 351
88cb7938 352
353 AliPHOSGetter * gime = AliPHOSGetter::Instance();
354 TClonesArray * digits = gime->Digits();
7b7c1533 355
9a1398dd 356 gMinuit->mncler(); // Reset Minuit's list of paramters
357 gMinuit->SetPrintLevel(-1) ; // No Printout
358 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
359 // To set the address of the minimization function
360
361 TList * toMinuit = new TList();
362 toMinuit->AddAt(emcRP,0) ;
7b7c1533 363 toMinuit->AddAt(digits,1) ;
9a1398dd 364
365 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
366
367 // filling initial values for fit parameters
368 AliPHOSDigit * digit ;
369
370 Int_t ierflg = 0;
371 Int_t index = 0 ;
372 Int_t nDigits = (Int_t) nPar / 3 ;
373
374 Int_t iDigit ;
375
7b7c1533 376 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
9a1398dd 377
378 for(iDigit = 0; iDigit < nDigits; iDigit++){
a0636361 379 digit = maxAt[iDigit];
9a1398dd 380
381 Int_t relid[4] ;
7b7c1533 382 Float_t x = 0.;
383 Float_t z = 0.;
384 geom->AbsToRelNumbering(digit->GetId(), relid) ;
385 geom->RelPosInModule(relid, x, z) ;
9a1398dd 386
387 Float_t energy = maxAtEnergy[iDigit] ;
388
389 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
390 index++ ;
391 if(ierflg != 0){
21cd0c07 392 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : x = %f\n", x ) ;
9a1398dd 393 return kFALSE;
394 }
395 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
396 index++ ;
397 if(ierflg != 0){
21cd0c07 398 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : z =%f\n", z ) ;
9a1398dd 399 return kFALSE;
400 }
401 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
402 index++ ;
403 if(ierflg != 0){
21cd0c07 404 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : energy = %f\n", energy ) ;
9a1398dd 405 return kFALSE;
406 }
407 }
408
409 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
410 // depends on it.
411 Double_t p1 = 1.0 ;
412 Double_t p2 = 0.0 ;
413
414 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
415 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
416 gMinuit->SetMaxIterations(5);
417 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
418
419 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
420
421 if(ierflg == 4){ // Minimum not found
21cd0c07 422 Warning("FindFit", "PHOS Unfolding fit not converged, cluster abandoned\n" );
9a1398dd 423 return kFALSE ;
424 }
425 for(index = 0; index < nPar; index++){
426 Double_t err ;
427 Double_t val ;
428 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
429 fitparameters[index] = val ;
430 }
431
432 delete toMinuit ;
433 return kTRUE;
c3c187e5 434
d15a28e7 435}
436
3758d9fc 437//____________________________________________________________________________
438void AliPHOSClusterizerv1::GetCalibrationParameters()
439{
5c149aeb 440 // Set calibration parameters:
eaa5012b 441 // if calibration database exists, they are read from database,
76d78b7a 442 // otherwise, reconstruction stops in the constructor of AliPHOSCalibData
5c149aeb 443 //
ef629623 444 // It is a user responsilibity to open CDB before reconstruction, for example:
445 // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
fbf811ec 446
28871337 447 fCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
f36dfa9a 448 if (fCalibData->GetCalibDataEmc() == 0)
449 AliFatal("Calibration parameters for PHOS EMC not found. Stop reconstruction.\n");
450 if (fCalibData->GetCalibDataCpv() == 0)
451 AliFatal("Calibration parameters for PHOS CPV not found. Stop reconstruction.\n");
452
3758d9fc 453}
548f0134 454
d15a28e7 455//____________________________________________________________________________
a4e98857 456void AliPHOSClusterizerv1::Init()
457{
2bd5457f 458 // Make all memory allocations which can not be done in default constructor.
459 // Attach the Clusterizer task to the list of PHOS tasks
ba898748 460
461 AliPHOSGetter* gime = AliPHOSGetter::Instance() ;
462 if(!gime)
463 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
fbf811ec 464
88cb7938 465 AliPHOSGeometry * geom = gime->PHOSGeometry();
ba898748 466
3758d9fc 467 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
468
7b7c1533 469 if(!gMinuit)
88cb7938 470 gMinuit = new TMinuit(100);
fbf811ec 471
ba898748 472 if ( !gime->Clusterizer() ) {
88cb7938 473 gime->PostClusterizer(this);
ba898748 474 }
9a1398dd 475}
7b7c1533 476
8d0f3f77 477//____________________________________________________________________________
478void AliPHOSClusterizerv1::InitParameters()
479{
480
481 fNumberOfCpvClusters = 0 ;
482 fNumberOfEmcClusters = 0 ;
483
484 fCpvClusteringThreshold = 0.0;
485 fEmcClusteringThreshold = 0.2;
486
487 fEmcLocMaxCut = 0.03 ;
488 fCpvLocMaxCut = 0.03 ;
ba898748 489
490 fEmcMinE = 0.01 ;
491 fCpvMinE = 0.0 ;
492
8d0f3f77 493 fW0 = 4.5 ;
494 fW0CPV = 4.0 ;
495
496 fEmcTimeGate = 1.e-8 ;
497
498 fToUnfold = kTRUE ;
88cb7938 499
8d0f3f77 500 fRecPointsInRun = 0 ;
c2cd5471 501
ba898748 502 fWrite = kTRUE ;
503
44ae287e 504 fCalibData = 0 ;
ae911977 505
eabde521 506 SetEventRange(0,-1) ;
f5eaa851 507
508 fIsOldRCUFormat = kFALSE;
8d0f3f77 509}
510
9a1398dd 511//____________________________________________________________________________
512Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
d15a28e7 513{
b2a60966 514 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
515 // = 1 are neighbour
516 // = 2 are not neighbour but do not continue searching
a4e98857 517 // neighbours are defined as digits having at least a common vertex
b2a60966 518 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
519 // which is compared to a digit (d2) not yet in a cluster
520
88cb7938 521 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
7b7c1533 522
d15a28e7 523 Int_t rv = 0 ;
524
d15a28e7 525 Int_t relid1[4] ;
7b7c1533 526 geom->AbsToRelNumbering(d1->GetId(), relid1) ;
d15a28e7 527
528 Int_t relid2[4] ;
7b7c1533 529 geom->AbsToRelNumbering(d2->GetId(), relid2) ;
d15a28e7 530
9688c1dd 531 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module
92862013 532 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
533 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
d15a28e7 534
92862013 535 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){
9688c1dd 536 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
d15a28e7 537 rv = 1 ;
538 }
539 else {
540 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
88cb7938 541 rv = 2; // Difference in row numbers is too large to look further
d15a28e7 542 }
543
544 }
545 else {
546
9688c1dd 547 if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )
d15a28e7 548 rv=2 ;
549
550 }
d72dfbc3 551
d15a28e7 552 return rv ;
553}
ba898748 554//____________________________________________________________________________
e95226ae 555void AliPHOSClusterizerv1::CleanDigits(TClonesArray * digits)
556{
557 // Remove digits with amplitudes below threshold
558
ba898748 559 for(Int_t i=0; i<digits->GetEntriesFast(); i++){
560 AliPHOSDigit * digit = static_cast<AliPHOSDigit*>(digits->At(i)) ;
e95226ae 561 if ( (IsInEmc(digit) && CalibrateEMC(digit->GetEnergy(),digit->GetId()) < fEmcMinE) ||
562 (IsInCpv(digit) && CalibrateCPV(digit->GetAmp() ,digit->GetId()) < fCpvMinE) )
ba898748 563 digits->RemoveAt(i) ;
564 }
565 digits->Compress() ;
566 for (Int_t i = 0 ; i < digits->GetEntriesFast() ; i++) {
567 AliPHOSDigit *digit = static_cast<AliPHOSDigit*>( digits->At(i) ) ;
568 digit->SetIndexInList(i) ;
569 }
d15a28e7 570
cbddd97f 571 //Overwrite digits tree
572 AliPHOSGetter* gime = AliPHOSGetter::Instance();
573 TTree * treeD = gime->TreeD();
574 treeD->Branch("PHOS", &digits);
575 treeD->Fill() ;
576 gime->WriteDigits("OVERWRITE");
577 gime->PhosLoader()->UnloadDigits() ;
ba898748 578}
d15a28e7 579//____________________________________________________________________________
9a1398dd 580Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
d15a28e7 581{
b2a60966 582 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
583
9f616d61 584 Bool_t rv = kFALSE ;
88cb7938 585 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
d15a28e7 586
2b629790 587 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
d15a28e7 588
2b629790 589 if(digit->GetId() <= nEMC ) rv = kTRUE;
ed4205d8 590
591 return rv ;
592}
593
ed4205d8 594//____________________________________________________________________________
9a1398dd 595Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
ed4205d8 596{
fad3e5b9 597 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
ed4205d8 598
599 Bool_t rv = kFALSE ;
88cb7938 600
601 AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
ed4205d8 602
2b629790 603 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
ed4205d8 604
2b629790 605 if(digit->GetId() > nEMC ) rv = kTRUE;
d15a28e7 606
607 return rv ;
608}
9a1398dd 609
610//____________________________________________________________________________
88cb7938 611void AliPHOSClusterizerv1::Unload()
612{
613 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
614 gime->PhosLoader()->UnloadDigits() ;
615 gime->PhosLoader()->UnloadRecPoints() ;
616}
617
618//____________________________________________________________________________
619void AliPHOSClusterizerv1::WriteRecPoints()
a4e98857 620{
7b7c1533 621
622 // Creates new branches with given title
a4e98857 623 // fills and writes into TreeR.
88cb7938 624
625 AliPHOSGetter * gime = AliPHOSGetter::Instance();
626
fbf811ec 627 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
628 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
629 TClonesArray * digits = gime->Digits() ;
ba898748 630
9a1398dd 631 Int_t index ;
ba898748 632 //Evaluate position, dispersion and other RecPoint properties..
633 Int_t nEmc = emcRecPoints->GetEntriesFast();
634 for(index = 0; index < nEmc; index++){
635 AliPHOSEmcRecPoint * rp = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) );
636 rp->Purify(fEmcMinE) ;
a978bcea 637 if(rp->GetMultiplicity()==0){
ba898748 638 emcRecPoints->RemoveAt(index) ;
639 delete rp ;
640 }
a978bcea 641
8c1fb709 642// No vertex is available now, calculate cirrections in PID
a978bcea 643 rp->EvalAll(fW0,digits) ;
8c1fb709 644 TVector3 fakeVtx(0.,0.,0.) ;
645 rp->EvalAll(fW0,fakeVtx,digits) ;
ba898748 646 }
647 emcRecPoints->Compress() ;
a978bcea 648// emcRecPoints->Sort() ; //Can not sort until position is calculated!
ba898748 649 // emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ;
650 for(index = 0; index < emcRecPoints->GetEntries(); index++){
88cb7938 651 dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) )->SetIndexInList(index) ;
ba898748 652 }
fbf811ec 653
af0b0570 654 //For each rec.point set the distance to the nearest bad crystal (BVP)
655 SetDistancesToBadChannels();
656
9a1398dd 657 //Now the same for CPV
ba898748 658 for(index = 0; index < cpvRecPoints->GetEntries(); index++){
659 AliPHOSCpvRecPoint * rp = dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(index) );
660 rp->EvalAll(fW0CPV,digits) ;
661 }
a978bcea 662// cpvRecPoints->Sort() ;
fbf811ec 663
7b7c1533 664 for(index = 0; index < cpvRecPoints->GetEntries(); index++)
2bb500e5 665 dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(index) )->SetIndexInList(index) ;
fbf811ec 666
7b7c1533 667 cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
9a1398dd 668
ba898748 669 if(fWrite){ //We write TreeR
670 TTree * treeR = gime->TreeR();
9a1398dd 671
ba898748 672 Int_t bufferSize = 32000 ;
673 Int_t splitlevel = 0 ;
9a1398dd 674
ba898748 675 //First EMC
676 TBranch * emcBranch = treeR->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
677 emcBranch->SetTitle(BranchName());
678
679 //Now CPV branch
680 TBranch * cpvBranch = treeR->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
681 cpvBranch->SetTitle(BranchName());
682
683 emcBranch ->Fill() ;
684 cpvBranch ->Fill() ;
685
686 gime->WriteRecPoints("OVERWRITE");
687 gime->WriteClusterizer("OVERWRITE");
688 }
9a1398dd 689}
690
691//____________________________________________________________________________
692void AliPHOSClusterizerv1::MakeClusters()
693{
694 // Steering method to construct the clusters stored in a list of Reconstructed Points
695 // A cluster is defined as a list of neighbour digits
d15a28e7 696
88cb7938 697
698 AliPHOSGetter * gime = AliPHOSGetter::Instance();
699
700 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
701 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
7b7c1533 702
7d493c2c 703 emcRecPoints->Delete() ;
704 cpvRecPoints->Delete() ;
7b7c1533 705
fbf811ec 706 TClonesArray * digits = gime->Digits() ;
88cb7938 707
ba898748 708 //Remove digits below threshold
709 CleanDigits(digits) ;
710
88cb7938 711 TClonesArray * digitsC = static_cast<TClonesArray*>( digits->Clone() ) ;
7b7c1533 712
d15a28e7 713 // Clusterization starts
9a1398dd 714
7b7c1533 715 TIter nextdigit(digitsC) ;
d15a28e7 716 AliPHOSDigit * digit ;
92862013 717 Bool_t notremoved = kTRUE ;
6ad0bfa0 718
88cb7938 719 while ( (digit = dynamic_cast<AliPHOSDigit *>( nextdigit()) ) ) { // scan over the list of digitsC
720
721
9a1398dd 722 AliPHOSRecPoint * clu = 0 ;
c161df70 723
d1de15f5 724 TArrayI clusterdigitslist(1500) ;
d15a28e7 725 Int_t index ;
f2bc1b87 726
e95226ae 727 if (( IsInEmc (digit) &&
728 CalibrateEMC(digit->GetEnergy(),digit->GetId()) > fEmcClusteringThreshold ) ||
729 ( IsInCpv (digit) &&
730 CalibrateCPV(digit->GetAmp() ,digit->GetId()) > fCpvClusteringThreshold ) ) {
d15a28e7 731 Int_t iDigitInCluster = 0 ;
7b7c1533 732
6ad0bfa0 733 if ( IsInEmc(digit) ) {
88cb7938 734 // start a new EMC RecPoint
735 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
736 emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
737
738 emcRecPoints->AddAt(new AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
739 clu = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(fNumberOfEmcClusters) ) ;
e95226ae 740 fNumberOfEmcClusters++ ;
741 clu->AddDigit(*digit, CalibrateEMC(digit->GetEnergy(),digit->GetId())) ;
742 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
743 iDigitInCluster++ ;
88cb7938 744 digitsC->Remove(digit) ;
ca77b032 745
d72dfbc3 746 } else {
88cb7938 747
748 // start a new CPV cluster
749 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
750 cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
751
752 cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
753
754 clu = dynamic_cast<AliPHOSCpvRecPoint *>( cpvRecPoints->At(fNumberOfCpvClusters) ) ;
755 fNumberOfCpvClusters++ ;
e95226ae 756 clu->AddDigit(*digit, CalibrateCPV(digit->GetAmp(),digit->GetId()) ) ;
88cb7938 757 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
758 iDigitInCluster++ ;
759 digitsC->Remove(digit) ;
d15a28e7 760 nextdigit.Reset() ;
88cb7938 761
762 // Here we remove remaining EMC digits, which cannot make a cluster
763
92862013 764 if( notremoved ) {
88cb7938 765 while( ( digit = dynamic_cast<AliPHOSDigit *>( nextdigit() ) ) ) {
9f616d61 766 if( IsInEmc(digit) )
88cb7938 767 digitsC->Remove(digit) ;
d15a28e7 768 else
88cb7938 769 break ;
770 }
771 notremoved = kFALSE ;
772 }
773
d15a28e7 774 } // else
775
776 nextdigit.Reset() ;
fad3e5b9 777
d15a28e7 778 AliPHOSDigit * digitN ;
779 index = 0 ;
9f616d61 780 while (index < iDigitInCluster){ // scan over digits already in cluster
88cb7938 781 digit = dynamic_cast<AliPHOSDigit*>( digits->At(clusterdigitslist[index]) ) ;
782 index++ ;
783 while ( (digitN = dynamic_cast<AliPHOSDigit *>( nextdigit() ) ) ) { // scan over the reduced list of digits
784 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
d15a28e7 785 switch (ineb ) {
9f616d61 786 case 0 : // not a neighbour
88cb7938 787 break ;
788 case 1 : // are neighbours
e95226ae 789 if (IsInEmc (digitN))
790 clu->AddDigit(*digitN, CalibrateEMC( digitN->GetEnergy(), digitN->GetId() ) );
791 else
792 clu->AddDigit(*digitN, CalibrateCPV( digitN->GetAmp() , digitN->GetId() ) );
793
88cb7938 794 clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ;
795 iDigitInCluster++ ;
796 digitsC->Remove(digitN) ;
797 break ;
9f616d61 798 case 2 : // too far from each other
88cb7938 799 goto endofloop;
800 } // switch
801
802 } // while digitN
803
d15a28e7 804 endofloop: ;
88cb7938 805 nextdigit.Reset() ;
806
d15a28e7 807 } // loop over cluster
ad8cfaf4 808
ad8cfaf4 809 } // energy theshold
ad8cfaf4 810
31aa6d6c 811
d15a28e7 812 } // while digit
813
7b7c1533 814 delete digitsC ;
9688c1dd 815
9a1398dd 816}
817
818//____________________________________________________________________________
a4e98857 819void AliPHOSClusterizerv1::MakeUnfolding()
820{
821 // Unfolds clusters using the shape of an ElectroMagnetic shower
9688c1dd 822 // Performs unfolding of all EMC/CPV clusters
9a1398dd 823
88cb7938 824 AliPHOSGetter * gime = AliPHOSGetter::Instance();
fc12304f 825
7b7c1533 826 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
88cb7938 827
828 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
829 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
fbf811ec 830 TClonesArray * digits = gime->Digits() ;
7b7c1533 831
a4e98857 832 // Unfold first EMC clusters
9a1398dd 833 if(fNumberOfEmcClusters > 0){
834
7b7c1533 835 Int_t nModulesToUnfold = geom->GetNModules() ;
9a1398dd 836
837 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
838 Int_t index ;
839 for(index = 0 ; index < numberofNotUnfolded ; index++){
840
88cb7938 841 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(index) ) ;
9a1398dd 842 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
88cb7938 843 break ;
9a1398dd 844
845 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
a0636361 846 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
9a1398dd 847 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
7b7c1533 848 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
9a1398dd 849
850 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
88cb7938 851 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
852 emcRecPoints->Remove(emcRecPoint);
853 emcRecPoints->Compress() ;
854 index-- ;
855 fNumberOfEmcClusters -- ;
856 numberofNotUnfolded-- ;
9a1398dd 857 }
c6bf27f2 858 else{
859 emcRecPoint->SetNExMax(1) ; //Only one local maximum
860 }
9a1398dd 861
862 delete[] maxAt ;
863 delete[] maxAtEnergy ;
864 }
865 }
a4e98857 866 // Unfolding of EMC clusters finished
9a1398dd 867
868
a4e98857 869 // Unfold now CPV clusters
9a1398dd 870 if(fNumberOfCpvClusters > 0){
871
9688c1dd 872 Int_t nModulesToUnfold = geom->GetNModules() ;
9a1398dd 873
874 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
875 Int_t index ;
876 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
877
88cb7938 878 AliPHOSRecPoint * recPoint = dynamic_cast<AliPHOSRecPoint *>( cpvRecPoints->At(index) ) ;
9a1398dd 879
880 if(recPoint->GetPHOSMod()> nModulesToUnfold)
88cb7938 881 break ;
9a1398dd 882
88cb7938 883 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint*>(recPoint) ;
9a1398dd 884
885 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
a0636361 886 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
9a1398dd 887 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
7b7c1533 888 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,digits) ;
9a1398dd 889
890 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
88cb7938 891 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
892 cpvRecPoints->Remove(emcRecPoint);
893 cpvRecPoints->Compress() ;
894 index-- ;
895 numberofCpvNotUnfolded-- ;
896 fNumberOfCpvClusters-- ;
9a1398dd 897 }
898
899 delete[] maxAt ;
900 delete[] maxAtEnergy ;
901 }
902 }
903 //Unfolding of Cpv clusters finished
904
905}
906
9a1398dd 907//____________________________________________________________________________
a978bcea 908Double_t AliPHOSClusterizerv1::ShowerShape(Double_t x, Double_t z)
9a1398dd 909{
910 // Shape of the shower (see PHOS TDR)
a4e98857 911 // If you change this function, change also the gradient evaluation in ChiSquare()
9a1398dd 912
a978bcea 913 //for the moment we neglect dependence on the incident angle.
914
915 Double_t r2 = x*x + z*z ;
916 Double_t r4 = r2*r2 ;
917 Double_t r295 = TMath::Power(r2, 2.95/2.) ;
9a1398dd 918 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
919 return shape ;
920}
921
922//____________________________________________________________________________
923void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
88cb7938 924 Int_t nMax,
925 AliPHOSDigit ** maxAt,
926 Float_t * maxAtEnergy)
9a1398dd 927{
928 // Performs the unfolding of a cluster with nMax overlapping showers
929
88cb7938 930 AliPHOSGetter * gime = AliPHOSGetter::Instance();
fc12304f 931
7b7c1533 932 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
88cb7938 933
fbf811ec 934 const TClonesArray * digits = gime->Digits() ;
935 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
936 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
7b7c1533 937
9a1398dd 938 Int_t nPar = 3 * nMax ;
939 Float_t * fitparameters = new Float_t[nPar] ;
940
941 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
942 if( !rv ) {
88cb7938 943 // Fit failed, return and remove cluster
c6bf27f2 944 iniEmc->SetNExMax(-1) ;
9a1398dd 945 delete[] fitparameters ;
946 return ;
947 }
948
949 // create ufolded rec points and fill them with new energy lists
950 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
951 // and later correct this number in acordance with actual energy deposition
952
953 Int_t nDigits = iniEmc->GetMultiplicity() ;
954 Float_t * efit = new Float_t[nDigits] ;
a978bcea 955 Float_t xDigit=0.,zDigit=0. ;
7b7c1533 956 Float_t xpar=0.,zpar=0.,epar=0. ;
9a1398dd 957 Int_t relid[4] ;
7b7c1533 958 AliPHOSDigit * digit = 0 ;
9a1398dd 959 Int_t * emcDigits = iniEmc->GetDigitsList() ;
960
a978bcea 961 TVector3 vIncid ;
962
9a1398dd 963 Int_t iparam ;
964 Int_t iDigit ;
965 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
88cb7938 966 digit = dynamic_cast<AliPHOSDigit*>( digits->At(emcDigits[iDigit] ) ) ;
7b7c1533 967 geom->AbsToRelNumbering(digit->GetId(), relid) ;
968 geom->RelPosInModule(relid, xDigit, zDigit) ;
9a1398dd 969 efit[iDigit] = 0;
970
971 iparam = 0 ;
972 while(iparam < nPar ){
973 xpar = fitparameters[iparam] ;
974 zpar = fitparameters[iparam+1] ;
975 epar = fitparameters[iparam+2] ;
976 iparam += 3 ;
a978bcea 977// geom->GetIncidentVector(fVtx,relid[0],xpar,zpar,vIncid) ;
978// efit[iDigit] += epar * ShowerShape(xDigit - xpar,zDigit - zpar,vIncid) ;
979 efit[iDigit] += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
9a1398dd 980 }
981 }
982
983
984 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
985 // so that energy deposited in each cell is distributed betwin new clusters proportionally
986 // to its contribution to efit
987
988 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
989 Float_t ratio ;
990
991 iparam = 0 ;
992 while(iparam < nPar ){
993 xpar = fitparameters[iparam] ;
994 zpar = fitparameters[iparam+1] ;
995 epar = fitparameters[iparam+2] ;
996 iparam += 3 ;
a978bcea 997// geom->GetIncidentVector(fVtx,iniEmc->GetPHOSMod(),xpar,zpar,vIncid) ;
9a1398dd 998
7b7c1533 999 AliPHOSEmcRecPoint * emcRP = 0 ;
9a1398dd 1000
1001 if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
1002
7b7c1533 1003 if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
88cb7938 1004 emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
9a1398dd 1005
73a68ccb 1006 (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
88cb7938 1007 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( emcRecPoints->At(fNumberOfEmcClusters) ) ;
9a1398dd 1008 fNumberOfEmcClusters++ ;
c6bf27f2 1009 emcRP->SetNExMax((Int_t)nPar/3) ;
9a1398dd 1010 }
1011 else{//create new entries in fCpvRecPoints
7b7c1533 1012 if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
88cb7938 1013 cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
9a1398dd 1014
73a68ccb 1015 (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
88cb7938 1016 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( cpvRecPoints->At(fNumberOfCpvClusters) ) ;
9a1398dd 1017 fNumberOfCpvClusters++ ;
1018 }
1019
1020 Float_t eDigit ;
1021 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
88cb7938 1022 digit = dynamic_cast<AliPHOSDigit*>( digits->At( emcDigits[iDigit] ) ) ;
7b7c1533 1023 geom->AbsToRelNumbering(digit->GetId(), relid) ;
1024 geom->RelPosInModule(relid, xDigit, zDigit) ;
a978bcea 1025// ratio = epar * ShowerShape(xDigit - xpar,zDigit - zpar,vIncid) / efit[iDigit] ;
1026 ratio = epar * ShowerShape(xDigit - xpar,zDigit - zpar) / efit[iDigit] ;
9a1398dd 1027 eDigit = emcEnergies[iDigit] * ratio ;
1028 emcRP->AddDigit( *digit, eDigit ) ;
88cb7938 1029 }
9a1398dd 1030 }
1031
1032 delete[] fitparameters ;
1033 delete[] efit ;
1034
1035}
1036
1037//_____________________________________________________________________________
1038void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
1039{
a4e98857 1040 // Calculates the Chi square for the cluster unfolding minimization
9a1398dd 1041 // Number of parameters, Gradient, Chi squared, parameters, what to do
1042
88cb7938 1043 TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;
9a1398dd 1044
88cb7938 1045 AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint*>( toMinuit->At(0) ) ;
1046 TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) ) ;
a978bcea 1047// TVector3 * vtx = dynamic_cast<TVector3*>(toMinuit->At(2)) ; //Vertex position
9a1398dd 1048
88cb7938 1049 // AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( gMinuit->GetObjectFit() ) ; // EmcRecPoint to fit
9a1398dd 1050
1051 Int_t * emcDigits = emcRP->GetDigitsList() ;
1052
7b7c1533 1053 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
9a1398dd 1054
1055 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
88cb7938 1056
1057 const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
a978bcea 1058// TVector3 vInc ;
9a1398dd 1059 fret = 0. ;
1060 Int_t iparam ;
1061
1062 if(iflag == 2)
1063 for(iparam = 0 ; iparam < nPar ; iparam++)
1064 Grad[iparam] = 0 ; // Will evaluate gradient
1065
1066 Double_t efit ;
1067
1068 AliPHOSDigit * digit ;
1069 Int_t iDigit ;
1070
7b7c1533 1071 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
9a1398dd 1072
88cb7938 1073 digit = dynamic_cast<AliPHOSDigit*>( digits->At( emcDigits[iDigit] ) );
9a1398dd 1074
1075 Int_t relid[4] ;
1076 Float_t xDigit ;
1077 Float_t zDigit ;
1078
1079 geom->AbsToRelNumbering(digit->GetId(), relid) ;
1080
1081 geom->RelPosInModule(relid, xDigit, zDigit) ;
1082
1083 if(iflag == 2){ // calculate gradient
1084 Int_t iParam = 0 ;
1085 efit = 0 ;
1086 while(iParam < nPar ){
a978bcea 1087 Double_t dx = (xDigit - x[iParam]) ;
88cb7938 1088 iParam++ ;
a978bcea 1089 Double_t dz = (zDigit - x[iParam]) ;
88cb7938 1090 iParam++ ;
a978bcea 1091// geom->GetIncidentVector(*vtx,emcRP->GetPHOSMod(),x[iParam-2],x[iParam-1],vInc) ;
1092// efit += x[iParam] * ShowerShape(dx,dz,vInc) ;
1093 efit += x[iParam] * ShowerShape(dx,dz) ;
88cb7938 1094 iParam++ ;
9a1398dd 1095 }
88cb7938 1096 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
9a1398dd 1097 iParam = 0 ;
1098 while(iParam < nPar ){
88cb7938 1099 Double_t xpar = x[iParam] ;
1100 Double_t zpar = x[iParam+1] ;
1101 Double_t epar = x[iParam+2] ;
a978bcea 1102// geom->GetIncidentVector(*vtx,emcRP->GetPHOSMod(),xpar,zpar,vInc) ;
88cb7938 1103 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
a978bcea 1104// Double_t shape = sum * ShowerShape(xDigit - xpar,zDigit - zpar,vInc) ;
1105 Double_t shape = sum * ShowerShape(xDigit - xpar,zDigit - zpar) ;
1106//DP: No incident angle dependence in gradient yet!!!!!!
88cb7938 1107 Double_t r4 = dr*dr*dr*dr ;
1108 Double_t r295 = TMath::Power(dr,2.95) ;
1109 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
1110 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
1111
1112 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
1113 iParam++ ;
1114 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
1115 iParam++ ;
1116 Grad[iParam] += shape ; // Derivative over energy
1117 iParam++ ;
9a1398dd 1118 }
1119 }
1120 efit = 0;
1121 iparam = 0 ;
1122
1123 while(iparam < nPar ){
1124 Double_t xpar = x[iparam] ;
1125 Double_t zpar = x[iparam+1] ;
1126 Double_t epar = x[iparam+2] ;
1127 iparam += 3 ;
a978bcea 1128// geom->GetIncidentVector(*vtx,emcRP->GetPHOSMod(),xpar,zpar,vInc) ;
1129// efit += epar * ShowerShape(xDigit - xpar,zDigit - zpar,vInc) ;
1130 efit += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
9a1398dd 1131 }
1132
88cb7938 1133 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
9a1398dd 1134 // Here we assume, that sigma = sqrt(E)
1135 }
83974468 1136
d15a28e7 1137}
1138
1139//____________________________________________________________________________
702ab87e 1140void AliPHOSClusterizerv1::Print(const Option_t *)const
d15a28e7 1141{
d1de15f5 1142 // Print clusterizer parameters
1143
21cd0c07 1144 TString message ;
1145 TString taskName(GetName()) ;
1146 taskName.ReplaceAll(Version(), "") ;
1147
1148 if( strcmp(GetName(), "") !=0 ) {
9a1398dd 1149 // Print parameters
21cd0c07 1150 message = "\n--------------- %s %s -----------\n" ;
1151 message += "Clusterizing digits from the file: %s\n" ;
1152 message += " Branch: %s\n" ;
1153 message += " EMC Clustering threshold = %f\n" ;
1154 message += " EMC Local Maximum cut = %f\n" ;
1155 message += " EMC Logarothmic weight = %f\n" ;
1156 message += " CPV Clustering threshold = %f\n" ;
1157 message += " CPV Local Maximum cut = %f\n" ;
1158 message += " CPV Logarothmic weight = %f\n" ;
9a1398dd 1159 if(fToUnfold)
21cd0c07 1160 message += " Unfolding on\n" ;
9a1398dd 1161 else
21cd0c07 1162 message += " Unfolding off\n" ;
9a1398dd 1163
21cd0c07 1164 message += "------------------------------------------------------------------" ;
9a1398dd 1165 }
21cd0c07 1166 else
1167 message = " AliPHOSClusterizerv1 not initialized " ;
1168
351dd634 1169 AliInfo(Form("%s, %s %s %s %s %s %s %s %s %s %s", message.Data(),
21cd0c07 1170 taskName.Data(),
1171 GetTitle(),
1172 taskName.Data(),
1173 GetName(),
1174 fEmcClusteringThreshold,
1175 fEmcLocMaxCut,
1176 fW0,
1177 fCpvClusteringThreshold,
1178 fCpvLocMaxCut,
351dd634 1179 fW0CPV )) ;
9a1398dd 1180}
a978bcea 1181//____________________________________________________________________________
1182//void AliPHOSClusterizerv1::GetVertex(void)
1183//{ //Extracts vertex posisition
1184//
1185 //ESD
1186//DP - todo if(){
1187//
1188// }
1189
1190// //MC Generator
1191// if(gAlice && gAlice->GetMCApp() && gAlice->Generator()){
1192// Float_t x,y,z ;
1193// gAlice->Generator()->GetOrigin(x,y,z) ;
1194// fVtx.SetXYZ(x,y,z) ;
1195// return ;
1196// }
1197//
1198// //No any source
1199// fVtx[0]=fVtx[1]=fVtx[2]=0. ;
1200//
1201//}
9a1398dd 1202//____________________________________________________________________________
a4e98857 1203void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1204{
1205 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
9a1398dd 1206
88cb7938 1207 AliPHOSGetter * gime = AliPHOSGetter::Instance();
1208
1209 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
1210 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
7b7c1533 1211
351dd634 1212 AliInfo(Form("\nevent %d \n Found %d EMC RecPoints and %d CPV RecPoints",
1213 gAlice->GetEvNumber(),
1214 emcRecPoints->GetEntriesFast(),
1215 cpvRecPoints->GetEntriesFast() )) ;
88cb7938 1216
1217 fRecPointsInRun += emcRecPoints->GetEntriesFast() ;
1218 fRecPointsInRun += cpvRecPoints->GetEntriesFast() ;
1219
11f9c5ff 1220
9a1398dd 1221 if(strstr(option,"all")) {
709e117a 1222 printf("\n EMC clusters \n") ;
1223 printf("Index Ene(MeV) Multi Module X Y Z Lambdas_1 Lambda_2 # of prim Primaries list\n") ;
9a1398dd 1224 Int_t index ;
7b7c1533 1225 for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1226 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ;
9a1398dd 1227 TVector3 locpos;
1228 rp->GetLocalPosition(locpos);
9a1398dd 1229 Float_t lambda[2];
1230 rp->GetElipsAxis(lambda);
9a1398dd 1231 Int_t * primaries;
1232 Int_t nprimaries;
1233 primaries = rp->GetPrimaries(nprimaries);
709e117a 1234 printf("\n%6d %8.2f %3d %2d %4.1f %4.1f %4.1f %4f %4f %2d : ",
11f9c5ff 1235 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(), rp->GetPHOSMod(),
1236 locpos.X(), locpos.Y(), locpos.Z(), lambda[0], lambda[1], nprimaries) ;
3bf72d32 1237
21cd0c07 1238 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
709e117a 1239 printf("%d ", primaries[iprimary] ) ;
21cd0c07 1240 }
709e117a 1241 printf("\n") ;
1242 }
1243
1244 //Now plot CPV recPoints
1245 printf("\n CPV clusters \n") ;
1246 printf("Index Ene(MeV) Module X Y Z \n") ;
1247 for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
2bb500e5 1248 AliPHOSCpvRecPoint * rp = (AliPHOSCpvRecPoint * )cpvRecPoints->At(index) ;
1249
1250 TVector3 locpos;
1251 rp->GetLocalPosition(locpos);
1252
1253 printf("\n%6d %8.2f %2d %4.1f %4.1f %4.1f \n",
1254 rp->GetIndexInList(), rp->GetEnergy(), rp->GetPHOSMod(),
1255 locpos.X(), locpos.Y(), locpos.Z()) ;
88cb7938 1256 }
9a1398dd 1257 }
d15a28e7 1258}
9a1398dd 1259
af0b0570 1260
1261//____________________________________________________________________________
1262void AliPHOSClusterizerv1::SetDistancesToBadChannels()
1263{
1264 //For each EMC rec. point set the distance to the nearest bad crystal.
1265 //Author: Boris Polichtchouk
1266
1267 if(!fCalibData->GetNumOfEmcBadChannels()) return;
1268 AliInfo(Form("%d bad channel(s) found.\n",fCalibData->GetNumOfEmcBadChannels()));
1269
1270 AliPHOSGetter* gime = AliPHOSGetter::Instance();
1271 AliPHOSGeometry* geom = gime->PHOSGeometry();
1272
1273 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
1274
1275 Int_t badIds[8000];
1276 fCalibData->EmcBadChannelIds(badIds);
1277
1278 AliPHOSEmcRecPoint* rp;
1279
1280 TMatrixF gmat;
1281 TVector3 gposRecPoint; // global (in ALICE frame) position of rec. point
1282 TVector3 gposBadChannel; // global position of bad crystal
1283 TVector3 dR;
1284
1285 Float_t dist,minDist;
1286
1287 for(Int_t iRP=0; iRP<emcRecPoints->GetEntries(); iRP++){
1288 rp = (AliPHOSEmcRecPoint*)emcRecPoints->At(iRP);
1289 minDist = 1.e+07;
1290
1291 for(Int_t iBad=0; iBad<fCalibData->GetNumOfEmcBadChannels(); iBad++) {
1292 rp->GetGlobalPosition(gposRecPoint,gmat);
1293 geom->RelPosInAlice(badIds[iBad],gposBadChannel);
1294 AliDebug(2,Form("BC position:[%.3f,%.3f,%.3f], RP position:[%.3f,%.3f,%.3f]. E=%.3f\n",
1295 gposBadChannel.X(),gposBadChannel.Y(),gposBadChannel.Z(),
1296 gposRecPoint.X(),gposRecPoint.Y(),gposRecPoint.Z(),rp->GetEnergy()));
1297 dR = gposBadChannel-gposRecPoint;
1298 dist = dR.Mag();
1299 if(dist<minDist) minDist = dist;
1300 }
1301
1302 rp->SetDistanceToBadCrystal(minDist);
1303 }
1304
1305}