]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/BASE/util/AliHLTRootFileWriterComponent.cxx
- implemented component registration via agents
[u/mrichter/AliRoot.git] / HLT / BASE / util / AliHLTRootFileWriterComponent.cxx
CommitLineData
4ddfc222 1// @(#) $Id$
2
9be2600f 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
4ddfc222 19/** @file AliHLTRootFileWriterComponent.cxx
20 @author Matthias Richter
21 @date
22 @brief Base class for writer components to store data in a ROOT file
23
24 */
25
26#include "AliHLTRootFileWriterComponent.h"
27#include "TFile.h"
28#include "TString.h"
29
4ddfc222 30/** ROOT macro for the implementation of ROOT specific class methods */
31ClassImp(AliHLTRootFileWriterComponent)
32
33AliHLTRootFileWriterComponent::AliHLTRootFileWriterComponent()
34 :
35 AliHLTFileWriter(),
36 fEventID(kAliHLTVoidEventID),
37 fCurrentFile(NULL)
38{
39 // see header file for class documentation
40 // or
41 // refer to README to build package
42 // or
43 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
44
45 // all blocks of one event go into the same file
46 SetMode(kConcatenateBlocks);
47}
48
49AliHLTRootFileWriterComponent::AliHLTRootFileWriterComponent(const AliHLTRootFileWriterComponent&)
50 :
51 AliHLTFileWriter(),
52 fEventID(kAliHLTVoidEventID),
53 fCurrentFile(NULL)
54{
55 // see header file for class documentation
56 HLTFatal("copy constructor untested");
57}
58
59AliHLTRootFileWriterComponent& AliHLTRootFileWriterComponent::operator=(const AliHLTRootFileWriterComponent&)
60{
61 // see header file for class documentation
62 HLTFatal("assignment operator untested");
63 return *this;
64}
65
66AliHLTRootFileWriterComponent::~AliHLTRootFileWriterComponent()
67{
68 // see header file for class documentation
69}
70
79c114b5 71int AliHLTRootFileWriterComponent::CloseWriter()
72{
73 // see header file for class documentation
74 if (fCurrentFile!=NULL) {
75 HLTDebug("close root file");
76 TFile* pFile=fCurrentFile; fCurrentFile=NULL;
77 pFile->Close(); delete pFile;
78 }
e83e889b 79 return 0;
79c114b5 80}
81
4ddfc222 82int AliHLTRootFileWriterComponent::DumpEvent( const AliHLTComponentEventData& evtData,
83 const AliHLTComponentBlockData* blocks,
84 AliHLTComponentTriggerData& trigData )
85{
86 // see header file for class documentation
87 int iResult=0;
88 if (evtData.fStructSize==0 && blocks==NULL && trigData.fStructSize==0) {
89 // this is just to get rid of the warning "unused parameter"
90 }
91 const TObject* pObj=GetFirstInputObject(kAliHLTAnyDataType);
92 HLTDebug("got first object %p", pObj);
93 int count=0;
94 while (pObj && iResult>=0) {
95 iResult=WriteObject(evtData.fEventID, pObj);
0a94ac22 96 if (iResult == 0) {
4ddfc222 97 count++;
98 HLTDebug("wrote object of class %s, data type %s", pObj->ClassName(), (DataType2Text(GetDataType(pObj)).c_str()));
99 }
100 pObj=GetNextInputObject();
101 }
79c114b5 102 HLTDebug("wrote %d object(s) from %d input blocks to file", count, GetNumberOfInputBlocks());
4ddfc222 103 return iResult;
104}
105
106int AliHLTRootFileWriterComponent::ScanArgument(int argc, const char** argv)
107{
108 // see header file for class documentation
109 // no other arguments known
110 if (argc==0 && argv==NULL) {
111 // this is just to get rid of the warning "unused parameter"
112 }
96bda103 113 int iResult=-EINVAL;
4ddfc222 114 return iResult;
115}
116
117int AliHLTRootFileWriterComponent::WriteObject(const AliHLTEventID_t eventID, const TObject *pOb)
118{
119 // see header file for class documentation
120 int iResult=0;
121 if (pOb) {
122 HLTDebug("write object %p (%s)", pOb, pOb->GetName());
123 if (CheckMode(kConcatenateEvents) && eventID != fEventID &&
124 fCurrentFile!=NULL && eventID!=kAliHLTVoidEventID) {
125 TFile* pFile=fCurrentFile; fCurrentFile=NULL;
126 pFile->Close(); delete pFile;
127 }
128 if (fCurrentFile==NULL) {
129 fCurrentFile=OpenFile(eventID, 0);
130 if (fCurrentFile) fEventID=eventID;
131 }
132 if (fCurrentFile) {
133 fCurrentFile->cd();
134 pOb->Write();
135 } else {
136 iResult=-EBADF;
137 }
138 }
139 return iResult;
140}
141
142TFile* AliHLTRootFileWriterComponent::OpenFile(const AliHLTEventID_t eventID, const int blockID, const char* option)
143{
144 // see header file for class documentation
145 TFile* pFile=NULL;
146 TString filename("");
5ab95e9a 147 if ((BuildFileName(eventID, blockID, kAliHLTVoidDataType, 0, filename))>=0 && filename.IsNull()==0) {
4ddfc222 148 pFile=new TFile(filename, option);
149 if (pFile) {
150 if (pFile->IsZombie()) {
151 delete pFile;
152 pFile=NULL;
153 HLTError("can not open ROOT file %s", filename.Data());
154 }
155 } else {
156 HLTFatal("memory allocation failed");
157 }
158 } else {
159 HLTError("failed to build a new file name for event %#8x", eventID);
160 }
161 return pFile;
162}