]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCEsdWriterComponent.cxx
- improvements in AliHLTFilePublisher/Writer
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCEsdWriterComponent.cxx
1 // @(#) $Id$
2
3 /** @file   AliHLTTPCEsdWriterComponent.cxx
4     @author Matthias Richter
5     @date   
6     @brief  Writer component to store tracks of the HLT TPC conformal
7             mapping tracker in the AliESD format
8
9                                                                           */
10 #include "AliHLTTPCEsdWriterComponent.h"
11 #include "AliESD.h"
12 #include "TTree.h"
13 #include "AliHLTTPCTrack.h"
14 #include "AliHLTTPCTrackArray.h"
15 #include "AliHLTTPCTrackletDataFormat.h"
16 #include "AliHLTTPCDefinitions.h"
17
18 /** global instance for component registration */
19 AliHLTTPCEsdWriterComponent gTPCEsdWriter;
20
21 /** ROOT macro for the implementation of ROOT specific class methods */
22 ClassImp(AliHLTTPCEsdWriterComponent)
23
24 AliHLTTPCEsdWriterComponent::AliHLTTPCEsdWriterComponent()
25   :
26   fTree(NULL),
27   fESD(NULL)
28 {
29 }
30
31 AliHLTTPCEsdWriterComponent::~AliHLTTPCEsdWriterComponent()
32 {
33 }
34
35 int AliHLTTPCEsdWriterComponent::InitWriter()
36 {
37   int iResult=0;
38   fESD = new AliESD;
39   if (fESD) {
40     fTree = new TTree("esdTree", "Tree with HLT ESD objects");
41     if (fTree) {
42       fTree->Branch("ESD", "AliESD", &fESD);
43     }
44     delete fESD;
45     fESD=NULL;
46   }
47   if (fTree==NULL) {
48     iResult=-ENOMEM;
49   }
50   return iResult;
51 }
52
53 int AliHLTTPCEsdWriterComponent::CloseWriter()
54 {
55   int iResult=0;
56   if (fTree) {
57     WriteObject(kAliHLTVoidEventID, fTree);
58     TTree* pTree=fTree;
59     fTree=NULL;
60     delete pTree;
61   } else {
62     HLTWarning("not initialized");
63   }
64   AliHLTRootFileWriterComponent::CloseWriter();
65 }
66
67 int AliHLTTPCEsdWriterComponent::DumpEvent( const AliHLTComponentEventData& evtData,
68                                             const AliHLTComponentBlockData* blocks, 
69                                             AliHLTComponentTriggerData& trigData )
70 {
71   int iResult=0;
72   TTree* pTree=fTree;
73   if (pTree) {
74     fESD = new AliESD;
75     if (fESD) {
76       AliESD* pESD=fESD;
77
78       const AliHLTComponentBlockData* iter = NULL;
79       AliHLTTPCTrackArray* tracks=NULL;
80       AliHLTTPCTrackletData* inPtr=NULL;
81  
82       for (int ndx=0; ndx<evtData.fBlockCnt && iResult>=0; ndx++) {
83         iter = blocks+ndx;
84         if ( iter->fDataType == AliHLTTPCDefinitions::gkTrackSegmentsDataType ) {
85           if (tracks==NULL) {
86             tracks=new AliHLTTPCTrackArray;
87             if (tracks) {
88               inPtr=(AliHLTTPCTrackletData*)iter->fPtr;
89               HLTDebug("reading block %d: %d tracklets", ndx, inPtr->fTrackletCnt);
90               tracks->FillTracks(inPtr->fTrackletCnt, inPtr->fTracklets);
91               if ((iResult=Tracks2ESD(tracks, pESD))>=0) {
92                 pTree->Fill();
93               }
94             } else {
95               iResult=-ENOMEM;
96             }
97           } else {
98             HLTWarning("can not process more than one track segment data block, "
99                        "please put a track merger in between");
100             break; // don't print the warning again
101           }
102         }
103       }
104
105
106       fESD=NULL;
107       delete pESD;
108     } else {
109       iResult=-ENOMEM;
110     }
111   }
112   return iResult;
113 }
114
115 int AliHLTTPCEsdWriterComponent::ScanArgument(int argc, const char** argv)
116 {
117   int iResult=AliHLTRootFileWriterComponent::ScanArgument(argc, argv);
118   return iResult;
119 }
120
121 int AliHLTTPCEsdWriterComponent::Tracks2ESD(AliHLTTPCTrackArray* pTracks, AliESD* pESD)
122 {
123   int iResult=0;
124   if (pTracks && pESD) {
125     HLTDebug("converting %d tracks from track array", pTracks->GetNTracks());
126     for (int i=0; i<pTracks->GetNTracks() && iResult>=0; i++) {
127       AliHLTTPCTrack* pTrack=(*pTracks)[i];
128       if (pTrack) {
129         int iLocal=pTrack->Convert2AliKalmanTrack();
130         if (iLocal>=0) {
131         AliESDtrack iotrack;
132         iotrack.UpdateTrackParams(pTrack,AliESDtrack::kTPCin);
133         iotrack.SetTPCPoints(pTrack->GetPoints());
134         pESD->AddTrack(&iotrack);
135         } else {
136           HLTError("conversion to AliKalmanTrack failed for track %d of %d", i, pTracks->GetNTracks());
137         }
138       } else {
139         HLTError("internal missmatch in array");
140         iResult=-EFAULT;
141       }
142     }
143     
144   } else {
145     iResult=-EINVAL;
146   }
147   return iResult;
148 }