]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/EMCAL/AliHLTEMCALTrackerComponent.cxx
Initial HLT-EMCAL components.
[u/mrichter/AliRoot.git] / HLT / EMCAL / AliHLTEMCALTrackerComponent.cxx
1 // $Id: AliHLTEMCALTrackerComponent.cxx 23618 2008-01-29 13:07:38Z hristov $
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   AliHLTEMCALTrackerComponent.cxx
20     @author Timm Steinbeck, Matthias Richter
21     @date   
22     @brief  A EMCALTracker processing component for the HLT. */
23
24 #if __GNUC__ >= 3
25 using namespace std;
26 #endif
27
28 #include "TFolder.h"
29 #include "TFile.h"
30
31 #include "AliHLTEMCALTrackerComponent.h"
32 #include "AliHLTEMCALDefinitions.h"
33 #include "AliCDBManager.h"
34 #include "AliEMCALTracker.h"
35 #include "AliEMCALReconstructor.h"
36 #include "AliESDEvent.h"
37 #include "AliMagFMaps.h"
38 #include "AliESDfriend.h"
39
40 #include <cstdlib>
41 #include <cerrno>
42 #include <string>
43
44 // this is a global object used for automatic component registration, do not use this
45 AliHLTEMCALTrackerComponent gAliHLTEMCALTrackerComponent;
46
47 ClassImp(AliHLTEMCALTrackerComponent);
48     
49 AliHLTEMCALTrackerComponent::AliHLTEMCALTrackerComponent()
50   : AliHLTProcessor()
51   , fOutputPercentage(100) // By default we copy to the output exactly what we got as input  
52   , fStrorageDBpath("local://$ALICE_ROOT")
53   , fCDB(NULL)
54   , fField(NULL)
55   , fGeometryFileName("")
56   , fGeometryFile(NULL)
57   , fGeoManager(NULL)
58   , fTracker(NULL)
59   , fInputFolder(new TFolder("EMCALtrackerFolder", "EMCALtrackerFolder"))
60 {
61   // Default constructor
62   fInputFolder->SetOwner(kTRUE);
63
64   fGeometryFileName = getenv("ALICE_ROOT");
65   fGeometryFileName += "/HLT/EMCAL/geometry.root";
66 }
67
68 AliHLTEMCALTrackerComponent::~AliHLTEMCALTrackerComponent()
69 {
70   // Destructor
71   delete fInputFolder;
72   fInputFolder = NULL;
73 }
74
75 const char* AliHLTEMCALTrackerComponent::GetComponentID()
76 {
77   // Return the component ID const char *
78   return "EMCALTracker"; // The ID of this component
79 }
80
81 void AliHLTEMCALTrackerComponent::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( AliHLTEMCALDefinitions::fgkClusterDataType );
86 }
87
88 AliHLTComponent_DataType AliHLTEMCALTrackerComponent::GetOutputDataType()
89 {
90   // Get the output data type
91   return AliHLTEMCALDefinitions::fgkEMCALESDDataType;
92 }
93
94 void AliHLTEMCALTrackerComponent::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* AliHLTEMCALTrackerComponent::Spawn()
103 {
104   // Spawn function, return new instance of this class
105   return new AliHLTEMCALTrackerComponent;
106 };
107
108 int AliHLTEMCALTrackerComponent::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::EMCALTracker::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::EMCALTracker::DoInit", "Missing Argument", "Missing output_percentage parameter");
122               return ENOTSUP;
123             }
124           Logging( kHLTLogDebug, "HLT::EMCALTracker::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::EMCALTracker::DoInit", "Wrong Argument", "Cannot convert output_percentage parameter '%s'", argv[i+1] );
129               return EINVAL;
130             }
131           Logging( kHLTLogInfo, "HLT::EMCALTracker::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::EMCALTracker::DoInit", "Missing Argument", "Missing -cdb argument");
141               return ENOTSUP;         
142             }
143           fStrorageDBpath = argv[i+1];
144           Logging( kHLTLogInfo, "HLT::EMCALTracker::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::EMCALTracker::DoInit", "Missing Argument", "Missing -geometry argument");
154               return ENOTSUP;         
155             }
156           fGeometryFileName = argv[i+1];
157           Logging( kHLTLogInfo, "HLT::EMCALTracker::DoInit", "GeomFile storage set", "GeomFile storage is %s", 
158                    fGeometryFileName.c_str() );   
159           i += 2;
160           continue;
161         }      
162
163       Logging(kHLTLogError, "HLT::EMCALTracker::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::EMCALTracker::DoInit", "Field", "Unable to init the field");
173
174   fCDB = AliCDBManager::Instance();
175   if (!fCDB)
176     {
177       Logging(kHLTLogError, "HLT::EMCALCalibration::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::EMCALCalibration::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 AliEMCALTracker;
191     }
192   else
193     {
194       Logging(kHLTLogError, "HLT::EMCALTracker::DoInit", "fGeometryFile", "Unable to open file. FATAL!");
195       return -1;
196     }
197
198   return 0;
199 }
200
201 int AliHLTEMCALTrackerComponent::DoDeinit()
202 {
203   // Deinitialization of the component
204
205   delete fField;
206   fField = 0;
207
208   delete fTracker;
209   fTracker = 0;
210   
211   if (fGeometryFile)
212     {
213       fGeometryFile->Close();
214       delete fGeometryFile;
215       fGeometryFile = 0;
216     }
217
218   fInputFolder->Clear();
219
220   return 0;
221 }
222
223 int AliHLTEMCALTrackerComponent::DoEvent( const AliHLTComponentEventData & evtData,
224                                         AliHLTComponentTriggerData & trigData )
225 {
226   // Process an event
227   
228   Logging( kHLTLogInfo, "HLT::EMCALTracker::DoEvent", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
229   Logging( kHLTLogInfo, "HLT::EMCALTracker::DoEvent", "BLOCKS", "NofBlocks %lu", evtData.fBlockCnt );
230
231   AliHLTUInt32_t dBlockSpecification = 0;
232
233   //implement a usage of the following
234 //   AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
235 //   AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
236 //   void *triggerData = trigData.fData;
237   Logging( kHLTLogDebug, "HLT::EMCALTracker::DoEvent", "Trigger data received", 
238            "Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
239
240   AliHLTComponentBlockData *dblock = (AliHLTComponentBlockData *)GetFirstInputBlock( AliHLTEMCALDefinitions::fgkClusterDataType );
241   if (dblock != 0)
242     {
243       dBlockSpecification = dblock->fSpecification;
244     }
245   else
246     {
247       Logging( kHLTLogWarning, "HLT::EMCALTracker::DoEvent", "DATAIN", "First Input Block not found! 0x%x", dblock);
248       return -1;
249     }
250
251   fInputFolder->Clear();
252   AliESDEvent *esd = 0; // we assume we receive this one from a global merger component
253   TTree *clusterTree = 0;
254   
255   int ibForce = 0;
256   TObject *tobjin = 0;
257
258   // here getfirstinput finds the first object with spec type..
259   // get the clusters tree
260   tobjin = (TObject *)GetFirstInputObject( AliHLTEMCALDefinitions::fgkClusterDataType, "TTree", ibForce);
261   while (tobjin != 0)
262     {
263       tobjin = (TObject *)GetNextInputObject( ibForce );
264       Logging( kHLTLogInfo, "HLT::EMCALTracker::DoEvent", "nextBLOCK", "Pointer = 0x%x", tobjin);
265       clusterTree = (TTree*)tobjin;
266       if (clusterTree != 0)
267         {
268           Int_t iNentries = clusterTree->GetEntries();
269           Logging( kHLTLogInfo, "HLT::EMCALTracker::DoEvent", "COUNT", "N of tree entries = %d", iNentries);
270           fTracker->LoadClusters(clusterTree);
271           fInputFolder->Add(clusterTree);      
272         }
273     }
274
275   // now get the ESD(s) - should in principle be only one...
276   tobjin = (TObject *)GetFirstInputObject( AliHLTEMCALDefinitions::fgkESDDataType, "AliESDevent", ibForce);
277   while (tobjin != 0)
278     {
279       esd = (AliESDEvent *)tobjin;
280       if (esd != 0)
281         {
282           HLTInfo("Got ESDevent");
283           fInputFolder->Add(esd);         
284         }
285       tobjin = (TObject *)GetNextInputObject( ibForce );
286       Logging( kHLTLogInfo, "HLT::EMCALTracker::DoEvent", "nextBLOCK", "Pointer = 0x%x", tobjin);
287       clusterTree = (TTree*)tobjin;
288     }
289
290   esd = (AliESDEvent*)fInputFolder->FindObject("AliESDevent");
291   if (esd != 0)
292     {
293       fTracker->PropagateBack(esd);
294       
295       //here transport the esd tracks further
296       Int_t nTracks = esd->GetNumberOfTracks();
297       HLTInfo("Number of tracks %d", nTracks);  
298       //esd->Print();
299       PushBack(esd, AliHLTEMCALDefinitions::fgkEMCALESDDataType);      
300     }
301   else
302     {
303       HLTError("No ESD events received!");
304     }
305
306   fTracker->UnloadClusters();  
307   fInputFolder->Clear();
308
309   HLTDebug("Event done.");
310   return 0;
311 }