]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDTrackerComponent.cxx
Correct initialization (freom rev. 31415)
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDTrackerComponent.cxx
1 // $Id$
2
3 /**************************************************************************
4  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Authors: Matthias Richter <Matthias.Richter@ift.uib.no>                *
7  *          Timm Steinbeck <timm@kip.uni-heidelberg.de>                   *
8  *          for The ALICE Off-line 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   AliHLTTRDTrackerComponent.cxx
20     @author Timm Steinbeck, Matthias Richter
21     @date   
22     @brief  A TRDTracker processing component for the HLT. */
23
24 #if __GNUC__ >= 3
25 using namespace std;
26 #endif
27
28 #include "TFile.h"
29 #include "TChain.h"
30
31 #include "AliHLTTRDTrackerComponent.h"
32 #include "AliHLTTRDDefinitions.h"
33 #include "AliCDBManager.h"
34
35 #include "AliTRDReconstructor.h"
36 #include "AliTRDrecoParam.h"
37 #include "AliESDEvent.h"
38 //#include "AliTRDtrackerHLT.h"
39 #include "AliTRDtracker.h"
40 #include "AliTRDCalibraFillHisto.h"
41 #include "AliMagF.h"
42 #include "AliTRDcluster.h"
43 #include "AliESDfriend.h"
44 #include <cstdlib>
45 #include <cerrno>
46 #include <string>
47
48 // this is a global object used for automatic component registration, do not use this
49 AliHLTTRDTrackerComponent gAliHLTTRDTrackerComponent;
50
51 ClassImp(AliHLTTRDTrackerComponent);
52     
53 AliHLTTRDTrackerComponent::AliHLTTRDTrackerComponent()
54   : AliHLTProcessor()
55   , fOutputPercentage(100) // By default we copy to the output exactly what we got as input  
56   , fStrorageDBpath("local://$ALICE_ROOT/OCDB")
57   , fCDB(NULL)
58   , fGeometryFileName("")
59   , fGeometryFile(NULL)
60   , fGeoManager(NULL)
61   , fTracker(NULL)
62 {
63   // Default constructor
64
65   fGeometryFileName = getenv("ALICE_ROOT");
66   fGeometryFileName += "/HLT/TRD/geometry.root";
67 }
68
69 AliHLTTRDTrackerComponent::~AliHLTTRDTrackerComponent()
70 {
71   // Destructor
72 }
73
74 const char* AliHLTTRDTrackerComponent::GetComponentID()
75 {
76   // Return the component ID const char *
77   return "TRDTracker"; // The ID of this component
78 }
79
80 void AliHLTTRDTrackerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
81 {
82   // Get the list of input data  
83   list.clear(); // We do not have any requirements for our input data type(s).
84   list.push_back( AliHLTTRDDefinitions::fgkClusterDataType );
85 }
86
87 AliHLTComponent_DataType AliHLTTRDTrackerComponent::GetOutputDataType()
88 {
89   // Get the output data type
90   return AliHLTTRDDefinitions::fgkClusterDataType;
91 }
92
93 void AliHLTTRDTrackerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
94 {
95   // Get the output data size
96   constBase = 0;
97   inputMultiplier = ((double)fOutputPercentage)/100.0;
98 }
99
100 // Spawn function, return new instance of this class
101 AliHLTComponent* AliHLTTRDTrackerComponent::Spawn()
102 {
103   // Spawn function, return new instance of this class
104   return new AliHLTTRDTrackerComponent;
105 };
106
107 int AliHLTTRDTrackerComponent::DoInit( int argc, const char** argv )
108 {
109   // perform initialization. We check whether our relative output size is specified in the arguments.
110   fOutputPercentage = 100;
111   int i = 0;
112   char* cpErr;
113   while ( i < argc )
114     {
115       HLTDebug("argv[%d] == %s", i, argv[i] );
116       if ( !strcmp( argv[i], "output_percentage" ) )
117         {
118           if ( i+1>=argc )
119             {
120               HLTError("Missing output_percentage parameter");
121               return ENOTSUP;
122             }
123           HLTDebug("Arguments", "argv[%d+1] == %s", i, argv[i+1] );
124           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
125           if ( *cpErr )
126             {
127               HLTError("Cannot convert output_percentage parameter '%s'", argv[i+1] );
128               return EINVAL;
129             }
130           HLTInfo("Output percentage set to %lu %%", fOutputPercentage );
131           i += 2;
132           continue;
133         }
134
135       if ( strcmp( argv[i], "-cdb" ) == 0)
136         {
137           if ( i+1 >= argc )
138             {
139               HLTError("Missing -cdb argument");
140               return ENOTSUP;         
141             }
142           fStrorageDBpath = argv[i+1];
143           HLTInfo("DB storage is %s", fStrorageDBpath.c_str() );          
144           i += 2;
145           continue;
146         }      
147
148       if ( strcmp( argv[i], "-geometry" ) == 0)
149         {
150           if ( i+1 >= argc )
151             {
152               HLTError("Missing -geometry argument");
153               return ENOTSUP;         
154             }
155           fGeometryFileName = argv[i+1];
156           HLTInfo("GeomFile storage is %s", fGeometryFileName.c_str() );          
157           i += 2;
158           continue;
159         }      
160
161       HLTError("Unknown option '%s'", argv[i] );
162       return EINVAL;
163     }
164
165
166   fCDB = AliCDBManager::Instance();
167   if (!fCDB)
168     {
169       HLTError("Could not get CDB instance", "fCDB 0x%x", fCDB);
170     }
171   else
172     {
173       fCDB->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
174       fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
175       HLTDebug("fCDB 0x%x", fCDB);
176     }
177     
178   fGeometryFile = TFile::Open(fGeometryFileName.c_str());
179   if (fGeometryFile)
180     {
181       fGeoManager = (TGeoManager *)fGeometryFile->Get("Geometry");
182       //fTracker = new AliTRDtrackerHLT(fGeometryFile);
183       AliTRDrecoParam *fPars = AliTRDrecoParam::GetLowFluxParam();
184       //fPars->SetSeeding(kTRUE);
185       //fPars->SetStreamLevel(0);
186       AliTRDReconstructor reconstructor; reconstructor.SetRecoParam(fPars);
187       // write clusters [cw] = true
188       // track seeding (stand alone tracking) [sa] = true
189       // PID method in reconstruction (NN) [nn] = true
190       // write online tracklets [tw] = false
191       // drift gas [ar] = false
192       reconstructor.SetOption("cw,sa");
193       fTracker = new AliTRDtracker(fGeometryFile);
194       //fTracker = new AliTRDtracker(fGeometryFile);
195     }
196   else
197     {
198       HLTError("Unable to open file. FATAL!");
199       return -1;
200     }
201
202   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
203   if (calibra == 0)
204     {
205       HLTError("Calibration Histos ::Instance failed");
206       return -1;      
207     }
208   else
209     {
210       calibra->Init2Dhistos();
211     }
212
213   return 0;
214 }
215
216 int AliHLTTRDTrackerComponent::DoDeinit()
217 {
218   // Deinitialization of the component
219
220   delete fTracker;
221   fTracker = 0;
222   
223   if (fGeometryFile)
224     {
225       fGeometryFile->Close();
226       delete fGeometryFile;
227       fGeometryFile = 0;
228     }
229
230   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
231   if (calibra)
232     {
233       // should not write in here!
234       calibra->Write2d();
235       calibra->Destroy();
236     }
237
238   return 0;
239 }
240
241 int AliHLTTRDTrackerComponent::DoEvent( const AliHLTComponentEventData & /*evtData*/,
242                                         AliHLTComponentTriggerData & /*trigData*/ )
243 {
244   // Process an event
245   
246   HLTInfo("Output percentage set to %lu %%", fOutputPercentage );
247   HLTInfo("NofBlocks %lu", GetNumberOfInputBlocks() );
248
249   AliHLTUInt32_t dBlockSpecification = 0;
250
251   //implement a usage of the following
252   //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
253   //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
254   //   void *triggerData = trigData.fData;
255   //HLTDebug("Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
256
257   AliHLTComponentBlockData *dblock = (AliHLTComponentBlockData *)GetFirstInputBlock( AliHLTTRDDefinitions::fgkClusterDataType );
258   if (dblock != 0)
259     {
260       dBlockSpecification = dblock->fSpecification;
261     }
262   else
263     {
264       HLTWarning("First Input Block not found! 0x%x", dblock);
265       return -1;
266     }
267
268   int ibForce = 0;
269   TObject *tobjin = (TObject *)GetFirstInputObject( AliHLTTRDDefinitions::fgkClusterDataType, "TTree", ibForce);
270   HLTInfo("Pointer = 0x%x", tobjin);
271
272   TTree *clusterTree = (TTree*)tobjin;
273   if (!clusterTree)
274     {
275       HLTWarning("First Input Block not a tree! 0x%x", tobjin);
276       return -1;
277     }
278
279   HLTInfo("Pointer = 0x%x Name = %s", clusterTree, clusterTree->GetName());
280
281   while (tobjin != 0)
282     {
283       if (clusterTree)
284         {
285           HLTInfo("Pointer = 0x%x Name = %s", clusterTree, clusterTree->GetName());
286           Int_t iNentries = clusterTree->GetEntries();
287           HLTInfo("N of tree entries = %d", iNentries);
288           fTracker->LoadClusters(clusterTree);
289         }
290       else
291         {
292           HLTError("Tree Pointer = 0x%x", clusterTree);
293         }
294
295       tobjin = (TObject *)GetNextInputObject( ibForce );
296       HLTInfo("Pointer = 0x%x", tobjin);
297       clusterTree = (TTree*)tobjin;
298     }
299
300   fTracker->SetAddTRDseeds();
301
302   AliESDfriend *esdFriend = new AliESDfriend();
303
304   AliESDEvent *esd = new AliESDEvent();
305   esd->CreateStdContent();
306   fTracker->PropagateBack(esd);
307   fTracker->RefitInward(esd);
308   fTracker->Clusters2Tracks(esd);
309
310   esd->GetESDfriend(esdFriend);
311   //here transport the esd tracks further
312
313   Int_t nTracks = esd->GetNumberOfTracks();
314   Int_t nTRDTracks = esd->GetNumberOfTrdTracks();
315   HLTInfo( "Number of tracks %d Number of TRD tracks %d", nTracks, nTRDTracks);
316
317
318   for (Int_t it = 0; it < nTracks; it++)
319     {
320       AliESDtrack* track = esd->GetTrack(it);
321       HLTInfo("Track %d 0x%x Pt %1.2f", it, track, track->Pt());
322       PushBack(track, AliHLTTRDDefinitions::fgkTRDSATracksDataType, ++dBlockSpecification);
323     }
324
325   delete esd;
326   delete esdFriend;
327
328   delete clusterTree;
329
330   return 0;
331 }