]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCEsdWriterComponent.cxx
added new helper components to libAliHLTUtil (EsdCollector and AliHLTOUTPublisher...
[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 "AliCDBEntry.h"
37 #include "AliCDBManager.h"
38 #include "TTree.h"
39 #include "TList.h"
40 #include "AliHLTTPCTrack.h"
41 #include "AliHLTTPCTrackArray.h"
42 #include "AliHLTTPCTrackletDataFormat.h"
43 #include "AliHLTTPCDefinitions.h"
44 #include "AliHLTTPCTransform.h"
45
46 /** ROOT macro for the implementation of ROOT specific class methods */
47 ClassImp(AliHLTTPCEsdWriterComponent)
48
49 AliHLTTPCEsdWriterComponent::AliHLTTPCEsdWriterComponent()
50   :
51   fSolenoidBz(0)
52 {
53   // see header file for class documentation
54   // or
55   // refer to README to build package
56   // or
57   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
58 }
59
60 AliHLTTPCEsdWriterComponent::~AliHLTTPCEsdWriterComponent()
61 {
62   // see header file for class documentation
63 }
64
65 AliHLTTPCEsdWriterComponent::AliWriter::AliWriter()
66   :
67   fTree(NULL),
68   fESD(NULL),
69   fBase(new AliHLTTPCEsdWriterComponent)
70 {
71   // see header file for class documentation
72   // or
73   // refer to README to build package
74   // or
75   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
76 }
77
78 AliHLTTPCEsdWriterComponent::AliWriter::~AliWriter()
79 {
80   // see header file for class documentation
81   if (fBase) delete fBase;
82   fBase=NULL;
83 }
84
85 void AliHLTTPCEsdWriterComponent::AliWriter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
86 {
87   // see header file for class documentation
88   list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
89   list.push_back(AliHLTTPCDefinitions::fgkTracksDataType);
90 }
91
92 int AliHLTTPCEsdWriterComponent::AliWriter::InitWriter()
93 {
94   // see header file for class documentation
95   int iResult=0;
96   fESD = new AliESDEvent;
97   if (fESD) {
98     fESD->CreateStdContent();
99     fTree = new TTree("esdTree", "Tree with HLT ESD objects");
100     if (fTree) {
101       fTree->SetDirectory(0);
102       fESD->WriteToTree(fTree);
103     }
104   }
105   if (fTree==NULL) {
106     iResult=-ENOMEM;
107   }
108
109   if (iResult>=0) {
110     iResult=fBase->Reconfigure(NULL, NULL);
111   }
112
113   return iResult;
114 }
115
116 int AliHLTTPCEsdWriterComponent::AliWriter::CloseWriter()
117 {
118   // see header file for class documentation
119   int iResult=0;
120   if (fTree) {
121     // the esd structure is written to the user info and is
122     // needed in te ReadFromTree method to read all objects correctly
123     if (fESD) fTree->GetUserInfo()->Add(fESD);
124     WriteObject(kAliHLTVoidEventID, fTree);
125     fTree->GetUserInfo()->Clear();
126     TTree* pTree=fTree;
127     fTree=NULL;
128     delete pTree;
129   } else {
130     HLTWarning("not initialized");
131   }
132
133   if (fESD) {
134     delete fESD;
135   }
136   iResult=AliHLTRootFileWriterComponent::CloseWriter();
137   return iResult;
138 }
139
140 int AliHLTTPCEsdWriterComponent::AliWriter::DumpEvent( const AliHLTComponentEventData& evtData,
141                                             const AliHLTComponentBlockData* blocks, 
142                                             AliHLTComponentTriggerData& /*trigData*/ )
143 {
144   // see header file for class documentation
145   int iResult=0;
146   TTree* pTree=fTree;
147   assert(fBase);
148   if (pTree && fBase) {
149     if (fESD) {
150       AliESDEvent* pESD=fESD;
151
152       iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt);
153
154     } else {
155       iResult=-ENOMEM;
156     }
157   }
158   return iResult;
159 }
160
161 int AliHLTTPCEsdWriterComponent::AliWriter::ScanArgument(int argc, const char** argv)
162 {
163   // see header file for class documentation
164   int iResult=AliHLTRootFileWriterComponent::ScanArgument(argc, argv);
165   return iResult;
166 }
167
168 int AliHLTTPCEsdWriterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD,
169                                                const AliHLTComponentBlockData* blocks,
170                                                int nBlocks, int* pMinSlice,
171                                                int* pMaxSlice)
172 {
173   // see header file for class documentation
174   int iResult=0;
175   if (pESD && blocks) {
176       pESD->Reset();
177     
178       pESD->SetMagneticField(fSolenoidBz);
179       const AliHLTComponentBlockData* iter = NULL;
180       AliHLTTPCTrackletData* inPtr=NULL;
181       int bIsTrackSegs=0;
182  
183       for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
184         iter = blocks+ndx;
185         if ( (bIsTrackSegs=(iter->fDataType == AliHLTTPCDefinitions::fgkTrackSegmentsDataType))==1 ||
186              iter->fDataType == AliHLTTPCDefinitions::fgkTracksDataType ) {
187           Int_t minslice=AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
188           Int_t maxslice=AliHLTTPCDefinitions::GetMaxSliceNr(iter->fSpecification);
189           if (bIsTrackSegs==0) {
190             // slice parameter and data specification ignored, tracks already in global coordinates
191             minslice=-1;
192             maxslice=-1;
193             if (pMinSlice) *pMinSlice=0;
194             if (pMaxSlice) *pMaxSlice=AliHLTTPCTransform::GetNSlice()-1;
195           } else {
196             if (pMinSlice && (*pMinSlice==-1 || *pMinSlice>minslice)) *pMinSlice=minslice;
197             if (pMaxSlice && (*pMaxSlice==-1 || *pMaxSlice<maxslice)) *pMaxSlice=maxslice;
198           }
199           //HLTDebug("dataspec %#x minslice %d", iter->fSpecification, minslice);
200           if (minslice >=-1 && minslice<AliHLTTPCTransform::GetNSlice()) {
201             if (minslice!=maxslice) {
202               HLTWarning("data from multiple sectors in one block: "
203                          "possible mismatch in treatment of local coordinate system");
204             }
205             AliHLTTPCTrackArray tracks;
206             inPtr=(AliHLTTPCTrackletData*)iter->fPtr;
207             HLTDebug("reading block %d (slice %d): %d tracklets", ndx, minslice, inPtr->fTrackletCnt);
208             tracks.FillTracks(inPtr->fTrackletCnt, inPtr->fTracklets, minslice, 0/*don't rotate*/);
209             if ((iResult=Tracks2ESD(&tracks, pESD))>=0) {
210             }
211           } else {
212             HLTError("invalid sector number");
213             iResult=-EBADF;
214           }
215         }
216       }
217       if (iResult>=0 && pTree) {
218         pTree->Fill();
219       }
220
221   } else {
222     iResult=-EINVAL;
223   }
224   return iResult;
225 }
226
227 int AliHLTTPCEsdWriterComponent::Tracks2ESD(AliHLTTPCTrackArray* pTracks, AliESDEvent* pESD)
228 {
229   // see header file for class documentation
230   int iResult=0;
231   if (pTracks && pESD) {
232     HLTDebug("converting %d tracks from track array", pTracks->GetNTracks());
233     for (int i=0; i<pTracks->GetNTracks() && iResult>=0; i++) {
234       AliHLTTPCTrack* pTrack=(*pTracks)[i];
235       if (pTrack) {
236         //HLTDebug("convert track %d", i);
237         //pTrack->Print();
238         int iLocal=pTrack->Convert2AliKalmanTrack();
239         if (iLocal>=0) {
240         AliESDtrack iotrack;
241         iotrack.UpdateTrackParams(pTrack,AliESDtrack::kTPCin);
242         iotrack.SetTPCPoints(pTrack->GetPoints());
243         pESD->AddTrack(&iotrack);
244         } else {
245           HLTError("conversion to AliKalmanTrack failed for track %d of %d", i, pTracks->GetNTracks());
246         }
247       } else {
248         HLTError("internal mismatch in array");
249         iResult=-EFAULT;
250       }
251     }
252     
253   } else {
254     iResult=-EINVAL;
255   }
256   return iResult;
257 }
258
259 int AliHLTTPCEsdWriterComponent::Configure(const char* arguments)
260 {
261   // see header file for class documentation
262   int iResult=0;
263   if (!arguments) return iResult;
264
265   TString allArgs=arguments;
266   TString argument;
267   int bMissingParam=0;
268
269   TObjArray* pTokens=allArgs.Tokenize(" ");
270   if (pTokens) {
271     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
272       argument=((TObjString*)pTokens->At(i))->GetString();
273       if (argument.IsNull()) continue;
274       
275       if (argument.CompareTo("-solenoidBz")==0) {
276         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
277         HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
278         fSolenoidBz=((TObjString*)pTokens->At(i))->GetString().Atof();
279         continue;
280       } else {
281         HLTError("unknown argument %s", argument.Data());
282         iResult=-EINVAL;
283         break;
284       }
285     }
286     delete pTokens;
287   }
288   if (bMissingParam) {
289     HLTError("missing parameter for argument %s", argument.Data());
290     iResult=-EINVAL;
291   }
292
293   return iResult;
294 }
295
296 int AliHLTTPCEsdWriterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
297 {
298   // see header file for class documentation
299   int iResult=0;
300   const char* path=kAliHLTCDBSolenoidBz;
301   const char* defaultNotify="";
302   if (cdbEntry) {
303     path=cdbEntry;
304     defaultNotify=" (default)";
305   }
306   if (path) {
307     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
308     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
309     if (pEntry) {
310       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
311       if (pString) {
312         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
313         iResult=Configure(pString->GetString().Data());
314       } else {
315         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
316       }
317     } else {
318       HLTError("can not fetch object \"%s\" from CDB", path);
319     }
320   }
321   
322   return iResult;
323 }
324
325 AliHLTTPCEsdWriterComponent::AliConverter::AliConverter()
326   :
327   fESD(NULL),
328   fBase(new AliHLTTPCEsdWriterComponent),
329   fWriteTree(1)
330 {
331   // see header file for class documentation
332   // or
333   // refer to README to build package
334   // or
335   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
336 }
337
338 AliHLTTPCEsdWriterComponent::AliConverter::~AliConverter()
339 {
340   // see header file for class documentation
341   if (fBase) delete fBase;
342   fBase=NULL;
343
344   if (fESD) delete fESD;
345   fESD=NULL;
346 }
347
348 void AliHLTTPCEsdWriterComponent::AliConverter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
349 {
350   // see header file for class documentation
351   list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
352   list.push_back(AliHLTTPCDefinitions::fgkTracksDataType);
353 }
354
355 AliHLTComponentDataType AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataType()
356 {
357   // see header file for class documentation
358   return kAliHLTDataTypeESDTree;
359 }
360
361 void AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
362 {
363   // see header file for class documentation
364   constBase=2000000;
365   inputMultiplier=10.0;
366 }
367
368 int AliHLTTPCEsdWriterComponent::AliConverter::DoInit(int argc, const char** argv)
369 {
370   // see header file for class documentation
371   int iResult=0;
372   TString argument="";
373   int bMissingParam=0;
374   for (int i=0; i<argc && iResult>=0; i++) {
375     argument=argv[i];
376     if (argument.IsNull()) continue;
377
378     // -notree
379     if (argument.CompareTo("-notree")==0) {
380       fWriteTree=0;
381
382       // -tree
383     } else if (argument.CompareTo("-tree")==0) {
384       fWriteTree=1;
385
386     } else {
387       HLTError("unknown argument %s", argument.Data());
388       break;
389     }
390   }
391   if (bMissingParam) {
392     HLTError("missing parameter for argument %s", argument.Data());
393     iResult=-EINVAL;
394   }
395
396   if (iResult>=0) {
397     iResult=fBase->Reconfigure(NULL, NULL);
398   }
399
400   return iResult;
401 }
402
403 int AliHLTTPCEsdWriterComponent::AliConverter::DoDeinit()
404 {
405   // see header file for class documentation
406   return 0;
407 }
408
409 int AliHLTTPCEsdWriterComponent::AliConverter::DoEvent(const AliHLTComponentEventData& evtData, 
410                                                        const AliHLTComponentBlockData* blocks, 
411                                                        AliHLTComponentTriggerData& /*trigData*/,
412                                                        AliHLTUInt8_t* /*outputPtr*/, 
413                                                        AliHLTUInt32_t& /*size*/,
414                                                        AliHLTComponentBlockDataList& /*outputBlocks*/ )
415 {
416   // see header file for class documentation
417   int iResult=0;
418   assert(fBase);
419   if (!fESD) {
420     fESD = new AliESDEvent;
421     if (fESD) {
422       fESD->CreateStdContent();
423     } else {
424       iResult=-ENOMEM;
425     }
426   }
427
428   AliESDEvent* pESD = fESD;
429   if (pESD && fBase) {
430     TTree* pTree = NULL;
431     // TODO: Matthias 06.12.2007
432     // Tried to write the ESD directly instead to a tree, but this did not work
433     // out. Information in the ESD is different, needs investigation.
434     if (fWriteTree)
435       pTree = new TTree("esdTree", "Tree with HLT ESD objects");
436     if (pTree) {
437       pTree->SetDirectory(0);
438       pESD->WriteToTree(pTree);
439     }
440
441     if ((iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt))>=0) {
442         // TODO: set the specification correctly
443       if (pTree) {
444         // the esd structure is written to the user info and is
445         // needed in te ReadFromTree method to read all objects correctly
446         pTree->GetUserInfo()->Add(pESD);
447         iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginTPC, 0);
448       } else {
449         iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginTPC, 0);
450       }
451     }
452     if (pTree) {
453       // clear user info list to prevent objects from being deleted
454       pTree->GetUserInfo()->Clear();
455       delete pTree;
456     }
457   }
458   return iResult;
459 }
460