]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTComponentHandler.cxx
get rid of compilation warnings
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTComponentHandler.cxx
1 // $Id$
2
3 /**************************************************************************
4  * This file is property of and copyright by the ALICE HLT Project        * 
5  * ALICE Experiment at CERN, All rights reserved.                         *
6  *                                                                        *
7  * Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8  *                  Timm Steinbeck <timm@kip.uni-heidelberg.de>           *
9  *                  for The ALICE HLT Project.                            *
10  *                                                                        *
11  * Permission to use, copy, modify and distribute this software and its   *
12  * documentation strictly for non-commercial purposes is hereby granted   *
13  * without fee, provided that the above copyright notice appears in all   *
14  * copies and that both the copyright notice and this permission notice   *
15  * appear in the supporting documentation. The authors make no claims     *
16  * about the suitability of this software for any purpose. It is          *
17  * provided "as is" without express or implied warranty.                  *
18  **************************************************************************/
19
20 /** @file   AliHLTComponentHandler.cxx
21     @author Matthias Richter, Timm Steinbeck
22     @date   
23     @brief  Implementation of HLT component handler. */
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 #if __GNUC__>= 3
32 using namespace std;
33 #endif
34 //#undef HAVE_DLFCN_H
35 #ifdef HAVE_DLFCN_H
36 #include <dlfcn.h>
37 #else
38 //#include <Riostream.h>
39 #include <TSystem.h>
40 #endif //HAVE_DLFCN_H
41 //#include "AliHLTStdIncludes.h"
42 #include "AliHLTComponentHandler.h"
43 #include "AliHLTComponent.h"
44 #include "AliHLTDataTypes.h"
45 #include "AliHLTModuleAgent.h"
46 #include "TString.h"
47
48 /** ROOT macro for the implementation of ROOT specific class methods */
49 ClassImp(AliHLTComponentHandler)
50
51 AliHLTComponentHandler::AliHLTComponentHandler()
52   :
53   fComponentList(),
54   fScheduleList(),
55   fLibraryList(),
56   fEnvironment(),
57   fOwnedComponents(),
58   fLibraryMode(kDynamic)
59 {
60   // see header file for class documentation
61   // or
62   // refer to README to build package
63   // or
64   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
65   memset(&fEnvironment, 0, sizeof(AliHLTComponentEnvironment));
66   AddStandardComponents();
67 }
68
69 AliHLTComponentHandler::AliHLTComponentHandler(AliHLTComponentEnvironment* pEnv)
70   :
71   AliHLTLogging(),
72   fComponentList(),
73   fScheduleList(),
74   fLibraryList(),
75   fEnvironment(),
76   fOwnedComponents(),
77   fLibraryMode(kDynamic)
78 {
79   // see header file for class documentation
80   if (pEnv) {
81     memcpy(&fEnvironment, pEnv, sizeof(AliHLTComponentEnvironment));
82     if (pEnv->fLoggingFunc) {
83       // the AliHLTLogging::Init method also sets the stream output
84       // and notification handler to AliLog. This should only be done
85       // if the logging environment contains a logging function
86       // for redirection
87       AliHLTLogging::Init(pEnv->fLoggingFunc);
88     }
89   }  else {
90     memset(&fEnvironment, 0, sizeof(AliHLTComponentEnvironment));
91   }
92   //#ifndef __DEBUG
93   //SetLocalLoggingLevel(kHLTLogError);
94   //#else
95   //SetLocalLoggingLevel(kHLTLogInfo);
96   //#endif
97
98   AddStandardComponents();
99 }
100
101 AliHLTComponentHandler::~AliHLTComponentHandler()
102 {
103   // see header file for class documentation
104   DeleteOwnedComponents();
105   UnloadLibraries();
106 }
107
108 int AliHLTComponentHandler::AnnounceVersion()
109 {
110   // see header file for class documentation
111   int iResult=0;
112 #ifdef PACKAGE_STRING
113   void HLTbaseCompileInfo( char*& date, char*& time);
114   char* date="";
115   char* time="";
116   HLTbaseCompileInfo(date, time);
117   if (!date) date="unknown";
118   if (!time) time="unknown";
119   HLTInfo("%s build on %s (%s)", PACKAGE_STRING, date, time);
120 #else
121   HLTInfo("ALICE High Level Trigger build on %s (%s) (embedded AliRoot build)", __DATE__, __TIME__);
122 #endif
123   return iResult;
124 }
125
126 Int_t AliHLTComponentHandler::AddComponent(AliHLTComponent* pSample)
127 {
128   // see header file for class documentation
129   Int_t iResult=0;
130   if (pSample==NULL) return -EINVAL;
131   if ((iResult=RegisterComponent(pSample))>=0) {
132     //HLTDebug("sample %s (%p) managed by handler", pSample->GetComponentID(), pSample);
133     fOwnedComponents.push_back(pSample);
134   }
135   return iResult;
136 }
137
138 Int_t AliHLTComponentHandler::RegisterComponent(AliHLTComponent* pSample)
139 {
140   // see header file for class documentation
141   Int_t iResult=0;
142   if (pSample) {
143     if (FindComponent(pSample->GetComponentID())==NULL) {
144       iResult=InsertComponent(pSample);
145       if (iResult>=0) {
146         HLTInfo("component %s registered", pSample->GetComponentID());
147       }
148     } else {
149       // component already registered
150       HLTDebug("component %s already registered, skipped", pSample->GetComponentID());
151       iResult=-EEXIST;
152     }
153   } else {
154     iResult=-EINVAL;
155   }
156   return iResult;
157 }
158
159 int AliHLTComponentHandler::DeregisterComponent( const char* componentID )
160 {
161   // see header file for class documentation
162
163   int iResult=0;
164   if (componentID) {
165     HLTWarning("not yet implemented, please notify the developers if you need this function");
166   } else {
167     iResult=-EINVAL;
168   }
169   return iResult;
170 }
171
172 Int_t AliHLTComponentHandler::ScheduleRegister(AliHLTComponent* pSample)
173 {
174   // see header file for class documentation
175   Int_t iResult=0;
176   if (pSample) {
177     fScheduleList.push_back(pSample);
178   } else {
179     iResult=-EINVAL;
180   }
181   return iResult;
182 }
183
184 int AliHLTComponentHandler::CreateComponent(const char* componentID, void* pEnvParam, int argc, const char** argv, AliHLTComponent*& component )
185 {
186   // see header file for class documentation
187   int iResult=0;
188   if (componentID) {
189     AliHLTComponent* pSample=FindComponent(componentID);
190     if (pSample!=NULL) {
191       component=pSample->Spawn();
192       if (component) {
193         HLTDebug("component \"%s\" created (%p)", componentID, component);
194         if ((iResult=component->Init(&fEnvironment, pEnvParam, argc, argv))!=0) {
195           HLTError("Initialization of component \"%s\" failed with error %d", componentID, iResult);
196           delete component;
197           component=NULL;
198         }
199       } else {
200         HLTError("can not spawn component \"%s\"", componentID);
201         iResult=-ENOENT;
202       }
203     } else {
204       HLTWarning("can not find component \"%s\"", componentID);
205       iResult=-ENOENT;
206     }
207   } else {
208     iResult=-EINVAL;
209   }
210   return iResult;
211 }
212
213 Int_t AliHLTComponentHandler::FindComponentIndex(const char* componentID)
214 {
215   // see header file for class documentation
216   Int_t iResult=0;
217   if (componentID) {
218     AliHLTComponentPList::iterator element=fComponentList.begin();
219     while (element!=fComponentList.end() && iResult>=0) {
220       if (strcmp(componentID, (*element)->GetComponentID())==0) {
221         break;
222       }
223       element++;
224       iResult++;
225     }
226     if (element==fComponentList.end()) iResult=-ENOENT;
227   } else {
228     iResult=-EINVAL;
229   }
230   return iResult;
231 }
232
233 AliHLTComponent* AliHLTComponentHandler::FindComponent(const char* componentID)
234 {
235   // see header file for class documentation
236   AliHLTComponent* pSample=NULL;
237   Int_t index=FindComponentIndex(componentID);
238   if (index>=0) {
239     pSample=(AliHLTComponent*)fComponentList.at(index);
240   }
241   return pSample;
242 }
243
244 Int_t AliHLTComponentHandler::InsertComponent(AliHLTComponent* pSample)
245 {
246   // see header file for class documentation
247   Int_t iResult=0;
248   if (pSample!=NULL) {
249     fComponentList.push_back(pSample);
250   } else {
251     iResult=-EINVAL;
252   }
253   return iResult;
254 }
255
256 void AliHLTComponentHandler::List() 
257 {
258   // see header file for class documentation
259   AliHLTComponentPList::iterator element=fComponentList.begin();
260   int index=0;
261   while (element!=fComponentList.end()) {
262     HLTInfo("%d. %s", index++, (*element++)->GetComponentID());
263   }
264 }
265
266 int AliHLTComponentHandler::HasOutputData( const char* componentID)
267 {
268   // see header file for class documentation
269   int iResult=0;
270   AliHLTComponent* pSample=FindComponent(componentID);
271   if (pSample) {
272     AliHLTComponent::TComponentType ct=AliHLTComponent::kUnknown;
273     ct=pSample->GetComponentType();
274     iResult=(ct==AliHLTComponent::kSource || ct==AliHLTComponent::kProcessor);
275   } else {
276     iResult=-ENOENT;
277   }
278   return iResult;
279 }
280
281 void AliHLTComponentHandler::SetEnvironment(AliHLTComponentEnvironment* pEnv) 
282 {
283   // see header file for class documentation
284   if (pEnv) {
285     memcpy(&fEnvironment, pEnv, sizeof(AliHLTComponentEnvironment));
286     if (fEnvironment.fLoggingFunc) {
287       // the AliHLTLogging::Init method also sets the stream output
288       // and notification handler to AliLog. This should only be done
289       // if the logging environment contains a logging function
290       // for redirection
291       AliHLTLogging::Init(fEnvironment.fLoggingFunc);
292     }
293   }
294 }
295
296 AliHLTComponentHandler::TLibraryMode AliHLTComponentHandler::SetLibraryMode(TLibraryMode mode)
297 {
298   // see header file for class documentation
299   TLibraryMode old=fLibraryMode;
300   fLibraryMode=mode;
301   return old;
302 }
303
304 int AliHLTComponentHandler::LoadLibrary( const char* libraryPath, int bActivateAgents)
305 {
306   // see header file for class documentation
307   int iResult=0;
308   if (libraryPath) {
309     // first activate all agents which are already loaded
310     if (bActivateAgents) ActivateAgents();
311
312     // set the global component handler for static component registration
313     AliHLTComponent::SetGlobalComponentHandler(this);
314
315     AliHLTLibHandle hLib;
316     const char* loadtype="";
317 #ifdef HAVE_DLFCN_H
318     // use interface to the dynamic linking loader
319     try {
320       hLib.fHandle=dlopen(libraryPath, RTLD_NOW);
321       loadtype="dlopen";
322     }
323     catch (...) {
324       // error message printed further down
325       loadtype="dlopen exeption";
326     }
327 #else
328     // use ROOT dynamic loader
329     // check if the library was already loaded, as Load returns
330     // 'failure' if the library was already loaded
331     try {
332     AliHLTLibHandle* pLib=FindLibrary(libraryPath);
333     if (pLib) {
334         int* pRootHandle=reinterpret_cast<int*>(pLib->fHandle);
335         (*pRootHandle)++;
336         HLTDebug("instance %d of library %s loaded", (*pRootHandle), libraryPath);
337         hLib.fHandle=pRootHandle;
338     }
339     
340     if (hLib.fHandle==NULL && gSystem->Load(libraryPath)==0) {
341       int* pRootHandle=new int;
342       if (pRootHandle) *pRootHandle=1;
343       hLib.fHandle=pRootHandle;
344       //HLTDebug("library %s loaded via gSystem", libraryPath);
345     }
346     loadtype="gSystem";
347     }
348     catch (...) {
349       // error message printed further down
350       loadtype="gSystem exeption";
351     }
352 #endif //HAVE_DLFCN_H
353     if (hLib.fHandle!=NULL) {
354       // create TString object to store library path and use pointer as handle 
355       hLib.fName=new TString(libraryPath);
356       hLib.fMode=fLibraryMode;
357       HLTInfo("library %s loaded (%s%s)", libraryPath, hLib.fMode==kStatic?"persistent, ":"", loadtype);
358       fLibraryList.insert(fLibraryList.begin(), hLib);
359       typedef void (*CompileInfo)( char*& date, char*& time);
360       CompileInfo fctInfo=(CompileInfo)FindSymbol(libraryPath, "CompileInfo");
361       if (fctInfo) {
362         char* date="";
363         char* time="";
364         (*fctInfo)(date, time);
365         if (!date) date="unknown";
366         if (!time) time="unknown";
367         HLTInfo("build on %s (%s)", date, time);
368       } else {
369         HLTInfo("no build info available (possible AliRoot embedded build)");
370       }
371
372       // static registration of components when library is loaded
373       iResult=RegisterScheduledComponents();
374
375     } else {
376       HLTError("can not load library %s (%s)", libraryPath, loadtype);
377 #ifdef HAVE_DLFCN_H
378       HLTError("dlopen error: %s", dlerror());
379 #endif //HAVE_DLFCN_H
380 #ifdef __APPLE__
381       iResult=-EFTYPE;
382 #else
383       iResult=-ELIBACC;
384 #endif
385     }
386     AliHLTComponent::UnsetGlobalComponentHandler();
387     
388     if (iResult>=0) {
389       // alternative dynamic registration by library agents
390       // !!! has to be done after UnsetGlobalComponentHandler
391       if (bActivateAgents) ActivateAgents();
392     }
393
394   } else {
395     iResult=-EINVAL;
396   }
397   return iResult;
398 }
399
400 int AliHLTComponentHandler::UnloadLibrary( const char* libraryPath )
401 {
402   // see header file for class documentation
403   int iResult=0;
404   if (libraryPath) {
405     vector<AliHLTLibHandle>::iterator element=fLibraryList.begin();
406     while (element!=fLibraryList.end()) {
407       TString* pName=reinterpret_cast<TString*>((*element).fName);
408       if (pName->CompareTo(libraryPath)==0) {
409         UnloadLibrary(*element);
410         fLibraryList.erase(element);
411         break;
412       }
413       element++;
414   }
415   } else {
416     iResult=-EINVAL;
417   }
418   return iResult;
419 }
420
421 int AliHLTComponentHandler::UnloadLibrary(AliHLTComponentHandler::AliHLTLibHandle &handle)
422 {
423   // see header file for class documentation
424   int iResult=0;
425   fgAliLoggingFunc=NULL;
426   TString* pName=reinterpret_cast<TString*>(handle.fName);
427   if (handle.fMode!=kStatic) {
428 #ifdef HAVE_DLFCN_H
429   try {
430     dlclose(handle.fHandle);
431   }
432   catch (...) {
433     HLTError("exeption caught during dlclose of library %s", pName!=NULL?pName->Data():"");
434   }
435 #else
436   int* pCount=reinterpret_cast<int*>(handle.fHandle);
437   if (--(*pCount)==0) {
438     if (pName) {
439       /** Matthias 26.04.2007
440        * I spent about a week to investigate a bug which seems to be in ROOT.
441        * Under certain circumstances, TSystem::Unload crashes. The crash occured
442        * for the first time, when libAliHLTUtil was loaded from AliHLTSystem right
443        * after the ComponentHandler was created. It does not occur when dlopen is
444        * used. 
445        * It has most likely to do with the garbage collection and automatic
446        * cleanup in ROOT. The crash occurs when ROOT is terminated and before
447        * an instance of AliHLTSystem was created.
448        *   root [0] AliHLTSystem gHLT
449        * It does not occur when the instance was created dynamically (but not even
450        * deleted)
451        *   root [0] AliHLTSystem* gHLT=new AliHLTSystem
452        *
453        * For that reason, the libraries are not unloaded here, even though there
454        * will be memory leaks.
455       gSystem->Unload(pName->Data());
456        */
457     }
458     else {
459       HLTError("missing library name, can not unload");
460     }
461     delete pCount;
462   }
463 #endif //HAVE_DLFCN_H
464   if (pName) {
465     HLTDebug("unload library %s", pName->Data());
466   } else {
467     HLTWarning("missing name for unloaded library");
468   }
469   }
470   handle.fName=NULL;
471   handle.fHandle=NULL;
472   if (pName) {
473     delete pName;
474   }
475   pName=NULL;
476   return iResult;
477 }
478
479 int AliHLTComponentHandler::UnloadLibraries()
480 {
481   // see header file for class documentation
482   int iResult=0;
483   vector<AliHLTLibHandle>::iterator element=fLibraryList.begin();
484   while (element!=fLibraryList.end()) {
485     UnloadLibrary(*element);
486     fLibraryList.erase(element);
487     element=fLibraryList.begin();
488   }
489   return iResult;
490 }
491
492 void* AliHLTComponentHandler::FindSymbol(const char* library, const char* symbol)
493 {
494   // see header file for class documentation
495   AliHLTLibHandle* hLib=FindLibrary(library);
496   if (hLib==NULL) return NULL;
497   void* pFunc=NULL;
498 #ifdef HAVE_DLFCN_H
499   pFunc=dlsym(hLib->fHandle, symbol);
500 #else
501   TString* name=reinterpret_cast<TString*>(hLib->fName);
502   pFunc=gSystem->DynFindSymbol(name->Data(), symbol);
503 #endif
504   return pFunc;
505 }
506
507 AliHLTComponentHandler::AliHLTLibHandle* AliHLTComponentHandler::FindLibrary(const char* library)
508 {
509   // see header file for class documentation
510   AliHLTLibHandle* hLib=NULL;
511   vector<AliHLTLibHandle>::iterator element=fLibraryList.begin();
512   while (element!=fLibraryList.end()) {
513     TString* name=reinterpret_cast<TString*>((*element).fName);
514     if (name->CompareTo(library)==0) {
515       hLib=&(*element);
516       break;
517     }
518     element++;
519   }
520   return hLib;
521 }
522
523 int AliHLTComponentHandler::AddStandardComponents()
524 {
525   // see header file for class documentation
526   int iResult=0;
527   AliHLTComponent::SetGlobalComponentHandler(this);
528   AliHLTComponent::UnsetGlobalComponentHandler();
529   iResult=RegisterScheduledComponents();
530   return iResult;
531 }
532
533 int AliHLTComponentHandler::RegisterScheduledComponents()
534 {
535   // see header file for class documentation
536   int iResult=0;
537   AliHLTComponentPList::iterator element=fScheduleList.begin();
538   int iLocalResult=0;
539   while (element!=fScheduleList.end()) {
540     iLocalResult=RegisterComponent(*element);
541     if (iResult==0) iResult=iLocalResult;
542     fScheduleList.erase(element);
543     element=fScheduleList.begin();
544   }
545   return iResult;
546 }
547
548 int AliHLTComponentHandler::ActivateAgents(const AliHLTModuleAgent** blackList, int size)
549 {
550   // see header file for class documentation
551   int iResult=0;
552   AliHLTModuleAgent* pAgent=AliHLTModuleAgent::GetFirstAgent();
553   while (pAgent && iResult>=0) {
554     if (blackList) {
555       int i=0;
556       for (; i<size; i++) {
557         if (blackList[i]==pAgent) break;
558       }
559       if (i<size) {
560         // this agent was in the list
561         pAgent=AliHLTModuleAgent::GetNextAgent();
562         continue;
563       }
564     }
565
566     pAgent->ActivateComponentHandler(this);
567     pAgent=AliHLTModuleAgent::GetNextAgent();
568   }
569   return iResult;
570 }
571
572 int AliHLTComponentHandler::DeleteOwnedComponents()
573 {
574   // see header file for class documentation
575   int iResult=0;
576   AliHLTComponentPList::iterator element=fOwnedComponents.begin();
577   while (element!=fOwnedComponents.end()) {
578     //DeregisterComponent((*element)->GetComponentID());
579     try {
580       delete *element;
581     }
582     catch (...) {
583       HLTError("delete managed sample %p", *element);
584     }
585     fOwnedComponents.erase(element);
586     element=fOwnedComponents.begin();
587   }
588   return iResult;
589 }