]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDTrackerComponent.cxx
Magnetic field settings fix for AliTracker::SetField
[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
30 #include "AliHLTTRDTrackerComponent.h"
31 #include "AliHLTTRDDefinitions.h"
32 #include "AliCDBManager.h"
33
34 #include "AliTRDclusterizerV1HLT.h"
35 #include "AliTRDReconstructor.h"
36 #include "AliESD.h"
37 #include "AliTRDtrackerHLT.h"
38 #include "AliTRDCalibraFillHisto.h"
39 #include "AliMagFMaps.h"
40 #include "AliTRDcluster.h"
41 #include "TObjArray.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 {
54   fOutputPercentage = 100; // By default we copy to the output exactly what we got as input
55
56   fStrorageDBpath = "local://$ALICE_ROOT";
57   fCDB = AliCDBManager::Instance();
58   //fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
59   fCDB->SetRun(0);
60
61   fGeometryFileName = getenv("ALICE_ROOT");
62   fGeometryFileName += "/HLT/TRD/geometry.root";
63 }
64
65 AliHLTTRDTrackerComponent::~AliHLTTRDTrackerComponent()
66 {
67 }
68
69 const char* AliHLTTRDTrackerComponent::GetComponentID()
70 {
71   return "TRDTracker"; // The ID of this component
72 }
73
74 void AliHLTTRDTrackerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
75 {
76   list.clear(); // We do not have any requirements for our input data type(s).
77   list.push_back( AliHLTTRDDefinitions::gkClusterDataType );
78 }
79
80 AliHLTComponent_DataType AliHLTTRDTrackerComponent::GetOutputDataType()
81 {
82   return AliHLTTRDDefinitions::gkClusterDataType;
83 }
84
85 void AliHLTTRDTrackerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
86 {
87   constBase = 0;
88   inputMultiplier = ((double)fOutputPercentage)/100.0;
89 }
90
91 // Spawn function, return new instance of this class
92 AliHLTComponent* AliHLTTRDTrackerComponent::Spawn()
93 {
94   return new AliHLTTRDTrackerComponent;
95 };
96
97 int AliHLTTRDTrackerComponent::DoInit( int argc, const char** argv )
98 {
99   // perform initialization. We check whether our relative output size is specified in the arguments.
100   fOutputPercentage = 100;
101   int i = 0;
102   char* cpErr;
103   while ( i < argc )
104     {
105       Logging( kHLTLogDebug, "HLT::TRDTracker::DoInit", "Arguments", "argv[%d] == %s", i, argv[i] );
106       if ( !strcmp( argv[i], "output_percentage" ) )
107         {
108           if ( i+1>=argc )
109             {
110               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing output_percentage parameter");
111               return ENOTSUP;
112             }
113           Logging( kHLTLogDebug, "HLT::TRDTracker::DoInit", "Arguments", "argv[%d+1] == %s", i, argv[i+1] );
114           fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
115           if ( *cpErr )
116             {
117               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Wrong Argument", "Cannot convert output_percentage parameter '%s'", argv[i+1] );
118               return EINVAL;
119             }
120           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
121           i += 2;
122           continue;
123         }
124
125       if ( strcmp( argv[i], "-cdb" ) == 0)
126         {
127           if ( i+1 >= argc )
128             {
129               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing -cdb argument");
130               return ENOTSUP;         
131             }
132           fStrorageDBpath = argv[i+1];
133           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "DB storage set", "DB storage is %s", fStrorageDBpath.c_str() );       
134           fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
135           i += 2;
136           continue;
137         }      
138
139       if ( strcmp( argv[i], "-geometry" ) == 0)
140         {
141           if ( i+1 >= argc )
142             {
143               Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Missing Argument", "Missing -geometry argument");
144               return ENOTSUP;         
145             }
146           fGeometryFileName = argv[i+1];
147           Logging( kHLTLogInfo, "HLT::TRDTracker::DoInit", "GeomFile storage set", "GeomFile storage is %s", 
148                    fGeometryFileName.c_str() );   
149           i += 2;
150           continue;
151         }      
152
153       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
154       return EINVAL;
155     }
156
157   fClusterizer = new AliTRDclusterizerV1HLT("TRCclusterizer", "TRCclusterizer");
158
159   //init alifield map - temporarly fixed - should come from a DB
160   fField = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
161   if (fField)
162     AliTracker::SetFieldMap(fField,1);
163   else
164     Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Field", "Unable to init the field");
165     
166   fGeometryFile = TFile::Open(fGeometryFileName.c_str());
167   if (fGeometryFile)
168     {
169       fGeoManager = (TGeoManager *)fGeometryFile->Get("Geometry");
170       fTracker = new AliTRDtrackerHLT(fGeometryFile);
171     }
172   else
173     {
174       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "fGeometryFile", "Unable to open file. FATAL!");
175       return -1;
176     }
177
178   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
179   if (calibra == 0)
180     {
181       Logging(kHLTLogError, "HLT::TRDTracker::DoInit", "Calibration Histos", "::Instance failed");
182       return -1;      
183     }
184   else
185     {
186       calibra->SetMITracking(kTRUE);
187       calibra->Init2Dhistos();
188     }
189
190   return 0;
191 }
192
193 int AliHLTTRDTrackerComponent::DoDeinit()
194 {
195   delete fClusterizer;
196   fClusterizer = 0;
197
198   delete fField;
199
200   delete fTracker;
201   fTracker = 0;
202   
203   if (fGeometryFile)
204     {
205       fGeometryFile->Close();
206       delete fGeometryFile;
207     }
208
209   AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
210   if (calibra)
211     {
212       calibra->Write2d();
213       calibra->Destroy();
214     }
215
216   return 0;
217 }
218
219 int AliHLTTRDTrackerComponent::DoEvent( const AliHLTComponentEventData & evtData,
220                                         AliHLTComponentTriggerData & trigData )
221 {
222   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "Output percentage set", "Output percentage set to %lu %%", fOutputPercentage );
223   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "BLOCKS", "NofBlocks %lu", evtData.fBlockCnt );
224
225   AliHLTUInt32_t fDblock_Specification = 0;
226
227   AliHLTComponentBlockData *dblock = (AliHLTComponentBlockData *)GetFirstInputBlock( AliHLTTRDDefinitions::gkClusterDataType );
228   if (dblock != 0)
229     {
230       fDblock_Specification = dblock->fSpecification;
231     }
232   else
233     {
234       Logging( kHLTLogWarning, "HLT::TRDTracker::DoEvent", "DATAIN", "First Input Block not found! 0x%x", dblock);
235       return -1;
236     }
237
238   int ibForce = 0;
239   TObject *tobjin = (TObject *)GetFirstInputObject( AliHLTTRDDefinitions::gkClusterDataType, "TObjArray", ibForce);
240   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "1stBLOCK", "Pointer = 0x%x", tobjin);
241
242 //   int iTotalClusterCounter = 0;
243   while (tobjin != 0)
244     {
245       TObjArray *clusters = (TObjArray *)tobjin;
246       if (clusters != 0)
247         {
248           //put back to the clusterizers tree
249           Int_t iSuggestedDet = -1; // take the det number from the first cluster
250           Bool_t kInsert = kFALSE;
251           if (clusters->GetEntries() > 0)
252             {
253               AliTRDcluster *cl = (AliTRDcluster*)clusters->At(0);
254               iSuggestedDet = cl->GetDetector();
255               kInsert = fClusterizer->InsertClusters(clusters, iSuggestedDet);      
256               Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "TOTREE", "Result = %d", kInsert);        
257             }
258             
259           Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "CLUSTERS", "Pointer = 0x%x", clusters);      
260 //        //dump the clusters to the log files
261 //        for (Int_t ic = 0; ic < clusters->GetEntries(); ic++)
262 //          {
263 //            AliTRDcluster *cl = (AliTRDcluster*)clusters->At(ic);
264 //            Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "CLUSTER", "%d : %d : Q = %1.2f", 
265 //                     iTotalClusterCounter, ic, cl->GetQ());         
266 //            iTotalClusterCounter++;
267 //          }
268           
269           //Pass the data further...
270           //PushBack(clusters, AliHLTTRDDefinitions::gkClusterDataType, fDblock_Specification);
271         }
272       else
273         {
274           Logging( kHLTLogError, "HLT::TRDTracker::DoEvent", "CLUSTERS", "Pointer = 0x%x", clusters);     
275         }
276
277       tobjin = (TObject *)GetNextInputObject( ibForce );
278       Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "nextBLOCK", "Pointer = 0x%x", tobjin);
279
280     }
281
282   Int_t iNclusters = fClusterizer->GetNclusters();
283   Logging( kHLTLogInfo, "HLT::TRDTracker::DoEvent", "COUNT", "N of Clusters = %d", iNclusters);
284   fTracker->LoadClusters(fClusterizer->GetClusterTree());
285
286   AliTRDReconstructor::SetSeedingOn(kTRUE);
287
288   AliESD *esd = new AliESD();
289
290   //fTracker->MakeSeedsMI(3, 5, esd);
291   fTracker->PropagateBack(esd);
292
293   //here transport the esd tracks further
294   //no receiver defined yet(!)
295   delete esd;
296   return 0;
297 }