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