]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCEsdWriterComponent.h
d7f4f9751a8462a9b53d8ffcd0dde37ff25d3ceb
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCEsdWriterComponent.h
1 // @(#) $Id$
2
3 #ifndef ALIHLTTPCESDWRITERCOMPONENT_H
4 #define ALIHLTTPCESDWRITERCOMPONENT_H
5 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
6  * See cxx source for full Copyright notice                               */
7
8 /** @file   AliHLTTPCEsdWriterComponent.h
9     @author Matthias Richter
10     @date   
11     @brief  Writer component to store tracks of the HLT TPC conformal
12             mapping tracker in the AliESD format
13
14                                                                           */
15 // see below for class documentation
16 // or
17 // refer to README to build package
18 // or
19 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
20
21 #include "AliHLTRootFileWriterComponent.h"
22 #include "AliHLTProcessor.h"
23
24 // forward declarations
25 class TTree;
26 class AliESDEvent;
27 class AliHLTTPCTrackArray;
28
29 /**
30  * @class AliHLTTPCEsdWriterComponent
31  * This class translates incoming track segments structures from the TPC
32  * conformal mapping tracker (datatype TRAKSEGS/TPC) or tracks in global 
33  * coordinates from the AliHLTTPCGlobalMergerComponent (TRACKS/TPC) into
34  * the ESD format.
35  *
36  * The \em TPCEsdWriter writes it to a ROOT file, the \em TPCEsdConverter
37  * to a TTree and sends the whole object to the output stream with data
38  * type @ref kAliHLTDataTypeESDTree and origin TPC.
39  *
40  * In case of TRAKSEGS, the component can only process data block containing
41  * data of one slice, but it can read an unlimeted number of data blocks.
42  *
43  * componentid: \b TPCEsdWriter <br>
44  * componentid: \b TPCEsdConverter <br>
45  * componentlibrary: \b libAliHLTTPC.so <br>
46  * Arguments TPCEsdWriter: <br>
47  * <!-- NOTE: ignore the \li. <i> and </i>: it's just doxygen formating -->
48  * \li -datafile     <i> filename   </i> <br>
49  *      file name base
50  * \li -directory    <i> directory  </i> <br>
51  *      target directory
52  *
53  * Arguments TPCEsdConverter: <br>
54  * <!-- NOTE: ignore the \li. <i> and </i>: it's just doxygen formating -->
55  * \li -notree
56  *      write ESD directly to output (@ref kAliHLTDataTypeESDObject)
57  * \li -tree
58  *      write ESD directly to TTree and to output (@ref kAliHLTDataTypeESDTree),
59  *      this is the default behavior.
60  *
61  * <pre>
62  * Example usage (HLT configuration file):
63  *         \<Proc ID="EsdWriter" type="prc">
64  *             \<Cmd>AliRootWrapperSubscriber -eventmodulo 1
65  *                 -componentid TPCEsdWriter
66  *                 -componentlibrary libAliHLTTPC.so
67  *                 -componentargs "-datafile AliESDs.root"
68  *            \</Cmd>
69  *
70  *            \<Parent>TR0-SC\</Parent>
71  *            \<Node>master\</Node>
72  *            \<Shm blocksize="1k" blockcount="1" type="sysv"/>
73  *        \</Proc>
74  * </pre>
75  *
76  * @see AliHLTFileWriter and AliHLTRootFileWriterComponent for more parameters
77  */
78 class AliHLTTPCEsdWriterComponent : public AliHLTLogging
79 {
80  public:
81   /** standard constructor */
82   AliHLTTPCEsdWriterComponent();
83   /** destructor */
84   ~AliHLTTPCEsdWriterComponent();
85
86   /**
87    * class AliHLTTPCEsdWriterComponent::AliWriter
88    * The writer component of the AliHLTTPCEsdWriterComponent.
89    */
90   class AliWriter : public AliHLTRootFileWriterComponent
91   {
92   public:
93   /** standard constructor */
94   AliWriter();
95   /** destructor */
96   virtual ~AliWriter();
97
98   /**
99    * The id of the component.
100    * @return component id (string)
101    */
102   const char* GetComponentID() {return "TPCEsdWriter";};
103
104   void GetInputDataTypes(AliHLTComponentDataTypeList& list);
105
106   /**
107    * Spawn function.
108    * @return new class instance
109    */
110   AliHLTComponent* Spawn() {return new AliHLTTPCEsdWriterComponent::AliWriter;}
111
112  protected:
113   /**
114    * Data processing method for the component.
115    * The function can be overloaded by specific ROOT file writer
116    * components.
117    * @param evtData       event data structure
118    * @param blocks        input data block descriptors
119    * @param trigData      trigger data structure
120    */
121   int DumpEvent( const AliHLTComponentEventData& evtData,
122                  const AliHLTComponentBlockData* blocks, 
123                  AliHLTComponentTriggerData& trigData );
124
125   using AliHLTRootFileWriterComponent::DumpEvent;
126
127   /**
128    * Scan one argument and adjacent parameters.
129    * @param argc           size of the argument array
130    * @param argv           agument array for component initialization
131    * @return number of processed members of the argv <br>
132    *         -EINVAL unknown argument <br>
133    *         -EPROTO parameter for argument missing <br>
134    */
135   int ScanArgument(int argc, const char** argv);
136  private:
137   /** copy constructor prohibited */
138   AliWriter(const AliWriter&);
139   /** assignment operator prohibited */
140   AliWriter& operator=(const AliWriter&);
141
142   /**
143    * Init the writer.
144    * The DoInit function is not available for this child class. InitWriter is the
145    * corresponding function for classes derived from AliHLTFileWriter.
146    */
147   int InitWriter();
148
149   /**
150    * Init the writer.
151    * The DoDeinit function is not available for this child class. CloseWriter is the
152    * corresponding function for classes derived from AliHLTFileWriter.
153    */
154   int CloseWriter();
155
156   /** the ESD tree */
157   TTree* fTree; //! transient value
158
159   /** the ESD */
160   AliESDEvent* fESD; //! transient value
161
162   /** pointer to the basic ESD conversion methods */
163   AliHLTTPCEsdWriterComponent* fBase; //! transient value
164   };
165
166   /**
167    * class AliHLTTPCEsdWriterComponent::AliConverter
168    * The converter component of the AliHLTTPCEsdWriterComponent.
169    * 
170    */
171   class AliConverter : public AliHLTProcessor
172   {
173   public:
174     /** standard constructor */
175     AliConverter();
176     /** destructor */
177     virtual ~AliConverter();
178
179     // interface methods of base class
180     const char* GetComponentID() {return "TPCEsdConverter";};
181     void GetInputDataTypes(AliHLTComponentDataTypeList& list);
182     AliHLTComponentDataType GetOutputDataType();
183     void GetOutputDataSize(unsigned long& constBase, double& inputMultiplier);
184     AliHLTComponent* Spawn() {return new AliHLTTPCEsdWriterComponent::AliConverter;}
185
186   protected:
187     // interface methods of base class
188     int DoInit(int argc, const char** argv);
189     int DoDeinit();
190     int DoEvent(const AliHLTComponentEventData& evtData,
191                 const AliHLTComponentBlockData* blocks, 
192                 AliHLTComponentTriggerData& trigData,
193                 AliHLTUInt8_t* outputPtr, 
194                 AliHLTUInt32_t& size,
195                 AliHLTComponentBlockDataList& outputBlocks );
196
197     using AliHLTProcessor::DoEvent;
198
199   private:
200     /** copy constructor prohibited */
201     AliConverter(const AliConverter&);
202     /** assignment operator prohibited */
203     AliConverter& operator=(const AliConverter&);
204
205     /** pointer to the basic ESD conversion methods */
206     AliHLTTPCEsdWriterComponent* fBase; //! transient value
207
208     /** write object to TTree or directly */
209     int fWriteTree; //!transient
210   };
211
212  protected:
213   /**
214    * Process the input data blocks.
215    * @param pTree    tree to be filled
216    * @param pESD     ESD to be filled
217    * @param blocks   data block descriptor array
218    * @param nBlocks  size of the array
219    * @param pMinSize [OUT] receives the minimum slice no
220    * @param pMaxSize [OUT] receives the maximum slice no
221    * @return neg. error code if failed
222    */
223   int ProcessBlocks(TTree* pTree, AliESDEvent* pESD, const AliHLTComponentBlockData* blocks,
224                     int nBlocks, int* pMinSlice=NULL, int* pMaxSlice=NULL);
225
226   /**
227    * Covert tracks to AliTPCtracks (AliKalmanTracks) and add them to ESD.
228    * @param pTracks  array of tracks
229    * @param pESD     pointer to ESD
230    * @return neg. error code if failed
231    */
232   int Tracks2ESD(AliHLTTPCTrackArray* pTracks, AliESDEvent* pESD);
233
234  private:
235   /** copy constructor prohibited */
236   AliHLTTPCEsdWriterComponent(const AliHLTTPCEsdWriterComponent&);
237   /** assignment operator prohibited */
238   AliHLTTPCEsdWriterComponent& operator=(const AliHLTTPCEsdWriterComponent&);
239
240   ClassDef(AliHLTTPCEsdWriterComponent, 1)
241 };
242 #endif