]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/CALO/AliHLTCaloClusterAnalyser.cxx
4d0f3ecaee80e3cdc679862ae777dfeacb5b0f02
[u/mrichter/AliRoot.git] / HLT / CALO / AliHLTCaloClusterAnalyser.cxx
1 // $Id: AliHLTCaloClusterAnalyser.cxx 35107 2009-09-30 01:45:06Z phille $
2
3 /**************************************************************************
4  * This file is property of and copyright by the ALICE HLT Project        * 
5  * All rights reserved.                                                   *
6  *                                                                        *
7  * Primary Authors: Oystein Djuvsland                                     *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          * 
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 /** 
19  * @file   AliHLTCaloClusterAnalyser.cxx
20  * @author Oystein Djuvsland
21  * @date 
22  * @brief  Cluster analyser for Calo HLT 
23  */
24
25 // see header file for class documentation
26 // or
27 // refer to README to build package
28 // or
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
30
31 #include "AliHLTCaloClusterAnalyser.h"
32 #include "AliHLTCaloRecPointHeaderStruct.h"
33 #include "AliHLTCaloRecPointDataStruct.h"
34 #include "AliHLTCaloClusterDataStruct.h"
35 //#include "AliHLTCaloPhysicsAnalyzer.h"
36 #include "AliHLTCaloGeometry.h"
37 #include "AliESDCaloCluster.h"
38 #include "TMath.h"
39 #include "TVector3.h"
40 #include "TH1F.h"
41 #include "TFile.h"
42 #include "AliHLTCaloClusterizer.h"
43
44 ClassImp(AliHLTCaloClusterAnalyser);
45
46 AliHLTCaloClusterAnalyser::AliHLTCaloClusterAnalyser() : 
47   //  AliHLTCaloBase(),
48   fLogWeight(4.5),
49   fRecPointArray(0),
50   fDigitDataArray(0),
51   fNRecPoints(0),
52   fCaloClusterDataPtr(0),
53   fCaloClusterHeaderPtr(0),
54   //fAnalyzerPtr(0),
55   fDoClusterFit(false),
56   fHaveCPVInfo(false),
57   fDoPID(false),
58   fHaveDistanceToBadChannel(false),
59   fGeometry(0),
60   fClusterType(AliESDCaloCluster::kPHOSCluster)
61 {
62   //See header file for documentation
63 }
64
65 AliHLTCaloClusterAnalyser::~AliHLTCaloClusterAnalyser() 
66 {
67   // See header file for class documentation
68 }
69
70 void 
71 AliHLTCaloClusterAnalyser::SetCaloClusterData(AliHLTCaloClusterDataStruct *caloClusterDataPtr)
72
73   //see header file for documentation
74   fCaloClusterDataPtr = caloClusterDataPtr; 
75 }
76
77 void
78 AliHLTCaloClusterAnalyser::SetRecPointArray(AliHLTCaloRecPointDataStruct **recPointDataPtr, Int_t nRecPoints)
79
80   fRecPointArray = recPointDataPtr; 
81   fNRecPoints = nRecPoints;
82 }
83
84 void 
85 AliHLTCaloClusterAnalyser::SetDigitDataArray(AliHLTCaloDigitDataStruct *digits) 
86
87 //   AliHLTCaloClusterizer cl("PHOS");
88   // cl.CheckDigits(fRecPointArray, digits, fNRecPoints);
89    fDigitDataArray = digits; 
90    //cl.CheckDigits(fRecPointArray, fDigitDataArray, fNRecPoints);
91 }
92
93 Int_t
94 AliHLTCaloClusterAnalyser::CalculateCenterOfGravity()
95 {
96   //see header file for documentation
97   Float_t wtot = 0.;
98   Float_t x = 0.;
99   Float_t z = 0.;
100   Float_t xi = 0.;
101   Float_t zi = 0.;
102
103   AliHLTCaloDigitDataStruct *digit = 0;
104
105   UInt_t iDigit = 0;
106
107   for(Int_t iRecPoint=0; iRecPoint < fNRecPoints; iRecPoint++) 
108     {
109       AliHLTCaloRecPointDataStruct *recPoint = fRecPointArray[iRecPoint];
110       //      digit = &(recPoint->fDigits);
111
112       Int_t *digitIndexPtr = &(recPoint->fDigits);
113
114       for(iDigit = 0; iDigit < recPoint->fMultiplicity; iDigit++)
115         {
116
117           digit = &(fDigitDataArray[*digitIndexPtr]);
118
119           xi = digit->fX+0.5;
120           zi = digit->fZ+0.5;
121
122           if (recPoint->fAmp > 0 && digit->fEnergy > 0) 
123             {
124               Float_t w = TMath::Max( 0., fLogWeight + TMath::Log( digit->fEnergy / recPoint->fAmp ) ) ;
125               x    += xi * w ;
126               z    += zi * w ;
127               wtot += w ;
128             }
129           digitIndexPtr++;
130         }
131
132       if (wtot>0) 
133         {
134           recPoint->fX = x/wtot ;
135           recPoint->fZ = z/wtot ;
136         }
137       else
138         {
139           recPoint->fAmp = 0;
140         }
141     }
142   return 0;
143 }
144
145
146 Int_t 
147 AliHLTCaloClusterAnalyser::CalculateRecPointMoments()
148 {       
149   //See header file for documentation
150   return 0;
151 }
152
153 Int_t 
154 AliHLTCaloClusterAnalyser::CalculateClusterMoments(AliHLTCaloRecPointDataStruct */*recPointPtr*/, AliHLTCaloClusterDataStruct* /*clusterPtr*/)
155 {
156   //See header file for documentation
157   return 0;
158 }
159
160
161 Int_t 
162 AliHLTCaloClusterAnalyser::DeconvoluteClusters()
163 {
164   //See header file for documentation
165   return 0;
166 }
167
168 Int_t 
169 AliHLTCaloClusterAnalyser::CreateClusters(Int_t nRecPoints, UInt_t availableSize, UInt_t& totSize)
170 {
171   //See header file for documentation
172
173    
174   totSize += sizeof(AliHLTCaloClusterDataStruct);
175   fNRecPoints = nRecPoints;
176
177   if(fGeometry == 0)
178   {
179      HLTError("No geometry object is initialised, creation of clusters stopped");
180   }
181
182   CalculateCenterOfGravity();
183
184   //  AliHLTCaloDigitDataStruct* digitPtr = &(recPointPtr->fDigits);  
185   AliHLTCaloDigitDataStruct* digitPtr = 0;
186
187   AliHLTCaloClusterDataStruct* caloClusterPtr = 0;
188   UShort_t* cellIDPtr = 0;
189   Float_t* cellAmpFracPtr = 0;;
190   
191 //  Int_t id = -1;
192   TVector3 globalPos;
193
194   for(Int_t i = 0; i < fNRecPoints; i++) //TODO needs fix when we start unfolding (number of clusters not necessarily same as number of recpoints gotten from the clusterizer
195     {
196       if((availableSize - totSize)  < sizeof(AliHLTCaloClusterDataStruct))
197       {
198          HLTError("Out of buffer");
199          return -ENOBUFS;
200       }
201       
202       caloClusterPtr = fCaloClusterDataPtr;
203      
204       cellIDPtr = &(caloClusterPtr->fCellsAbsId);
205       cellAmpFracPtr = &(caloClusterPtr->fCellsAmpFraction);
206      
207       AliHLTCaloRecPointDataStruct *recPointPtr = fRecPointArray[i];
208       
209       AliHLTCaloGlobalCoordinate globalCoord;
210       fGeometry->GetGlobalCoordinates(*recPointPtr, globalCoord);
211
212       caloClusterPtr->fGlobalPos[0] =  globalCoord.fX;
213       caloClusterPtr->fGlobalPos[1] =  globalCoord.fY;
214       caloClusterPtr->fGlobalPos[2] =  globalCoord.fZ;
215
216       HLTDebug("Cluster local position: x = %f, z = %f, module = %d", recPointPtr->fX, recPointPtr->fZ, recPointPtr->fModule);
217       HLTDebug("Cluster global position: x = %f, y = %f, z = %f", globalCoord.fX, globalCoord.fY, globalCoord.fZ);
218       
219       //caloClusterPtr->fNCells = 0;//recPointPtr->fMultiplicity;
220       caloClusterPtr->fNCells = recPointPtr->fMultiplicity;
221
222       caloClusterPtr->fClusterType = fClusterType;
223 //      Int_t tmpSize = 0;//totSize + (caloClusterPtr->fNCells-1)*(sizeof(Short_t) + sizeof(Float_t));
224
225       //TODO remove hardcoded 10; 
226       memset(caloClusterPtr->fTracksMatched, 0xff, sizeof(Int_t)*10);
227
228       Int_t tmpSize = totSize + (caloClusterPtr->fNCells-1)*(sizeof(Short_t) + sizeof(Float_t));
229
230       if((availableSize - totSize)  < tmpSize)
231       {
232          HLTError("Out of buffer");
233          return -ENOBUFS;
234       }
235       Int_t *digitIndexPtr = &(recPointPtr->fDigits);
236       Int_t id = 0;
237
238       for(UInt_t j = 0; j < caloClusterPtr->fNCells; j++)
239         {
240            digitPtr = &(fDigitDataArray[*digitIndexPtr]);
241            id++;
242            fGeometry->GetCellAbsId(recPointPtr->fModule, digitPtr->fX, digitPtr->fZ, id);
243           *cellIDPtr = id;
244           *cellAmpFracPtr = digitPtr->fEnergy/recPointPtr->fAmp;
245           //printf("Cell ID pointer: %x\n", cellIDPtr);
246          //printf("Cell Amp Pointer: %x\n", cellAmpFracPtr);
247          //printf("Digit ID: %d\n", *cellIDPtr);
248           digitPtr++;
249           //cellIDPtr = reinterpret_cast<UShort_t*>(reinterpret_cast<char*>(cellAmpFracPtr) + sizeof(Float_t)); 
250           cellIDPtr = reinterpret_cast<UShort_t*>(reinterpret_cast<char*>(cellIDPtr) + sizeof(Float_t) + sizeof(Short_t)); 
251           cellAmpFracPtr = reinterpret_cast<Float_t*>(reinterpret_cast<char*>(cellAmpFracPtr) + sizeof(Short_t) + sizeof(Float_t));
252           
253           digitIndexPtr++;
254           
255         }
256
257       totSize += tmpSize;
258
259       caloClusterPtr->fEnergy = recPointPtr->fAmp;
260
261       if(fDoClusterFit)
262         {
263           FitCluster(recPointPtr);
264         }
265       else
266         {
267           caloClusterPtr->fDispersion = 0;
268           caloClusterPtr->fFitQuality = 0;
269           caloClusterPtr->fM20 = 0;
270           caloClusterPtr->fM02 = 0;
271
272         }
273       if(fHaveCPVInfo)
274         {
275           caloClusterPtr->fEmcCpvDistance = GetCPVDistance(recPointPtr);
276         }
277       else
278         {
279           caloClusterPtr->fEmcCpvDistance = -1;
280         }
281       if(fDoPID)
282         {
283           DoParticleIdentification(caloClusterPtr);
284         }
285       else
286         {
287           for(Int_t k = 0; k < AliPID::kSPECIESN; k++)
288             {
289               caloClusterPtr->fPID[k] = 0;
290             }
291         }
292       if(fHaveDistanceToBadChannel)
293         {
294           caloClusterPtr->fDistanceToBadChannel = GetDistanceToBadChannel(caloClusterPtr);
295         }
296       else
297         {
298           caloClusterPtr->fDistanceToBadChannel = -1;
299         }
300
301       caloClusterPtr->fClusterType = fClusterType;
302           //totSize += sizeof(AliHLTCaloClusterDataStruct) + (caloClusterPtr->fNCells)*(sizeof(Short_t) +sizeof(Float_t)-1);   
303       //totSize += sizeof(AliHLTCaloClusterDataStruct) + (caloClusterPtr->fNCells-1)*(sizeof(Short_t) + sizeof(Float_t));   
304       //printf("CaloClusterPtr: %x, energy: %f\n", caloClusterPtr, caloClusterPtr->fEnergy);
305       
306       //      caloClusterPtr = reinterpret_cast<AliHLTCaloClusterDataStruct*>(cellAmpFracPtr);
307       //caloClusterPtr = reinterpret_cast<AliHLTCaloClusterDataStruct*>(cellIDPtr);
308       fCaloClusterDataPtr = reinterpret_cast<AliHLTCaloClusterDataStruct*>(cellIDPtr);
309       recPointPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(digitPtr);
310       //digitPtr = &(recPointPtr->fDigits);  
311     }
312
313 return fNRecPoints;
314
315 }
316