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