]>
Commit | Line | Data |
---|---|---|
ab48128d | 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 | /* $Id$ */ | |
16 | //_________________________________________________________________________ | |
70a93198 | 17 | // Reconstructed Points for the EMCAL |
18 | // A RecPoint is a cluster of digits | |
40164976 | 19 | // |
20 | // | |
d64c959b | 21 | //*-- Author: Yves Schutz (SUBATECH) |
70a93198 | 22 | //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH) |
23 | //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04 | |
ab48128d | 24 | |
25 | // --- ROOT system --- | |
9aa6a5f6 | 26 | #include "TPad.h" |
27 | #include "TGraph.h" | |
28 | #include "TPaveText.h" | |
29 | #include "TClonesArray.h" | |
30 | #include "TMath.h" | |
31 | #include "TGeoMatrix.h" | |
32 | #include "TGeoManager.h" | |
33 | #include "TGeoPhysicalNode.h" | |
25bb3dcb | 34 | #include "TRandom.h" |
ab48128d | 35 | |
36 | // --- Standard library --- | |
9aa6a5f6 | 37 | #include <Riostream.h> |
ab48128d | 38 | |
39 | // --- AliRoot header files --- | |
1d59832c | 40 | //#include "AliGenerator.h" |
41 | class AliGenerator; | |
1d59832c | 42 | class AliEMCAL; |
9aa6a5f6 | 43 | #include "AliLog.h" |
44 | #include "AliGeomManager.h" | |
ab48128d | 45 | #include "AliEMCALGeometry.h" |
4635df1f | 46 | #include "AliEMCALHit.h" |
ab48128d | 47 | #include "AliEMCALDigit.h" |
48 | #include "AliEMCALRecPoint.h" | |
40164976 | 49 | #include "AliCaloCalibPedestal.h" |
50 | #include "AliEMCALGeoParams.h" | |
ab48128d | 51 | |
52 | ClassImp(AliEMCALRecPoint) | |
53 | ||
ab48128d | 54 | //____________________________________________________________________________ |
55 | AliEMCALRecPoint::AliEMCALRecPoint() | |
9aa6a5f6 | 56 | : AliCluster(), fGeomPtr(0), |
57 | fAmp(0), fIndexInList(-1), //to be set when the point is already stored | |
25bb3dcb | 58 | fGlobPos(0,0,0),fLocPos(0,0,0), |
9aa6a5f6 | 59 | fMaxDigit(100), fMulDigit(0), fMaxTrack(200), |
60 | fMulTrack(0), fDigitsList(0), fTracksList(0), | |
61 | fClusterType(-1), fCoreEnergy(0), fDispersion(0), | |
62 | fEnergyList(0), fTimeList(0), fAbsIdList(0), | |
0561e246 | 63 | fTime(0.), fNExMax(0), fCoreRadius(10), //HG check this |
9aa6a5f6 | 64 | fDETracksList(0), fMulParent(0), fMaxParent(0), |
65 | fParentsList(0), fDEParentsList(0), fSuperModuleNumber(0), | |
25bb3dcb | 66 | fDigitIndMax(-1), fDistToBadTower(-1), fSharedCluster(kFALSE) |
ab48128d | 67 | { |
68 | // ctor | |
9aa6a5f6 | 69 | fGeomPtr = AliEMCALGeometry::GetInstance(); |
70 | ||
706863b6 | 71 | fLambda[0] = 0; |
72 | fLambda[1] = 0; | |
9aa6a5f6 | 73 | |
ab48128d | 74 | } |
75 | ||
76 | //____________________________________________________________________________ | |
9aa6a5f6 | 77 | AliEMCALRecPoint::AliEMCALRecPoint(const char *) |
78 | : AliCluster(), fGeomPtr(0), | |
79 | fAmp(0), fIndexInList(-1), //to be set when the point is already stored | |
25bb3dcb | 80 | fGlobPos(0,0,0), fLocPos(0,0,0), |
9aa6a5f6 | 81 | fMaxDigit(100), fMulDigit(0), fMaxTrack(1000), fMulTrack(0), |
82 | fDigitsList(new Int_t[fMaxDigit]), fTracksList(new Int_t[fMaxTrack]), | |
83 | fClusterType(-1), fCoreEnergy(0), fDispersion(0), | |
84 | fEnergyList(new Float_t[fMaxDigit]), fTimeList(new Float_t[fMaxDigit]), | |
0561e246 | 85 | fAbsIdList(new Int_t[fMaxDigit]), fTime(-1.), fNExMax(0), fCoreRadius(10), |
9aa6a5f6 | 86 | fDETracksList(new Float_t[fMaxTrack]), fMulParent(0), fMaxParent(1000), |
87 | fParentsList(new Int_t[fMaxParent]), fDEParentsList(new Float_t[fMaxParent]), | |
25bb3dcb | 88 | fSuperModuleNumber(0), fDigitIndMax(-1), fDistToBadTower(-1),fSharedCluster(kFALSE) |
ab48128d | 89 | { |
90 | // ctor | |
af5bdd85 | 91 | for (Int_t i = 0; i < fMaxTrack; i++) |
92 | fDETracksList[i] = 0; | |
94478418 | 93 | for (Int_t i = 0; i < fMaxParent; i++) { |
94 | fParentsList[i] = -1; | |
af5bdd85 | 95 | fDEParentsList[i] = 0; |
94478418 | 96 | } |
18a21c7c | 97 | |
9aa6a5f6 | 98 | fGeomPtr = AliEMCALGeometry::GetInstance(); |
706863b6 | 99 | fLambda[0] = 0; |
100 | fLambda[1] = 0; | |
70a93198 | 101 | } |
18a21c7c | 102 | |
70a93198 | 103 | //____________________________________________________________________________ |
18a21c7c | 104 | AliEMCALRecPoint::AliEMCALRecPoint(const AliEMCALRecPoint & rp) |
9aa6a5f6 | 105 | : AliCluster(rp), fGeomPtr(rp.fGeomPtr), |
106 | fAmp(rp.fAmp), fIndexInList(rp.fIndexInList), | |
25bb3dcb | 107 | fGlobPos(rp.fGlobPos),fLocPos(rp.fLocPos), |
9aa6a5f6 | 108 | fMaxDigit(rp.fMaxDigit), fMulDigit(rp.fMulDigit), |
109 | fMaxTrack(rp.fMaxTrack), fMulTrack(rp.fMaxTrack), | |
110 | fDigitsList(new Int_t[rp.fMaxDigit]), fTracksList(new Int_t[rp.fMaxTrack]), | |
111 | fClusterType(rp.fClusterType), fCoreEnergy(rp.fCoreEnergy), | |
18a21c7c | 112 | fDispersion(rp.fDispersion), |
9aa6a5f6 | 113 | fEnergyList(new Float_t[rp.fMaxDigit]), fTimeList(new Float_t[rp.fMaxDigit]), |
0561e246 | 114 | fAbsIdList(new Int_t[rp.fMaxDigit]), fTime(rp.fTime), fNExMax(rp.fNExMax),fCoreRadius(rp.fCoreRadius), |
9aa6a5f6 | 115 | fDETracksList(new Float_t[rp.fMaxTrack]), fMulParent(rp.fMulParent), |
116 | fMaxParent(rp.fMaxParent), fParentsList(new Int_t[rp.fMaxParent]), | |
117 | fDEParentsList(new Float_t[rp.fMaxParent]), | |
40164976 | 118 | fSuperModuleNumber(rp.fSuperModuleNumber), fDigitIndMax(rp.fDigitIndMax), |
25bb3dcb | 119 | fDistToBadTower(rp.fDistToBadTower), fSharedCluster(rp.fSharedCluster) |
0a4cb131 | 120 | { |
121 | //copy ctor | |
0a4cb131 | 122 | fLambda[0] = rp.fLambda[0]; |
123 | fLambda[1] = rp.fLambda[1]; | |
18a21c7c | 124 | |
0a4cb131 | 125 | for(Int_t i = 0; i < rp.fMulDigit; i++) { |
126 | fEnergyList[i] = rp.fEnergyList[i]; | |
127 | fTimeList[i] = rp.fTimeList[i]; | |
128 | fAbsIdList[i] = rp.fAbsIdList[i]; | |
129 | } | |
9aa6a5f6 | 130 | |
af5bdd85 | 131 | for(Int_t i = 0; i < rp.fMulTrack; i++) fDETracksList[i] = rp.fDETracksList[i]; |
9aa6a5f6 | 132 | |
133 | for(Int_t i = 0; i < rp.fMulParent; i++) { | |
134 | fParentsList[i] = rp.fParentsList[i]; | |
135 | fDEParentsList[i] = rp.fDEParentsList[i]; | |
136 | } | |
0a4cb131 | 137 | |
138 | } | |
139 | //____________________________________________________________________________ | |
70a93198 | 140 | AliEMCALRecPoint::~AliEMCALRecPoint() |
141 | { | |
142 | // dtor | |
143 | if ( fEnergyList ) | |
144 | delete[] fEnergyList ; | |
85c60a8e | 145 | if ( fTimeList ) |
146 | delete[] fTimeList ; | |
e52475ed | 147 | if ( fAbsIdList ) |
148 | delete[] fAbsIdList ; | |
af5bdd85 | 149 | if ( fDETracksList) |
150 | delete[] fDETracksList; | |
87cdc3be | 151 | if ( fParentsList) |
152 | delete[] fParentsList; | |
af5bdd85 | 153 | if ( fDEParentsList) |
154 | delete[] fDEParentsList; | |
25bb3dcb | 155 | |
9aa6a5f6 | 156 | delete [] fDigitsList ; |
157 | delete [] fTracksList ; | |
158 | } | |
159 | ||
160 | //____________________________________________________________________________ | |
161 | AliEMCALRecPoint& AliEMCALRecPoint::operator= (const AliEMCALRecPoint &rp) | |
162 | { | |
37890aaf | 163 | // assignment operator |
164 | ||
9aa6a5f6 | 165 | if(&rp == this) return *this; |
166 | ||
167 | fGeomPtr = rp.fGeomPtr; | |
168 | fAmp = rp.fAmp; | |
169 | fIndexInList = rp.fIndexInList; | |
25bb3dcb | 170 | fGlobPos = rp.fGlobPos; |
171 | fLocPos = rp.fLocPos; | |
9aa6a5f6 | 172 | fMaxDigit = rp.fMaxDigit; |
173 | fMulDigit = rp.fMulDigit; | |
174 | fMaxTrack = rp.fMaxTrack; | |
175 | fMulTrack = rp.fMaxTrack; | |
176 | for(Int_t i = 0; i<fMaxDigit; i++) fDigitsList[i] = rp.fDigitsList[i]; | |
177 | for(Int_t i = 0; i<fMaxTrack; i++) fTracksList[i] = rp.fTracksList[i]; | |
178 | fClusterType = rp.fClusterType; | |
6ce4e2ed | 179 | fCoreEnergy = rp.fCoreEnergy; |
180 | fDispersion = rp.fDispersion; | |
9aa6a5f6 | 181 | for(Int_t i = 0; i<fMaxDigit; i++) { |
182 | fEnergyList[i] = rp.fEnergyList[i]; | |
6ce4e2ed | 183 | fTimeList[i] = rp.fTimeList[i]; |
184 | fAbsIdList[i] = rp.fAbsIdList[i]; | |
9aa6a5f6 | 185 | } |
186 | fTime = rp.fTime; | |
0561e246 | 187 | fNExMax = rp.fNExMax; |
9aa6a5f6 | 188 | fCoreRadius = rp.fCoreRadius; |
189 | for(Int_t i = 0; i < fMaxTrack; i++) fDETracksList[i] = rp.fDETracksList[i]; | |
190 | fMulParent = rp.fMulParent; | |
191 | fMaxParent = rp.fMaxParent; | |
192 | for(Int_t i = 0; i < fMaxParent; i++) { | |
193 | fParentsList[i] = rp.fParentsList[i]; | |
194 | fDEParentsList[i] = rp.fDEParentsList[i]; | |
195 | } | |
196 | fSuperModuleNumber = rp.fSuperModuleNumber; | |
197 | fDigitIndMax = rp.fDigitIndMax; | |
198 | ||
199 | fLambda[0] = rp.fLambda[0]; | |
200 | fLambda[1] = rp.fLambda[1]; | |
40164976 | 201 | |
202 | fDistToBadTower = rp.fDistToBadTower; | |
25bb3dcb | 203 | fSharedCluster = rp.fSharedCluster; |
40164976 | 204 | |
9aa6a5f6 | 205 | return *this; |
206 | ||
70a93198 | 207 | } |
208 | ||
209 | //____________________________________________________________________________ | |
6ce4e2ed | 210 | void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t energy, Bool_t shared) |
70a93198 | 211 | { |
212 | // Adds a digit to the RecPoint | |
213 | // and accumulates the total amplitude and the multiplicity | |
214 | ||
215 | if(fEnergyList == 0) | |
216 | fEnergyList = new Float_t[fMaxDigit]; | |
85c60a8e | 217 | if(fTimeList == 0) |
218 | fTimeList = new Float_t[fMaxDigit]; | |
e52475ed | 219 | if(fAbsIdList == 0) { |
220 | fAbsIdList = new Int_t[fMaxDigit]; | |
e52475ed | 221 | } |
70a93198 | 222 | |
223 | if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists | |
224 | fMaxDigit*=2 ; | |
6ce4e2ed | 225 | Int_t * tempo = new Int_t[fMaxDigit]; |
226 | Float_t * tempoE = new Float_t[fMaxDigit]; | |
227 | Float_t * tempoT = new Float_t[fMaxDigit]; | |
228 | Int_t * tempoId = new Int_t[fMaxDigit]; | |
70a93198 | 229 | |
230 | Int_t index ; | |
231 | for ( index = 0 ; index < fMulDigit ; index++ ){ | |
e52475ed | 232 | tempo[index] = fDigitsList[index] ; |
233 | tempoE[index] = fEnergyList[index] ; | |
85c60a8e | 234 | tempoT[index] = fTimeList[index] ; |
e52475ed | 235 | tempoId[index] = fAbsIdList[index] ; |
70a93198 | 236 | } |
237 | ||
94478418 | 238 | delete [] fDigitsList ; |
70a93198 | 239 | delete [] fEnergyList ; |
85c60a8e | 240 | delete [] fTimeList ; |
e52475ed | 241 | delete [] fAbsIdList ; |
e52475ed | 242 | |
94478418 | 243 | fDigitsList = tempo; |
244 | fEnergyList = tempoE; | |
6ce4e2ed | 245 | fTimeList = tempoT; |
246 | fAbsIdList = tempoId; | |
70a93198 | 247 | } // if |
248 | ||
249 | fDigitsList[fMulDigit] = digit.GetIndexInList() ; | |
6ce4e2ed | 250 | fEnergyList[fMulDigit] = energy ; |
94478418 | 251 | fTimeList[fMulDigit] = digit.GetTimeR() ; |
e52475ed | 252 | fAbsIdList[fMulDigit] = digit.GetId(); |
70a93198 | 253 | fMulDigit++ ; |
6ce4e2ed | 254 | fAmp += energy ; |
25bb3dcb | 255 | |
256 | if(shared) fSharedCluster = kTRUE; | |
257 | ||
258 | //GCB, May-2010, setting moved to EvalAll method, set the super module number for the largest energy digit position. | |
9aa6a5f6 | 259 | //JLK 10-Oct-2007 this hasn't been filled before because it was in |
260 | //the wrong place in previous versions. | |
261 | //Now we evaluate it only if the supermodulenumber for this recpoint | |
262 | //has not yet been set (or is the 0th one) | |
25bb3dcb | 263 | //if(fSuperModuleNumber == 0) |
264 | //fSuperModuleNumber = fGeomPtr->GetSuperModuleNumber(digit.GetId()); | |
9aa6a5f6 | 265 | |
70a93198 | 266 | } |
267 | //____________________________________________________________________________ | |
268 | Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const | |
269 | { | |
270 | // Tells if (true) or not (false) two digits are neighbours | |
271 | // A neighbour is defined as being two digits which share a corner | |
25bb3dcb | 272 | // ONLY USED IN CASE OF UNFOLDING |
273 | ||
47583f1b | 274 | Bool_t areNeighbours = kFALSE ; |
275 | Int_t nSupMod=0, nModule=0, nIphi=0, nIeta=0; | |
276 | Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0; | |
277 | Int_t relid1[2] , relid2[2] ; // ieta, iphi | |
278 | Int_t rowdiff=0, coldiff=0; | |
70a93198 | 279 | |
e52475ed | 280 | areNeighbours = kFALSE ; |
281 | ||
2bb3725c | 282 | fGeomPtr->GetCellIndex(digit1->GetId(), nSupMod,nModule,nIphi,nIeta); |
283 | fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, relid1[0],relid1[1]); | |
e52475ed | 284 | |
2bb3725c | 285 | fGeomPtr->GetCellIndex(digit2->GetId(), nSupMod1,nModule1,nIphi1,nIeta1); |
286 | fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, relid2[0],relid2[1]); | |
70a93198 | 287 | |
25bb3dcb | 288 | // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2-1 |
289 | // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0 | |
290 | if(fSharedCluster){ | |
291 | if(nSupMod1%2) relid1[1]+=AliEMCALGeoParams::fgkEMCALCols; | |
292 | else relid2[1]+=AliEMCALGeoParams::fgkEMCALCols; | |
293 | } | |
294 | ||
e52475ed | 295 | rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ; |
296 | coldiff = TMath::Abs( relid1[1] - relid2[1] ) ; | |
70a93198 | 297 | |
298 | if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) | |
e52475ed | 299 | areNeighbours = kTRUE ; |
ab48128d | 300 | |
70a93198 | 301 | return areNeighbours; |
302 | } | |
303 | ||
304 | //____________________________________________________________________________ | |
305 | Int_t AliEMCALRecPoint::Compare(const TObject * obj) const | |
306 | { | |
307 | // Compares two RecPoints according to their position in the EMCAL modules | |
308 | ||
25bb3dcb | 309 | Float_t delta = 1 ; //Width of "Sorting row". |
310 | ||
53e430a3 | 311 | Int_t rv = 2 ; |
70a93198 | 312 | |
313 | AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ; | |
314 | ||
315 | TVector3 locpos1; | |
316 | GetLocalPosition(locpos1); | |
317 | TVector3 locpos2; | |
318 | clu->GetLocalPosition(locpos2); | |
319 | ||
9848d950 | 320 | Int_t rowdif = (Int_t)(TMath::Ceil(locpos1.X()/delta)-TMath::Ceil(locpos2.X()/delta)) ; |
70a93198 | 321 | if (rowdif> 0) |
322 | rv = 1 ; | |
323 | else if(rowdif < 0) | |
324 | rv = -1 ; | |
325 | else if(locpos1.Y()>locpos2.Y()) | |
326 | rv = -1 ; | |
327 | else | |
328 | rv = 1 ; | |
329 | ||
330 | return rv ; | |
ab48128d | 331 | } |
332 | ||
25bb3dcb | 333 | // GCB, May-2010, Method not used, just comment it but remove? |
ab48128d | 334 | //____________________________________________________________________________ |
25bb3dcb | 335 | //Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py) |
336 | //{ | |
337 | // // Compute distance from point px,py to a AliEMCALRecPoint considered as a Tmarker | |
338 | // // Compute the closest distance of approach from point px,py to this marker. | |
339 | // // The distance is computed in pixels units. | |
340 | // // HG Still need to update -> Not sure what this should achieve | |
341 | // | |
342 | // TVector3 pos(0.,0.,0.) ; | |
343 | // GetLocalPosition(pos) ; | |
344 | // Float_t x = pos.X() ; | |
345 | // Float_t y = pos.Y() ; | |
346 | // const Int_t kMaxDiff = 10; | |
347 | // Int_t pxm = gPad->XtoAbsPixel(x); | |
348 | // Int_t pym = gPad->YtoAbsPixel(y); | |
349 | // Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym); | |
350 | // | |
351 | // if (dist > kMaxDiff) return 9999; | |
352 | // return dist; | |
353 | //} | |
ab48128d | 354 | |
355 | //___________________________________________________________________________ | |
356 | void AliEMCALRecPoint::Draw(Option_t *option) | |
357 | { | |
358 | // Draw this AliEMCALRecPoint with its current attributes | |
359 | ||
360 | AppendPad(option); | |
361 | } | |
362 | ||
25bb3dcb | 363 | // GCB, May-2010, Method not used, just comment it but remove? |
ab48128d | 364 | //______________________________________________________________________________ |
25bb3dcb | 365 | //void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t) |
366 | //{ | |
367 | // // Execute action corresponding to one event | |
368 | // // This member function is called when a AliEMCALRecPoint is clicked with the locator | |
369 | // // | |
370 | // // If Left button is clicked on AliEMCALRecPoint, the digits are switched on | |
371 | // // and switched off when the mouse button is released. | |
372 | // | |
373 | // // static Int_t pxold, pyold; | |
374 | // | |
375 | // /* static TGraph * digitgraph = 0 ; | |
376 | // static TPaveText* clustertext = 0 ; | |
377 | // | |
378 | // if (!gPad->IsEditable()) return; | |
379 | // | |
380 | // switch (event) { | |
381 | // | |
382 | // | |
383 | // case kButton1Down:{ | |
384 | // AliEMCALDigit * digit ; | |
385 | // | |
386 | // Int_t iDigit; | |
387 | // Int_t relid[2] ; | |
388 | // | |
389 | // const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ; | |
390 | // Float_t * xi = new Float_t [kMulDigit] ; | |
391 | // Float_t * zi = new Float_t [kMulDigit] ; | |
392 | // | |
393 | // for(iDigit = 0; iDigit < kMulDigit; iDigit++) { | |
394 | // Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code"); | |
395 | // digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]); | |
396 | // fGeomPtr->AbsToRelNumbering(digit->GetId(), relid) ; | |
397 | // fGeomPtr->PosInAlice(relid, xi[iDigit], zi[iDigit]) ; | |
398 | // } | |
399 | // | |
400 | // if (!digitgraph) { | |
401 | // digitgraph = new TGraph(fMulDigit,xi,zi); | |
402 | // digitgraph-> SetMarkerStyle(5) ; | |
403 | // digitgraph-> SetMarkerSize(1.) ; | |
404 | // digitgraph-> SetMarkerColor(1) ; | |
405 | // digitgraph-> Draw("P") ; | |
406 | // } | |
407 | // if (!clustertext) { | |
408 | // | |
409 | // TVector3 pos(0.,0.,0.) ; | |
410 | // GetLocalPosition(pos) ; | |
411 | // clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ; | |
412 | // Text_t line1[40] ; | |
413 | // Text_t line2[40] ; | |
414 | // sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ; | |
415 | // sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ; | |
416 | // clustertext ->AddText(line1) ; | |
417 | // clustertext ->AddText(line2) ; | |
418 | // clustertext ->Draw(""); | |
419 | // } | |
420 | // gPad->Update() ; | |
421 | // Print("") ; | |
422 | // delete[] xi ; | |
423 | // delete[] zi ; | |
424 | // } | |
425 | // | |
426 | //break; | |
427 | // | |
428 | // case kButton1Up: | |
429 | // if (digitgraph) { | |
430 | // delete digitgraph ; | |
431 | // digitgraph = 0 ; | |
432 | // } | |
433 | // if (clustertext) { | |
434 | // delete clustertext ; | |
435 | // clustertext = 0 ; | |
436 | // } | |
437 | // | |
438 | // break; | |
439 | // | |
440 | // }*/ | |
441 | //} | |
ab48128d | 442 | |
70a93198 | 443 | //____________________________________________________________________________ |
444 | void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits) | |
445 | { | |
25bb3dcb | 446 | // Evaluates cluster parameters |
447 | ||
448 | // First calculate the index of digit with maximum amplitude and get | |
449 | // the supermodule number where it sits. | |
6ce4e2ed | 450 | |
25bb3dcb | 451 | fDigitIndMax = GetMaximalEnergyIndex(); |
452 | fSuperModuleNumber = fGeomPtr->GetSuperModuleNumber(GetAbsIdMaxDigit()); | |
453 | ||
454 | //Evaluate global and local position | |
455 | EvalGlobalPosition(logWeight, digits) ; | |
70a93198 | 456 | EvalLocalPosition(logWeight, digits) ; |
25bb3dcb | 457 | |
458 | //Evaluate shower parameters | |
70a93198 | 459 | EvalElipsAxis(logWeight, digits) ; |
460 | EvalDispersion(logWeight, digits) ; | |
25bb3dcb | 461 | |
4635df1f | 462 | //EvalCoreEnergy(logWeight, digits); |
70a93198 | 463 | EvalTime(digits) ; |
87cdc3be | 464 | EvalPrimaries(digits) ; |
465 | EvalParents(digits); | |
40164976 | 466 | |
9aa6a5f6 | 467 | //Called last because it sets the global position of the cluster? |
468 | EvalLocal2TrackingCSTransform(); | |
469 | ||
70a93198 | 470 | } |
471 | ||
472 | //____________________________________________________________________________ | |
473 | void AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits) | |
474 | { | |
475 | // Calculates the dispersion of the shower at the origin of the RecPoint | |
1d46d1f6 | 476 | // in cell units - Nov 16,2006 |
70a93198 | 477 | |
1d46d1f6 | 478 | Double_t d = 0., wtot = 0., w = 0.; |
af5bdd85 | 479 | Int_t iDigit=0, nstat=0; |
53e430a3 | 480 | AliEMCALDigit * digit=0; |
25bb3dcb | 481 | |
1d46d1f6 | 482 | // Calculates the dispersion in cell units |
483 | Double_t etai, phii, etaMean=0.0, phiMean=0.0; | |
2bb3725c | 484 | int nSupMod=0, nModule=0, nIphi=0, nIeta=0; |
1d46d1f6 | 485 | int iphi=0, ieta=0; |
486 | // Calculate mean values | |
70a93198 | 487 | for(iDigit=0; iDigit < fMulDigit; iDigit++) { |
488 | digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ; | |
1963b290 | 489 | |
1d46d1f6 | 490 | if (fAmp>0 && fEnergyList[iDigit]>0) { |
2bb3725c | 491 | fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta); |
492 | fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta); | |
25bb3dcb | 493 | |
494 | // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2 | |
495 | // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0 | |
496 | if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols; | |
497 | ||
1d46d1f6 | 498 | etai=(Double_t)ieta; |
499 | phii=(Double_t)iphi; | |
500 | w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ; | |
501 | ||
502 | if(w>0.0) { | |
503 | phiMean += phii*w; | |
504 | etaMean += etai*w; | |
505 | wtot += w; | |
506 | } | |
507 | } | |
508 | } | |
509 | if (wtot>0) { | |
510 | phiMean /= wtot ; | |
511 | etaMean /= wtot ; | |
512 | } else AliError(Form("Wrong weight %f\n", wtot)); | |
70a93198 | 513 | |
1d46d1f6 | 514 | // Calculate dispersion |
515 | for(iDigit=0; iDigit < fMulDigit; iDigit++) { | |
516 | digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ; | |
517 | ||
518 | if (fAmp>0 && fEnergyList[iDigit]>0) { | |
2bb3725c | 519 | fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta); |
520 | fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta); | |
25bb3dcb | 521 | |
522 | // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2 | |
523 | // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0 | |
524 | if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols; | |
525 | ||
1d46d1f6 | 526 | etai=(Double_t)ieta; |
527 | phii=(Double_t)iphi; | |
528 | w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ; | |
529 | ||
530 | if(w>0.0) { | |
531 | nstat++; | |
532 | d += w*((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean)); | |
e52475ed | 533 | } |
534 | } | |
ab48128d | 535 | } |
70a93198 | 536 | |
e52475ed | 537 | if ( wtot > 0 && nstat>1) d /= wtot ; |
538 | else d = 0. ; | |
70a93198 | 539 | |
540 | fDispersion = TMath::Sqrt(d) ; | |
25bb3dcb | 541 | //printf("AliEMCALRecPoint::EvalDispersion() : Dispersion %f \n",fDispersion); |
ab48128d | 542 | } |
70a93198 | 543 | |
40164976 | 544 | //____________________________________________________________________________ |
545 | void AliEMCALRecPoint::EvalDistanceToBadChannels(AliCaloCalibPedestal* caloped) | |
546 | { | |
547 | //For each EMC rec. point set the distance to the nearest bad channel. | |
548 | //AliInfo(Form("%d bad channel(s) found.\n", caloped->GetDeadTowerCount())); | |
2ebdefe6 | 549 | //It is done in cell units and not in global or local position as before (Sept 2010) |
40164976 | 550 | |
551 | if(!caloped->GetDeadTowerCount()) return; | |
25bb3dcb | 552 | |
553 | //Get channels map of the supermodule where the cluster is. | |
554 | TH2D* hMap = caloped->GetDeadMap(fSuperModuleNumber); | |
40164976 | 555 | |
2ebdefe6 | 556 | Int_t dRrow, dReta; |
557 | Float_t minDist = 10000.; | |
558 | Float_t dist = 0.; | |
559 | Int_t nSupMod, nModule; | |
560 | Int_t nIphi, nIeta; | |
561 | Int_t iphi, ieta; | |
562 | fDigitIndMax = GetMaximalEnergyIndex(); | |
563 | fGeomPtr->GetCellIndex(fAbsIdList[fDigitIndMax], nSupMod,nModule,nIphi,nIeta); | |
564 | fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta); | |
565 | ||
566 | // TVector3 dR; | |
567 | // TVector3 cellpos; | |
568 | // Float_t minDist = 100000; | |
569 | // Float_t dist = 0; | |
570 | // Int_t absId = -1; | |
571 | ||
572 | //Loop on tower status map | |
40164976 | 573 | for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){ |
574 | for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){ | |
575 | //Check if tower is bad. | |
576 | if(hMap->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue; | |
2ebdefe6 | 577 | //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow); |
578 | ||
579 | dRrow=TMath::Abs(irow-iphi); | |
580 | dReta=TMath::Abs(icol-ieta); | |
581 | dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta); | |
25bb3dcb | 582 | if(dist < minDist) minDist = dist; |
2ebdefe6 | 583 | |
584 | // //Tower is bad, get the absId of the index. | |
585 | // absId = fGeomPtr->GetAbsCellIdFromCellIndexes(fSuperModuleNumber, irow, icol); | |
586 | // | |
587 | // //Get the position of this tower. | |
588 | // | |
589 | // //Calculate the distance in local coordinates | |
590 | // //fGeomPtr->RelPosCellInSModule(absId,cellpos); | |
591 | // //Calculate distance between this tower and cluster, set if is smaller than previous. | |
592 | // //dR = cellpos-fLocPos; | |
593 | // | |
594 | // //Calculate the distance in global coordinates | |
595 | // fGeomPtr->GetGlobal(absId,cellpos); | |
596 | // //Calculate distance between this tower and cluster, set if it is smaller than previous. | |
597 | // dR = cellpos-fGlobPos; | |
598 | // | |
599 | // dist = dR.Mag(); | |
600 | // if(dist < minDist) minDist = dist; | |
601 | ||
40164976 | 602 | } |
89dc5acd | 603 | } |
2ebdefe6 | 604 | |
25bb3dcb | 605 | //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module |
606 | if (fSharedCluster) { | |
607 | TH2D* hMap2 = 0; | |
608 | Int_t nSupMod2 = -1; | |
609 | ||
610 | //The only possible combinations are (0,1), (2,3) ... (10,11) | |
611 | if(fSuperModuleNumber%2) nSupMod2 = fSuperModuleNumber-1; | |
612 | else nSupMod2 = fSuperModuleNumber+1; | |
613 | hMap2 = caloped->GetDeadMap(nSupMod2); | |
2ebdefe6 | 614 | |
25bb3dcb | 615 | //Loop on tower status map of second super module |
616 | for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){ | |
617 | for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){ | |
618 | //Check if tower is bad. | |
619 | if(hMap2->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue; | |
620 | //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow); | |
2ebdefe6 | 621 | |
622 | dRrow=TMath::Abs(irow-iphi); | |
623 | ||
624 | if(fSuperModuleNumber%2) { | |
625 | dReta=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+ieta)); | |
626 | } | |
627 | else { | |
628 | dReta=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-ieta); | |
629 | } | |
630 | ||
631 | dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta); | |
632 | if(dist < minDist) minDist = dist; | |
633 | ||
634 | // | |
635 | // //Tower is bad, get the absId of the index. | |
636 | // absId = fGeomPtr->GetAbsCellIdFromCellIndexes(nSupMod2, irow, icol); | |
637 | // | |
638 | // //Get the position of this tower. | |
639 | // | |
640 | // //Calculate the distance in global coordinates | |
641 | // fGeomPtr->GetGlobal(absId,cellpos); | |
642 | // //Calculate distance between this tower and cluster, set if it is smaller than previous. | |
643 | // dR = cellpos-fGlobPos; | |
644 | // | |
645 | // dist = dR.Mag(); | |
646 | // if(dist < minDist) minDist = dist; | |
25bb3dcb | 647 | } |
648 | } | |
649 | ||
650 | }// shared cluster in 2 SuperModules | |
2ebdefe6 | 651 | |
40164976 | 652 | fDistToBadTower = minDist; |
25bb3dcb | 653 | //printf("AliEMCALRecPoint::EvalDistanceToBadChannel() - Distance to Bad is %f cm, shared cluster? %d \n",fDistToBadTower,fSharedCluster); |
40164976 | 654 | } |
655 | ||
656 | ||
ab48128d | 657 | //____________________________________________________________________________ |
70a93198 | 658 | void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits) |
88cb7938 | 659 | { |
25bb3dcb | 660 | // Calculates the center of gravity in the local EMCAL-module coordinates |
661 | // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing | |
662 | ||
53e430a3 | 663 | AliEMCALDigit * digit=0; |
25bb3dcb | 664 | Int_t i=0, nstat=0; |
665 | ||
47583f1b | 666 | Double_t dist = TmaxInCm(Double_t(fAmp)); |
25bb3dcb | 667 | //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it? |
668 | ||
669 | Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.; | |
670 | ||
671 | //printf(" dist : %f e : %f \n", dist, fAmp); | |
672 | for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) { | |
673 | digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ; | |
674 | ||
2ebdefe6 | 675 | if(!digit) { |
676 | AliError("No Digit!!"); | |
677 | continue; | |
678 | } | |
679 | ||
25bb3dcb | 680 | //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]); |
681 | fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]); | |
682 | ||
683 | //Temporal patch, due to mapping problem, need to swap "y" in one of the 2 SM, although no effect in position calculation. GCB 05/2010 | |
684 | if(fSharedCluster && fSuperModuleNumber != fGeomPtr->GetSuperModuleNumber(digit->GetId())) xyzi[1]*=-1; | |
685 | ||
686 | //printf("EvalLocalPosition Cell: Id %i, SM %i : dist %f Local x,y,z %f %f %f \n", | |
687 | // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()), dist, xyzi[0], xyzi[1], xyzi[2]); | |
688 | ||
689 | if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp )); | |
690 | else w = fEnergyList[iDigit]; // just energy | |
691 | ||
692 | if(w>0.0) { | |
693 | wtot += w ; | |
694 | nstat++; | |
695 | for(i=0; i<3; i++ ) { | |
696 | clXYZ[i] += (w*xyzi[i]); | |
697 | clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]); | |
698 | } | |
699 | } | |
700 | } | |
701 | // cout << " wtot " << wtot << endl; | |
702 | if ( wtot > 0 ) { | |
703 | // xRMS = TMath::Sqrt(x2m - xMean*xMean); | |
704 | for(i=0; i<3; i++ ) { | |
705 | clXYZ[i] /= wtot; | |
706 | if(nstat>1) { | |
707 | clRmsXYZ[i] /= (wtot*wtot); | |
708 | clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i]; | |
709 | if(clRmsXYZ[i] > 0.0) { | |
710 | clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]); | |
711 | } else clRmsXYZ[i] = 0; | |
712 | } else clRmsXYZ[i] = 0; | |
713 | } | |
714 | } else { | |
715 | for(i=0; i<3; i++ ) { | |
716 | clXYZ[i] = clRmsXYZ[i] = -1.; | |
717 | } | |
718 | } | |
719 | // clRmsXYZ[i] ?? | |
720 | ||
721 | // // Cluster of one single digit, smear the position to avoid discrete position | |
722 | // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm. | |
723 | // // Rndm generates a number in ]0,1] | |
724 | // if (fMulDigit==1) { | |
725 | // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm()); | |
726 | // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm()); | |
727 | // } | |
728 | ||
729 | //Set position in local vector | |
730 | fLocPos.SetX(clXYZ[0]); | |
731 | fLocPos.SetY(clXYZ[1]); | |
732 | fLocPos.SetZ(clXYZ[2]); | |
733 | ||
734 | if (gDebug==2) | |
735 | printf("EvalLocalPosition Cluster: Local (x,y,z) = (%f,%f,%f) \n", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; | |
736 | ||
737 | } | |
738 | ||
739 | ||
740 | //____________________________________________________________________________ | |
741 | void AliEMCALRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits) | |
742 | { | |
743 | // Calculates the center of gravity in the global ALICE coordinates | |
e52475ed | 744 | // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing |
70a93198 | 745 | |
53e430a3 | 746 | AliEMCALDigit * digit=0; |
25bb3dcb | 747 | Int_t i=0, nstat=0; |
748 | ||
47583f1b | 749 | Double_t dist = TmaxInCm(Double_t(fAmp)); |
25bb3dcb | 750 | //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it? |
751 | ||
752 | Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, lxyzi[3], xyzi[3], wtot=0., w=0.; | |
70a93198 | 753 | |
1ae500a2 | 754 | //printf(" dist : %f e : %f \n", dist, fAmp); |
e52475ed | 755 | for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) { |
70a93198 | 756 | digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ; |
e52475ed | 757 | |
2ebdefe6 | 758 | if(!digit) { |
759 | AliError("No Digit!!"); | |
760 | continue; | |
761 | } | |
762 | ||
25bb3dcb | 763 | //Get the local coordinates of the cell |
764 | //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, lxyzi[0], lxyzi[1], lxyzi[2]); | |
765 | fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, lxyzi[0], lxyzi[1], lxyzi[2]); | |
766 | ||
767 | //Now get the global coordinate | |
768 | fGeomPtr->GetGlobal(lxyzi,xyzi, fGeomPtr->GetSuperModuleNumber(digit->GetId())); | |
769 | //TVector3 pos(xyzi[0], xyzi[1], xyzi[2]); | |
770 | //printf("EvalGlobalPosition Cell: Id %i, SM %i : dist %f Local (x,y,z) = (%f %f %f), eta %f, phi%f \n", | |
771 | // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()),dist, xyzi[0], xyzi[1], xyzi[2],pos.Eta(),pos.Phi()*TMath::RadToDeg()); | |
772 | ||
af5bdd85 | 773 | if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp )); |
774 | else w = fEnergyList[iDigit]; // just energy | |
ab48128d | 775 | |
e52475ed | 776 | if(w>0.0) { |
777 | wtot += w ; | |
778 | nstat++; | |
779 | for(i=0; i<3; i++ ) { | |
780 | clXYZ[i] += (w*xyzi[i]); | |
781 | clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]); | |
782 | } | |
783 | } | |
784 | } | |
785 | // cout << " wtot " << wtot << endl; | |
70a93198 | 786 | if ( wtot > 0 ) { |
e52475ed | 787 | // xRMS = TMath::Sqrt(x2m - xMean*xMean); |
788 | for(i=0; i<3; i++ ) { | |
789 | clXYZ[i] /= wtot; | |
790 | if(nstat>1) { | |
791 | clRmsXYZ[i] /= (wtot*wtot); | |
792 | clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i]; | |
793 | if(clRmsXYZ[i] > 0.0) { | |
794 | clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]); | |
795 | } else clRmsXYZ[i] = 0; | |
796 | } else clRmsXYZ[i] = 0; | |
797 | } | |
70a93198 | 798 | } else { |
e52475ed | 799 | for(i=0; i<3; i++ ) { |
800 | clXYZ[i] = clRmsXYZ[i] = -1.; | |
801 | } | |
70a93198 | 802 | } |
e52475ed | 803 | // clRmsXYZ[i] ?? |
25bb3dcb | 804 | |
805 | // // Cluster of one single digit, smear the position to avoid discrete position | |
806 | // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm. | |
807 | // // Rndm generates a number in ]0,1] | |
808 | // if (fMulDigit==1) { | |
809 | // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm()); | |
810 | // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm()); | |
811 | // } | |
812 | ||
813 | //Set position in global vector | |
814 | fGlobPos.SetX(clXYZ[0]); | |
815 | fGlobPos.SetY(clXYZ[1]); | |
816 | fGlobPos.SetZ(clXYZ[2]); | |
817 | ||
818 | if (gDebug==2) | |
819 | printf("EvalGlobalPosition Cluster: (x ,y ,z) = (%f,%f,%f), eta %f,phi %f\n", | |
820 | fGlobPos.X(), fGlobPos.Y(), fGlobPos.Z(),fGlobPos.Eta(),fGlobPos.Phi()*TMath::RadToDeg()) ; | |
ab48128d | 821 | } |
822 | ||
1ae500a2 | 823 | //____________________________________________________________________________ |
824 | void AliEMCALRecPoint::EvalLocalPositionFit(Double_t deff, Double_t logWeight, | |
825 | Double_t phiSlope, TClonesArray * digits) | |
826 | { | |
47583f1b | 827 | // Evaluates local position of clusters in SM |
828 | ||
829 | Double_t ycorr=0; | |
830 | AliEMCALDigit *digit=0; | |
25bb3dcb | 831 | Int_t i=0, nstat=0; |
1ae500a2 | 832 | Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.; |
833 | ||
47583f1b | 834 | Double_t dist = TmaxInCm(Double_t(fAmp)); |
25bb3dcb | 835 | //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it? |
836 | ||
1ae500a2 | 837 | for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) { |
838 | digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ; | |
a51e676d | 839 | if(digit){ |
840 | dist = deff; | |
841 | //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]); | |
842 | fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]); | |
843 | ||
844 | if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp )); | |
845 | else w = fEnergyList[iDigit]; // just energy | |
846 | ||
847 | if(w>0.0) { | |
848 | wtot += w ; | |
849 | nstat++; | |
850 | for(i=0; i<3; i++ ) { | |
851 | clXYZ[i] += (w*xyzi[i]); | |
852 | clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]); | |
853 | } | |
1ae500a2 | 854 | } |
a51e676d | 855 | }else AliError("Digit null"); |
856 | }//loop | |
1ae500a2 | 857 | // cout << " wtot " << wtot << endl; |
858 | if ( wtot > 0 ) { | |
859 | // xRMS = TMath::Sqrt(x2m - xMean*xMean); | |
860 | for(i=0; i<3; i++ ) { | |
861 | clXYZ[i] /= wtot; | |
862 | if(nstat>1) { | |
863 | clRmsXYZ[i] /= (wtot*wtot); | |
864 | clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i]; | |
865 | if(clRmsXYZ[i] > 0.0) { | |
866 | clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]); | |
867 | } else clRmsXYZ[i] = 0; | |
868 | } else clRmsXYZ[i] = 0; | |
869 | } | |
870 | } else { | |
871 | for(i=0; i<3; i++ ) { | |
872 | clXYZ[i] = clRmsXYZ[i] = -1.; | |
873 | } | |
874 | } | |
875 | // clRmsXYZ[i] ?? | |
876 | if(phiSlope != 0.0 && logWeight > 0.0 && wtot) { | |
877 | // Correction in phi direction (y - coords here); Aug 16; | |
878 | // May be put to global level or seperate method | |
879 | ycorr = clXYZ[1] * (1. + phiSlope); | |
880 | //printf(" y %f : ycorr %f : slope %f \n", clXYZ[1], ycorr, phiSlope); | |
881 | clXYZ[1] = ycorr; | |
882 | } | |
25bb3dcb | 883 | |
1ae500a2 | 884 | fLocPos.SetX(clXYZ[0]); |
885 | fLocPos.SetY(clXYZ[1]); | |
886 | fLocPos.SetZ(clXYZ[2]); | |
887 | ||
888 | // if (gDebug==2) | |
889 | // printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; | |
1ae500a2 | 890 | } |
891 | ||
9aa6a5f6 | 892 | //_____________________________________________________________________________ |
1ae500a2 | 893 | Bool_t AliEMCALRecPoint::EvalLocalPosition2(TClonesArray * digits, TArrayD &ed) |
894 | { | |
895 | // Evaluated local position of rec.point using digits | |
896 | // and parametrisation of w0 and deff | |
897 | //printf(" <I> AliEMCALRecPoint::EvalLocalPosition2() \n"); | |
898 | return AliEMCALRecPoint::EvalLocalPositionFromDigits(digits, ed, fLocPos); | |
899 | } | |
900 | ||
9aa6a5f6 | 901 | //_____________________________________________________________________________ |
1ae500a2 | 902 | Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(TClonesArray *digits, TArrayD &ed, TVector3 &locPos) |
903 | { | |
904 | // Used when digits should be recalibrated | |
47583f1b | 905 | Double_t deff=0, w0=0, esum=0; |
906 | Int_t iDigit=0; | |
907 | // AliEMCALDigit *digit; | |
1ae500a2 | 908 | |
909 | if(ed.GetSize() && (digits->GetEntries()!=ed.GetSize())) return kFALSE; | |
910 | ||
911 | // Calculate sum energy of digits | |
912 | esum = 0.0; | |
913 | for(iDigit=0; iDigit<ed.GetSize(); iDigit++) esum += ed[iDigit]; | |
914 | ||
915 | GetDeffW0(esum, deff, w0); | |
916 | ||
917 | return EvalLocalPositionFromDigits(esum, deff, w0, digits, ed, locPos); | |
918 | } | |
919 | ||
9aa6a5f6 | 920 | //_____________________________________________________________________________ |
921 | Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(const Double_t esum, const Double_t deff, const Double_t w0, TClonesArray *digits, TArrayD &ed, TVector3 &locPos) | |
1ae500a2 | 922 | { |
40164976 | 923 | //Evaluate position of digits in supermodule. |
47583f1b | 924 | AliEMCALDigit *digit=0; |
1ae500a2 | 925 | |
25bb3dcb | 926 | Int_t i=0, nstat=0; |
1ae500a2 | 927 | Double_t clXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.; |
25bb3dcb | 928 | //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it? |
929 | ||
9aa6a5f6 | 930 | // Get pointer to EMCAL geometry |
931 | // (can't use fGeomPtr in static method) | |
932 | AliEMCALGeometry* geo = AliEMCALGeometry::GetInstance(); | |
933 | ||
1ae500a2 | 934 | for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) { |
935 | digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)); | |
a51e676d | 936 | if(digit){ |
937 | //geo->RelPosCellInSModule(digit->GetId(), idMax, deff, xyzi[0], xyzi[1], xyzi[2]); | |
938 | geo->RelPosCellInSModule(digit->GetId(), deff, xyzi[0], xyzi[1], xyzi[2]); | |
939 | ||
940 | if(w0 > 0.0) w = TMath::Max( 0., w0 + TMath::Log(ed[iDigit] / esum)); | |
941 | else w = ed[iDigit]; // just energy | |
942 | ||
943 | if(w>0.0) { | |
944 | wtot += w ; | |
945 | nstat++; | |
946 | for(i=0; i<3; i++ ) { | |
947 | clXYZ[i] += (w*xyzi[i]); | |
948 | } | |
1ae500a2 | 949 | } |
a51e676d | 950 | }else AliError("Digit null"); |
951 | }//loop | |
1ae500a2 | 952 | // cout << " wtot " << wtot << endl; |
953 | if (wtot > 0) { | |
954 | for(i=0; i<3; i++ ) { | |
955 | clXYZ[i] /= wtot; | |
956 | } | |
957 | locPos.SetX(clXYZ[0]); | |
958 | locPos.SetY(clXYZ[1]); | |
959 | locPos.SetZ(clXYZ[2]); | |
960 | return kTRUE; | |
961 | } else { | |
962 | return kFALSE; | |
963 | } | |
964 | ||
965 | } | |
966 | ||
9aa6a5f6 | 967 | //_____________________________________________________________________________ |
1ae500a2 | 968 | void AliEMCALRecPoint::GetDeffW0(const Double_t esum , Double_t &deff, Double_t &w0) |
969 | { | |
970 | // | |
971 | // Aug 31, 2001 | |
972 | // Applied for simulation data with threshold 3 adc | |
973 | // Calculate efective distance (deff) and weigh parameter (w0) | |
974 | // for coordinate calculation; 0.5 GeV < esum <100 GeV. | |
975 | // Look to: http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/CALIB/GEOMCORR/deffandW0VaEgamma_2.gif | |
976 | // | |
47583f1b | 977 | Double_t e=0.0; |
37890aaf | 978 | const Double_t kdp0=9.25147, kdp1=1.16700; // Hard coded now |
979 | const Double_t kwp0=4.83713, kwp1=-2.77970e-01, kwp2 = 4.41116; | |
1ae500a2 | 980 | |
981 | // No extrapolation here | |
982 | e = esum<0.5?0.5:esum; | |
983 | e = e>100.?100.:e; | |
984 | ||
37890aaf | 985 | deff = kdp0 + kdp1*TMath::Log(e); |
986 | w0 = kwp0 / (1. + TMath::Exp(kwp1*(e+kwp2))); | |
1ae500a2 | 987 | //printf("<I> AliEMCALRecPoint::GetDeffW0 esum %5.2f : deff %5.2f : w0 %5.2f \n", esum, deff, w0); |
988 | } | |
e52475ed | 989 | |
70a93198 | 990 | //______________________________________________________________________________ |
991 | void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits) | |
992 | { | |
993 | // This function calculates energy in the core, | |
4635df1f | 994 | // i.e. within a radius rad = fCoreEnergy around the center. Beyond this radius |
70a93198 | 995 | // in accordance with shower profile the energy deposition |
996 | // should be less than 2% | |
1d46d1f6 | 997 | // Unfinished - Nov 15,2006 |
998 | // Distance is calculate in (phi,eta) units | |
70a93198 | 999 | |
53e430a3 | 1000 | AliEMCALDigit * digit = 0 ; |
5dee926e | 1001 | |
53e430a3 | 1002 | Int_t iDigit=0; |
70a93198 | 1003 | |
e52475ed | 1004 | if (!fLocPos.Mag()) { |
70a93198 | 1005 | EvalLocalPosition(logWeight, digits); |
1006 | } | |
1007 | ||
1d46d1f6 | 1008 | Double_t phiPoint = fLocPos.Phi(), etaPoint = fLocPos.Eta(); |
1009 | Double_t eta, phi, distance; | |
70a93198 | 1010 | for(iDigit=0; iDigit < fMulDigit; iDigit++) { |
1011 | digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ; | |
4635df1f | 1012 | |
1d46d1f6 | 1013 | eta = phi = 0.0; |
1014 | fGeomPtr->EtaPhiFromIndex(digit->GetId(),eta, phi) ; | |
1015 | phi = phi * TMath::DegToRad(); | |
70a93198 | 1016 | |
1d46d1f6 | 1017 | distance = TMath::Sqrt((eta-etaPoint)*(eta-etaPoint)+(phi-phiPoint)*(phi-phiPoint)); |
70a93198 | 1018 | if(distance < fCoreRadius) |
1019 | fCoreEnergy += fEnergyList[iDigit] ; | |
1020 | } | |
1021 | ||
1022 | } | |
ab48128d | 1023 | //____________________________________________________________________________ |
70a93198 | 1024 | void AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits) |
ab48128d | 1025 | { |
70a93198 | 1026 | // Calculates the axis of the shower ellipsoid in eta and phi |
1d46d1f6 | 1027 | // in cell units |
ab48128d | 1028 | |
47583f1b | 1029 | TString gn(fGeomPtr->GetName()); |
1d46d1f6 | 1030 | |
1031 | Double_t wtot = 0.; | |
70a93198 | 1032 | Double_t x = 0.; |
1033 | Double_t z = 0.; | |
1034 | Double_t dxx = 0.; | |
1035 | Double_t dzz = 0.; | |
1036 | Double_t dxz = 0.; | |
ab48128d | 1037 | |
1d46d1f6 | 1038 | AliEMCALDigit * digit = 0; |
25bb3dcb | 1039 | |
53e430a3 | 1040 | Double_t etai =0, phii=0, w=0; |
2bb3725c | 1041 | int nSupMod=0, nModule=0, nIphi=0, nIeta=0; |
1d46d1f6 | 1042 | int iphi=0, ieta=0; |
1043 | for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) { | |
70a93198 | 1044 | digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ; |
1d46d1f6 | 1045 | etai = phii = 0.; |
1d46d1f6 | 1046 | // Nov 15,2006 - use cell numbers as coordinates |
1047 | // Copied for shish-kebab geometry, ieta,iphi is cast as double as eta,phi | |
1048 | // We can use the eta,phi(or coordinates) of cell | |
dc77cc84 | 1049 | nSupMod = nModule = nIphi = nIeta = iphi = ieta = 0; |
1d46d1f6 | 1050 | |
dc77cc84 | 1051 | fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta); |
1052 | fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta); | |
25bb3dcb | 1053 | |
1054 | // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2 | |
1055 | // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0 | |
1056 | if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols; | |
1057 | ||
dc77cc84 | 1058 | etai=(Double_t)ieta; |
1059 | phii=(Double_t)iphi; | |
25bb3dcb | 1060 | |
1d46d1f6 | 1061 | w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ; |
1062 | // fAmp summed amplitude of digits, i.e. energy of recpoint | |
1063 | // Gives smaller value of lambda than log weight | |
1064 | // w = fEnergyList[iDigit] / fAmp; // Nov 16, 2006 - try just energy | |
ff1e7e2f | 1065 | |
70a93198 | 1066 | dxx += w * etai * etai ; |
1067 | x += w * etai ; | |
1068 | dzz += w * phii * phii ; | |
1069 | z += w * phii ; | |
1963b290 | 1070 | |
ff1e7e2f | 1071 | dxz += w * etai * phii ; |
1963b290 | 1072 | |
70a93198 | 1073 | wtot += w ; |
1074 | } | |
ff1e7e2f | 1075 | |
70a93198 | 1076 | if ( wtot > 0 ) { |
1077 | dxx /= wtot ; | |
1078 | x /= wtot ; | |
1079 | dxx -= x * x ; | |
1080 | dzz /= wtot ; | |
1081 | z /= wtot ; | |
1082 | dzz -= z * z ; | |
1083 | dxz /= wtot ; | |
1084 | dxz -= x * z ; | |
ab48128d | 1085 | |
70a93198 | 1086 | fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ; |
1087 | if(fLambda[0] > 0) | |
1088 | fLambda[0] = TMath::Sqrt(fLambda[0]) ; | |
1089 | else | |
1090 | fLambda[0] = 0; | |
1091 | ||
1092 | fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ; | |
ff1e7e2f | 1093 | |
70a93198 | 1094 | if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda. |
1095 | fLambda[1] = TMath::Sqrt(fLambda[1]) ; | |
1096 | else | |
1097 | fLambda[1]= 0. ; | |
1098 | } else { | |
1099 | fLambda[0]= 0. ; | |
1100 | fLambda[1]= 0. ; | |
ab48128d | 1101 | } |
ff1e7e2f | 1102 | |
25bb3dcb | 1103 | //printf("AliEMCALRecPoint::EvalElipsAxis() lambdas = %f,%f \n", fLambda[0],fLambda[1]) ; |
ff1e7e2f | 1104 | |
ab48128d | 1105 | } |
1106 | ||
1107 | //______________________________________________________________________________ | |
1108 | void AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits) | |
1109 | { | |
af5bdd85 | 1110 | // Constructs the list of primary particles (tracks) which |
1111 | // have contributed to this RecPoint and calculate deposited energy | |
1112 | // for each track | |
ab48128d | 1113 | |
53e430a3 | 1114 | AliEMCALDigit * digit =0; |
af5bdd85 | 1115 | Int_t * primArray = new Int_t[fMaxTrack] ; |
7f3091e7 | 1116 | memset(primArray,-1,sizeof(Int_t)*fMaxTrack); |
af5bdd85 | 1117 | Float_t * dEPrimArray = new Float_t[fMaxTrack] ; |
7f3091e7 | 1118 | memset(dEPrimArray,-1,sizeof(Int_t)*fMaxTrack); |
2ebdefe6 | 1119 | |
ab48128d | 1120 | Int_t index ; |
1121 | for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits | |
1122 | digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; | |
2ebdefe6 | 1123 | if(!digit) { |
1124 | AliError("No Digit!!"); | |
1125 | continue; | |
1126 | } | |
1127 | ||
ab48128d | 1128 | Int_t nprimaries = digit->GetNprimary() ; |
5c0368b8 | 1129 | if ( nprimaries == 0 ) continue ; |
ab48128d | 1130 | Int_t jndex ; |
1131 | for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit | |
1132 | if ( fMulTrack > fMaxTrack ) { | |
2ebdefe6 | 1133 | fMulTrack = fMaxTrack ; |
1134 | Error("EvalPrimaries", "increase fMaxTrack ") ; | |
1135 | break ; | |
ab48128d | 1136 | } |
af5bdd85 | 1137 | Int_t newPrimary = digit->GetPrimary(jndex+1); |
1138 | Float_t dEPrimary = digit->GetDEPrimary(jndex+1); | |
ab48128d | 1139 | Int_t kndex ; |
1140 | Bool_t already = kFALSE ; | |
1141 | for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored | |
2ebdefe6 | 1142 | if ( newPrimary == primArray[kndex] ){ |
1143 | already = kTRUE ; | |
1144 | dEPrimArray[kndex] += dEPrimary; | |
1145 | break ; | |
1146 | } | |
ab48128d | 1147 | } // end of check |
5c0368b8 | 1148 | if ( !already && (fMulTrack < fMaxTrack)) { // store it |
2ebdefe6 | 1149 | primArray[fMulTrack] = newPrimary ; |
1150 | dEPrimArray[fMulTrack] = dEPrimary ; | |
1151 | fMulTrack++ ; | |
ab48128d | 1152 | } // store it |
1153 | } // all primaries in digit | |
ab48128d | 1154 | } // all digits |
2ebdefe6 | 1155 | |
af5bdd85 | 1156 | Int_t *sortIdx = new Int_t[fMulTrack]; |
1157 | TMath::Sort(fMulTrack,dEPrimArray,sortIdx); | |
1158 | for(index = 0; index < fMulTrack; index++) { | |
1159 | fTracksList[index] = primArray[sortIdx[index]] ; | |
1160 | fDETracksList[index] = dEPrimArray[sortIdx[index]] ; | |
1161 | } | |
1162 | delete [] sortIdx; | |
1163 | delete [] primArray ; | |
1164 | delete [] dEPrimArray ; | |
2ebdefe6 | 1165 | |
ab48128d | 1166 | } |
7ee5c5be | 1167 | |
87cdc3be | 1168 | //______________________________________________________________________________ |
1169 | void AliEMCALRecPoint::EvalParents(TClonesArray * digits) | |
1170 | { | |
1171 | // Constructs the list of parent particles (tracks) which have contributed to this RecPoint | |
2ebdefe6 | 1172 | |
53e430a3 | 1173 | AliEMCALDigit * digit=0 ; |
af5bdd85 | 1174 | Int_t * parentArray = new Int_t[fMaxTrack] ; |
7f3091e7 | 1175 | memset(parentArray,-1,sizeof(Int_t)*fMaxTrack); |
af5bdd85 | 1176 | Float_t * dEParentArray = new Float_t[fMaxTrack] ; |
7f3091e7 | 1177 | memset(dEParentArray,-1,sizeof(Int_t)*fMaxTrack); |
2ebdefe6 | 1178 | |
87cdc3be | 1179 | Int_t index ; |
1180 | for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits | |
94478418 | 1181 | if (fDigitsList[index] >= digits->GetEntries() || fDigitsList[index] < 0) |
2ebdefe6 | 1182 | AliError(Form("Trying to get invalid digit %d (idx in WriteRecPoint %d)",fDigitsList[index],index)); |
87cdc3be | 1183 | digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; |
2ebdefe6 | 1184 | if(!digit) { |
1185 | AliError("No Digit!!"); | |
1186 | continue; | |
1187 | } | |
1188 | ||
87cdc3be | 1189 | Int_t nparents = digit->GetNiparent() ; |
5c0368b8 | 1190 | if ( nparents == 0 ) continue ; |
2ebdefe6 | 1191 | |
87cdc3be | 1192 | Int_t jndex ; |
1193 | for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit | |
1194 | if ( fMulParent > fMaxParent ) { | |
2ebdefe6 | 1195 | fMulTrack = - 1 ; |
1196 | Error("EvalParents", "increase fMaxParent") ; | |
1197 | break ; | |
87cdc3be | 1198 | } |
af5bdd85 | 1199 | Int_t newParent = digit->GetIparent(jndex+1) ; |
1200 | Float_t newdEParent = digit->GetDEParent(jndex+1) ; | |
87cdc3be | 1201 | Int_t kndex ; |
1202 | Bool_t already = kFALSE ; | |
f1d429fd | 1203 | for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored |
2ebdefe6 | 1204 | if ( newParent == parentArray[kndex] ){ |
1205 | dEParentArray[kndex] += newdEParent; | |
1206 | already = kTRUE ; | |
1207 | break ; | |
1208 | } | |
87cdc3be | 1209 | } // end of check |
94478418 | 1210 | if ( !already && (fMulParent < fMaxParent)) { // store it |
2ebdefe6 | 1211 | parentArray[fMulParent] = newParent ; |
1212 | dEParentArray[fMulParent] = newdEParent ; | |
1213 | fMulParent++ ; | |
87cdc3be | 1214 | } // store it |
1215 | } // all parents in digit | |
87cdc3be | 1216 | } // all digits |
2ebdefe6 | 1217 | |
27e2a47c | 1218 | if (fMulParent>0) { |
af5bdd85 | 1219 | Int_t *sortIdx = new Int_t[fMulParent]; |
1220 | TMath::Sort(fMulParent,dEParentArray,sortIdx); | |
1221 | for(index = 0; index < fMulParent; index++) { | |
1222 | fParentsList[index] = parentArray[sortIdx[index]] ; | |
1223 | fDEParentsList[index] = dEParentArray[sortIdx[index]] ; | |
1224 | } | |
1225 | delete [] sortIdx; | |
27e2a47c | 1226 | } |
2ebdefe6 | 1227 | |
af5bdd85 | 1228 | delete [] parentArray; |
1229 | delete [] dEParentArray; | |
87cdc3be | 1230 | } |
1231 | ||
70a93198 | 1232 | //____________________________________________________________________________ |
1233 | void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const | |
1234 | { | |
9aa6a5f6 | 1235 | // returns the position of the cluster in the local reference system |
1236 | // of the sub-detector | |
70a93198 | 1237 | |
9aa6a5f6 | 1238 | lpos = fLocPos; |
70a93198 | 1239 | } |
1240 | ||
ab48128d | 1241 | //____________________________________________________________________________ |
1242 | void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const | |
1243 | { | |
1244 | // returns the position of the cluster in the global reference system of ALICE | |
70a93198 | 1245 | // These are now the Cartesian X, Y and Z |
e52475ed | 1246 | // cout<<" geom "<<geom<<endl; |
25bb3dcb | 1247 | // fGeomPtr->GetGlobal(fLocPos, gpos, fSuperModuleNumber); |
1248 | gpos = fGlobPos; | |
1249 | ||
9aa6a5f6 | 1250 | } |
1251 | ||
1252 | //____________________________________________________________________________ | |
7cfcebd3 | 1253 | //void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos, TMatrixF & gmat) const |
1254 | //{ | |
1255 | // // returns the position of the cluster in the global reference system of ALICE | |
1256 | // // These are now the Cartesian X, Y and Z | |
1257 | // // cout<<" geom "<<geom<<endl; | |
1258 | // | |
1259 | // //To be implemented | |
1260 | // fGeomPtr->GetGlobalEMCAL(this, gpos, gmat); | |
1261 | // | |
1262 | //} | |
9aa6a5f6 | 1263 | |
1264 | //_____________________________________________________________________________ | |
1265 | void AliEMCALRecPoint::EvalLocal2TrackingCSTransform() | |
1266 | { | |
1267 | //Evaluates local to "tracking" c.s. transformation (B.P.). | |
1268 | //All evaluations should be completed before calling for this | |
1269 | //function. | |
1270 | //See ALICE PPR Chapter 5 p.18 for "tracking" c.s. definition, | |
1271 | //or just ask Jouri Belikov. :) | |
1272 | ||
1273 | SetVolumeId(AliGeomManager::LayerToVolUID(AliGeomManager::kEMCAL,GetSuperModuleNumber())); | |
1274 | ||
1275 | const TGeoHMatrix* tr2loc = GetTracking2LocalMatrix(); | |
1276 | if(!tr2loc) AliFatal(Form("No Tracking2LocalMatrix found.")); | |
1277 | ||
1278 | Double_t lxyz[3] = {fLocPos.X(),fLocPos.Y(),fLocPos.Z()}; | |
1279 | Double_t txyz[3] = {0,0,0}; | |
1280 | ||
1281 | tr2loc->MasterToLocal(lxyz,txyz); | |
1282 | SetX(txyz[0]); SetY(txyz[1]); SetZ(txyz[2]); | |
1283 | ||
1284 | if(AliLog::GetGlobalDebugLevel()>0) { | |
7cfcebd3 | 1285 | TVector3 gpos; //TMatrixF gmat; |
1286 | //GetGlobalPosition(gpos,gmat); //Not doing anythin special, replace by next line. | |
1287 | fGeomPtr->GetGlobal(fLocPos, gpos, GetSuperModuleNumber()); | |
1288 | ||
9aa6a5f6 | 1289 | Float_t gxyz[3]; |
1290 | GetGlobalXYZ(gxyz); | |
1291 | AliInfo(Form("lCS-->(%.3f,%.3f,%.3f), tCS-->(%.3f,%.3f,%.3f), gCS-->(%.3f,%.3f,%.3f), gCScalc-\ | |
1292 | ->(%.3f,%.3f,%.3f), supermodule %d", | |
1293 | fLocPos.X(),fLocPos.Y(),fLocPos.Z(), | |
1294 | GetX(),GetY(),GetZ(), | |
1295 | gpos.X(),gpos.Y(),gpos.Z(), | |
1296 | gxyz[0],gxyz[1],gxyz[2],GetSuperModuleNumber())); | |
1297 | } | |
1298 | ||
70a93198 | 1299 | } |
1300 | ||
1301 | //____________________________________________________________________________ | |
1302 | Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const | |
1303 | { | |
1304 | // Finds the maximum energy in the cluster | |
ab48128d | 1305 | |
70a93198 | 1306 | Float_t menergy = 0. ; |
1307 | ||
1308 | Int_t iDigit; | |
70a93198 | 1309 | for(iDigit=0; iDigit<fMulDigit; iDigit++) { |
1310 | ||
1311 | if(fEnergyList[iDigit] > menergy) | |
1312 | menergy = fEnergyList[iDigit] ; | |
1313 | } | |
1314 | return menergy ; | |
ab48128d | 1315 | } |
1316 | ||
25bb3dcb | 1317 | //____________________________________________________________________________ |
1318 | Int_t AliEMCALRecPoint::GetMaximalEnergyIndex(void) const | |
1319 | { | |
1320 | // Finds the maximum energy in the cluster | |
1321 | ||
1322 | Float_t menergy = 0. ; | |
44fe3410 | 1323 | Int_t mid = 0 ; |
25bb3dcb | 1324 | Int_t iDigit; |
1325 | ||
1326 | for(iDigit=0; iDigit<fMulDigit; iDigit++) { | |
1327 | ||
1328 | if(fEnergyList[iDigit] > menergy){ | |
1329 | menergy = fEnergyList[iDigit] ; | |
1330 | mid = iDigit ; | |
1331 | } | |
1332 | }//loop on cluster digits | |
1333 | ||
1334 | return mid ; | |
1335 | } | |
1336 | ||
1337 | ||
aad8e277 | 1338 | //____________________________________________________________________________ |
70a93198 | 1339 | Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const |
aad8e277 | 1340 | { |
70a93198 | 1341 | // Calculates the multiplicity of digits with energy larger than H*energy |
1342 | ||
1343 | Int_t multipl = 0 ; | |
1344 | Int_t iDigit ; | |
1345 | for(iDigit=0; iDigit<fMulDigit; iDigit++) { | |
1346 | ||
1347 | if(fEnergyList[iDigit] > H * fAmp) | |
1348 | multipl++ ; | |
1349 | } | |
1350 | return multipl ; | |
1351 | } | |
1352 | ||
1353 | //____________________________________________________________________________ | |
1354 | Int_t AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy, | |
1355 | Float_t locMaxCut,TClonesArray * digits) const | |
1356 | { | |
1357 | // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum | |
1358 | // energy difference between two local maxima | |
1359 | ||
53e430a3 | 1360 | AliEMCALDigit * digit = 0; |
1361 | AliEMCALDigit * digitN = 0; | |
70a93198 | 1362 | |
53e430a3 | 1363 | Int_t iDigitN = 0 ; |
1364 | Int_t iDigit = 0 ; | |
70a93198 | 1365 | |
1366 | for(iDigit = 0; iDigit < fMulDigit; iDigit++) | |
1367 | maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ; | |
1368 | ||
1369 | for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) { | |
1370 | if(maxAt[iDigit]) { | |
1371 | digit = maxAt[iDigit] ; | |
1372 | ||
1373 | for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) { | |
53e430a3 | 1374 | digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; |
70a93198 | 1375 | |
1376 | if ( AreNeighbours(digit, digitN) ) { | |
1377 | if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) { | |
1378 | maxAt[iDigitN] = 0 ; | |
1379 | // but may be digit too is not local max ? | |
1380 | if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) | |
1381 | maxAt[iDigit] = 0 ; | |
1382 | } | |
1383 | else { | |
1384 | maxAt[iDigit] = 0 ; | |
1385 | // but may be digitN too is not local max ? | |
1386 | if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) | |
1387 | maxAt[iDigitN] = 0 ; | |
1388 | } | |
1389 | } // if Areneighbours | |
1390 | } // while digitN | |
1391 | } // slot not empty | |
1392 | } // while digit | |
1393 | ||
1394 | iDigitN = 0 ; | |
1395 | for(iDigit = 0; iDigit < fMulDigit; iDigit++) { | |
1396 | if(maxAt[iDigit] ){ | |
1397 | maxAt[iDigitN] = maxAt[iDigit] ; | |
1398 | maxAtEnergy[iDigitN] = fEnergyList[iDigit] ; | |
1399 | iDigitN++ ; | |
1400 | } | |
1401 | } | |
1402 | return iDigitN ; | |
1403 | } | |
4635df1f | 1404 | |
1405 | //____________________________________________________________________________ | |
1406 | Int_t AliEMCALRecPoint::GetPrimaryIndex() const | |
1407 | { | |
1408 | // Get the primary track index in TreeK which deposits the most energy | |
af5bdd85 | 1409 | // in Digits which forms RecPoint. |
4635df1f | 1410 | |
af5bdd85 | 1411 | if (fMulTrack) |
1412 | return fTracksList[0]; | |
1413 | return -12345; | |
4635df1f | 1414 | } |
1415 | ||
70a93198 | 1416 | //____________________________________________________________________________ |
1417 | void AliEMCALRecPoint::EvalTime(TClonesArray * digits){ | |
1418 | // time is set to the time of the digit with the maximum energy | |
1419 | ||
53e430a3 | 1420 | Float_t maxE = 0; |
1421 | Int_t maxAt = 0; | |
70a93198 | 1422 | for(Int_t idig=0; idig < fMulDigit; idig++){ |
1423 | if(fEnergyList[idig] > maxE){ | |
53e430a3 | 1424 | maxE = fEnergyList[idig] ; |
70a93198 | 1425 | maxAt = idig; |
1426 | } | |
1427 | } | |
1428 | fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ; | |
aad8e277 | 1429 | |
aad8e277 | 1430 | } |
ab48128d | 1431 | |
1432 | //______________________________________________________________________________ | |
1433 | void AliEMCALRecPoint::Paint(Option_t *) | |
1434 | { | |
1435 | // Paint this ALiRecPoint as a TMarker with its current attributes | |
1436 | ||
1437 | TVector3 pos(0.,0.,0.) ; | |
1438 | GetLocalPosition(pos) ; | |
1439 | Coord_t x = pos.X() ; | |
1440 | Coord_t y = pos.Z() ; | |
1441 | Color_t markercolor = 1 ; | |
53e430a3 | 1442 | Size_t markersize = 1.; |
ab48128d | 1443 | Style_t markerstyle = 5 ; |
1444 | ||
1445 | if (!gPad->IsBatch()) { | |
1446 | gVirtualX->SetMarkerColor(markercolor) ; | |
1447 | gVirtualX->SetMarkerSize (markersize) ; | |
1448 | gVirtualX->SetMarkerStyle(markerstyle) ; | |
1449 | } | |
1450 | gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ; | |
1451 | gPad->PaintPolyMarker(1,&x,&y,"") ; | |
1452 | } | |
70a93198 | 1453 | |
9aa6a5f6 | 1454 | //_____________________________________________________________________ |
1ae500a2 | 1455 | Double_t AliEMCALRecPoint::TmaxInCm(const Double_t e , const Int_t key) |
1456 | { | |
9aa6a5f6 | 1457 | // e energy in GeV) |
1ae500a2 | 1458 | // key = 0(gamma, default) |
1459 | // != 0(electron) | |
47583f1b | 1460 | const Double_t ca = 4.82; // shower max parameter - first guess; ca=TMath::Log(1000./8.07) |
1461 | const Double_t x0 = 1.23; // radiation lenght (cm) | |
1462 | Double_t tmax = 0.; // position of electromagnetic shower max in cm | |
1ae500a2 | 1463 | |
1ae500a2 | 1464 | if(e>0.1) { |
1465 | tmax = TMath::Log(e) + ca; | |
1466 | if (key==0) tmax += 0.5; | |
1467 | else tmax -= 0.5; | |
37890aaf | 1468 | tmax *= x0; // convert to cm |
1ae500a2 | 1469 | } |
1470 | return tmax; | |
1471 | } | |
1472 | ||
70a93198 | 1473 | //______________________________________________________________________________ |
1474 | Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const | |
1475 | { | |
1476 | //Converts Theta (Radians) to Eta(Radians) | |
1477 | return (2.*TMath::ATan(TMath::Exp(-arg))); | |
1478 | } | |
1479 | ||
1480 | //______________________________________________________________________________ | |
1481 | Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const | |
1482 | { | |
1483 | //Converts Eta (Radians) to Theta(Radians) | |
1484 | return (-1 * TMath::Log(TMath::Tan(0.5 * arg))); | |
1485 | } | |
261b1065 | 1486 | |
1487 | //____________________________________________________________________________ | |
e1a51e6e | 1488 | void AliEMCALRecPoint::Print(Option_t *opt) const |
261b1065 | 1489 | { |
1490 | // Print the list of digits belonging to the cluster | |
e1a51e6e | 1491 | if(strlen(opt)==0) return; |
261b1065 | 1492 | TString message ; |
4800667c | 1493 | message = "AliEMCALRecPoint:\n" ; |
261b1065 | 1494 | message += " digits # = " ; |
1495 | Info("Print", message.Data()) ; | |
1496 | ||
1497 | Int_t iDigit; | |
1498 | for(iDigit=0; iDigit<fMulDigit; iDigit++) | |
1499 | printf(" %d ", fDigitsList[iDigit] ) ; | |
e52475ed | 1500 | printf("\n"); |
1501 | ||
261b1065 | 1502 | Info("Print", " Energies = ") ; |
1503 | for(iDigit=0; iDigit<fMulDigit; iDigit++) | |
1504 | printf(" %f ", fEnergyList[iDigit] ) ; | |
e52475ed | 1505 | printf("\n"); |
1506 | ||
1507 | Info("Print", "\n Abs Ids = ") ; | |
1508 | for(iDigit=0; iDigit<fMulDigit; iDigit++) | |
1509 | printf(" %i ", fAbsIdList[iDigit] ) ; | |
1510 | printf("\n"); | |
1511 | ||
1512 | Info("Print", " Primaries ") ; | |
261b1065 | 1513 | for(iDigit = 0;iDigit < fMulTrack; iDigit++) |
1514 | printf(" %d ", fTracksList[iDigit]) ; | |
e52475ed | 1515 | |
1516 | printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]); | |
1517 | ||
85c60a8e | 1518 | message = " ClusterType = %d" ; |
1519 | message += " Multiplicity = %d" ; | |
261b1065 | 1520 | message += " Cluster Energy = %f" ; |
1521 | message += " Core energy = %f" ; | |
1522 | message += " Core radius = %f" ; | |
1523 | message += " Number of primaries %d" ; | |
1524 | message += " Stored at position %d" ; | |
85c60a8e | 1525 | Info("Print", message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList() ) ; |
261b1065 | 1526 | } |
1d46d1f6 | 1527 | |
9aa6a5f6 | 1528 | //___________________________________________________________ |
1d46d1f6 | 1529 | Double_t AliEMCALRecPoint::GetPointEnergy() const |
1530 | { | |
37890aaf | 1531 | //Returns energy .... |
47583f1b | 1532 | Double_t e=0.0; |
1d46d1f6 | 1533 | for(int ic=0; ic<GetMultiplicity(); ic++) e += double(fEnergyList[ic]); |
1534 | return e; | |
1535 | } |