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