]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSClusterizerv1.cxx
pedestal values subtracted
[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 *
a28333c5 20 * $Log: AliPHOSClusterizerv1.cxx,v $
21 * Revision 1.118 2007/12/11 21:23:26 kharlov
22 * Added possibility to swith off unfolding
23 *
a42611a9 24 * Revision 1.117 2007/10/18 08:42:05 kharlov
25 * Bad channels cleaned before clusterization
26 *
454c62a4 27 * Revision 1.116 2007/10/01 20:24:08 kharlov
28 * Memory leaks fixed
29 *
e68222ce 30 * Revision 1.115 2007/09/26 14:22:17 cvetan
31 * Important changes to the reconstructor classes. Complete elimination of the run-loaders, which are now steered only from AliReconstruction. Removal of the corresponding Reconstruct() and FillESD() methods.
32 *
d76c31f4 33 * Revision 1.114 2007/09/06 16:06:44 kharlov
34 * Absence of sorting results in loose of all unfolded clusters
35 *
6d0a5c2e 36 * Revision 1.113 2007/08/28 12:55:07 policheh
37 * Loaders removed from the reconstruction code (C.Cheshkov)
38 *
9a2cdbdf 39 * Revision 1.112 2007/08/22 09:20:50 hristov
40 * Updated QA classes (Yves)
41 *
a5fa6165 42 * Revision 1.111 2007/08/08 12:11:28 kharlov
43 * Protection against uninitialized fQADM
44 *
7ae07b77 45 * Revision 1.110 2007/08/07 14:16:00 kharlov
46 * Quality assurance added (Yves Schutz)
47 *
ddd1a39c 48 * Revision 1.109 2007/07/24 17:20:35 policheh
49 * Usage of RecoParam objects instead of hardcoded parameters in reconstruction.
50 * (See $ALICE_ROOT/PHOS/macros/BeamTest2006/RawReconstruction.C).
51 *
3799bcb5 52 * Revision 1.108 2007/06/18 07:00:51 kharlov
53 * Bug fix for attempt to use AliPHOSEmcRecPoint after its deletion
54 *
cd79ec76 55 * Revision 1.107 2007/05/25 14:12:26 policheh
56 * Local to tracking CS transformation added for CPV rec. points
57 *
76e42795 58 * Revision 1.106 2007/05/24 13:01:22 policheh
59 * Local to tracking CS transformation invoked for each EMC rec.point
60 *
788ab371 61 * Revision 1.105 2007/05/02 13:41:22 kharlov
62 * Mode protection against absence of calib.data from AliPHOSCalibData to AliPHOSClusterizerv1::GetCalibrationParameters()
63 *
f36dfa9a 64 * Revision 1.104 2007/04/27 16:55:53 kharlov
65 * Calibration stops if PHOS CDB objects do not exist
66 *
76d78b7a 67 * Revision 1.103 2007/04/11 11:55:45 policheh
68 * SetDistancesToBadChannels() added.
69 *
af0b0570 70 * Revision 1.102 2007/03/28 19:18:15 kharlov
71 * RecPoints recalculation in TSM removed
72 *
8c1fb709 73 * Revision 1.101 2007/03/06 06:51:27 kharlov
74 * Calculation of cluster properties dep. on vertex posponed to TrackSegmentMaker
75 *
a978bcea 76 * Revision 1.100 2007/01/10 11:07:26 kharlov
77 * Raw digits writing to file (B.Polichtchouk)
78 *
cbddd97f 79 * Revision 1.99 2006/11/07 16:49:51 kharlov
80 * Corrections for next event switch in case of raw data (B.Polichtchouk)
81 *
e1c29d81 82 * Revision 1.98 2006/10/27 17:14:27 kharlov
83 * Introduce AliDebug and AliLog (B.Polichtchouk)
84 *
e1afea68 85 * Revision 1.97 2006/08/29 11:41:19 kharlov
86 * Missing implementation of ctors and = operator are added
87 *
f1133801 88 * Revision 1.96 2006/08/25 16:56:30 kharlov
89 * Compliance with Effective C++
90 *
0378398c 91 * Revision 1.95 2006/08/11 12:36:26 cvetan
92 * Update of the PHOS code needed in order to read and reconstruct the beam test raw data (i.e. without an existing galice.root)
93 *
52c5f046 94 * Revision 1.94 2006/08/07 12:27:49 hristov
95 * Removing obsolete code which affected the event numbering scheme
96 *
4d4750bd 97 * Revision 1.93 2006/08/01 12:20:17 cvetan
98 * 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
99 *
f5eaa851 100 * Revision 1.92 2006/04/29 20:26:46 hristov
101 * Separate EMC and CPV calibration (Yu.Kharlov)
102 *
e95226ae 103 * Revision 1.91 2006/04/22 10:30:17 hristov
104 * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
105 *
27a73a5d 106 * Revision 1.90 2006/04/11 15:22:59 hristov
107 * run number in query set to -1: forces AliCDBManager to use its run number (A.Colla)
108 *
28871337 109 * Revision 1.89 2006/03/13 14:05:42 kharlov
110 * Calibration objects for EMC and CPV
111 *
fc6706cb 112 * Revision 1.88 2006/01/11 08:54:52 hristov
113 * Additional protection in case no calibration entry was found
114 *
e521c5a0 115 * Revision 1.87 2005/11/22 08:46:43 kharlov
116 * Updated to new CDB (Boris Polichtchouk)
117 *
ef629623 118 * Revision 1.86 2005/11/14 21:52:43 hristov
119 * Coding conventions
120 *
fdf65bb5 121 * Revision 1.85 2005/09/27 16:08:08 hristov
122 * New version of CDB storage framework (A.Colla)
123 *
9e1ceb13 124 * Revision 1.84 2005/09/21 10:02:47 kharlov
125 * Reading calibration from CDB (Boris Polichtchouk)
126 *
8d6d4fb5 127 * Revision 1.82 2005/09/02 15:43:13 kharlov
128 * Add comments in GetCalibrationParameters and Calibrate
129 *
5c149aeb 130 * Revision 1.81 2005/09/02 14:32:07 kharlov
131 * Calibration of raw data
132 *
44ae287e 133 * Revision 1.80 2005/08/24 15:31:36 kharlov
134 * Setting raw digits flag
135 *
9735c615 136 * Revision 1.79 2005/07/25 15:53:53 kharlov
137 * Read raw data
138 *
3cf4f75f 139 * Revision 1.78 2005/05/28 14:19:04 schutz
140 * Compilation warnings fixed by T.P.
141 *
702ab87e 142 */
143
9a1398dd 144//*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
d15a28e7 145//////////////////////////////////////////////////////////////////////////////
9a1398dd 146// Clusterization class. Performs clusterization (collects neighbouring active cells) and
a4e98857 147// unfolds the clusters having several local maxima.
148// Results are stored in TreeR#, branches PHOSEmcRP (EMC recPoints),
9a1398dd 149// PHOSCpvRP (CPV RecPoints) and AliPHOSClusterizer (Clusterizer with all
f035f6ce 150// parameters including input digits branch title, thresholds etc.)
a4e98857 151// This TTask is normally called from Reconstructioner, but can as well be used in
152// standalone mode.
153// Use Case:
9a2cdbdf 154// root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1(<pointer_to_phos_geometry_onject>)
155// root [1] cl->Digits2Clusters(digitsTree,clusterTree)
156// //finds RecPoints in the current event
a4e98857 157// root [2] cl->SetDigitsBranch("digits2")
f035f6ce 158// //sets another title for Digitis (input) branch
a4e98857 159// root [3] cl->SetRecPointsBranch("recp2")
f035f6ce 160// //sets another title four output branches
a4e98857 161// root [4] cl->SetEmcLocalMaxCut(0.03)
9a1398dd 162// //set clusterization parameters
ed4205d8 163
d15a28e7 164// --- ROOT system ---
165
166#include "TMath.h"
9a1398dd 167#include "TMinuit.h"
168#include "TTree.h"
9a1398dd 169#include "TBenchmark.h"
9a2cdbdf 170#include "TClonesArray.h"
d15a28e7 171
172// --- Standard library ---
173
d15a28e7 174// --- AliRoot header files ---
351dd634 175#include "AliLog.h"
9a2cdbdf 176#include "AliConfig.h"
e957fea8 177#include "AliPHOSGeometry.h"
d15a28e7 178#include "AliPHOSClusterizerv1.h"
88cb7938 179#include "AliPHOSEmcRecPoint.h"
9a1398dd 180#include "AliPHOSCpvRecPoint.h"
d15a28e7 181#include "AliPHOSDigit.h"
9a1398dd 182#include "AliPHOSDigitizer.h"
ba898748 183#include "AliPHOSCalibrationDB.h"
9e1ceb13 184#include "AliCDBManager.h"
eaa5012b 185#include "AliCDBStorage.h"
fdf65bb5 186#include "AliCDBEntry.h"
3799bcb5 187#include "AliPHOSRecoParam.h"
9a2cdbdf 188#include "AliPHOSReconstructor.h"
a28333c5 189#include "AliPHOSCalibData.h"
d15a28e7 190
191ClassImp(AliPHOSClusterizerv1)
f035f6ce 192
d15a28e7 193//____________________________________________________________________________
0378398c 194AliPHOSClusterizerv1::AliPHOSClusterizerv1() :
195 AliPHOSClusterizer(),
196 fDefaultInit(0), fEmcCrystals(0), fToUnfold(0),
d07a44a8 197 fWrite(0),
198 fNumberOfEmcClusters(0), fNumberOfCpvClusters(0),
a28333c5 199 fEmcClusteringThreshold(0), fCpvClusteringThreshold(0),
0378398c 200 fEmcLocMaxCut(0), fW0(0), fCpvLocMaxCut(0),
d55c9739 201 fW0CPV(0), fEmcTimeGate(0), fEcoreRadius(0)
d15a28e7 202{
f035f6ce 203 // default ctor (to be used mainly by Streamer)
f035f6ce 204
8d0f3f77 205 InitParameters() ;
92f521a9 206 fDefaultInit = kTRUE ;
9a1398dd 207}
7b7c1533 208
9a1398dd 209//____________________________________________________________________________
9a2cdbdf 210AliPHOSClusterizerv1::AliPHOSClusterizerv1(AliPHOSGeometry *geom) :
211 AliPHOSClusterizer(geom),
0378398c 212 fDefaultInit(0), fEmcCrystals(0), fToUnfold(0),
d07a44a8 213 fWrite(0),
214 fNumberOfEmcClusters(0), fNumberOfCpvClusters(0),
a28333c5 215 fEmcClusteringThreshold(0), fCpvClusteringThreshold(0),
0378398c 216 fEmcLocMaxCut(0), fW0(0), fCpvLocMaxCut(0),
d55c9739 217 fW0CPV(0), fEmcTimeGate(0), fEcoreRadius(0)
9a1398dd 218{
a4e98857 219 // ctor with the indication of the file where header Tree and digits Tree are stored
9a1398dd 220
8d0f3f77 221 InitParameters() ;
2bd5457f 222 Init() ;
92f521a9 223 fDefaultInit = kFALSE ;
9a1398dd 224}
fc12304f 225
9688c1dd 226//____________________________________________________________________________
227 AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
228{
0bc3b8ed 229 // dtor
230
9688c1dd 231}
3758d9fc 232//____________________________________________________________________________
9a2cdbdf 233void AliPHOSClusterizerv1::Digits2Clusters(Option_t *option)
a4e98857 234{
9a2cdbdf 235 // Steering method to perform clusterization for one event
236 // The input is the tree with digits.
237 // The output is the tree with clusters.
9a1398dd 238
239 if(strstr(option,"tim"))
7d493c2c 240 gBenchmark->Start("PHOSClusterizer");
9a1398dd 241
eabde521 242 if(strstr(option,"print")) {
88cb7938 243 Print() ;
eabde521 244 return ;
245 }
bed9e3fb 246
9a2cdbdf 247 MakeClusters() ;
9a1398dd 248
9a2cdbdf 249 AliDebug(2,Form(" ---- Printing clusters (%d)\n",
250 fEMCRecPoints->GetEntries()));
251 if(AliLog::GetGlobalDebugLevel()>1)
252 fEMCRecPoints->Print();
88cb7938 253
9a2cdbdf 254 if(fToUnfold)
255 MakeUnfolding();
ddd1a39c 256
9a2cdbdf 257 WriteRecPoints();
7b7c1533 258
9a2cdbdf 259 if(strstr(option,"deb"))
260 PrintRecPoints(option) ;
94de8339 261
9a1398dd 262 if(strstr(option,"tim")){
263 gBenchmark->Stop("PHOSClusterizer");
9a2cdbdf 264 AliInfo(Form("took %f seconds for Clusterizing\n",
265 gBenchmark->GetCpuTime("PHOSClusterizer")));
8d8258f6 266 }
6483babc 267 fEMCRecPoints->Delete();
268 fCPVRecPoints->Delete();
9a1398dd 269}
270
271//____________________________________________________________________________
a0636361 272Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
d1de15f5 273 Int_t nPar, Float_t * fitparameters) const
9a1398dd 274{
275 // Calls TMinuit to fit the energy distribution of a cluster with several maxima
a4e98857 276 // The initial values for fitting procedure are set equal to the positions of local maxima.
f035f6ce 277 // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
9a1398dd 278
88cb7938 279
9a1398dd 280 gMinuit->mncler(); // Reset Minuit's list of paramters
281 gMinuit->SetPrintLevel(-1) ; // No Printout
282 gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;
283 // To set the address of the minimization function
284
285 TList * toMinuit = new TList();
286 toMinuit->AddAt(emcRP,0) ;
9a2cdbdf 287 toMinuit->AddAt(fDigitsArr,1) ;
288 toMinuit->AddAt(fGeom,2) ;
9a1398dd 289
290 gMinuit->SetObjectFit(toMinuit) ; // To tranfer pointer to UnfoldingChiSquare
291
292 // filling initial values for fit parameters
293 AliPHOSDigit * digit ;
294
295 Int_t ierflg = 0;
296 Int_t index = 0 ;
297 Int_t nDigits = (Int_t) nPar / 3 ;
298
299 Int_t iDigit ;
300
9a1398dd 301 for(iDigit = 0; iDigit < nDigits; iDigit++){
a0636361 302 digit = maxAt[iDigit];
9a1398dd 303
304 Int_t relid[4] ;
7b7c1533 305 Float_t x = 0.;
306 Float_t z = 0.;
9a2cdbdf 307 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
308 fGeom->RelPosInModule(relid, x, z) ;
9a1398dd 309
310 Float_t energy = maxAtEnergy[iDigit] ;
311
312 gMinuit->mnparm(index, "x", x, 0.1, 0, 0, ierflg) ;
313 index++ ;
314 if(ierflg != 0){
21cd0c07 315 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : x = %f\n", x ) ;
9a1398dd 316 return kFALSE;
317 }
318 gMinuit->mnparm(index, "z", z, 0.1, 0, 0, ierflg) ;
319 index++ ;
320 if(ierflg != 0){
21cd0c07 321 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : z =%f\n", z ) ;
9a1398dd 322 return kFALSE;
323 }
324 gMinuit->mnparm(index, "Energy", energy , 0.05*energy, 0., 4.*energy, ierflg) ;
325 index++ ;
326 if(ierflg != 0){
21cd0c07 327 Warning("FindFit", "PHOS Unfolding unable to set initial value for fit procedure : energy = %f\n", energy ) ;
9a1398dd 328 return kFALSE;
329 }
330 }
331
332 Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
333 // depends on it.
334 Double_t p1 = 1.0 ;
335 Double_t p2 = 0.0 ;
336
337 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
338 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
339 gMinuit->SetMaxIterations(5);
340 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
341
342 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
343
344 if(ierflg == 4){ // Minimum not found
21cd0c07 345 Warning("FindFit", "PHOS Unfolding fit not converged, cluster abandoned\n" );
9a1398dd 346 return kFALSE ;
347 }
348 for(index = 0; index < nPar; index++){
349 Double_t err ;
350 Double_t val ;
351 gMinuit->GetParameter(index, val, err) ; // Returns value and error of parameter index
352 fitparameters[index] = val ;
353 }
354
355 delete toMinuit ;
356 return kTRUE;
c3c187e5 357
d15a28e7 358}
359
548f0134 360
d15a28e7 361//____________________________________________________________________________
a4e98857 362void AliPHOSClusterizerv1::Init()
363{
2bd5457f 364 // Make all memory allocations which can not be done in default constructor.
365 // Attach the Clusterizer task to the list of PHOS tasks
ba898748 366
9a2cdbdf 367 fEmcCrystals = fGeom->GetNModules() * fGeom->GetNCristalsInModule() ;
3758d9fc 368
7b7c1533 369 if(!gMinuit)
88cb7938 370 gMinuit = new TMinuit(100);
fbf811ec 371
a28333c5 372 if (!fgCalibData)
373 fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
374 if (fgCalibData->GetCalibDataEmc() == 0)
375 AliFatal("Calibration parameters for PHOS EMC not found. Stop reconstruction.\n");
e68222ce 376
9a1398dd 377}
7b7c1533 378
8d0f3f77 379//____________________________________________________________________________
380void AliPHOSClusterizerv1::InitParameters()
381{
382
383 fNumberOfCpvClusters = 0 ;
384 fNumberOfEmcClusters = 0 ;
3799bcb5 385
386 const AliPHOSRecoParam* parEmc = AliPHOSReconstructor::GetRecoParamEmc();
387 if(!parEmc) AliFatal("Reconstruction parameters for EMC not set!");
388
389 const AliPHOSRecoParam* parCpv = AliPHOSReconstructor::GetRecoParamCpv();
390 if(!parCpv) AliFatal("Reconstruction parameters for CPV not set!");
391
392 fCpvClusteringThreshold = parCpv->GetClusteringThreshold();
393 fEmcClusteringThreshold = parEmc->GetClusteringThreshold();
8d0f3f77 394
3799bcb5 395 fEmcLocMaxCut = parEmc->GetLocalMaxCut();
396 fCpvLocMaxCut = parCpv->GetLocalMaxCut();
ba898748 397
3799bcb5 398 fW0 = parEmc->GetLogWeight();
399 fW0CPV = parCpv->GetLogWeight();
8d0f3f77 400
454c62a4 401 fEmcTimeGate = 1.e-6 ;
d55c9739 402 fEcoreRadius = parEmc->GetEcoreRadius();
8d0f3f77 403
a42611a9 404 fToUnfold = parEmc->ToUnfold() ;
88cb7938 405
ba898748 406 fWrite = kTRUE ;
8d0f3f77 407}
408
9a1398dd 409//____________________________________________________________________________
410Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
d15a28e7 411{
b2a60966 412 // Gives the neighbourness of two digits = 0 are not neighbour but continue searching
413 // = 1 are neighbour
414 // = 2 are not neighbour but do not continue searching
d07a44a8 415 // =-1 are not neighbour, continue searching, but do not look before d2 next time
a4e98857 416 // neighbours are defined as digits having at least a common vertex
b2a60966 417 // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster
418 // which is compared to a digit (d2) not yet in a cluster
419
d15a28e7 420 Int_t relid1[4] ;
9a2cdbdf 421 fGeom->AbsToRelNumbering(d1->GetId(), relid1) ;
d15a28e7 422
423 Int_t relid2[4] ;
9a2cdbdf 424 fGeom->AbsToRelNumbering(d2->GetId(), relid2) ;
d15a28e7 425
9688c1dd 426 if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module
92862013 427 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
428 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
d15a28e7 429
454c62a4 430 if (( coldiff <= 1 ) && ( rowdiff <= 1 )){ //At least common vertex
431 // if (( relid1[2]==relid2[2] && coldiff <= 1 ) || ( relid1[3]==relid2[3] && rowdiff <= 1 )){ //common side
9688c1dd 432 if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
d07a44a8 433 return 1 ;
d15a28e7 434 }
435 else {
436 if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1))
d07a44a8 437 return 2; // Difference in row numbers is too large to look further
d15a28e7 438 }
d07a44a8 439 return 0 ;
d15a28e7 440
441 }
442 else {
d07a44a8 443 if(relid1[0] > relid2[0] && relid1[1]==relid2[1] ) //we switched to the next module
444 return -1 ;
445 if(relid1[1] < relid2[1]) //we switched from EMC(0) to CPV(-1)
446 return -1 ;
d15a28e7 447
d07a44a8 448 return 2 ;
d15a28e7 449
450 }
d72dfbc3 451
d07a44a8 452 return 0 ;
d15a28e7 453}
d15a28e7 454//____________________________________________________________________________
9a1398dd 455Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
d15a28e7 456{
b2a60966 457 // Tells if (true) or not (false) the digit is in a PHOS-EMC module
458
9f616d61 459 Bool_t rv = kFALSE ;
d15a28e7 460
9a2cdbdf 461 Int_t nEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
d15a28e7 462
2b629790 463 if(digit->GetId() <= nEMC ) rv = kTRUE;
ed4205d8 464
465 return rv ;
466}
467
ed4205d8 468//____________________________________________________________________________
9a1398dd 469Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
ed4205d8 470{
fad3e5b9 471 // Tells if (true) or not (false) the digit is in a PHOS-CPV module
ed4205d8 472
473 Bool_t rv = kFALSE ;
88cb7938 474
9a2cdbdf 475 Int_t nEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
ed4205d8 476
2b629790 477 if(digit->GetId() > nEMC ) rv = kTRUE;
d15a28e7 478
479 return rv ;
480}
9a1398dd 481
88cb7938 482//____________________________________________________________________________
483void AliPHOSClusterizerv1::WriteRecPoints()
a4e98857 484{
7b7c1533 485
486 // Creates new branches with given title
a4e98857 487 // fills and writes into TreeR.
88cb7938 488
9a1398dd 489 Int_t index ;
ba898748 490 //Evaluate position, dispersion and other RecPoint properties..
9a2cdbdf 491 Int_t nEmc = fEMCRecPoints->GetEntriesFast();
a28333c5 492 Float_t emcMinE= AliPHOSReconstructor::GetRecoParamEmc()->GetMinE(); //Minimal digit energy
e347da7b 493 TVector3 fakeVtx(0.,0.,0.) ;
ba898748 494 for(index = 0; index < nEmc; index++){
cd79ec76 495 AliPHOSEmcRecPoint * rp =
9a2cdbdf 496 dynamic_cast<AliPHOSEmcRecPoint *>( fEMCRecPoints->At(index) );
a28333c5 497 rp->Purify(emcMinE) ;
a978bcea 498 if(rp->GetMultiplicity()==0){
9a2cdbdf 499 fEMCRecPoints->RemoveAt(index) ;
ba898748 500 delete rp ;
cd79ec76 501 continue;
ba898748 502 }
a978bcea 503
cd79ec76 504 // No vertex is available now, calculate corrections in PID
e347da7b 505 rp->EvalAll(fDigitsArr) ;
d55c9739 506 rp->EvalCoreEnergy(fW0,fEcoreRadius,fDigitsArr) ;
9a2cdbdf 507 rp->EvalAll(fW0,fakeVtx,fDigitsArr) ;
cd79ec76 508 rp->EvalLocal2TrackingCSTransform();
ba898748 509 }
9a2cdbdf 510 fEMCRecPoints->Compress() ;
6d0a5c2e 511 fEMCRecPoints->Sort() ;
9a2cdbdf 512 // fEMCRecPoints->Expand(fEMCRecPoints->GetEntriesFast()) ;
513 for(index = 0; index < fEMCRecPoints->GetEntries(); index++){
514 dynamic_cast<AliPHOSEmcRecPoint *>( fEMCRecPoints->At(index) )->SetIndexInList(index) ;
ba898748 515 }
fbf811ec 516
af0b0570 517 //For each rec.point set the distance to the nearest bad crystal (BVP)
518 SetDistancesToBadChannels();
519
9a1398dd 520 //Now the same for CPV
9a2cdbdf 521 for(index = 0; index < fCPVRecPoints->GetEntries(); index++){
522 AliPHOSCpvRecPoint * rp = dynamic_cast<AliPHOSCpvRecPoint *>( fCPVRecPoints->At(index) );
e347da7b 523 rp->EvalAll(fDigitsArr) ;
524 rp->EvalAll(fW0CPV,fakeVtx,fDigitsArr) ;
76e42795 525 rp->EvalLocal2TrackingCSTransform();
ba898748 526 }
6d0a5c2e 527 fCPVRecPoints->Sort() ;
fbf811ec 528
9a2cdbdf 529 for(index = 0; index < fCPVRecPoints->GetEntries(); index++)
530 dynamic_cast<AliPHOSCpvRecPoint *>( fCPVRecPoints->At(index) )->SetIndexInList(index) ;
fbf811ec 531
9a2cdbdf 532 fCPVRecPoints->Expand(fCPVRecPoints->GetEntriesFast()) ;
9a1398dd 533
ba898748 534 if(fWrite){ //We write TreeR
9a2cdbdf 535 fTreeR->Fill();
ba898748 536 }
9a1398dd 537}
538
539//____________________________________________________________________________
540void AliPHOSClusterizerv1::MakeClusters()
541{
542 // Steering method to construct the clusters stored in a list of Reconstructed Points
543 // A cluster is defined as a list of neighbour digits
d15a28e7 544
8d8258f6 545 fNumberOfCpvClusters = 0 ;
546 fNumberOfEmcClusters = 0 ;
547
d07a44a8 548 //Mark all digits as unused yet
549 Int_t nDigits=fDigitsArr->GetEntriesFast() ;
550 for(Int_t i=0; i<nDigits; i++){
551 fDigitsUsed[i]=0 ;
552 }
553 Int_t iFirst = 0 ; //first index of digit which potentially can be a part of cluster
554 //e.g. first digit in this module, first CPV digit etc.
d15a28e7 555 AliPHOSDigit * digit ;
d07a44a8 556 TArrayI clusterdigitslist(1500) ;
557 AliPHOSRecPoint * clu = 0 ;
558 for(Int_t i=0; i<nDigits; i++){
559 if(fDigitsUsed[i])
560 continue ;
6ad0bfa0 561
d07a44a8 562 digit=static_cast<AliPHOSDigit*>(fDigitsArr->At(i)) ;
88cb7938 563
d07a44a8 564 clu=0 ;
c161df70 565
d15a28e7 566 Int_t index ;
f2bc1b87 567
d07a44a8 568 //is this digit so energetic that start cluster?
a28333c5 569 if (( IsInEmc(digit) && digit->GetEnergy() > fEmcClusteringThreshold ) ||
570 ( IsInCpv(digit) && digit->GetEnergy() > fCpvClusteringThreshold ) ) {
d15a28e7 571 Int_t iDigitInCluster = 0 ;
7b7c1533 572
6ad0bfa0 573 if ( IsInEmc(digit) ) {
88cb7938 574 // start a new EMC RecPoint
9a2cdbdf 575 if(fNumberOfEmcClusters >= fEMCRecPoints->GetSize())
576 fEMCRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
88cb7938 577
9a2cdbdf 578 fEMCRecPoints->AddAt(new AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
d07a44a8 579 clu = static_cast<AliPHOSEmcRecPoint *>( fEMCRecPoints->At(fNumberOfEmcClusters) ) ;
e95226ae 580 fNumberOfEmcClusters++ ;
a28333c5 581 clu->AddDigit(*digit, digit->GetEnergy()) ;
e95226ae 582 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
583 iDigitInCluster++ ;
d07a44a8 584 fDigitsUsed[i]=kTRUE ;
d72dfbc3 585 } else {
88cb7938 586 // start a new CPV cluster
9a2cdbdf 587 if(fNumberOfCpvClusters >= fCPVRecPoints->GetSize())
588 fCPVRecPoints->Expand(2*fNumberOfCpvClusters+1);
88cb7938 589
9a2cdbdf 590 fCPVRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
d07a44a8 591 clu = static_cast<AliPHOSCpvRecPoint *>( fCPVRecPoints->At(fNumberOfCpvClusters) ) ;
88cb7938 592 fNumberOfCpvClusters++ ;
a28333c5 593 clu->AddDigit(*digit, digit->GetEnergy()) ;
88cb7938 594 clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;
595 iDigitInCluster++ ;
d07a44a8 596 fDigitsUsed[i]=kTRUE ;
d15a28e7 597 } // else
598
d07a44a8 599 //Now scan remaining digits in list to find neigbours of our seed
600
d15a28e7 601 AliPHOSDigit * digitN ;
602 index = 0 ;
9f616d61 603 while (index < iDigitInCluster){ // scan over digits already in cluster
d07a44a8 604 digit = static_cast<AliPHOSDigit*>( fDigitsArr->At(clusterdigitslist[index]) ) ;
88cb7938 605 index++ ;
d07a44a8 606 for(Int_t j=iFirst; j<nDigits; j++){
607 if(fDigitsUsed[j])
608 continue ; //look through remaining digits
609 digitN = static_cast<AliPHOSDigit*>( fDigitsArr->At(j) ) ;
88cb7938 610 Int_t ineb = AreNeighbours(digit, digitN); // call (digit,digitN) in THAT oder !!!!!
d15a28e7 611 switch (ineb ) {
d07a44a8 612 case -1: //too early (e.g. previous module), do not look before j at subsequent passes
613 iFirst=j ;
614 break ;
9f616d61 615 case 0 : // not a neighbour
88cb7938 616 break ;
617 case 1 : // are neighbours
a28333c5 618 clu->AddDigit(*digitN, digitN->GetEnergy());
d07a44a8 619 clusterdigitslist[iDigitInCluster] = j ;
88cb7938 620 iDigitInCluster++ ;
d07a44a8 621 fDigitsUsed[j]=kTRUE ;
88cb7938 622 break ;
9f616d61 623 case 2 : // too far from each other
88cb7938 624 goto endofloop;
625 } // switch
626
d07a44a8 627 }
88cb7938 628
d07a44a8 629 endofloop: ; //scanned all possible neighbours for this digit
88cb7938 630
d15a28e7 631 } // loop over cluster
ad8cfaf4 632 } // energy theshold
d07a44a8 633 }
9688c1dd 634
9a1398dd 635}
636
637//____________________________________________________________________________
a4e98857 638void AliPHOSClusterizerv1::MakeUnfolding()
639{
640 // Unfolds clusters using the shape of an ElectroMagnetic shower
9688c1dd 641 // Performs unfolding of all EMC/CPV clusters
9a1398dd 642
a4e98857 643 // Unfold first EMC clusters
9a1398dd 644 if(fNumberOfEmcClusters > 0){
645
9a2cdbdf 646 Int_t nModulesToUnfold = fGeom->GetNModules() ;
9a1398dd 647
648 Int_t numberofNotUnfolded = fNumberOfEmcClusters ;
649 Int_t index ;
650 for(index = 0 ; index < numberofNotUnfolded ; index++){
651
9a2cdbdf 652 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint *>( fEMCRecPoints->At(index) ) ;
9a1398dd 653 if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
88cb7938 654 break ;
9a1398dd 655
656 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
a0636361 657 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
9a1398dd 658 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
9a2cdbdf 659 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,fDigitsArr) ;
9a1398dd 660
661 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
88cb7938 662 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
a28333c5 663
9a2cdbdf 664 fEMCRecPoints->Remove(emcRecPoint);
665 fEMCRecPoints->Compress() ;
88cb7938 666 index-- ;
667 fNumberOfEmcClusters -- ;
668 numberofNotUnfolded-- ;
9a1398dd 669 }
c6bf27f2 670 else{
671 emcRecPoint->SetNExMax(1) ; //Only one local maximum
672 }
9a1398dd 673
674 delete[] maxAt ;
675 delete[] maxAtEnergy ;
676 }
677 }
a4e98857 678 // Unfolding of EMC clusters finished
9a1398dd 679
680
a4e98857 681 // Unfold now CPV clusters
9a1398dd 682 if(fNumberOfCpvClusters > 0){
683
9a2cdbdf 684 Int_t nModulesToUnfold = fGeom->GetNModules() ;
9a1398dd 685
686 Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;
687 Int_t index ;
688 for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
689
9a2cdbdf 690 AliPHOSRecPoint * recPoint = dynamic_cast<AliPHOSRecPoint *>( fCPVRecPoints->At(index) ) ;
9a1398dd 691
692 if(recPoint->GetPHOSMod()> nModulesToUnfold)
88cb7938 693 break ;
9a1398dd 694
88cb7938 695 AliPHOSEmcRecPoint * emcRecPoint = dynamic_cast<AliPHOSEmcRecPoint*>(recPoint) ;
9a1398dd 696
697 Int_t nMultipl = emcRecPoint->GetMultiplicity() ;
a0636361 698 AliPHOSDigit ** maxAt = new AliPHOSDigit*[nMultipl] ;
9a1398dd 699 Float_t * maxAtEnergy = new Float_t[nMultipl] ;
9a2cdbdf 700 Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,fDigitsArr) ;
9a1398dd 701
702 if( nMax > 1 ) { // if cluster is very flat (no pronounced maximum) then nMax = 0
88cb7938 703 UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
9a2cdbdf 704 fCPVRecPoints->Remove(emcRecPoint);
705 fCPVRecPoints->Compress() ;
88cb7938 706 index-- ;
707 numberofCpvNotUnfolded-- ;
708 fNumberOfCpvClusters-- ;
9a1398dd 709 }
710
711 delete[] maxAt ;
712 delete[] maxAtEnergy ;
713 }
714 }
715 //Unfolding of Cpv clusters finished
716
717}
718
9a1398dd 719//____________________________________________________________________________
a978bcea 720Double_t AliPHOSClusterizerv1::ShowerShape(Double_t x, Double_t z)
9a1398dd 721{
722 // Shape of the shower (see PHOS TDR)
a4e98857 723 // If you change this function, change also the gradient evaluation in ChiSquare()
9a1398dd 724
a978bcea 725 //for the moment we neglect dependence on the incident angle.
726
727 Double_t r2 = x*x + z*z ;
728 Double_t r4 = r2*r2 ;
729 Double_t r295 = TMath::Power(r2, 2.95/2.) ;
9a1398dd 730 Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
731 return shape ;
732}
733
734//____________________________________________________________________________
735void AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc,
88cb7938 736 Int_t nMax,
737 AliPHOSDigit ** maxAt,
738 Float_t * maxAtEnergy)
9a1398dd 739{
740 // Performs the unfolding of a cluster with nMax overlapping showers
741
742 Int_t nPar = 3 * nMax ;
743 Float_t * fitparameters = new Float_t[nPar] ;
744
745 Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
a28333c5 746
9a1398dd 747 if( !rv ) {
88cb7938 748 // Fit failed, return and remove cluster
c6bf27f2 749 iniEmc->SetNExMax(-1) ;
9a1398dd 750 delete[] fitparameters ;
751 return ;
752 }
753
754 // create ufolded rec points and fill them with new energy lists
755 // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
756 // and later correct this number in acordance with actual energy deposition
757
758 Int_t nDigits = iniEmc->GetMultiplicity() ;
759 Float_t * efit = new Float_t[nDigits] ;
a978bcea 760 Float_t xDigit=0.,zDigit=0. ;
7b7c1533 761 Float_t xpar=0.,zpar=0.,epar=0. ;
9a1398dd 762 Int_t relid[4] ;
7b7c1533 763 AliPHOSDigit * digit = 0 ;
9a1398dd 764 Int_t * emcDigits = iniEmc->GetDigitsList() ;
765
a978bcea 766 TVector3 vIncid ;
767
9a1398dd 768 Int_t iparam ;
769 Int_t iDigit ;
770 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
9a2cdbdf 771 digit = dynamic_cast<AliPHOSDigit*>( fDigitsArr->At(emcDigits[iDigit] ) ) ;
772 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
773 fGeom->RelPosInModule(relid, xDigit, zDigit) ;
9a1398dd 774 efit[iDigit] = 0;
775
776 iparam = 0 ;
777 while(iparam < nPar ){
778 xpar = fitparameters[iparam] ;
779 zpar = fitparameters[iparam+1] ;
780 epar = fitparameters[iparam+2] ;
781 iparam += 3 ;
9a2cdbdf 782// fGeom->GetIncidentVector(fVtx,relid[0],xpar,zpar,vIncid) ;
a978bcea 783// efit[iDigit] += epar * ShowerShape(xDigit - xpar,zDigit - zpar,vIncid) ;
784 efit[iDigit] += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
9a1398dd 785 }
786 }
787
9a1398dd 788 // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
789 // so that energy deposited in each cell is distributed betwin new clusters proportionally
790 // to its contribution to efit
791
792 Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
793 Float_t ratio ;
794
795 iparam = 0 ;
796 while(iparam < nPar ){
797 xpar = fitparameters[iparam] ;
798 zpar = fitparameters[iparam+1] ;
799 epar = fitparameters[iparam+2] ;
800 iparam += 3 ;
9a2cdbdf 801// fGeom->GetIncidentVector(fVtx,iniEmc->GetPHOSMod(),xpar,zpar,vIncid) ;
a28333c5 802
7b7c1533 803 AliPHOSEmcRecPoint * emcRP = 0 ;
9a1398dd 804
9a2cdbdf 805 if(iniEmc->IsEmc()){ //create new entries in fEMCRecPoints...
9a1398dd 806
9a2cdbdf 807 if(fNumberOfEmcClusters >= fEMCRecPoints->GetSize())
808 fEMCRecPoints->Expand(2*fNumberOfEmcClusters) ;
9a1398dd 809
9a2cdbdf 810 (*fEMCRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
811 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( fEMCRecPoints->At(fNumberOfEmcClusters) ) ;
9a1398dd 812 fNumberOfEmcClusters++ ;
c6bf27f2 813 emcRP->SetNExMax((Int_t)nPar/3) ;
9a1398dd 814 }
9a2cdbdf 815 else{//create new entries in fCPVRecPoints
816 if(fNumberOfCpvClusters >= fCPVRecPoints->GetSize())
817 fCPVRecPoints->Expand(2*fNumberOfCpvClusters) ;
9a1398dd 818
9a2cdbdf 819 (*fCPVRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
820 emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( fCPVRecPoints->At(fNumberOfCpvClusters) ) ;
9a1398dd 821 fNumberOfCpvClusters++ ;
822 }
823
824 Float_t eDigit ;
825 for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
9a2cdbdf 826 digit = dynamic_cast<AliPHOSDigit*>( fDigitsArr->At( emcDigits[iDigit] ) ) ;
827 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
828 fGeom->RelPosInModule(relid, xDigit, zDigit) ;
a978bcea 829// ratio = epar * ShowerShape(xDigit - xpar,zDigit - zpar,vIncid) / efit[iDigit] ;
830 ratio = epar * ShowerShape(xDigit - xpar,zDigit - zpar) / efit[iDigit] ;
9a1398dd 831 eDigit = emcEnergies[iDigit] * ratio ;
832 emcRP->AddDigit( *digit, eDigit ) ;
88cb7938 833 }
9a1398dd 834 }
835
836 delete[] fitparameters ;
837 delete[] efit ;
838
839}
840
841//_____________________________________________________________________________
842void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
843{
a4e98857 844 // Calculates the Chi square for the cluster unfolding minimization
9a1398dd 845 // Number of parameters, Gradient, Chi squared, parameters, what to do
846
88cb7938 847 TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;
9a1398dd 848
88cb7938 849 AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint*>( toMinuit->At(0) ) ;
850 TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) ) ;
9a2cdbdf 851 // A bit buggy way to get an access to the geometry
852 // To be revised!
853 AliPHOSGeometry *geom = dynamic_cast<AliPHOSGeometry *>(toMinuit->At(2));
854
855// TVector3 * vtx = dynamic_cast<TVector3*>(toMinuit->At(3)) ; //Vertex position
9a1398dd 856
88cb7938 857 // AliPHOSEmcRecPoint * emcRP = dynamic_cast<AliPHOSEmcRecPoint *>( gMinuit->GetObjectFit() ) ; // EmcRecPoint to fit
9a1398dd 858
859 Int_t * emcDigits = emcRP->GetDigitsList() ;
860
7b7c1533 861 Int_t nOdigits = emcRP->GetDigitsMultiplicity() ;
9a1398dd 862
863 Float_t * emcEnergies = emcRP->GetEnergiesList() ;
88cb7938 864
a978bcea 865// TVector3 vInc ;
9a1398dd 866 fret = 0. ;
867 Int_t iparam ;
868
869 if(iflag == 2)
870 for(iparam = 0 ; iparam < nPar ; iparam++)
871 Grad[iparam] = 0 ; // Will evaluate gradient
872
873 Double_t efit ;
874
875 AliPHOSDigit * digit ;
876 Int_t iDigit ;
877
7b7c1533 878 for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
9a1398dd 879
88cb7938 880 digit = dynamic_cast<AliPHOSDigit*>( digits->At( emcDigits[iDigit] ) );
9a1398dd 881
882 Int_t relid[4] ;
883 Float_t xDigit ;
884 Float_t zDigit ;
885
886 geom->AbsToRelNumbering(digit->GetId(), relid) ;
887
888 geom->RelPosInModule(relid, xDigit, zDigit) ;
889
890 if(iflag == 2){ // calculate gradient
891 Int_t iParam = 0 ;
892 efit = 0 ;
893 while(iParam < nPar ){
a978bcea 894 Double_t dx = (xDigit - x[iParam]) ;
88cb7938 895 iParam++ ;
a978bcea 896 Double_t dz = (zDigit - x[iParam]) ;
88cb7938 897 iParam++ ;
9a2cdbdf 898// fGeom->GetIncidentVector(*vtx,emcRP->GetPHOSMod(),x[iParam-2],x[iParam-1],vInc) ;
a978bcea 899// efit += x[iParam] * ShowerShape(dx,dz,vInc) ;
900 efit += x[iParam] * ShowerShape(dx,dz) ;
88cb7938 901 iParam++ ;
9a1398dd 902 }
88cb7938 903 Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E)
9a1398dd 904 iParam = 0 ;
905 while(iParam < nPar ){
88cb7938 906 Double_t xpar = x[iParam] ;
907 Double_t zpar = x[iParam+1] ;
908 Double_t epar = x[iParam+2] ;
9a2cdbdf 909// fGeom->GetIncidentVector(*vtx,emcRP->GetPHOSMod(),xpar,zpar,vInc) ;
88cb7938 910 Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
a978bcea 911// Double_t shape = sum * ShowerShape(xDigit - xpar,zDigit - zpar,vInc) ;
912 Double_t shape = sum * ShowerShape(xDigit - xpar,zDigit - zpar) ;
913//DP: No incident angle dependence in gradient yet!!!!!!
88cb7938 914 Double_t r4 = dr*dr*dr*dr ;
915 Double_t r295 = TMath::Power(dr,2.95) ;
916 Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
917 0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
918
919 Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ; // Derivative over x
920 iParam++ ;
921 Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ; // Derivative over z
922 iParam++ ;
923 Grad[iParam] += shape ; // Derivative over energy
924 iParam++ ;
9a1398dd 925 }
926 }
927 efit = 0;
928 iparam = 0 ;
929
930 while(iparam < nPar ){
931 Double_t xpar = x[iparam] ;
932 Double_t zpar = x[iparam+1] ;
933 Double_t epar = x[iparam+2] ;
934 iparam += 3 ;
9a2cdbdf 935// fGeom->GetIncidentVector(*vtx,emcRP->GetPHOSMod(),xpar,zpar,vInc) ;
a978bcea 936// efit += epar * ShowerShape(xDigit - xpar,zDigit - zpar,vInc) ;
937 efit += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
9a1398dd 938 }
939
88cb7938 940 fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ;
9a1398dd 941 // Here we assume, that sigma = sqrt(E)
942 }
83974468 943
d15a28e7 944}
945
946//____________________________________________________________________________
702ab87e 947void AliPHOSClusterizerv1::Print(const Option_t *)const
d15a28e7 948{
d1de15f5 949 // Print clusterizer parameters
950
21cd0c07 951 TString message ;
952 TString taskName(GetName()) ;
953 taskName.ReplaceAll(Version(), "") ;
954
955 if( strcmp(GetName(), "") !=0 ) {
9a1398dd 956 // Print parameters
21cd0c07 957 message = "\n--------------- %s %s -----------\n" ;
958 message += "Clusterizing digits from the file: %s\n" ;
959 message += " Branch: %s\n" ;
960 message += " EMC Clustering threshold = %f\n" ;
961 message += " EMC Local Maximum cut = %f\n" ;
962 message += " EMC Logarothmic weight = %f\n" ;
963 message += " CPV Clustering threshold = %f\n" ;
964 message += " CPV Local Maximum cut = %f\n" ;
965 message += " CPV Logarothmic weight = %f\n" ;
9a1398dd 966 if(fToUnfold)
21cd0c07 967 message += " Unfolding on\n" ;
9a1398dd 968 else
21cd0c07 969 message += " Unfolding off\n" ;
9a1398dd 970
21cd0c07 971 message += "------------------------------------------------------------------" ;
9a1398dd 972 }
21cd0c07 973 else
974 message = " AliPHOSClusterizerv1 not initialized " ;
975
351dd634 976 AliInfo(Form("%s, %s %s %s %s %s %s %s %s %s %s", message.Data(),
21cd0c07 977 taskName.Data(),
978 GetTitle(),
979 taskName.Data(),
980 GetName(),
981 fEmcClusteringThreshold,
982 fEmcLocMaxCut,
983 fW0,
984 fCpvClusteringThreshold,
985 fCpvLocMaxCut,
351dd634 986 fW0CPV )) ;
9a1398dd 987}
a978bcea 988//____________________________________________________________________________
a4e98857 989void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
990{
991 // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
9a1398dd 992
9a2cdbdf 993 AliInfo(Form("\nFound %d EMC RecPoints and %d CPV RecPoints",
994 fEMCRecPoints->GetEntriesFast(),
995 fCPVRecPoints->GetEntriesFast() )) ;
88cb7938 996
9a1398dd 997 if(strstr(option,"all")) {
709e117a 998 printf("\n EMC clusters \n") ;
999 printf("Index Ene(MeV) Multi Module X Y Z Lambdas_1 Lambda_2 # of prim Primaries list\n") ;
9a1398dd 1000 Int_t index ;
9a2cdbdf 1001 for (index = 0 ; index < fEMCRecPoints->GetEntries() ; index++) {
1002 AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )fEMCRecPoints->At(index) ;
9a1398dd 1003 TVector3 locpos;
1004 rp->GetLocalPosition(locpos);
9a1398dd 1005 Float_t lambda[2];
1006 rp->GetElipsAxis(lambda);
9a1398dd 1007 Int_t * primaries;
1008 Int_t nprimaries;
1009 primaries = rp->GetPrimaries(nprimaries);
709e117a 1010 printf("\n%6d %8.2f %3d %2d %4.1f %4.1f %4.1f %4f %4f %2d : ",
11f9c5ff 1011 rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(), rp->GetPHOSMod(),
1012 locpos.X(), locpos.Y(), locpos.Z(), lambda[0], lambda[1], nprimaries) ;
3bf72d32 1013
21cd0c07 1014 for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
709e117a 1015 printf("%d ", primaries[iprimary] ) ;
21cd0c07 1016 }
709e117a 1017 printf("\n") ;
1018 }
1019
1020 //Now plot CPV recPoints
1021 printf("\n CPV clusters \n") ;
1022 printf("Index Ene(MeV) Module X Y Z \n") ;
9a2cdbdf 1023 for (index = 0 ; index < fCPVRecPoints->GetEntries() ; index++) {
1024 AliPHOSCpvRecPoint * rp = (AliPHOSCpvRecPoint * )fCPVRecPoints->At(index) ;
2bb500e5 1025
1026 TVector3 locpos;
1027 rp->GetLocalPosition(locpos);
1028
1029 printf("\n%6d %8.2f %2d %4.1f %4.1f %4.1f \n",
1030 rp->GetIndexInList(), rp->GetEnergy(), rp->GetPHOSMod(),
1031 locpos.X(), locpos.Y(), locpos.Z()) ;
88cb7938 1032 }
9a1398dd 1033 }
d15a28e7 1034}
9a1398dd 1035
af0b0570 1036
1037//____________________________________________________________________________
1038void AliPHOSClusterizerv1::SetDistancesToBadChannels()
1039{
1040 //For each EMC rec. point set the distance to the nearest bad crystal.
1041 //Author: Boris Polichtchouk
1042
e68222ce 1043 if(!fgCalibData->GetNumOfEmcBadChannels()) return;
1044 AliInfo(Form("%d bad channel(s) found.\n",fgCalibData->GetNumOfEmcBadChannels()));
af0b0570 1045
af0b0570 1046 Int_t badIds[8000];
e68222ce 1047 fgCalibData->EmcBadChannelIds(badIds);
af0b0570 1048
1049 AliPHOSEmcRecPoint* rp;
1050
1051 TMatrixF gmat;
1052 TVector3 gposRecPoint; // global (in ALICE frame) position of rec. point
1053 TVector3 gposBadChannel; // global position of bad crystal
1054 TVector3 dR;
1055
1056 Float_t dist,minDist;
79987dcd 1057 Int_t relid[4]={0,0,0,0} ;
1058 TVector3 lpos ;
9a2cdbdf 1059 for(Int_t iRP=0; iRP<fEMCRecPoints->GetEntries(); iRP++){
1060 rp = (AliPHOSEmcRecPoint*)fEMCRecPoints->At(iRP);
79987dcd 1061 //evaluate distance to border
1062 relid[0]=rp->GetPHOSMod() ;
1063 relid[2]=1 ;
1064 relid[3]=1 ;
1065 Float_t xcorner,zcorner;
1066 fGeom->RelPosInModule(relid,xcorner,zcorner) ; //coordinate of the corner cell
1067 rp->GetLocalPosition(lpos) ;
1068 minDist = 2.2+TMath::Min(-xcorner-TMath::Abs(lpos.X()),-zcorner-TMath::Abs(lpos.Z())); //2.2 - crystal size
e68222ce 1069 for(Int_t iBad=0; iBad<fgCalibData->GetNumOfEmcBadChannels(); iBad++) {
e250de8f 1070 fGeom->AbsToRelNumbering(badIds[iBad],relid) ;
79987dcd 1071 if(relid[0]!=rp->GetPHOSMod()) //We can not evaluate global position directly since
1072 continue ; //bad channels can be in the module which does not exist in simulations.
af0b0570 1073 rp->GetGlobalPosition(gposRecPoint,gmat);
9a2cdbdf 1074 fGeom->RelPosInAlice(badIds[iBad],gposBadChannel);
af0b0570 1075 AliDebug(2,Form("BC position:[%.3f,%.3f,%.3f], RP position:[%.3f,%.3f,%.3f]. E=%.3f\n",
1076 gposBadChannel.X(),gposBadChannel.Y(),gposBadChannel.Z(),
1077 gposRecPoint.X(),gposRecPoint.Y(),gposRecPoint.Z(),rp->GetEnergy()));
1078 dR = gposBadChannel-gposRecPoint;
1079 dist = dR.Mag();
1080 if(dist<minDist) minDist = dist;
1081 }
1082
1083 rp->SetDistanceToBadCrystal(minDist);
1084 }
1085
1086}