Adding custom linkdef file so that dictionaries are correctly generated for AliHLTTri...
[u/mrichter/AliRoot.git] / HLT / BASE / util / AliHLTRootFilePublisherComponent.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  *                  for The ALICE HLT Project.                            *
9  *                                                                        *
10  * Permission to use, copy, modify and distribute this software and its   *
11  * documentation strictly for non-commercial purposes is hereby granted   *
12  * without fee, provided that the above copyright notice appears in all   *
13  * copies and that both the copyright notice and this permission notice   *
14  * appear in the supporting documentation. The authors make no claims     *
15  * about the suitability of this software for any purpose. It is          *
16  * provided "as is" without express or implied warranty.                  *
17  **************************************************************************/
18
19 /** @file   AliHLTRootFilePublisherComponent.cxx
20     @author Matthias Richter, Jochen Thaeder
21     @date   
22     @brief  HLT file publisher component implementation. */
23
24 #include "AliHLTRootFilePublisherComponent.h"
25 //#include <TObjString.h>
26 //#include <TMath.h>
27 //#include <TFile.h>
28
29 #include "TList.h"
30 #include "TTree.h"
31 #include "TKey.h"
32 #include "TFile.h"
33
34
35 /** ROOT macro for the implementation of ROOT specific class methods */
36 ClassImp(AliHLTRootFilePublisherComponent)
37
38 /*
39  * ---------------------------------------------------------------------------------
40  *                            Constructor / Destructor
41  * ---------------------------------------------------------------------------------
42  */
43
44 // #################################################################################
45 AliHLTRootFilePublisherComponent::AliHLTRootFilePublisherComponent()
46   :
47   AliHLTFilePublisher(),
48   fpCurrentEvent(NULL),
49   fObjectName("") {
50   // see header file for class documentation
51   // or
52   // refer to README to build package
53   // or
54   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
55
56   // Set file to ROOT-File
57   SetIsRawFile( kFALSE );
58 }
59
60 // #################################################################################
61 AliHLTRootFilePublisherComponent::~AliHLTRootFilePublisherComponent() {
62   // see header file for class documentation
63
64   // file list and file name list are owner of their objects and
65   // delete all the objects
66 }
67
68 /*
69  * ---------------------------------------------------------------------------------
70  * Public functions to implement AliHLTComponent's interface.
71  * These functions are required for the registration process
72  * ---------------------------------------------------------------------------------
73  */
74
75 // #################################################################################
76 const char* AliHLTRootFilePublisherComponent::GetComponentID() {
77   // see header file for class documentation
78   return "ROOTFilePublisher";
79 }
80
81 // #################################################################################
82 AliHLTComponent* AliHLTRootFilePublisherComponent::Spawn() {
83   // see header file for class documentation
84   return new AliHLTRootFilePublisherComponent;
85 }
86
87 /*
88  * ---------------------------------------------------------------------------------
89  * Protected functions to implement AliHLTComponent's interface.
90  * These functions provide initialization as well as the actual processing
91  * capabilities of the component. 
92  * ---------------------------------------------------------------------------------
93  */
94
95 // #################################################################################
96 Int_t AliHLTRootFilePublisherComponent::ScanArgument(Int_t argc, const char** argv) {
97   // see header file for class documentation
98
99   Int_t iResult = 0;
100
101   TString argument = "";
102   TString parameter = "";
103   Int_t bMissingParam = 0;
104   fObjectName = "";  // Reset this to the default: read all objects in the file.
105   
106   argument=argv[iResult];
107   if (argument.IsNull()) return -EINVAL;
108
109   // -objectname
110   if ( !argument.CompareTo("-objectname") ) {
111     if ( ! (bMissingParam=(++iResult>=argc)) ) {
112       parameter = argv[iResult];
113       parameter.Remove(TString::kLeading, ' '); // remove all blanks
114       fObjectName = parameter;
115     } 
116   }
117   else {
118     HLTError("unknown argument %s", argument.Data());
119     iResult = -EINVAL;    
120   }
121   
122   if ( bMissingParam ) {
123     HLTError("missing parameter for argument %s", argument.Data());
124     iResult = -EPROTO;
125   }
126
127   return iResult;
128 }
129
130  // #################################################################################
131 Int_t AliHLTRootFilePublisherComponent::GetEvent( const AliHLTComponentEventData& /*evtData*/,
132                                                 AliHLTComponentTriggerData& /*trigData*/,
133                                                 AliHLTUInt8_t* /*outputPtr*/, 
134                                                 AliHLTUInt32_t& size,
135                                                 vector<AliHLTComponentBlockData>& /*outputBlocks*/ ) {
136   // see header file for class documentation
137
138   if ( !IsDataEvent() ) return 0;
139
140   Int_t iResult=0;
141   size=0;
142
143   // -- Ptr to current event
144   TObjLink *lnk = fpCurrentEvent;
145   if ( lnk == NULL) {
146     lnk = GetEventList()->FirstLink();
147     fpCurrentEvent = lnk;
148   }
149
150   if ( lnk ) {
151     EventFiles* pEventDesc = dynamic_cast<EventFiles*>( lnk->GetObject() );
152     if (pEventDesc) {
153     
154       HLTDebug("publishing files for event %p", pEventDesc);
155       TList& files=*pEventDesc; // type conversion operator defined
156       TObjLink *flnk=files.FirstLink();
157
158       while (flnk && iResult>=0) {
159
160         FileDesc* pFileDesc=dynamic_cast<FileDesc*>(flnk->GetObject());
161
162         if (not fOpenFilesAtStart) pFileDesc->OpenFile();
163         TFile* pFile=NULL;
164
165         if (pFileDesc && (pFile=*pFileDesc)!=NULL) {
166
167           for ( Int_t i = 0; i < pFile->GetListOfKeys()->GetEntries(); i++  ){
168             TKey * key= dynamic_cast<TKey*>( pFile->GetListOfKeys()->At(i) );
169
170             if ( fObjectName != "" ) {
171               if ( !( ((TString) key->GetName()).CompareTo(fObjectName) ) )
172                 PushBack( key->ReadObj(), *pFileDesc, *pFileDesc ); 
173             }
174             else 
175               PushBack( key->ReadObj(), *pFileDesc, *pFileDesc );             
176               // above : type conversion operator defined for DataType and Spec
177           }
178
179           if (not fOpenFilesAtStart) pFileDesc->CloseFile();
180         } else {
181           HLTError("no file available");
182           iResult=-EFAULT;
183         }
184         flnk = flnk->Next();
185       }
186     } else {
187       HLTError("can not get event descriptor from list link");
188       iResult=-EFAULT;
189     }
190   } else {
191     iResult=-ENOENT;
192   }
193   if (iResult>=0 && fpCurrentEvent) fpCurrentEvent=fpCurrentEvent->Next();
194
195   return iResult;
196 }