]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCEsdWriterComponent.cxx
minor corrections, remove obsolete include file
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCEsdWriterComponent.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   AliHLTTPCEsdWriterComponent.cxx
20     @author Matthias Richter
21     @date   
22     @brief  Writer component to store tracks of the HLT TPC conformal
23             mapping tracker in the AliESD format
24
25                                                                           */
26 // see header file for class documentation
27 // or
28 // refer to README to build package
29 // or
30 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
31
32 #include <cassert>
33 #include "AliHLTTPCEsdWriterComponent.h"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "TTree.h"
37 #include "AliHLTTPCTrack.h"
38 #include "AliHLTTPCTrackArray.h"
39 #include "AliHLTTPCTrackletDataFormat.h"
40 #include "AliHLTTPCDefinitions.h"
41 #include "AliHLTTPCTransform.h"
42
43 /** ROOT macro for the implementation of ROOT specific class methods */
44 ClassImp(AliHLTTPCEsdWriterComponent)
45
46 AliHLTTPCEsdWriterComponent::AliHLTTPCEsdWriterComponent()
47 {
48   // see header file for class documentation
49   // or
50   // refer to README to build package
51   // or
52   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
53 }
54
55 AliHLTTPCEsdWriterComponent::~AliHLTTPCEsdWriterComponent()
56 {
57   // see header file for class documentation
58 }
59
60 AliHLTTPCEsdWriterComponent::AliWriter::AliWriter()
61   :
62   fTree(NULL),
63   fESD(NULL),
64   fBase(new AliHLTTPCEsdWriterComponent)
65 {
66   // see header file for class documentation
67   // or
68   // refer to README to build package
69   // or
70   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
71 }
72
73 AliHLTTPCEsdWriterComponent::AliWriter::~AliWriter()
74 {
75   // see header file for class documentation
76   if (fBase) delete fBase;
77   fBase=NULL;
78 }
79
80 void AliHLTTPCEsdWriterComponent::AliWriter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
81 {
82   // see header file for class documentation
83   list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
84 }
85
86 int AliHLTTPCEsdWriterComponent::AliWriter::InitWriter()
87 {
88   // see header file for class documentation
89   int iResult=0;
90   fESD = new AliESDEvent;
91   if (fESD) {
92     fESD->CreateStdContent();
93     fTree = new TTree("esdTree", "Tree with HLT ESD objects");
94     if (fTree) {
95       fTree->SetDirectory(0);
96       fESD->WriteToTree(fTree);
97     }
98   }
99   if (fTree==NULL) {
100     iResult=-ENOMEM;
101   }
102   return iResult;
103 }
104
105 int AliHLTTPCEsdWriterComponent::AliWriter::CloseWriter()
106 {
107   // see header file for class documentation
108   int iResult=0;
109   if (fTree) {
110     WriteObject(kAliHLTVoidEventID, fTree);
111     TTree* pTree=fTree;
112     fTree=NULL;
113     delete pTree;
114   } else {
115     HLTWarning("not initialized");
116   }
117   iResult=AliHLTRootFileWriterComponent::CloseWriter();
118   return iResult;
119 }
120
121 int AliHLTTPCEsdWriterComponent::AliWriter::DumpEvent( const AliHLTComponentEventData& evtData,
122                                             const AliHLTComponentBlockData* blocks, 
123                                             AliHLTComponentTriggerData& /*trigData*/ )
124 {
125   // see header file for class documentation
126   int iResult=0;
127   TTree* pTree=fTree;
128   assert(fBase);
129   if (pTree && fBase) {
130     if (fESD) {
131       AliESDEvent* pESD=fESD;
132
133       iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt);
134
135     } else {
136       iResult=-ENOMEM;
137     }
138   }
139   return iResult;
140 }
141
142 int AliHLTTPCEsdWriterComponent::AliWriter::ScanArgument(int argc, const char** argv)
143 {
144   // see header file for class documentation
145   int iResult=AliHLTRootFileWriterComponent::ScanArgument(argc, argv);
146   return iResult;
147 }
148
149 int AliHLTTPCEsdWriterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD,
150                                                const AliHLTComponentBlockData* blocks,
151                                                int nBlocks, int* pMinSlice,
152                                                int* pMaxSlice)
153 {
154   // see header file for class documentation
155   int iResult=0;
156   if (pESD && blocks) {
157       const AliHLTComponentBlockData* iter = NULL;
158       AliHLTTPCTrackletData* inPtr=NULL;
159       int bIsTrackSegs=0;
160  
161       for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
162         iter = blocks+ndx;
163         if ( (bIsTrackSegs=(iter->fDataType == AliHLTTPCDefinitions::fgkTrackSegmentsDataType))==1 ||
164              iter->fDataType == AliHLTTPCDefinitions::fgkTracksDataType ) {
165           Int_t minslice=AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
166           Int_t maxslice=AliHLTTPCDefinitions::GetMaxSliceNr(iter->fSpecification);
167           if (bIsTrackSegs==0) {
168             // slice parameter and data specification ignored, tracks already in global coordinates
169             minslice=-1;
170             maxslice=-1;
171             if (pMinSlice) *pMinSlice=0;
172             if (pMaxSlice) *pMaxSlice=AliHLTTPCTransform::GetNSlice()-1;
173           } else {
174             if (pMinSlice && (*pMinSlice==-1 || *pMinSlice>minslice)) *pMinSlice=minslice;
175             if (pMaxSlice && (*pMaxSlice==-1 || *pMaxSlice<maxslice)) *pMaxSlice=maxslice;
176           }
177           //HLTDebug("dataspec %#x minslice %d", iter->fSpecification, minslice);
178           if (minslice >=-1 && minslice<AliHLTTPCTransform::GetNSlice()) {
179             if (minslice!=maxslice) {
180               HLTWarning("data from multiple sectors in one block: "
181                          "possible missmatch in treatment of local coordinate system");
182             }
183             AliHLTTPCTrackArray tracks;
184             inPtr=(AliHLTTPCTrackletData*)iter->fPtr;
185             HLTDebug("reading block %d (slice %d): %d tracklets", ndx, minslice, inPtr->fTrackletCnt);
186             tracks.FillTracks(inPtr->fTrackletCnt, inPtr->fTracklets, minslice, 0/*don't rotate*/);
187             if ((iResult=Tracks2ESD(&tracks, pESD))>=0) {
188             }
189           } else {
190             HLTError("invalid sector number");
191             iResult=-EBADF;
192           }
193         }
194       }
195       if (iResult>=0 && pTree) {
196         pTree->Fill();
197       }
198
199       pESD->Reset();
200     
201   } else {
202     iResult=-EINVAL;
203   }
204   return iResult;
205 }
206
207 int AliHLTTPCEsdWriterComponent::Tracks2ESD(AliHLTTPCTrackArray* pTracks, AliESDEvent* pESD)
208 {
209   // see header file for class documentation
210   int iResult=0;
211   if (pTracks && pESD) {
212     HLTDebug("converting %d tracks from track array", pTracks->GetNTracks());
213     for (int i=0; i<pTracks->GetNTracks() && iResult>=0; i++) {
214       AliHLTTPCTrack* pTrack=(*pTracks)[i];
215       if (pTrack) {
216         //HLTDebug("convert track %d", i);
217         //pTrack->Print();
218         int iLocal=pTrack->Convert2AliKalmanTrack();
219         if (iLocal>=0) {
220         AliESDtrack iotrack;
221         iotrack.UpdateTrackParams(pTrack,AliESDtrack::kTPCin);
222         iotrack.SetTPCPoints(pTrack->GetPoints());
223         pESD->AddTrack(&iotrack);
224         } else {
225           HLTError("conversion to AliKalmanTrack failed for track %d of %d", i, pTracks->GetNTracks());
226         }
227       } else {
228         HLTError("internal missmatch in array");
229         iResult=-EFAULT;
230       }
231     }
232     
233   } else {
234     iResult=-EINVAL;
235   }
236   return iResult;
237 }
238
239 AliHLTTPCEsdWriterComponent::AliConverter::AliConverter()
240   :
241   fBase(new AliHLTTPCEsdWriterComponent),
242   fWriteTree(1)
243 {
244   // see header file for class documentation
245   // or
246   // refer to README to build package
247   // or
248   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
249 }
250
251 AliHLTTPCEsdWriterComponent::AliConverter::~AliConverter()
252 {
253   // see header file for class documentation
254   if (fBase) delete fBase;
255   fBase=NULL;
256 }
257
258 void AliHLTTPCEsdWriterComponent::AliConverter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
259 {
260   // see header file for class documentation
261   list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
262 }
263
264 AliHLTComponentDataType AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataType()
265 {
266   // see header file for class documentation
267   return kAliHLTDataTypeESDTree;
268 }
269
270 void AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
271 {
272   // see header file for class documentation
273   constBase=1000000;
274   inputMultiplier=1.0;
275 }
276
277 int AliHLTTPCEsdWriterComponent::AliConverter::DoInit(int argc, const char** argv)
278 {
279   // see header file for class documentation
280   int iResult=0;
281   TString argument="";
282   int bMissingParam=0;
283   for (int i=0; i<argc && iResult>=0; i++) {
284     argument=argv[i];
285     if (argument.IsNull()) continue;
286
287     // -notree
288     if (argument.CompareTo("-notree")==0) {
289       fWriteTree=0;
290
291       // -tree
292     } else if (argument.CompareTo("-tree")==0) {
293       fWriteTree=1;
294
295     } else {
296       HLTError("unknown argument %s", argument.Data());
297       break;
298     }
299   }
300   if (bMissingParam) {
301     HLTError("missing parameter for argument %s", argument.Data());
302     iResult=-EINVAL;
303   }
304
305   return iResult;
306 }
307
308 int AliHLTTPCEsdWriterComponent::AliConverter::DoDeinit()
309 {
310   // see header file for class documentation
311   return 0;
312 }
313
314 int AliHLTTPCEsdWriterComponent::AliConverter::DoEvent(const AliHLTComponentEventData& evtData, 
315                                                        const AliHLTComponentBlockData* blocks, 
316                                                        AliHLTComponentTriggerData& /*trigData*/,
317                                                        AliHLTUInt8_t* /*outputPtr*/, 
318                                                        AliHLTUInt32_t& /*size*/,
319                                                        AliHLTComponentBlockDataList& /*outputBlocks*/ )
320 {
321   // see header file for class documentation
322   int iResult=0;
323   assert(fBase);
324   AliESDEvent* pESD = new AliESDEvent;
325   if (pESD && fBase) {
326     pESD->CreateStdContent();
327     TTree* pTree = NULL;
328     // TODO: Matthias 06.12.2007
329     // Tried to write the ESD directly instead to a tree, but this did not work
330     // out. Information in the ESD is different, needs investigation.
331     if (fWriteTree)
332       pTree = new TTree("esdTree", "Tree with HLT ESD objects");
333     if (pTree) {
334       pTree->SetDirectory(0);
335       pESD->WriteToTree(pTree);
336     }
337
338     if ((iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt))>=0) {
339         // TODO: set the specification correctly
340       if (pTree)
341         iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginTPC, 0);
342       else
343         iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginTPC, 0);
344     }
345     if (pTree)
346       delete pTree;
347
348     delete pESD;
349   }
350   return iResult;
351 }
352