]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCEsdWriterComponent.cxx
documentation
[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 #include "AliHLTTPCEsdWriterComponent.h"
27 #include "AliESD.h"
28 #include "TTree.h"
29 #include "AliHLTTPCTrack.h"
30 #include "AliHLTTPCTrackArray.h"
31 #include "AliHLTTPCTrackletDataFormat.h"
32 #include "AliHLTTPCDefinitions.h"
33
34 /** global instance for component registration */
35 AliHLTTPCEsdWriterComponent gTPCEsdWriter;
36
37 /** ROOT macro for the implementation of ROOT specific class methods */
38 ClassImp(AliHLTTPCEsdWriterComponent)
39
40 AliHLTTPCEsdWriterComponent::AliHLTTPCEsdWriterComponent()
41   :
42   fTree(NULL),
43   fESD(NULL)
44 {
45   // see header file for class documentation
46   // or
47   // refer to README to build package
48   // or
49   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
50 }
51
52 AliHLTTPCEsdWriterComponent::~AliHLTTPCEsdWriterComponent()
53 {
54   // see header file for class documentation
55 }
56
57 int AliHLTTPCEsdWriterComponent::InitWriter()
58 {
59   // see header file for class documentation
60   int iResult=0;
61   fESD = new AliESD;
62   if (fESD) {
63     fTree = new TTree("esdTree", "Tree with HLT ESD objects");
64     if (fTree) {
65       fTree->Branch("ESD", "AliESD", &fESD);
66     }
67     delete fESD;
68     fESD=NULL;
69   }
70   if (fTree==NULL) {
71     iResult=-ENOMEM;
72   }
73   return iResult;
74 }
75
76 int AliHLTTPCEsdWriterComponent::CloseWriter()
77 {
78   // see header file for class documentation
79   int iResult=0;
80   if (fTree) {
81     WriteObject(kAliHLTVoidEventID, fTree);
82     TTree* pTree=fTree;
83     fTree=NULL;
84     delete pTree;
85   } else {
86     HLTWarning("not initialized");
87   }
88   iResult=AliHLTRootFileWriterComponent::CloseWriter();
89   return iResult;
90 }
91
92 int AliHLTTPCEsdWriterComponent::DumpEvent( const AliHLTComponentEventData& evtData,
93                                             const AliHLTComponentBlockData* blocks, 
94                                             AliHLTComponentTriggerData& trigData )
95 {
96   // see header file for class documentation
97   int iResult=0;
98   TTree* pTree=fTree;
99   if (pTree) {
100     fESD = new AliESD;
101     if (fESD) {
102       AliESD* pESD=fESD;
103
104       const AliHLTComponentBlockData* iter = NULL;
105       AliHLTTPCTrackletData* inPtr=NULL;
106       int bIsTrackSegs=0;
107  
108       for (int ndx=0; ndx<(int)evtData.fBlockCnt && iResult>=0; ndx++) {
109         iter = blocks+ndx;
110         if ( (bIsTrackSegs=(iter->fDataType == AliHLTTPCDefinitions::fgkTrackSegmentsDataType))==1 ||
111              iter->fDataType == AliHLTTPCDefinitions::fgkTracksDataType ) {
112           Int_t minslice=AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
113           Int_t maxslice=AliHLTTPCDefinitions::GetMaxSliceNr(iter->fSpecification);
114           if (bIsTrackSegs==0) {
115             // slice parameter and data specification ignored, tracks already in global coordinates
116             minslice=0;
117             maxslice=0;
118           }
119           //HLTDebug("dataspec %#x minslice %d", iter->fSpecification, minslice);
120           if (minslice >=0 && minslice<36) {
121             if (minslice!=maxslice) {
122               HLTWarning("data from multiple sectors in one block: "
123                          "possible missmatch in treatment of local coordinate system");
124             }
125             AliHLTTPCTrackArray tracks;
126             inPtr=(AliHLTTPCTrackletData*)iter->fPtr;
127             HLTDebug("reading block %d (slice %d): %d tracklets", ndx, minslice, inPtr->fTrackletCnt);
128             tracks.FillTracks(inPtr->fTrackletCnt, inPtr->fTracklets, minslice, 0/*don't rotate*/);
129             if ((iResult=Tracks2ESD(&tracks, pESD))>=0) {
130             }
131           } else {
132             HLTError("invalid sector number");
133             iResult=-EBADF;
134           }
135         }
136       }
137       if (iResult>=0) {
138         pTree->Fill();
139       }
140
141       fESD=NULL;
142       delete pESD;
143     } else {
144       iResult=-ENOMEM;
145     }
146   }
147   return iResult;
148 }
149
150 int AliHLTTPCEsdWriterComponent::ScanArgument(int argc, const char** argv)
151 {
152   // see header file for class documentation
153   int iResult=AliHLTRootFileWriterComponent::ScanArgument(argc, argv);
154   return iResult;
155 }
156
157 int AliHLTTPCEsdWriterComponent::Tracks2ESD(AliHLTTPCTrackArray* pTracks, AliESD* pESD)
158 {
159   // see header file for class documentation
160   int iResult=0;
161   if (pTracks && pESD) {
162     HLTDebug("converting %d tracks from track array", pTracks->GetNTracks());
163     for (int i=0; i<pTracks->GetNTracks() && iResult>=0; i++) {
164       AliHLTTPCTrack* pTrack=(*pTracks)[i];
165       if (pTrack) {
166         //HLTDebug("convert track %d", i);
167         //pTrack->Print();
168         int iLocal=pTrack->Convert2AliKalmanTrack();
169         if (iLocal>=0) {
170         AliESDtrack iotrack;
171         iotrack.UpdateTrackParams(pTrack,AliESDtrack::kTPCin);
172         iotrack.SetTPCPoints(pTrack->GetPoints());
173         pESD->AddTrack(&iotrack);
174         } else {
175           HLTError("conversion to AliKalmanTrack failed for track %d of %d", i, pTracks->GetNTracks());
176         }
177       } else {
178         HLTError("internal missmatch in array");
179         iResult=-EFAULT;
180       }
181     }
182     
183   } else {
184     iResult=-EINVAL;
185   }
186   return iResult;
187 }