]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDTrackerComponent.cxx
Update and cleanup triggered by the offline developments.
[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
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   fCDB = AliCDBManager::Instance();
65   //fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
66   fCDB->SetRun(0);
67
68   fGeometryFileName = getenv("ALICE_ROOT");
69   fGeometryFileName += "/HLT/TRD/geometry.root";
70 }
71
72 AliHLTTRDTrackerComponent::~AliHLTTRDTrackerComponent()
73 {
74   // Destructor
75 }
76
77 const char* AliHLTTRDTrackerComponent::GetComponentID()
78 {
79   // Return the component ID const char *
80   return "TRDTracker"; // The ID of this component
81 }
82
83 void AliHLTTRDTrackerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
84 {
85   // Get the list of input data  
86   list.clear(); // We do not have any requirements for our input data type(s).
87   list.push_back( AliHLTTRDDefinitions::fgkClusterDataType );
88 }
89
90 AliHLTComponent_DataType AliHLTTRDTrackerComponent::GetOutputDataType()
91 {
92   // Get the output data type
93   return AliHLTTRDDefinitions::fgkClusterDataType;
94 }
95
96 void AliHLTTRDTrackerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
97 {
98   // Get the output data size
99   constBase = 0;
100   inputMultiplier = ((double)fOutputPercentage)/100.0;
101 }
102
103 // Spawn function, return new instance of this class
104 AliHLTComponent* AliHLTTRDTrackerComponent::Spawn()
105 {
106   // Spawn function, return new instance of this class
107   return new AliHLTTRDTrackerComponent;
108 };
109
110 int AliHLTTRDTrackerComponent::DoInit( int argc, const char** argv )
111 {
112   // perform initialization. We check whether our relative output size is specified in the arguments.
113   fOutputPercentage = 100;
114   int i = 0;
115   char* cpErr;
116   while ( i < argc )
117     {
118       Logging( kHLTLogDebug, "HLT::TRDTracker::DoInit", "Arguments", "argv[%d] == %s", i, argv[i] );
119       if ( !strcmp( argv[i], "output_percentage" ) )
120         {
121           if ( i+1>=argc )
122             {
123               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing output_percentage parameter");
124               return ENOTSUP;
125             }
126           Logging( kHLTLogDebug, "HLT::TRDTracker::DoInit", "Arguments", "argv[%d+1] == %s", i, argv[i+1] );
127           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
128           if ( *cpErr )
129             {
130               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Wrong Argument", "Cannot convert output_percentage parameter '%s'", argv[i+1] );
131               return EINVAL;
132             }
133           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
134           i += 2;
135           continue;
136         }
137
138       if ( strcmp( argv[i], "-cdb" ) == 0)
139         {
140           if ( i+1 >= argc )
141             {
142               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing -cdb argument");
143               return ENOTSUP;         
144             }
145           fStrorageDBpath = argv[i+1];
146           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "DB storage set", "DB storage is %s", fStrorageDBpath.c_str() );       
147           fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
148           i += 2;
149           continue;
150         }      
151
152       if ( strcmp( argv[i], "-geometry" ) == 0)
153         {
154           if ( i+1 >= argc )
155             {
156               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing -geometry argument");
157               return ENOTSUP;         
158             }
159           fGeometryFileName = argv[i+1];
160           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "GeomFile storage set", "GeomFile storage is %s", 
161                    fGeometryFileName.c_str() );   
162           i += 2;
163           continue;
164         }      
165
166       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
167       return EINVAL;
168     }
169
170   //init alifield map - temporarly fixed - should come from a DB
171   fField = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
172   if (fField)
173     AliTracker::SetFieldMap(fField,1);
174   else
175     Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Field", "Unable to init the field");
176     
177   fGeometryFile = TFile::Open(fGeometryFileName.c_str());
178   if (fGeometryFile)
179     {
180       fGeoManager = (TGeoManager *)fGeometryFile->Get("Geometry");
181       //fTracker = new AliTRDtrackerHLT(fGeometryFile);
182       fTracker = new AliTRDtracker(fGeometryFile);
183       //fTracker = new AliTRDtracker(fGeometryFile);
184     }
185   else
186     {
187       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "fGeometryFile", "Unable to open file. FATAL!");
188       return -1;
189     }
190
191   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
192   if (calibra == 0)
193     {
194       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Calibration Histos", "::Instance failed");
195       return -1;      
196     }
197   else
198     {
199       calibra->SetMITracking(kTRUE);
200       calibra->Init2Dhistos();
201     }
202
203   return 0;
204 }
205
206 int AliHLTTRDTrackerComponent::DoDeinit()
207 {
208   // Deinitialization of the component
209
210   delete fField;
211   fField = 0;
212
213   delete fTracker;
214   fTracker = 0;
215   
216   if (fGeometryFile)
217     {
218       fGeometryFile->Close();
219       delete fGeometryFile;
220       fGeometryFile = 0;
221     }
222
223   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
224   if (calibra)
225     {
226       // should not write in here!
227       calibra->Write2d();
228       Logging( kHLTLogInfo, "HLT::TRDTracker::DoDeinit", "CALIBRA", "before destroy");
229       calibra->Destroy();
230       Logging( kHLTLogInfo, "HLT::TRDTracker::DoDeinit", "CALIBRA", "after destroy");
231     }
232
233   return 0;
234 }
235
236 int AliHLTTRDTrackerComponent::DoEvent( const AliHLTComponentEventData & evtData,
237                                         AliHLTComponentTriggerData & trigData )
238 {
239   // Process an event
240   
241   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
242   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "BLOCKS", "NofBlocks %lu", evtData.fBlockCnt );
243
244   AliHLTUInt32_t fDblock_Specification = 0;
245
246   //implement a usage of the following
247 //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
248 //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
249 //   void *triggerData = trigData.fData;
250   Logging( kHLTLogDebug, "HLT::TRDClusterizer::DoEvent", "Trigger data received", 
251            "Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
252
253   AliHLTComponentBlockData *dblock = (AliHLTComponentBlockData *)GetFirstInputBlock( AliHLTTRDDefinitions::fgkClusterDataType );
254   if (dblock != 0)
255     {
256       fDblock_Specification = dblock->fSpecification;
257     }
258   else
259     {
260       Logging( kHLTLogWarning, "HLT::TRDTracker::DoEvent", "DATAIN", "First Input Block not found! 0x%x", dblock);
261       return -1;
262     }
263
264   int ibForce = 0;
265   TObject *tobjin = (TObject *)GetFirstInputObject( AliHLTTRDDefinitions::fgkClusterDataType, "TTree", ibForce);
266   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "1stBLOCK", "Pointer = 0x%x", tobjin);
267
268   TTree *clusterTree = (TTree*)tobjin;
269   if (!clusterTree)
270     {
271       Logging( kHLTLogWarning, "HLT::TRDTracker::DoEvent", "DATAIN", "First Input Block not a tree! 0x%x", tobjin);
272       return -1;
273     }
274
275   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "1stBLOCK", "Pointer = 0x%x Name = %s", clusterTree, clusterTree->GetName());
276
277   while (tobjin != 0)
278     {
279       if (clusterTree)
280         {
281           Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "CLUSTERS", "Pointer = 0x%x Name = %s", clusterTree, clusterTree->GetName());
282           Int_t iNentries = clusterTree->GetEntries();
283           Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "COUNT", "N of tree entries = %d", iNentries);
284           fTracker->LoadClusters(clusterTree);
285         }
286       else
287         {
288           Logging( kHLTLogError, "HLT::TRDTracker::DoEvent", "CLUSTERS", "Tree Pointer = 0x%x", clusterTree);
289         }
290
291       tobjin = (TObject *)GetNextInputObject( ibForce );
292       Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "nextBLOCK", "Pointer = 0x%x", tobjin);
293       clusterTree = (TTree*)tobjin;
294     }
295
296   AliTRDReconstructor::SetSeedingOn(kTRUE);
297
298   AliESDEvent *esd = new AliESDEvent();
299   esd->CreateStdContent();
300   fTracker->PropagateBack(esd);
301   fTracker->RefitInward(esd);
302
303   //here transport the esd tracks further
304
305   Int_t nTracks = esd->GetNumberOfTracks();
306   Int_t nTRDTracks = esd->GetNumberOfTrdTracks();
307   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "Number of tracks %d Number of TRD tracks %d", nTracks, nTRDTracks);
308
309 //   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
310 //   calibra->Init2Dhistostrack();
311
312   for (Int_t it = 0; it < nTracks; it++)
313     {
314       AliESDtrack* track = esd->GetTrack(it);
315       Int_t nCalibObjects = 0;
316       Int_t idx = 0;
317       while (track->GetCalibObject(idx) != 0)
318         {
319           nCalibObjects++;
320           idx++;
321         }
322       Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "Track 0x%x NcalibObjects %d", track, nCalibObjects);
323
324       PushBack(track, AliHLTTRDDefinitions::fgkTRDSATracksDataType, fDblock_Specification);
325 //       if (calibra->GetMItracking())
326 //      {
327 //        calibra->UpdateHistograms(track);
328 //      }
329     }
330
331   //no receiver defined yet(!)
332   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "now deleting");
333   delete esd;
334   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "after delete esd");
335   delete clusterTree;
336
337   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "DONE", "after delete clusterTree");
338   return 0;
339 }