]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalHistoCollector.cxx
- new gain calibb
[u/mrichter/AliRoot.git] / HLT / global / AliHLTGlobalHistoCollector.cxx
1
2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project        *
4 //* ALICE Experiment at CERN, All rights reserved.                         *
5 //*                                                                        *
6 //* Primary Authors: Kalliopi Kanaki <Kalliopi.Kanaki@ift.uib.no>          *
7 //*                  for The ALICE HLT Project.                            *
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 /** @file   AliHLTGlobalHistoCollector.cxx
19     @author Kalliopi Kanaki
20     @date   
21     @brief  The Histogram Handler component
22 */
23
24 #if __GNUC__>= 3
25 using namespace std;
26 #endif
27 #include "AliHLTGlobalHistoCollector.h"
28 #include "AliCDBEntry.h"
29 #include "AliCDBManager.h"
30 #include "TString.h"
31 #include "TObjArray.h"
32 #include "TObjString.h"
33 #include "TH1.h"
34 #include "TTimeStamp.h"
35 #include "TSystem.h"
36
37 ClassImp(AliHLTGlobalHistoCollector) //ROOT macro for the implementation of ROOT specific class methods
38
39   AliHLTGlobalHistoCollector::AliHLTGlobalHistoCollector()
40     :    
41     fUID(0),
42     fStore(),
43     fBenchmark("GlobalHistoCollector")
44 {
45   // see header file for class documentation
46   // or
47   // refer to README to build package
48   // or
49   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
50 }
51
52 AliHLTGlobalHistoCollector::~AliHLTGlobalHistoCollector() { 
53   // see header file for class documentation
54   Clear();
55 }
56
57
58
59 // Public functions to implement AliHLTComponent's interface.
60 // These functions are required for the registration process
61
62 const char* AliHLTGlobalHistoCollector::GetComponentID() 
63
64   // see header file for class documentation
65   return "GlobalHistoCollector";
66 }
67
68 void AliHLTGlobalHistoCollector::GetInputDataTypes( vector<AliHLTComponentDataType>& list) 
69
70   // see header file for class documentation
71
72   list.clear(); 
73   list.push_back( kAliHLTAllDataTypes );
74 }
75
76 AliHLTComponentDataType AliHLTGlobalHistoCollector::GetOutputDataType() 
77
78   // see header file for class documentation
79   return kAliHLTAllDataTypes;
80 }
81
82
83 void AliHLTGlobalHistoCollector::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) 
84
85   // see header file for class documentation
86   constBase=0;
87   inputMultiplier=1.0;
88 }
89
90 AliHLTComponent* AliHLTGlobalHistoCollector::Spawn() 
91
92   // see header file for class documentation
93   return new AliHLTGlobalHistoCollector();
94 }
95         
96 int AliHLTGlobalHistoCollector::DoInit( int argc, const char** argv ) 
97
98   // see header file for class documentation
99
100   Clear(); 
101   fBenchmark.Reset();
102   fBenchmark.SetTimer(0,"total");
103   fBenchmark.SetTimer(1,"merging");
104
105   int iResult=0;
106   
107   TString configuration="";
108   TString argument="";
109   for (int j=0; j<argc && iResult>=0; j++) {
110     
111     argument=argv[j];
112     if (!configuration.IsNull()) configuration+=" ";
113     configuration+=argument;    
114   }
115    
116   if (!configuration.IsNull()) {
117     iResult=Configure(configuration.Data());
118   } else {
119     iResult=Reconfigure(NULL, NULL);
120   }
121   fUID = 0;
122   return iResult;
123 }
124
125
126 int AliHLTGlobalHistoCollector::DoDeinit() 
127
128   // see header file for class documentation 
129     
130   Clear();
131   fUID = 0;
132   return 0;
133 }
134
135 int AliHLTGlobalHistoCollector::Configure(const char* arguments) 
136
137   // see header file for class documentation
138   
139   int iResult=0;
140   if (!arguments) return iResult;
141   HLTInfo("parsing configuration string \'%s\'", arguments);
142
143   TString allArgs=arguments;
144   TString argument;
145   int bMissingParam=0;
146
147   TObjArray* pTokens=allArgs.Tokenize(" ");
148   if (pTokens) {
149     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
150       argument=((TObjString*)pTokens->At(i))->GetString();
151       if (argument.IsNull()) continue;
152      
153       //if (argument.CompareTo("-sum-noise-histograms")==0) {
154       //fNoiseHistograms = kTRUE;
155       //HLTInfo("got \'-sum-noise-histograms\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
156       //else 
157       {
158         HLTError("unknown argument %s", argument.Data());
159         iResult=-EINVAL;
160         break;
161       }
162     } // end for
163     
164     delete pTokens;
165   
166   } // end if pTokens
167   
168   if (bMissingParam) {
169     HLTError("missing parameter for argument %s", argument.Data());
170     iResult=-EINVAL;
171   }
172   return iResult;
173 }
174
175
176 int AliHLTGlobalHistoCollector::Reconfigure(const char* /*cdbEntry*/, const char* /*chainId*/) { 
177   // see header file for class documentation
178   
179   return 0; // no CDB entry exist
180   /*
181   int iResult=0;  
182   const char* path="HLT/ConfigGlobal/GlobalHistoCollector";
183   const char* defaultNotify="";
184   if (cdbEntry) {
185     path=cdbEntry;
186     defaultNotify=" (default)";
187   }
188   
189   if (path) {
190     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
191     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
192     if (pEntry) {
193       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
194       if (pString) {
195         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
196         iResult=Configure(pString->GetString().Data());
197       } else {
198         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
199       }
200     } else {
201       HLTError("cannot fetch object \"%s\" from CDB", path);
202     }
203   }
204   return iResult;
205 */
206 }
207
208
209 void AliHLTGlobalHistoCollector::Clear() 
210
211   // reset the store
212
213   for( unsigned int i=0; i<fStore.size(); i++ ){
214     for( unsigned int j=0; j<fStore[i].fInstances.size(); j++ ){
215       delete fStore[i].fInstances[j].fObject;
216     }
217     delete fStore[i].fMergedObject;
218   }
219   fStore.clear();
220 }
221
222
223
224
225 int AliHLTGlobalHistoCollector::DoEvent(const AliHLTComponentEventData & evtData, AliHLTComponentTriggerData& /*trigData*/)
226 {
227   // see header file for class documentation
228   //cout<<"\n\nDoEvent called"<<endl;
229
230   if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
231
232   fBenchmark.StartNewEvent();
233   fBenchmark.Start(0);
234
235   if( fUID == 0 ){
236     TTimeStamp t;
237     fUID = ( gSystem->GetPid() + t.GetNanoSec())*10 + evtData.fEventID;
238   }
239
240   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(); i!=NULL; i=GetNextInputBlock() ){
241     fBenchmark.AddInput(i->fSize);
242   }
243
244   const TObject *iter = NULL;
245   for(iter = GetFirstInputObject(); iter != NULL; iter = GetNextInputObject()){
246             
247     if( !iter->InheritsFrom(TH1::Class()) 
248         && !iter->InheritsFrom(TSeqCollection::Class()) ) continue;
249
250     //cout<<"received object "<<iter->GetName()<<" with id="<<GetSpecification(iter)<<endl;
251
252     //search for the base entry, if not exist then create a new entry   
253
254     int iColl = -1;
255     for( unsigned int i=0; i<fStore.size(); i++ ){
256       if( fStore[i].fHLTDataType != GetDataType(iter) ) continue;
257       if( fStore[i].fInstances.size()<1 ) continue; 
258       TObject * obj = fStore[i].fInstances[0].fObject;
259       if( !obj ) continue;
260       if( TString(obj->GetName()).CompareTo(iter->GetName())==0){
261         iColl = i;
262         break;
263       }
264     }
265     //cout<<"Collection found: "<<iColl<<endl;
266     if( iColl<0 ){
267       AliHLTGlobalHCCollection c;
268       c.fHLTDataType = GetDataType(iter);
269       c.fMergedObject = 0;
270       c.fNeedToMerge = 1;
271       fStore.push_back(c);
272       iColl = fStore.size()-1;
273     }else{
274       fStore[iColl].fNeedToMerge = 1;
275     }
276
277     // search for the specific entry, if not exist then create a new one
278     
279     AliHLTGlobalHCCollection &c = fStore[iColl];
280    
281     int iSpec=-1;
282     for( unsigned int i=0; i<c.fInstances.size(); i++ ){
283       AliHLTGlobalHCInstance &inst = c.fInstances[i];
284       if( inst.fHLTSpecification == GetSpecification(iter) ){
285         iSpec = i;
286         break;
287       }
288     }
289     //cout<<"Instance found:"<<iSpec<<endl;
290     if( iSpec<0 ){
291       AliHLTGlobalHCInstance inst;
292       inst.fHLTSpecification = GetSpecification(iter);
293       inst.fObject = 0;
294       c.fInstances.push_back(inst);
295       iSpec = c.fInstances.size()-1;      
296     }else{
297       delete c.fInstances[iSpec].fObject;
298     }
299
300     c.fInstances[iSpec].fObject = iter->Clone();
301     
302     //cout<<"index = "<<iColl<<","<<iSpec<<endl;
303
304   } // end for loop over input blocks
305
306   fBenchmark.Start(1);
307  
308   
309   // merge histos 
310
311   for( unsigned int iColl = 0; iColl<fStore.size(); iColl++){
312     AliHLTGlobalHCCollection &c = fStore[iColl];
313     if( !c.fNeedToMerge && c.fMergedObject ) continue;
314     if( c.fInstances.size() <1 ) continue;
315     delete c.fMergedObject;
316     c.fMergedObject = c.fInstances[0].fObject->Clone();
317     TList l;
318     for( unsigned int i=1; i<c.fInstances.size(); i++ ){
319       l.Add(c.fInstances[i].fObject);
320     }
321
322     if( c.fMergedObject->InheritsFrom(TH1::Class()) ){
323       TH1 *histo = dynamic_cast<TH1*>(c.fMergedObject);
324       if( histo ) histo->Merge(&l);
325     }
326     else if( c.fMergedObject->InheritsFrom(TSeqCollection::Class()) ){
327       TSeqCollection *list = dynamic_cast<TSeqCollection*>(c.fMergedObject);
328       if( list ) list->Merge(&l);
329     }    
330     c.fNeedToMerge = 0;
331   }
332   fBenchmark.Stop(1);
333  
334   // send output 
335
336   for( unsigned int i=0; i<fStore.size(); i++ ){
337     if( fStore[i].fMergedObject ){
338       PushBack((TObject*) fStore[i].fMergedObject, fStore[i].fHLTDataType, fUID );
339       fBenchmark.AddOutput(GetLastObjectSize());
340     }
341   }
342
343   fBenchmark.Stop(0);
344   HLTInfo(fBenchmark.GetStatistics());
345
346   return 0;
347 }
348