]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/PHOS/AliHLTPHOSRcuHistogramProducerComponent.cxx
New classes that updates histograms contniously during run.
[u/mrichter/AliRoot.git] / HLT / PHOS / AliHLTPHOSRcuHistogramProducerComponent.cxx
1 /**************************************************************************
2  * Copyright(c) 2006, ALICE Experiment at CERN, All rights reserved.      *
3  *                                                                        *
4  * Authors: Boris Polichtchouk & Per Thomas Hille for the ALICE           *
5  * offline/HLT Project. Contributors are mentioned in the code where      *
6  * appropriate.                                                           *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 #include <iostream>
18 #include "stdio.h"
19 #include <cstdlib>
20 #include "AliHLTPHOSRcuCellEnergyDataStruct.h"
21 //#include "AliHLTPHOSModuleCellAccumulatedEnergyDataStruct.h"
22 #include "AliHLTPHOSRcuHistogramProducer.h"
23 #include "AliHLTPHOSRcuHistogramProducerComponent.h"
24
25
26
27 const AliHLTComponentDataType  AliHLTPHOSRcuHistogramProducerComponent::inputDataTypes[]={kAliHLTVoidDataType,{0,"",""}}; //'zero' terminated array
28 const AliHLTComponentDataType  AliHLTPHOSRcuHistogramProducerComponent::outputDataType=kAliHLTVoidDataType;
29 AliHLTPHOSRcuHistogramProducerComponent gAliHLTPHOSRcuHistogramProducerComponent;
30
31 //AliHLTPHOSHistogramProducerComponent gAliHLTPHOSHistogramProducerComponent;
32 /*************************************************************************
33 * Class AliHLTPHOSRcuHistogramProducerComponent accumulating histograms  *
34 * with amplitudes per PHOS channel                                       *
35 * It is intended to run at the HLT farm                                  *
36 * and it fills the histograms with amplitudes per channel.               * 
37 * Usage example see in PHOS/macros/Shuttle/AliPHOSCalibHistoProducer.C   *
38 **************************************************************************/
39 AliHLTPHOSRcuHistogramProducerComponent:: AliHLTPHOSRcuHistogramProducerComponent():AliHLTProcessor(), fEventCount(0), fRcuHistoProducerPtr(0)
40 {
41   //  Reset();
42
43
44
45
46 AliHLTPHOSRcuHistogramProducerComponent::~ AliHLTPHOSRcuHistogramProducerComponent()
47 {
48
49 }
50
51
52 AliHLTPHOSRcuHistogramProducerComponent::AliHLTPHOSRcuHistogramProducerComponent(const  AliHLTPHOSRcuHistogramProducerComponent & ) : AliHLTProcessor(), fEventCount(0), fRcuHistoProducerPtr(0)
53 {
54
55 }
56
57
58 int 
59 AliHLTPHOSRcuHistogramProducerComponent::Deinit()
60 {
61   cout << "AliHLTPHOSRcuHistogramProducerComponent::Deinit()" << endl;
62   fRcuHistoProducerPtr->WriteEnergyHistograms();
63   return 0;
64 }
65
66
67 int 
68 AliHLTPHOSRcuHistogramProducerComponent::DoDeinit()
69 {
70   Logging(kHLTLogInfo, "HLT", "PHOS", ",AliHLTPHOSRcuHistogramProducer DoDeinit");
71   return 0;
72 }
73
74
75 const char* 
76 AliHLTPHOSRcuHistogramProducerComponent::GetComponentID()
77 {
78   return "RcuHistogramProducer";
79 }
80
81
82 void
83  AliHLTPHOSRcuHistogramProducerComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
84 {
85   const AliHLTComponentDataType* pType=inputDataTypes;
86   while (pType->fID!=0) 
87     {
88       list.push_back(*pType);
89       pType++;
90     }
91 }
92
93
94 AliHLTComponentDataType 
95 AliHLTPHOSRcuHistogramProducerComponent::GetOutputDataType()
96 {
97   return AliHLTPHOSDefinitions::gkCellEnergyDataType;
98 }
99
100
101 void
102 AliHLTPHOSRcuHistogramProducerComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier )
103 {
104   constBase = 30;
105   inputMultiplier = 1;
106 }
107
108
109 int  AliHLTPHOSRcuHistogramProducerComponent::DoEvent( const AliHLTComponentEventData& evtData, const AliHLTComponentBlockData* blocks, 
110                                               AliHLTComponentTriggerData& trigData, AliHLTUInt8_t* outputPtr, 
111                                               AliHLTUInt32_t& size, vector<AliHLTComponentBlockData>& outputBlocks )
112 {
113   unsigned long ndx       = 0;
114   UInt_t offset           = 0; 
115   UInt_t mysize           = 0;
116   UInt_t tSize            = 0;
117   const AliHLTComponentBlockData* iter = NULL;   
118   AliHLTPHOSRcuCellEnergyDataStruct *cellDataPtr;
119   AliHLTUInt8_t* outBPtr;
120  
121   // outBPtr = outputPtr;
122   // fOutPtr =  (AliHLTPHOSRcuCellAccumulatedEnergyDataStruct*)outBPtr;
123   int tmpCnt;
124
125   for( ndx = 0; ndx < evtData.fBlockCnt; ndx++ )
126     {
127       iter = blocks+ndx;
128       cellDataPtr = (AliHLTPHOSRcuCellEnergyDataStruct*)( iter->fPtr);
129       tmpCnt =  cellDataPtr->fCnt;
130
131       for(int i= 0; i <= tmpCnt; i ++)
132         {
133           fRcuHistoProducerPtr->FillEnergy(cellDataPtr->fValidData[i].fX,
134                                            cellDataPtr->fValidData[i].fZ, 
135                                            cellDataPtr->fValidData[i].fGain, 
136                                            cellDataPtr->fValidData[i].fEnergy);
137         }
138     }
139   
140   outBPtr = outputPtr;
141   fOutPtr =  (AliHLTPHOSRcuCellAccumulatedEnergyDataStruct*)outBPtr;
142   const AliHLTPHOSRcuCellAccumulatedEnergyDataStruct  &innPtr = fRcuHistoProducerPtr->GetCellAccumulatedEnergies();
143
144   fOutPtr->fModuleID = fModuleID;
145   fOutPtr->fRcuX     = fRcuX;
146   fOutPtr->fRcuZ     = fRcuZ;
147
148
149   for(int x=0; x < N_XCOLUMNS_RCU; x ++)
150     {
151       for(int z=0; z < N_XCOLUMNS_RCU; z ++)
152         {
153           for(int gain =0;  gain < N_GAINS; gain ++)
154             {
155               fOutPtr->fAccumulatedEnergies[x][z][gain] = innPtr.fAccumulatedEnergies[x][z][gain];
156               fOutPtr->fHits[x][z][gain] = innPtr.fHits[x][z][gain];
157             }
158         }
159     }
160
161
162   //pushing data to shared output memory
163   mysize += sizeof(AliHLTPHOSRcuCellAccumulatedEnergyDataStruct);
164   AliHLTComponentBlockData bd;
165   FillBlockData( bd );
166   bd.fOffset = offset;
167   bd.fSize = mysize;
168   bd.fDataType = AliHLTPHOSDefinitions::gkCellAccumulatedEnergyDataType;
169   bd.fSpecification = 0xFFFFFFFF;
170   outputBlocks.push_back( bd );
171   tSize += mysize;
172   outBPtr += mysize;
173
174   if( tSize > size )
175     {
176       Logging( kHLTLogFatal, "HLT::AliHLTRcuHistogramProducerComponent::DoEvent", "Too much data",
177                "Data written over allowed buffer. Amount written: %lu, allowed amount: %lu."
178                , tSize, size );
179       return EMSGSIZE;
180     }
181
182   fEventCount++; 
183   return 0;
184 }//end DoEvent
185
186
187 int
188 AliHLTPHOSRcuHistogramProducerComponent::DoInit( int argc, const char** argv )
189 {
190   int iResult=0;
191   TString argument="";
192   //  fRcuHistoProducerPtr = new AliHLTPHOSRcuHistogramProducer();
193   AliHLTUInt8_t tmpRcuX;
194   AliHLTUInt8_t tmpRcuZ; 
195   AliHLTUInt8_t tmpModuleID;
196   Bool_t isSetEquippmentID = kFALSE;
197
198   for(int i=0; i<argc && iResult>=0; i++) 
199     {
200       argument=argv[i];
201       
202       if (argument.IsNull()) 
203         {
204           continue;
205         }
206       if (argument.CompareTo("-equipmentID") == 0) 
207         {
208           if(i+1 <= argc)
209             {
210               fEquippmentID = atoi(argv[i+1]);
211               isSetEquippmentID = kTRUE;
212             }
213           else
214             {
215                iResult= -1;
216                Logging( kHLTLogFatal, "HLT::AliHLTPHOSRcuHistogramProducerComponent::DoInt( int argc, const char** argv )", "Missing argument",
217                         "The argument -equippmentID expects a number");
218                return  iResult;   
219             }
220         }
221
222       int rcuIndex =  (fEquippmentID - 1792)%N_RCUS_PER_MODULE;
223       //     fModuleID = (fEquippmentID  -1792 -rcuIndex)/N_RCUS_PER_MODULE;
224       tmpModuleID = ((fEquippmentID -1792 -rcuIndex)/N_RCUS_PER_MODULE);
225       SetModuleID(tmpModuleID);
226
227       if(rcuIndex == 0)
228         {
229           tmpRcuX = 0; 
230           tmpRcuZ = 0;
231         }
232       
233       if(rcuIndex == 1)
234         {
235          tmpRcuX = 0; 
236          tmpRcuZ = 1;
237         }
238       
239       if(rcuIndex == 2)
240         {
241           tmpRcuX = 1; 
242           tmpRcuZ = 0;
243         }
244       
245       if(rcuIndex == 3)
246         {
247           tmpRcuX = 1; 
248           tmpRcuZ = 1;
249         }
250
251       SetRcuX(tmpRcuX);
252       SetRcuZ(tmpRcuZ); 
253       cout <<"********InitInfo************"<< endl;
254       cout <<"AliHLTPHOSRcuHistogramProducerComponent::SetCoordinate"<< endl;
255       cout <<"Equpippment ID =\t"<< fEquippmentID <<endl;
256       cout <<"Module ID =\t"<<  (int) tmpModuleID<<endl;
257       cout <<"RCUX =\t\t" << (int)tmpRcuX << endl;
258       cout <<"RCUZ =\t\t" << (int)tmpRcuZ << endl;
259      }
260
261   if(isSetEquippmentID == kFALSE)
262     {
263        Logging( kHLTLogFatal, "HLT::AliHLTPHOSRcuHistogramProducerComponent::DoInt( int argc, const char** argv )", "Missing argument",
264                         "The argument equippmentID is not set: set it with a component argumet like this: -equippmentID  <number>");
265        iResult = -2; 
266     }
267
268   
269   fRcuHistoProducerPtr = new AliHLTPHOSRcuHistogramProducer( tmpModuleID, tmpRcuX, tmpRcuZ);
270
271
272   return  iResult;
273 }
274
275
276 void 
277 AliHLTPHOSRcuHistogramProducerComponent::SetRcuX(AliHLTUInt8_t X)
278 {
279   fRcuX = X;
280 }
281
282
283 void 
284 AliHLTPHOSRcuHistogramProducerComponent::SetRcuZ(AliHLTUInt8_t Z)
285 {
286   fRcuZ = Z;
287 }
288
289
290 void 
291 AliHLTPHOSRcuHistogramProducerComponent::SetModuleID(AliHLTUInt8_t moduleID)
292 {
293   fModuleID = moduleID;
294 }
295
296
297 void 
298 AliHLTPHOSRcuHistogramProducerComponent::SetEquippmentId(int id)
299 {
300   fEquippmentID = id;
301   fRcuHistoProducerPtr->SetEquippmentId(id);
302 }
303
304
305 int 
306 AliHLTPHOSRcuHistogramProducerComponent::GetEquippmentId()
307 {
308   return  fEquippmentID;
309 }
310
311
312 AliHLTComponent*
313 AliHLTPHOSRcuHistogramProducerComponent::Spawn()
314 {
315   return new AliHLTPHOSRcuHistogramProducerComponent;
316 }
317
318