]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDTrackerComponent.cxx
Fix needed by the new TRD reconstruction (Christoph)
[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 "AliESDEvent.h"
37 //#include "AliTRDtrackerHLT.h"
38 #include "AliTRDtracker.h"
39 #include "AliTRDCalibraFillHisto.h"
40 #include "AliMagFMaps.h"
41 #include "AliTRDcluster.h"
42 #include "AliESDfriend.h"
43 #include <cstdlib>
44 #include <cerrno>
45 #include <string>
46
47 // this is a global object used for automatic component registration, do not use this
48 AliHLTTRDTrackerComponent gAliHLTTRDTrackerComponent;
49
50 ClassImp(AliHLTTRDTrackerComponent);
51     
52 AliHLTTRDTrackerComponent::AliHLTTRDTrackerComponent()
53   : AliHLTProcessor()
54   , fOutputPercentage(100) // By default we copy to the output exactly what we got as input  
55   , fStrorageDBpath("local://$ALICE_ROOT")
56   , fCDB(NULL)
57   , fField(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       Logging( kHLTLogDebug, "HLT::TRDTracker::DoInit", "Arguments", "argv[%d] == %s", i, argv[i] );
116       if ( !strcmp( argv[i], "output_percentage" ) )
117         {
118           if ( i+1>=argc )
119             {
120               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing output_percentage parameter");
121               return ENOTSUP;
122             }
123           Logging( kHLTLogDebug, "HLT::TRDTracker::DoInit", "Arguments", "argv[%d+1] == %s", i, argv[i+1] );
124           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
125           if ( *cpErr )
126             {
127               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Wrong Argument", "Cannot convert output_percentage parameter '%s'", argv[i+1] );
128               return EINVAL;
129             }
130           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "Output percentage set", "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               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing -cdb argument");
140               return ENOTSUP;         
141             }
142           fStrorageDBpath = argv[i+1];
143           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "DB storage set", "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               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing -geometry argument");
153               return ENOTSUP;         
154             }
155           fGeometryFileName = argv[i+1];
156           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "GeomFile storage set", "GeomFile storage is %s", 
157                    fGeometryFileName.c_str() );   
158           i += 2;
159           continue;
160         }      
161
162       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
163       return EINVAL;
164     }
165
166   //init alifield map - temporarly fixed - should come from a DB
167   fField = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
168   if (fField)
169     AliTracker::SetFieldMap(fField,1);
170   else
171     Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Field", "Unable to init the field");
172
173   fCDB = AliCDBManager::Instance();
174   if (!fCDB)
175     {
176       Logging(kHLTLogError, "HLT::TRDCalibration::DoInit", "Could not get CDB instance", "fCDB 0x%x", fCDB);
177     }
178   else
179     {
180       fCDB->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
181       fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
182       Logging(kHLTLogDebug, "HLT::TRDCalibration::DoInit", "CDB instance", "fCDB 0x%x", fCDB);
183     }
184     
185   fGeometryFile = TFile::Open(fGeometryFileName.c_str());
186   if (fGeometryFile)
187     {
188       fGeoManager = (TGeoManager *)fGeometryFile->Get("Geometry");
189       //fTracker = new AliTRDtrackerHLT(fGeometryFile);
190       fTracker = new AliTRDtracker(fGeometryFile);
191       //fTracker = new AliTRDtracker(fGeometryFile);
192     }
193   else
194     {
195       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "fGeometryFile", "Unable to open file. FATAL!");
196       return -1;
197     }
198
199   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
200   if (calibra == 0)
201     {
202       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Calibration Histos", "::Instance failed");
203       return -1;      
204     }
205   else
206     {
207       calibra->Init2Dhistos();
208     }
209
210   return 0;
211 }
212
213 int AliHLTTRDTrackerComponent::DoDeinit()
214 {
215   // Deinitialization of the component
216
217   delete fField;
218   fField = 0;
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       Logging( kHLTLogInfo, "HLT::TRDTracker::DoDeinit", "CALIBRA", "before destroy");
236       calibra->Destroy();
237       Logging( kHLTLogInfo, "HLT::TRDTracker::DoDeinit", "CALIBRA", "after destroy");
238     }
239
240   return 0;
241 }
242
243 int AliHLTTRDTrackerComponent::DoEvent( const AliHLTComponentEventData & evtData,
244                                         AliHLTComponentTriggerData & trigData )
245 {
246   // Process an event
247   
248   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
249   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "BLOCKS", "NofBlocks %lu", evtData.fBlockCnt );
250
251   AliHLTUInt32_t dBlockSpecification = 0;
252
253   //implement a usage of the following
254 //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
255 //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
256 //   void *triggerData = trigData.fData;
257   Logging( kHLTLogDebug, "HLT::TRDTracker::DoEvent", "Trigger data received", 
258            "Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
259
260   AliHLTComponentBlockData *dblock = (AliHLTComponentBlockData *)GetFirstInputBlock( AliHLTTRDDefinitions::fgkClusterDataType );
261   if (dblock != 0)
262     {
263       dBlockSpecification = dblock->fSpecification;
264     }
265   else
266     {
267       Logging( kHLTLogWarning, "HLT::TRDTracker::DoEvent", "DATAIN", "First Input Block not found! 0x%x", dblock);
268       return -1;
269     }
270
271   int ibForce = 0;
272   TObject *tobjin = (TObject *)GetFirstInputObject( AliHLTTRDDefinitions::fgkClusterDataType, "TTree", ibForce);
273   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "1stBLOCK", "Pointer = 0x%x", tobjin);
274
275   TTree *clusterTree = (TTree*)tobjin;
276   if (!clusterTree)
277     {
278       Logging( kHLTLogWarning, "HLT::TRDTracker::DoEvent", "DATAIN", "First Input Block not a tree! 0x%x", tobjin);
279       return -1;
280     }
281
282   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "1stBLOCK", "Pointer = 0x%x Name = %s", clusterTree, clusterTree->GetName());
283
284   while (tobjin != 0)
285     {
286       if (clusterTree)
287         {
288           Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "CLUSTERS", "Pointer = 0x%x Name = %s", clusterTree, clusterTree->GetName());
289           Int_t iNentries = clusterTree->GetEntries();
290           Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "COUNT", "N of tree entries = %d", iNentries);
291           fTracker->LoadClusters(clusterTree);
292         }
293       else
294         {
295           Logging( kHLTLogError, "HLT::TRDTracker::DoEvent", "CLUSTERS", "Tree Pointer = 0x%x", clusterTree);
296         }
297
298       tobjin = (TObject *)GetNextInputObject( ibForce );
299       Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "nextBLOCK", "Pointer = 0x%x", tobjin);
300       clusterTree = (TTree*)tobjin;
301     }
302
303   AliTRDReconstructor::SetSeedingOn(kTRUE);
304
305   fTracker->SetAddTRDseeds();
306
307   AliESDfriend *esdFriend = new AliESDfriend();
308
309   AliESDEvent *esd = new AliESDEvent();
310   esd->CreateStdContent();
311   fTracker->PropagateBack(esd);
312   fTracker->RefitInward(esd);
313   fTracker->Clusters2Tracks(esd);
314
315   esd->GetESDfriend(esdFriend);
316   //here transport the esd tracks further
317
318   Int_t nTracks = esd->GetNumberOfTracks();
319   Int_t nTRDTracks = esd->GetNumberOfTrdTracks();
320   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "Number of tracks %d Number of TRD tracks %d", nTracks, nTRDTracks);
321
322 //   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
323 //   calibra->Init2Dhistostrack();
324
325   for (Int_t it = 0; it < nTracks; it++)
326     {
327       AliESDtrack* track = esd->GetTrack(it);
328
329 //       Int_t nCalibObjects = 0;
330 //       Int_t idx = 0;
331 //       while (track->GetCalibObject(idx) != 0)
332 //      {
333 //        nCalibObjects++;
334 //        idx++;
335 //      }
336 //       Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "Track 0x%x NcalibObjects %d", track, nCalibObjects);
337
338       Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "Track %d 0x%x Pt %1.2f", it, track, track->Pt());
339       PushBack(track, AliHLTTRDDefinitions::fgkTRDSATracksDataType, ++dBlockSpecification);
340 //       if (calibra->GetMItracking())
341 //      {
342 //        calibra->UpdateHistograms(track);
343 //      }
344     }
345
346   //PushBack(esd, AliHLTTRDDefinitions::fgkTRDSAEsdDataType, dBlockSpecification);
347   //PushBack(esdFriend, AliHLTTRDDefinitions::fgkTRDSAEsdDataType, dBlockSpecification);
348
349   //no receiver defined yet(!)
350   //Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "now deleting");
351   delete esd;
352   delete esdFriend;
353
354   //Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "after delete esd");
355   delete clusterTree;
356
357   //Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "after delete clusterTree");
358   return 0;
359 }