]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/offline/AliHLTTPCOfflineClustererComponent.cxx
further work on the TPC offline wrappers, test macro added (Jacek & Matthias)
[u/mrichter/AliRoot.git] / HLT / TPCLib / offline / AliHLTTPCOfflineClustererComponent.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                                                                        *
9 //* Permission to use, copy, modify and distribute this software and its   *
10 //* documentation strictly for non-commercial purposes is hereby granted   *
11 //* without fee, provided that the above copyright notice appears in all   *
12 //* copies and that both the copyright notice and this permission notice   *
13 //* appear in the supporting documentation. The authors make no claims     *
14 //* about the suitability of this software for any purpose. It is          *
15 //* provided "as is" without express or implied warranty.                  *
16 //**************************************************************************
17
18 /** @file   AliHLTTPCOfflineClustererComponent.cxx
19     @author Matthias Richter
20     @date   
21     @brief  Wrapper component to the TPC offline cluster finder
22 */
23
24 #include "AliHLTTPCOfflineClustererComponent.h"
25 #include "AliHLTTPCDefinitions.h"
26 #include "AliCDBManager.h"
27 #include "AliGeomManager.h"
28 #include "AliTPCRecoParam.h"
29 #include "AliTPCParam.h"
30 #include "AliTPCParamSR.h"
31 #include "AliRawReaderMemory.h"
32 #include "AliTPCclustererMI.h"
33 #include "AliTPCClustersRow.h"
34 #include "AliMagFMaps.h"
35 #include "AliTracker.h"
36 #include "AliDAQ.h"
37 #include "TString.h"
38 #include "TObjArray.h"
39 #include "TObjString.h"
40 #include "TTree.h"
41
42 /** ROOT macro for the implementation of ROOT specific class methods */
43 ClassImp(AliHLTTPCOfflineClustererComponent)
44
45 AliHLTTPCOfflineClustererComponent::AliHLTTPCOfflineClustererComponent() 
46 : AliHLTProcessor(),
47 fOutputPercentage(100),
48 fGeometryFileName(""),
49 fTPCRecoParam(0),
50 fTPCGeomParam(0),
51 fRawReader(0),
52 fClusterer(0),
53 fMagField(0)
54 {
55   // Default constructor
56   fGeometryFileName = getenv("ALICE_ROOT");
57   fGeometryFileName += "/HLT/TPCLib/offline/geometry.root";
58 }
59
60 AliHLTTPCOfflineClustererComponent::~AliHLTTPCOfflineClustererComponent()
61 {
62   // Destructor
63 }
64
65 const char* AliHLTTPCOfflineClustererComponent::GetComponentID()
66 {
67   // Return component ID
68   return "TPCOfflineClusterer";
69 }
70
71 void AliHLTTPCOfflineClustererComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
72 {
73   // Get the list of input data types
74   list.push_back(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC);
75 }
76
77 AliHLTComponentDataType AliHLTTPCOfflineClustererComponent::GetOutputDataType()
78 {
79   // Return output data type
80   // Tree or TObjArray of clusters
81   return  kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC/*AliHLTTPCDefinitions::fgkOfflineClustersDataType*/;
82
83 }
84
85 void AliHLTTPCOfflineClustererComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
86 {
87   // Get output data size
88   constBase = 0;
89   inputMultiplier = ((double)fOutputPercentage)/100.0;
90 }
91
92 AliHLTComponent* AliHLTTPCOfflineClustererComponent::Spawn()
93 {
94   // Return a new instance of the class 
95   return new AliHLTTPCOfflineClustererComponent;
96 }
97
98 int AliHLTTPCOfflineClustererComponent::DoInit( int argc, const char** argv )
99 {
100   // Perfom initialisation. Checks arguments (argc,argv)  
101   // and make initialisation. Returns 0 if success.  
102   //
103   int iResult=0;
104   HLTInfo("parsing %d arguments", argc);
105
106   TString argument="";
107   TString configuration=""; 
108   int bMissingParam=0;
109
110   // loop over input parameters
111   for (int i=0; i<argc && iResult>=0; i++) {
112     argument=argv[i];
113     if (argument.IsNull()) continue;
114
115     if (argument.CompareTo("-geometry")==0) {
116       if ((bMissingParam=(++i>=argc))) break;
117
118       HLTInfo("got \'-geometry\' argument: %s", argv[i]);
119       fGeometryFileName = argv[i];
120       HLTInfo("Geometry file is: %s", fGeometryFileName.c_str());
121
122       // the remaining arguments are treated as configuration
123     } else {
124       if (!configuration.IsNull()) configuration+=" ";
125       configuration+=argument;
126     }
127   } // end loop
128
129   if (bMissingParam) {
130     HLTError("missing parameter for argument %s", argument.Data());
131     iResult=-EINVAL;
132   }
133
134   if (iResult>=0 && !configuration.IsNull()) {
135     iResult=Configure(configuration.Data());
136   } else {
137     iResult=Reconfigure(NULL, NULL);
138   }
139
140   //
141   // initialisation
142   //
143
144   /*
145    
146   // Load geometry
147   HLTInfo("Geometry file %s",fGeometryFileName.c_str());
148   AliGeomManager::LoadGeometry(fGeometryFileName.c_str());
149   if((AliGeomManager::GetGeometry()) == 0) {
150     HLTError("Cannot load geometry from file %s",fGeometryFileName.c_str());
151     iResult=-EINVAL;
152   }
153  
154   // Magnetic field
155   fMagField = new AliMagFMaps("Maps","Maps", 2, 1.0, 10., 2);
156   AliTracker::SetFieldMap(fMagField,kFALSE);
157
158   */
159
160   // Raw Reader
161   fRawReader = new AliRawReaderMemory;
162
163   // TPC reconstruction parameters
164   fTPCRecoParam = AliTPCRecoParam::GetLowFluxParam();
165   if (fTPCRecoParam) {
166     fTPCRecoParam->SetClusterSharing(kTRUE);
167   }
168
169   // TPC geometry parameters
170   fTPCGeomParam = new AliTPCParamSR;
171   if (fTPCGeomParam) {
172     fTPCGeomParam->ReadGeoMatrices();
173   }
174
175   // Init clusterer
176   fClusterer = new AliTPCclustererMI(fTPCGeomParam,fTPCRecoParam);
177
178   if (!fRawReader || !fClusterer || !fTPCRecoParam || !fTPCGeomParam) {
179     HLTError("failed creating internal objects");
180     iResult=-ENOMEM;
181   }
182   return iResult;
183 }
184
185 int AliHLTTPCOfflineClustererComponent::DoDeinit()
186 {
187   // Deinitialisation of the component
188   if (fTPCRecoParam) delete fTPCRecoParam; fTPCRecoParam=0;
189   if (fTPCGeomParam) delete fTPCGeomParam; fTPCGeomParam=0;
190   if (fRawReader) delete fRawReader; fRawReader=0;
191   if (fClusterer) delete fClusterer; fClusterer=0;
192   if (fMagField) delete fMagField; fMagField=0;
193
194   return 0;
195 }
196
197 int AliHLTTPCOfflineClustererComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
198 {
199   // see header file for class documentation
200   HLTInfo("processing data");
201
202   int iResult=0;
203   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC);
204        pBlock!=NULL && iResult>=0;
205        pBlock=GetNextInputBlock()) {
206     int slice=AliHLTTPCDefinitions::GetMinSliceNr(pBlock->fSpecification);
207     int patch=AliHLTTPCDefinitions::GetMinPatchNr(pBlock->fSpecification);
208
209     if (slice!=AliHLTTPCDefinitions::GetMaxSliceNr(pBlock->fSpecification) ||
210             patch!=AliHLTTPCDefinitions::GetMaxPatchNr(pBlock->fSpecification)) {
211           HLTError("ambiguous readout partition (specification 0x%08x), skipping input block", pBlock->fSpecification);
212           break;
213     }
214     if (slice<0 || slice>35 || patch<0 || patch>5) {
215       HLTError("invalid readout partition %d/%d (specification 0x%08x, skipping input block", slice, patch,  pBlock->fSpecification);
216       break;
217     }
218
219     if (fRawReader && fClusterer) {
220       // setup raw reader and cluster finder
221       fRawReader->SetMemory( reinterpret_cast<UChar_t*>( pBlock->fPtr ), pBlock->fSize );
222       int ddlId=AliDAQ::DdlIDOffset("TPC");
223       if (patch<2) {
224         ddlId+=2*slice+patch;
225       } else {
226         ddlId+=72;
227         ddlId+=4*slice+patch-2;   
228       }
229       fRawReader->SetEquipmentID(ddlId);
230
231       // run the cluster finder
232       fClusterer->Digits2Clusters(fRawReader);
233
234       AliTPCClustersRow *clrow = 0x0;
235       Int_t nbClusters = 0;
236       Int_t lower   = fClusterer->GetOutputArray()->LowerBound();
237       Int_t entries = fClusterer->GetOutputArray()->GetEntriesFast();
238
239       for (Int_t i=lower; i<entries; i++) {
240         clrow = (AliTPCClustersRow*) fClusterer->GetOutputArray()->At(i);
241         if(!clrow) continue;
242         if(!clrow->GetArray()) continue;
243
244         nbClusters += clrow->GetArray()->GetEntries() ;
245       }
246
247       // insert TObjArray of clusters into output stream
248       PushBack(fClusterer->GetOutputArray(), kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC/*AliHLTTPCDefinitions::fgkOfflineClustersDataType*/, pBlock->fSpecification);
249       HLTInfo("processing done: DDL %d Number of clusters %d",ddlId, nbClusters);
250
251     } else {
252       HLTError("component not initialized");
253       iResult=-EFAULT;
254     }
255   }
256
257   return iResult;
258 }
259
260 int AliHLTTPCOfflineClustererComponent::Configure(const char* arguments)
261 {
262   // see header file for class documentation
263   int iResult=0;
264   if (!arguments) return iResult;
265
266   TString allArgs=arguments;
267   TString argument;
268   int bMissingParam=0;
269
270   TObjArray* pTokens=allArgs.Tokenize(" ");
271   if (pTokens) {
272     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
273       argument=((TObjString*)pTokens->At(i))->GetString();
274       if (argument.IsNull()) continue;
275
276       if (argument.CompareTo("-something")==0) {
277         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
278
279       } else {
280         HLTError("unknown argument %s", argument.Data());
281         iResult=-EINVAL;
282         break;
283       }
284     }
285     delete pTokens;
286   }
287   if (bMissingParam) {
288     HLTError("missing parameter for argument %s", argument.Data());
289     iResult=-EINVAL;
290   }
291   return iResult;
292 }
293
294 int AliHLTTPCOfflineClustererComponent::Reconfigure(const char* /*cdbEntry*/, const char* /*chainId*/)
295 {
296   // see header file for class documentation
297   int iResult=0;
298   // CDB stuff needs to be implemented
299   return iResult;
300 }