]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/CALO/AliHLTCaloClusterAnalyser.cxx
- fixing compilation warnings
[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        wtot = 0;
114        x = 0; 
115        z = 0;
116
117        /*Float_t maxAmp = 0;
118        Int_t maxX = 0;
119        Int_t maxZ = 0;
120       */for(iDigit = 0; iDigit < recPoint->fMultiplicity; iDigit++)
121         {
122
123           digit = &(fDigitDataArray[*digitIndexPtr]);
124
125           xi = digit->fX+0.5;
126           zi = digit->fZ+0.5;
127
128           if (recPoint->fAmp > 0 && digit->fEnergy > 0) 
129             {
130               Float_t w = TMath::Max( 0., fLogWeight + TMath::Log( digit->fEnergy / recPoint->fAmp ) ) ;
131               x    += xi * w ;
132               z    += zi * w ;
133               wtot += w ;
134 /*            if(digit->fEnergy > maxAmp)
135               {
136                  maxAmp = digit->fEnergy;
137                  maxX = digit->fX + 0.5;
138                  maxZ = digit->fZ + 0.5;
139               }*/
140             }
141           digitIndexPtr++;
142         }
143
144       if (wtot>0) 
145         {
146           recPoint->fX = x/wtot ;
147           recPoint->fZ = z/wtot ;
148         }
149       else
150         {
151           recPoint->fAmp = 0;
152         }
153 //      printf("Max digit: E = %f, x = %d, z= %d, cluster: E = %f, x = %f, z = %f\n" , maxAmp, maxX, maxZ, recPoint->fAmp, recPoint->fX, recPoint->fZ);
154     }
155   return 0;
156 }
157
158
159 Int_t 
160 AliHLTCaloClusterAnalyser::CalculateRecPointMoments()
161 {       
162   //See header file for documentation
163   return 0;
164 }
165
166 Int_t 
167 AliHLTCaloClusterAnalyser::CalculateClusterMoments(AliHLTCaloRecPointDataStruct */*recPointPtr*/, AliHLTCaloClusterDataStruct* /*clusterPtr*/)
168 {
169   //See header file for documentation
170   return 0;
171 }
172
173
174 Int_t 
175 AliHLTCaloClusterAnalyser::DeconvoluteClusters()
176 {
177   //See header file for documentation
178   return 0;
179 }
180
181 Int_t 
182 AliHLTCaloClusterAnalyser::CreateClusters(Int_t nRecPoints, UInt_t availableSize, UInt_t& totSize)
183 {
184   //See header file for documentation
185
186    
187   fNRecPoints = nRecPoints;
188
189   if(fGeometry == 0)
190   {
191      HLTError("No geometry object is initialised, creation of clusters stopped");
192   }
193
194   CalculateCenterOfGravity();
195
196   //  AliHLTCaloDigitDataStruct* digitPtr = &(recPointPtr->fDigits);  
197   AliHLTCaloDigitDataStruct* digitPtr = 0;
198
199   AliHLTCaloClusterDataStruct* caloClusterPtr = 0;
200
201   //  Int_t id = -1;
202   TVector3 globalPos;
203
204   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
205     {
206       if((availableSize - totSize)  < sizeof(AliHLTCaloClusterDataStruct))
207       {
208          HLTError("Out of buffer: available size is: %d, total size used: %d", availableSize, totSize);
209          return -ENOBUFS;
210       }
211       totSize += sizeof(AliHLTCaloClusterDataStruct);
212       
213       caloClusterPtr = fCaloClusterDataPtr;
214      
215       AliHLTCaloRecPointDataStruct *recPointPtr = fRecPointArray[i];
216       
217       AliHLTCaloGlobalCoordinate globalCoord;
218       fGeometry->GetGlobalCoordinates(*recPointPtr, globalCoord);
219
220       caloClusterPtr->fGlobalPos[0] =  globalCoord.fX;
221       caloClusterPtr->fGlobalPos[1] =  globalCoord.fY;
222       caloClusterPtr->fGlobalPos[2] =  globalCoord.fZ;
223
224       HLTDebug("Cluster local position: x = %f, z = %f, module = %d", recPointPtr->fX, recPointPtr->fZ, recPointPtr->fModule);
225       HLTDebug("Cluster global position: x = %f, y = %f, z = %f", globalCoord.fX, globalCoord.fY, globalCoord.fZ);
226       
227       //caloClusterPtr->fNCells = 0;//recPointPtr->fMultiplicity;
228       caloClusterPtr->fNCells = recPointPtr->fMultiplicity;
229
230       caloClusterPtr->fClusterType = fClusterType;
231 //      Int_t tmpSize = 0;//totSize + (caloClusterPtr->fNCells-1)*(sizeof(Short_t) + sizeof(Float_t));
232
233       //TODO remove hardcoded 10; 
234       memset(caloClusterPtr->fTracksMatched, 0xff, sizeof(Int_t)*10);
235
236       //Int_t tmpSize = totSize + (caloClusterPtr->fNCells-1)*(sizeof(Short_t) + sizeof(Float_t));
237       UInt_t tmpSize = (caloClusterPtr->fNCells-1)*sizeof(AliHLTCaloCellDataStruct);
238
239       if((availableSize - totSize)  < tmpSize)
240       {
241          HLTError("Out of buffer, available size is: %d, total size used: %d, extra size needed: %d", availableSize, totSize, tmpSize);
242          return -ENOBUFS;
243       }
244       Int_t *digitIndexPtr = &(recPointPtr->fDigits);
245       Int_t id = 0;
246
247        AliHLTCaloCellDataStruct *cellPtr = &(caloClusterPtr->fCaloCells);
248        
249       for(UInt_t j = 0; j < caloClusterPtr->fNCells; j++)
250         {
251            digitPtr = &(fDigitDataArray[*digitIndexPtr]);
252            id++;
253            fGeometry->GetCellAbsId(recPointPtr->fModule, digitPtr->fX, digitPtr->fZ, id);
254            
255           cellPtr->fCellsAbsId= id;
256           cellPtr->fCellsAmpFraction = digitPtr->fEnergy/recPointPtr->fAmp;
257           //printf("Cell ID pointer: %x\n", cellIDPtr);
258          //printf("Cell Amp Pointer: %x\n", cellAmpFracPtr);
259          //printf("Cell pos: x = %d, z = %d\n", digitPtr->fX, digitPtr->fZ);
260          //printf("Cell ID: %d, pointer: %x\n", *cellIDPtr, cellIDPtr);
261          //printf("Cell Amp: %f, pointer: %x\n", *cellAmpFracPtr, cellAmpFracPtr);
262           cellPtr++;
263           digitIndexPtr++;
264           
265         }
266
267       totSize += tmpSize;
268
269       caloClusterPtr->fEnergy = recPointPtr->fAmp;
270
271       HLTDebug("Cluster global position: x = %f, y = %f, z = %f, energy: %f, number of cells: %d, cluster pointer: %x", globalCoord.fX, globalCoord.fY, globalCoord.fZ, caloClusterPtr->fEnergy, caloClusterPtr->fNCells,  caloClusterPtr);
272
273       if(fDoClusterFit)
274         {
275           FitCluster(recPointPtr);
276         }
277       else
278         {
279           caloClusterPtr->fDispersion = 0;
280           caloClusterPtr->fFitQuality = 0;
281           caloClusterPtr->fM20 = 0;
282           caloClusterPtr->fM02 = 0;
283
284         }
285       if(fHaveCPVInfo)
286         {
287           caloClusterPtr->fEmcCpvDistance = GetCPVDistance(recPointPtr);
288         }
289       else
290         {
291           caloClusterPtr->fEmcCpvDistance = -1;
292         }
293       if(fDoPID)
294         {
295           DoParticleIdentification(caloClusterPtr);
296         }
297       else
298         {
299           for(Int_t k = 0; k < AliPID::kSPECIESN; k++)
300             {
301               caloClusterPtr->fPID[k] = 0;
302             }
303         }
304       if(fHaveDistanceToBadChannel)
305         {
306           caloClusterPtr->fDistanceToBadChannel = GetDistanceToBadChannel(caloClusterPtr);
307         }
308       else
309         {
310           caloClusterPtr->fDistanceToBadChannel = -1;
311         }
312
313       caloClusterPtr->fClusterType = fClusterType;
314           //totSize += sizeof(AliHLTCaloClusterDataStruct) + (caloClusterPtr->fNCells)*(sizeof(Short_t) +sizeof(Float_t)-1);   
315       //totSize += sizeof(AliHLTCaloClusterDataStruct) + (caloClusterPtr->fNCells-1)*(sizeof(Short_t) + sizeof(Float_t));   
316       //printf("CaloClusterPtr: %x, energy: %f\n", caloClusterPtr, caloClusterPtr->fEnergy);
317       
318       //      caloClusterPtr = reinterpret_cast<AliHLTCaloClusterDataStruct*>(cellAmpFracPtr);
319       //caloClusterPtr = reinterpret_cast<AliHLTCaloClusterDataStruct*>(cellIDPtr);
320       fCaloClusterDataPtr = reinterpret_cast<AliHLTCaloClusterDataStruct*>(cellPtr);
321       recPointPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(digitPtr);
322       //digitPtr = &(recPointPtr->fDigits);  
323     }
324
325 return fNRecPoints;
326
327 }
328