c07f268c7e7e0d34ca565d2351b816b49949f525
[u/mrichter/AliRoot.git] / HLT / CALO / AliHLTCaloClusterizer.cxx
1 // $Id$\r
2 \r
3 /**************************************************************************\r
4  * This file is property of and copyright by the ALICE HLT Project        * \r
5  * All rights reserved.                                                   *\r
6  *                                                                        *\r
7  * Primary Authors: Oystein Djuvsland                                     *\r
8  *                                                                        *\r
9  * Permission to use, copy, modify and distribute this software and its   *\r
10  * documentation strictly for non-commercial purposes is hereby granted   *\r
11  * without fee, provided that the above copyright notice appears in all   *\r
12  * copies and that both the copyright notice and this permission notice   *\r
13  * appear in the supporting documentation. The authors make no claims     *\r
14  * about the suitability of this software for any purpose. It is          * \r
15  * provided "as is" without express or implied warranty.                  *\r
16  **************************************************************************/\r
17 \r
18 /** \r
19  * @file   AliHLTCaloClusterizer.cxx\r
20  * @author Oystein Djuvsland\r
21  * @date \r
22  * @brief  Clusterizer for PHOS HLT \r
23  */\r
24 \r
25 // see header file for class documentation\r
26 // or\r
27 // refer to README to build package\r
28 // or\r
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt\r
30 \r
31 #include "AliHLTCaloClusterizer.h"\r
32 #include "AliHLTLogging.h"\r
33 #include "TMath.h"\r
34 #include "AliHLTCaloRecPointDataStruct.h"\r
35 #include "AliHLTCaloDigitDataStruct.h"\r
36 #include "AliHLTCaloDigitContainerDataStruct.h"\r
37 #include "AliHLTCaloConstantsHandler.h"\r
38 \r
39 ClassImp(AliHLTCaloClusterizer);\r
40 \r
41 AliHLTCaloClusterizer::AliHLTCaloClusterizer(TString det):\r
42   AliHLTCaloConstantsHandler(det),\r
43   fRecPointArray(0),\r
44   fRecPointDataPtr(0),\r
45   fFirstRecPointPtr(0),\r
46   fArraySize(0),\r
47   fAvailableSize(0),\r
48   fUsedSize(0),\r
49   fNRecPoints(0),\r
50   fDigitIndexPtr(0),\r
51   fEmcClusteringThreshold(0),\r
52   fEmcMinEnergyThreshold(0),\r
53   fEmcTimeGate(0),\r
54   fDigitsInCluster(0),\r
55   fDigitsPointerArray(0),\r
56   fDigitContainerPtr(0),\r
57   fMaxDigitIndexDiff(0),\r
58   fNDigits(0)\r
59 {\r
60   //See header file for documentation\r
61   //fEmcClusteringThreshold = 0.2;\r
62   //fEmcMinEnergyThreshold = 0.03;\r
63 \r
64   fEmcClusteringThreshold = 0.1;\r
65   fEmcMinEnergyThreshold = 0.01;\r
66   fEmcTimeGate = 1.e-6 ;\r
67   \r
68   fMaxDigitIndexDiff = 2*fCaloConstants->GetNZROWSMOD();\r
69   \r
70   \r
71   fArraySize = 10;\r
72   fRecPointArray = new AliHLTCaloRecPointDataStruct*[fArraySize];\r
73   \r
74   fAvailableSize = sizeof(AliHLTCaloRecPointDataStruct) * 20;\r
75   fFirstRecPointPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(new UChar_t[fAvailableSize]);\r
76   fRecPointDataPtr = fFirstRecPointPtr;  \r
77 \r
78 }//end\r
79 \r
80 AliHLTCaloClusterizer::~AliHLTCaloClusterizer()  \r
81 {\r
82   //See header file for documentation\r
83 }\r
84 \r
85 void \r
86 AliHLTCaloClusterizer::SetRecPointDataPtr(AliHLTCaloRecPointDataStruct* recPointDataPtr)\r
87 {\r
88   // See header file for documentation\r
89   fRecPointDataPtr = recPointDataPtr;\r
90 }\r
91 \r
92 Int_t \r
93 AliHLTCaloClusterizer::ClusterizeEvent(Int_t nDigits)\r
94 {\r
95   //see header file for documentation\r
96   Int_t nRecPoints = 0;\r
97   fNRecPoints = 0;\r
98   fUsedSize = 0;\r
99   fNDigits = nDigits;\r
100   fRecPointDataPtr = fFirstRecPointPtr;\r
101 \r
102   //Clusterization starts\r
103   for(Int_t i = 0; i < nDigits; i++)\r
104     { \r
105       fDigitsInCluster = 0;\r
106 \r
107       if(fDigitsPointerArray[i]->fEnergy < fEmcClusteringThreshold)\r
108         {\r
109           continue;\r
110         }\r
111       CheckArray();\r
112       CheckBuffer();\r
113 \r
114       // First digit is placed at the fDigits member variable in the recpoint\r
115       fDigitIndexPtr = &(fRecPointDataPtr->fDigits);\r
116       //fUsedSize += sizeof(AliHLTCaloRecPointDataStruct);\r
117       \r
118       fRecPointDataPtr->fAmp = 0;\r
119       fRecPointDataPtr->fModule = fDigitsPointerArray[i]->fModule;\r
120 \r
121       // Assigning the digit to this rec point\r
122       fRecPointDataPtr->fDigits = i;\r
123       fUsedSize += sizeof(AliHLTCaloRecPointDataStruct);\r
124       \r
125       // Incrementing the pointer to be ready for new entry\r
126       fDigitIndexPtr++;\r
127 \r
128       fRecPointDataPtr->fAmp += fDigitsPointerArray[i]->fEnergy;\r
129       fDigitsPointerArray[i]->fEnergy = 0;\r
130       fDigitsInCluster++;\r
131       nRecPoints++;\r
132 \r
133       // Scanning for the neighbours\r
134       if(ScanForNeighbourDigits(i, fRecPointDataPtr) != 0)\r
135       {\r
136          return -1;\r
137       }\r
138 \r
139       //fUsedSize += sizeof(AliHLTCaloRecPointDataStruct) + (fDigitsInCluster-1)*sizeof(AliHLTCaloDigitDataStruct);   \r
140       \r
141       fRecPointDataPtr->fMultiplicity = fDigitsInCluster;     \r
142       fRecPointArray[fNRecPoints] = fRecPointDataPtr; \r
143       \r
144       fRecPointDataPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(fDigitIndexPtr);\r
145 \r
146       fNRecPoints++;\r
147       \r
148     }//end of clusterization\r
149 \r
150     return nRecPoints;\r
151 }\r
152 \r
153 Int_t\r
154 AliHLTCaloClusterizer::ScanForNeighbourDigits(Int_t index, AliHLTCaloRecPointDataStruct* recPoint)\r
155 {\r
156   //see header file for documentation\r
157   Int_t max = TMath::Min(fNDigits, (Int_t)fMaxDigitIndexDiff+index);\r
158   Int_t min = TMath::Max(0, (Int_t)(index - (Int_t)fMaxDigitIndexDiff));\r
159 \r
160   max = fNDigits;\r
161   min = 0;\r
162   for(Int_t j = min; j < max; j++)\r
163     {\r
164       if(fDigitsPointerArray[j]->fEnergy > fEmcMinEnergyThreshold)\r
165         {\r
166           if(j != index)\r
167             {\r
168               if(AreNeighbours(fDigitsPointerArray[index],\r
169                                fDigitsPointerArray[j]))\r
170                 {\r
171 //                if((fAvailableSize - fUsedSize) < sizeof(Int_t))\r
172 //                   {\r
173 //                      UChar_t *tmp = new UChar_t[fAvailableSize*2];\r
174 //                      memcpy(tmp, fRecPointDataPtr, fAvailableSize);\r
175 //                      for(Int_t n = 0; n < fNRecPoints; n++)\r
176 //                      {\r
177 //                         fRecPointArray[n] = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(reinterpret_cast<UChar_t*>(fRecPointArray[n]) - reinterpret_cast<UChar_t*>(fFirstRecPointPtr) + reinterpret_cast<UChar_t*>(tmp));\r
178 //                      }\r
179 //                      fRecPointDataPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(tmp);\r
180 //                      fFirstRecPointPtr = fRecPointDataPtr;\r
181 //                      fUsedSize = 0;\r
182 //                   }  \r
183                   CheckBuffer();\r
184                   \r
185                   // Assigning index to digit\r
186                   *fDigitIndexPtr = j;\r
187                   fUsedSize += sizeof(Int_t);\r
188                   \r
189                   // Incrementing digit pointer to be ready for new entry\r
190                   fDigitIndexPtr++;\r
191                   \r
192                   fRecPointDataPtr->fAmp += fDigitsPointerArray[j]->fEnergy;\r
193                   fDigitsPointerArray[j]->fEnergy = 0;        \r
194                   fDigitsInCluster++;\r
195                   ScanForNeighbourDigits(j, recPoint);\r
196                 }\r
197             }\r
198         }\r
199     }\r
200   return 0;\r
201 }\r
202 \r
203 Int_t \r
204 AliHLTCaloClusterizer::AreNeighbours(AliHLTCaloDigitDataStruct* digit1, \r
205                                             AliHLTCaloDigitDataStruct* digit2)\r
206 {\r
207   //see header file for documentation\r
208   if ( (digit1->fModule == digit2->fModule) /*&& (coord1[1]==coord2[1])*/ ) // inside the same PHOS module\r
209     { \r
210       Int_t rowdiff = TMath::Abs( digit1->fZ - digit2->fZ );  \r
211       Int_t coldiff = TMath::Abs( digit1->fX - digit2->fX ); \r
212       if (( coldiff <= 1   &&  rowdiff == 0 ) || ( coldiff == 0 &&  rowdiff <= 1 ))\r
213         {\r
214 //        cout << "Are neighbours: digit (E = "  << digit1->fEnergy << ") with x = " << digit1->fX << " and z = " << digit1->fZ << \r
215 //          " is neighbour with digit (E = " << digit2->fEnergy << ") with x = " << digit2->fX << " and z = " << digit2->fZ << endl;\r
216 \r
217           if(TMath::Abs(digit1->fTime - digit2->fTime ) < fEmcTimeGate)\r
218             {\r
219               return 1; \r
220             }\r
221         }\r
222 \r
223    /*   Float_t rowdiff = TMath::Abs( digit1->fZ - digit2->fZ );  \r
224       Float_t coldiff = TMath::Abs( digit1->fX - digit2->fX ); \r
225       if (( coldiff <= 2.4   &&  rowdiff < 0.4 ) || ( coldiff < 0.4 &&  rowdiff <= 2.4 ))\r
226         {\r
227           //      cout << "Are neighbours: digit (E = "  << digit1->fEnergy << ") with x = " << digit1->fX << " and z = " << digit1->fZ << \r
228           //        " is neighbour with digit (E = " << digit2->fEnergy << ") with x = " << digit2->fX << " and z = " << digit2->fZ << endl;\r
229 \r
230           if(TMath::Abs(digit1->fTime - digit2->fTime ) < fEmcTimeGate)\r
231             {\r
232               return 1; \r
233             }\r
234         }\r
235    */\r
236       else\r
237         {\r
238 //          cout << "Not neighbours: digit (E = "  << digit1->fEnergy << ") with x = " << digit1->fX << " and z = " << digit1->fZ << \r
239 //          " is not neighbour with digit (E = " << digit2->fEnergy << ") with x = " << digit2->fX << " and z = " << digit2->fZ << endl;\r
240         }\r
241     }\r
242   return 0;\r
243 }\r
244 \r
245 \r
246 \r
247 Int_t AliHLTCaloClusterizer::CheckArray()\r
248 {\r
249       if(fArraySize == fNRecPoints)\r
250         {\r
251            fArraySize *= 2;\r
252            AliHLTCaloRecPointDataStruct **tmp = new AliHLTCaloRecPointDataStruct*[fArraySize];\r
253            memcpy(tmp, fRecPointArray, fArraySize/2 * sizeof(AliHLTCaloRecPointDataStruct*));\r
254            delete [] fRecPointArray;\r
255            fRecPointArray = tmp;\r
256            //fRecPointArray[fNRecPoints-1] = fRecPointDataPtr;\r
257            //Int_t recPointOffset = reinterpret_cast<UChar_t*>(fRecPointDataPtr) - reinterpret_cast<UChar_t*>(fFirstRecPointPtr);\r
258            //fRecPointDataPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(reinterpret_cast<UChar_t*>(tmp) + recPointOffset);\r
259            \r
260         }\r
261    return 0;\r
262 }\r
263 \r
264 Int_t AliHLTCaloClusterizer::CheckBuffer()\r
265 {\r
266    // See header file for class documentation \r
267         if((fAvailableSize - fUsedSize) < (Int_t)sizeof(AliHLTCaloRecPointDataStruct))\r
268         {\r
269             Int_t recPointOffset = reinterpret_cast<UChar_t*>(fRecPointDataPtr) - reinterpret_cast<UChar_t*>(fFirstRecPointPtr);\r
270             Int_t digitIndexOffset = reinterpret_cast<UChar_t*>(fDigitIndexPtr) - reinterpret_cast<UChar_t*>(fRecPointDataPtr);\r
271             UChar_t *tmp = new UChar_t[fAvailableSize*2];\r
272 \r
273             memcpy(tmp, fFirstRecPointPtr, fUsedSize);\r
274             fAvailableSize *= 2;\r
275             for(Int_t n = 0; n < fNRecPoints; n++)\r
276             {\r
277                fRecPointArray[n] = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(reinterpret_cast<UChar_t*>(fRecPointArray[n]) - reinterpret_cast<UChar_t*>(fFirstRecPointPtr) + reinterpret_cast<UChar_t*>(tmp));\r
278             }\r
279             delete [] fFirstRecPointPtr;\r
280             fFirstRecPointPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(tmp);\r
281             fRecPointDataPtr = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(tmp + recPointOffset);\r
282             fDigitIndexPtr = reinterpret_cast<Int_t*>(reinterpret_cast<UChar_t*>(fRecPointDataPtr) + digitIndexOffset);\r
283             //fUsedSize = 0;\r
284         }\r
285    return 0;\r
286 }\r
287 \r
288 Int_t AliHLTCaloClusterizer::CheckDigits(AliHLTCaloRecPointDataStruct** recArray, AliHLTCaloDigitDataStruct** digitArray, Int_t nRP)\r
289 {\r
290   AliHLTCaloRecPointDataStruct **recpoints = recArray;\r
291   AliHLTCaloDigitDataStruct **digits = digitArray;\r
292   Int_t nRecPoints = nRP;\r
293   \r
294   if(recArray == 0)\r
295   {\r
296      recpoints = fRecPointArray;\r
297   }\r
298   if(digitArray == 0)\r
299   {\r
300      digits = fDigitsPointerArray;\r
301   }\r
302   if(nRP == 0)\r
303   {\r
304      nRecPoints = fNRecPoints;\r
305   }\r
306   for(Int_t i = 0; i < nRecPoints; i++)\r
307   {\r
308           \r
309      AliHLTCaloRecPointDataStruct *recPoint = recpoints[i];\r
310 \r
311      //AliHLTCaloRecPointDataStruct *recPoint = fRecPointArray[0];\r
312      Int_t multiplicity = recPoint->fMultiplicity;\r
313      Int_t *digitIndexPtr = &(recPoint->fDigits);\r
314      for(Int_t j = 0; j < multiplicity; j++)\r
315      {\r
316         //AliHLTCaloRecPointDataStruct *recPoint = fRecPointArray[j];\r
317 //      AliHLTCaloDigitDataStruct *digit = digits[*digitIndexPtr];\r
318         digitIndexPtr++;\r
319         //recPoint = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(digitIndexPtr);\r
320      }\r
321   }\r
322      \r
323      return 0;\r
324 }\r
325 \r
326 Int_t AliHLTCaloClusterizer::CheckDigits(AliHLTCaloRecPointDataStruct** recArray, AliHLTCaloDigitDataStruct* digitArray, Int_t nRP)\r
327 {\r
328   AliHLTCaloRecPointDataStruct **recpoints = recArray;\r
329   AliHLTCaloDigitDataStruct *digits = digitArray;\r
330   Int_t nRecPoints = nRP;\r
331   \r
332   if(recArray == 0)\r
333   {\r
334      recpoints = fRecPointArray;\r
335   }\r
336     if(nRP == 0)\r
337   {\r
338      nRecPoints = fNRecPoints;\r
339   }\r
340   for(Int_t i = 0; i < nRecPoints; i++)\r
341   {\r
342           \r
343      AliHLTCaloRecPointDataStruct *recPoint = recpoints[i];\r
344 \r
345      //AliHLTCaloRecPointDataStruct *recPoint = fRecPointArray[0];\r
346      Int_t multiplicity = recPoint->fMultiplicity;\r
347      Int_t *digitIndexPtr = &(recPoint->fDigits);\r
348      for(Int_t j = 0; j < multiplicity; j++)\r
349      {\r
350         //AliHLTCaloRecPointDataStruct *recPoint = fRecPointArray[j];\r
351         AliHLTCaloDigitDataStruct digit = digits[*digitIndexPtr];\r
352         digitIndexPtr++;\r
353         //recPoint = reinterpret_cast<AliHLTCaloRecPointDataStruct*>(digitIndexPtr);\r
354      }\r
355   }\r
356    return 0;\r
357      \r
358 }\r