Line data Source code
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: *
8 : //* for The ALICE HLT 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 AliHLTTRDClusterizerComponent.cxx
20 : @author Theodor Rascanu
21 : @date
22 : @brief A TRDClusterizer processing component for the HLT.
23 : */
24 :
25 : // see header file for class documentation //
26 : // or //
27 : // refer to README to build package //
28 : // or //
29 : // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt //
30 :
31 : #include "TTree.h"
32 : #include "TFile.h"
33 : #include "TBranch.h"
34 :
35 : #include "AliHLTTRDClusterizerComponent.h"
36 : #include "AliHLTTRDDefinitions.h"
37 : #include "AliHLTTRDClusterizer.h"
38 : #include "AliHLTTRDUtils.h"
39 :
40 : #include "AliGeomManager.h"
41 : #include "AliTRDReconstructor.h"
42 : #include "AliCDBManager.h"
43 : #include "AliCDBStorage.h"
44 : #include "AliCDBEntry.h"
45 : #include "AliTRDrecoParam.h"
46 : #include "AliTRDrawStream.h"
47 : #include "AliTRDcluster.h"
48 :
49 : #include "AliRawReaderMemory.h"
50 :
51 : #ifdef HAVE_VALGRIND_CALLGRIND_H
52 : #include <valgrind/callgrind.h>
53 : #else
54 : #define CALLGRIND_START_INSTRUMENTATION (void)0
55 : #define CALLGRIND_STOP_INSTRUMENTATION (void)0
56 : #endif
57 :
58 : #include <cstdlib>
59 : #include <cerrno>
60 : #include <string>
61 :
62 : using namespace std;
63 :
64 6 : ClassImp(AliHLTTRDClusterizerComponent)
65 :
66 : AliHLTTRDClusterizerComponent::AliHLTTRDClusterizerComponent()
67 6 : : AliHLTProcessor(),
68 6 : fOutputPercentage(100),
69 6 : fOutputConst(0),
70 6 : fClusterizer(NULL),
71 6 : fRecoParam(NULL),
72 6 : fMemReader(NULL),
73 6 : fReconstructor(NULL),
74 6 : fRecoParamType(-1),
75 6 : fRecoDataType(-1),
76 6 : fRawDataVersion(2),
77 6 : fyPosMethod(1),
78 6 : fgeometryFileName(""),
79 6 : fProcessTracklets(kFALSE),
80 6 : fHLTstreamer(kTRUE),
81 6 : fTC(kFALSE),
82 6 : fHLTflag(kTRUE),
83 6 : fHighLevelOutput(kFALSE),
84 6 : fEmulateHLTClusters(kFALSE)
85 24 : {
86 : // Default constructor
87 :
88 9 : }
89 :
90 : AliHLTTRDClusterizerComponent::~AliHLTTRDClusterizerComponent()
91 24 : {
92 : // Destructor
93 : // Work is Done in DoDeInit()
94 12 : }
95 :
96 :
97 : const char* AliHLTTRDClusterizerComponent::GetComponentID()
98 : {
99 : // Return the component ID const char *
100 162 : return "TRDClusterizer"; // The ID of this component
101 : }
102 :
103 : void AliHLTTRDClusterizerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
104 : {
105 : // Get the list of input data
106 0 : list.clear(); // We do not have any requirements for our input data type(s).
107 0 : list.push_back(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
108 0 : }
109 :
110 : AliHLTComponentDataType AliHLTTRDClusterizerComponent::GetOutputDataType()
111 : {
112 : // Get the output data type
113 0 : return kAliHLTMultipleDataType;
114 : }
115 :
116 : int AliHLTTRDClusterizerComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
117 : {
118 : // Get the output data type
119 0 : tgtList.clear();
120 0 : tgtList.push_back(AliHLTTRDDefinitions::fgkClusterDataType);
121 0 : tgtList.push_back(AliHLTTRDDefinitions::fgkMCMtrackletDataType);
122 0 : return tgtList.size();
123 : }
124 :
125 :
126 : void AliHLTTRDClusterizerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
127 : {
128 : // Get the output data size
129 0 : constBase = fOutputConst;
130 0 : inputMultiplier = ((double)fOutputPercentage)*4/100.0;
131 0 : }
132 :
133 : AliHLTComponent* AliHLTTRDClusterizerComponent::Spawn()
134 : {
135 : // Spawn function, return new instance of this class
136 0 : return new AliHLTTRDClusterizerComponent;
137 0 : };
138 :
139 : int AliHLTTRDClusterizerComponent::DoInit( int argc, const char** argv )
140 : {
141 : // perform initialization. We check whether our relative output size is specified in the arguments.
142 : int iResult=0;
143 :
144 0 : fReconstructor = new AliTRDReconstructor();
145 : HLTDebug("TRDReconstructor at 0x%x", fReconstructor);
146 :
147 0 : TString configuration="";
148 0 : TString argument="";
149 0 : for (int i=0; i<argc && iResult>=0; i++) {
150 0 : argument=argv[i];
151 0 : if (!configuration.IsNull()) configuration+=" ";
152 0 : configuration+=argument;
153 : }
154 :
155 0 : if (!configuration.IsNull()) {
156 0 : iResult=Configure(configuration.Data());
157 0 : } else {
158 0 : iResult=Reconfigure(NULL, NULL);
159 : }
160 :
161 0 : if(!fClusterizer){
162 0 : HLTFatal("Clusterizer was not initialized!");
163 0 : return -1;
164 : }
165 :
166 0 : if(iResult<0) return iResult;
167 :
168 0 : fMemReader = new AliRawReaderMemory;
169 0 : fClusterizer->SetReconstructor(fReconstructor);
170 0 : fClusterizer->SetUseLabels(kFALSE);
171 :
172 0 : if(fReconstructor->IsProcessingTracklets())
173 0 : fOutputConst = fClusterizer->GetTrMemBlockSize();
174 :
175 0 : return iResult;
176 0 : }
177 :
178 : int AliHLTTRDClusterizerComponent::DoDeinit()
179 : {
180 : // Deinitialization of the component
181 0 : delete fMemReader;
182 0 : fMemReader = 0;
183 0 : delete fClusterizer;
184 0 : fClusterizer = 0;
185 :
186 : //fReconstructor->SetClusters(0x0);
187 0 : delete fReconstructor;
188 0 : fReconstructor = 0x0;
189 :
190 0 : if (fRecoParam)
191 : {
192 : HLTDebug("Deleting fRecoParam");
193 0 : delete fRecoParam;
194 0 : fRecoParam = 0;
195 0 : }
196 :
197 0 : return 0;
198 : }
199 :
200 : int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponentEventData& evtData,
201 : const AliHLTComponentBlockData* blocks,
202 : AliHLTComponent_TriggerData& /*trigData*/,
203 : AliHLTUInt8_t* outputPtr,
204 : AliHLTUInt32_t& size,
205 : vector<AliHLTComponent_BlockData>& outputBlocks )
206 : {
207 : // Process an event
208 :
209 : #ifdef HAVE_VALGRIND_CALLGRIND_H
210 : if (evtData.fEventID == 10)
211 : CALLGRIND_START_INSTRUMENTATION;
212 :
213 : if(GetFirstInputBlock(kAliHLTDataTypeEOR))
214 : CALLGRIND_STOP_INSTRUMENTATION;
215 : #endif
216 :
217 0 : if(!IsDataEvent())return 0;
218 :
219 : HLTDebug( "NofBlocks %i", evtData.fBlockCnt );
220 : // Process an event
221 : AliHLTUInt32_t totalSize = 0, offset = 0;
222 :
223 : //implement a usage of the following
224 : // AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
225 : // AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
226 : // void *triggerData = trigData.fData;
227 : //HLTDebug( "Trigger data received. Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
228 :
229 : // Loop over all input blocks in the event
230 0 : AliHLTComponentDataType expectedDataType = (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
231 0 : for ( UInt_t iBlock = 0; iBlock < evtData.fBlockCnt; iBlock++ )
232 : {
233 0 : const AliHLTComponentBlockData &block = blocks[iBlock];
234 : // lets not use the internal TRD data types here : AliHLTTRDDefinitions::fgkDDLRawDataType
235 : // which is depreciated - we use HLT global defs instead
236 : // if ( block.fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD) )
237 0 : AliHLTComponentDataType inputDataType = block.fDataType;
238 0 : if ( inputDataType != expectedDataType)
239 : {
240 : HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) Wrong received datatype: %s - required datatype: %s; Skipping",
241 : iBlock, evtData.fBlockCnt,
242 : evtData.fEventID, evtData.fEventID,
243 : DataType2Text(inputDataType).c_str(),
244 : DataType2Text(expectedDataType).c_str());
245 0 : continue;
246 : }
247 : else
248 : {
249 : HLTDebug("We get the right data type: Block # %i/%i; Event 0x%08LX (%Lu) Received datatype: %s; Block Size: %i",
250 : iBlock, evtData.fBlockCnt,
251 : evtData.fEventID, evtData.fEventID,
252 : DataType2Text(inputDataType).c_str(),
253 : block.fSize);
254 : }
255 :
256 : #ifndef NDEBUG
257 0 : unsigned long constBase;
258 0 : double inputMultiplier;
259 0 : GetOutputDataSize(constBase,inputMultiplier);
260 0 : if(size<(constBase+block.fSize*inputMultiplier)){
261 0 : HLTWarning("Memory Block given might be too small: %i < %i; Event %Lu", size, constBase+block.fSize*inputMultiplier, evtData.fEventID);
262 : }
263 : #endif
264 :
265 : // fMemReader->Reset();
266 0 : fMemReader->SetMemory((UChar_t*) block.fPtr, block.fSize);
267 :
268 0 : AliHLTUInt32_t spec = block.fSpecification;
269 :
270 0 : Int_t id = AliHLTTRDUtils::GetSM(spec) + 1024;
271 :
272 0 : fMemReader->SetEquipmentID(id);
273 :
274 0 : fClusterizer->SetMemBlock(outputPtr+offset);
275 0 : Bool_t bclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
276 0 : if(bclustered)
277 : {
278 : HLTDebug("Clustered successfully");
279 : }
280 : else
281 : {
282 0 : HLTError("Clustering ERROR");
283 0 : return -1;
284 : }
285 :
286 : AliHLTUInt32_t addedSize;
287 0 : if(fReconstructor->IsProcessingTracklets()){
288 0 : addedSize = fClusterizer->GetAddedTrSize();
289 0 : totalSize += fClusterizer->GetTrMemBlockSize(); //if IsProcessingTracklets() is enabled we always reserve a data block of size GetTrMemBlockSize() for the tracklets
290 0 : if (addedSize > 0){
291 : // Using low-level interface
292 : // with interface classes
293 0 : if ( totalSize > size )
294 : {
295 0 : HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
296 : totalSize, size );
297 0 : return EMSGSIZE;
298 : }
299 :
300 : // Fill block
301 0 : AliHLTComponentBlockData bd;
302 0 : FillBlockData( bd );
303 0 : bd.fOffset = offset;
304 0 : bd.fSize = addedSize;
305 0 : bd.fSpecification = block.fSpecification;
306 0 : bd.fDataType = AliHLTTRDDefinitions::fgkMCMtrackletDataType;
307 0 : outputBlocks.push_back( bd );
308 : HLTDebug( "BD ptr 0x%x, offset %i, size %i, dataType %s, spec 0x%x ", bd.fPtr, bd.fOffset, bd.fSize, DataType2Text(bd.fDataType).c_str(), spec);
309 0 : }
310 : offset = totalSize;
311 0 : }
312 :
313 0 : addedSize = fClusterizer->GetAddedClSize();
314 0 : if (addedSize > 0){
315 :
316 0 : Int_t* nTimeBins = (Int_t*)(outputPtr+offset+fClusterizer->GetAddedClSize());
317 0 : *nTimeBins = fClusterizer->GetNTimeBins();
318 0 : addedSize += sizeof(*nTimeBins);
319 :
320 0 : totalSize += addedSize;
321 0 : if ( totalSize > size )
322 : {
323 0 : HLTError("Too much data; Data written over allowed buffer. Amount written: %lu, allowed amount: %lu.",
324 : totalSize, size );
325 0 : return EMSGSIZE;
326 : }
327 :
328 : // Fill block
329 0 : AliHLTComponentBlockData bd;
330 0 : FillBlockData( bd );
331 0 : bd.fOffset = offset;
332 0 : bd.fSize = addedSize;
333 0 : bd.fSpecification = block.fSpecification;
334 0 : bd.fDataType = AliHLTTRDDefinitions::fgkClusterDataType;
335 0 : outputBlocks.push_back( bd );
336 : HLTDebug( "BD ptr 0x%x, offset %i, size %i, dataType %s, spec 0x%x ", bd.fPtr, bd.fOffset, bd.fSize, DataType2Text(bd.fDataType).c_str(), spec);
337 : offset = totalSize;
338 0 : }
339 : else{
340 : HLTDebug("Array of clusters is empty!");
341 : }
342 0 : }
343 : //fReconstructor->SetClusters(0x0);
344 :
345 0 : size = totalSize;
346 : HLTDebug("Event is done. size written to the output is %i", size);
347 0 : return 0;
348 0 : }
349 :
350 : void AliHLTTRDClusterizerComponent::PrintObject(
351 : #ifdef __DEBUG
352 : TClonesArray* inClustersArray
353 : #else
354 : TClonesArray*
355 : #endif
356 : )
357 : {
358 : #ifdef __DEBUG
359 : AliTRDcluster* cluster=0x0;
360 :
361 : for (Int_t i=0; i < inClustersArray->GetEntriesFast(); i++){
362 : cluster = dynamic_cast<AliTRDcluster*>(inClustersArray->At(i));
363 : HLTDebug("cluster[%i]",i);
364 : HLTDebug(" PadCol = %i; PadRow = %i; PadTime = %i", cluster->GetPadCol(), cluster->GetPadRow(), cluster->GetPadTime());
365 : HLTDebug(" Detector = %i, Amplitude = %f, Center = %f", cluster->GetDetector(), cluster->GetQ(), cluster->GetCenter());
366 : HLTDebug(" LocalTimeBin = %i; NPads = %i; maskedPosition: %s, status: %s", cluster->GetLocalTimeBin(), cluster->GetNPads(),cluster->GetPadMaskedPosition(),cluster->GetPadMaskedPosition());
367 : }
368 : #endif
369 0 : }
370 :
371 : int AliHLTTRDClusterizerComponent::Configure(const char* arguments){
372 : int iResult=0;
373 0 : if (!arguments) return iResult;
374 :
375 0 : TString allArgs=arguments;
376 0 : TString argument;
377 : int bMissingParam=0;
378 :
379 0 : TObjArray* pTokens=allArgs.Tokenize(" ");
380 0 : if (pTokens) {
381 0 : for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
382 0 : argument=((TObjString*)pTokens->At(i))->GetString();
383 0 : if (argument.IsNull()) continue;
384 :
385 0 : if (argument.CompareTo("output_percentage")==0) {
386 0 : if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
387 0 : HLTInfo("Setting output percentage to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
388 0 : fOutputPercentage=((TObjString*)pTokens->At(i))->GetString().Atoi();
389 0 : continue;
390 : }
391 0 : else if (argument.CompareTo("-geometry")==0) {
392 0 : if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
393 0 : HLTInfo("Setting geometry to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
394 0 : fgeometryFileName=((TObjString*)pTokens->At(i))->GetString();
395 0 : continue;
396 : }
397 0 : else if (argument.CompareTo("-lowflux")==0) {
398 0 : fRecoParamType = 0;
399 0 : HLTInfo("Low flux reconstruction selected");
400 : continue;
401 : }
402 0 : else if (argument.CompareTo("-highflux")==0) {
403 0 : fRecoParamType = 1;
404 0 : HLTInfo("High flux reconstruction selected");
405 : continue;
406 : }
407 0 : else if (argument.CompareTo("-cosmics")==0) {
408 0 : fRecoParamType = 2;
409 0 : HLTInfo("Cosmics reconstruction selected");
410 : continue;
411 : }
412 0 : else if (argument.CompareTo("-simulation")==0) {
413 0 : fRecoDataType = 0;
414 0 : HLTInfo("Awaiting simulated data");
415 : continue;
416 : }
417 0 : else if (argument.CompareTo("-experiment")==0) {
418 0 : fRecoDataType = 1;
419 0 : HLTInfo("Awaiting real data");
420 : continue;
421 : }
422 0 : else if (argument.CompareTo("-processTracklets")==0) {
423 0 : fProcessTracklets = kTRUE;
424 0 : HLTInfo("Writing L1 tracklets to output");
425 : continue;
426 : }
427 0 : else if (argument.CompareTo("-noZS")==0) {
428 0 : fOutputPercentage = 10;
429 0 : HLTInfo("Awaiting non zero surpressed data");
430 : continue;
431 : }
432 0 : else if (argument.CompareTo("-HLTflag")==0) {
433 0 : if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
434 0 : TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
435 0 : if (toCompareTo.CompareTo("yes")==0){
436 0 : HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
437 0 : fHLTflag=kTRUE;
438 0 : }
439 0 : else if (toCompareTo.CompareTo("no")==0){
440 0 : HLTInfo("Setting HLTflag to: %s", toCompareTo.Data());
441 0 : fHLTflag=kFALSE;
442 : }
443 : else {
444 0 : HLTError("unknown argument for HLTflag: %s", toCompareTo.Data());
445 : iResult=-EINVAL;
446 0 : break;
447 : }
448 0 : continue;
449 0 : }
450 0 : else if (argument.CompareTo("-faststreamer")==0) {
451 0 : fHLTstreamer = kTRUE;
452 0 : HLTInfo("Useing fast raw streamer");
453 : continue;
454 : }
455 0 : else if (argument.CompareTo("-nofaststreamer")==0) {
456 0 : fHLTstreamer = kFALSE;
457 0 : HLTInfo("Don't use fast raw streamer");
458 : continue;
459 : }
460 0 : else if (argument.CompareTo("-tailcancellation")==0) {
461 0 : fTC = kTRUE;
462 0 : HLTInfo("Useing tailcancellation");
463 : continue;
464 : }
465 0 : else if (argument.CompareTo("-rawver")==0) {
466 0 : if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
467 0 : HLTInfo("Raw data version is: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
468 0 : fRawDataVersion=((TObjString*)pTokens->At(i))->GetString().Atoi();
469 0 : continue;
470 : }
471 0 : else if (argument.CompareTo("-highLevelOutput")==0) {
472 0 : if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
473 0 : TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
474 0 : if (toCompareTo.CompareTo("yes")==0){
475 0 : HLTWarning("Setting highLevelOutput to: %s", toCompareTo.Data());
476 0 : fHighLevelOutput=kTRUE;
477 0 : }
478 0 : else if (toCompareTo.CompareTo("no")==0){
479 0 : HLTInfo("Setting highLevelOutput to: %s", toCompareTo.Data());
480 0 : fHighLevelOutput=kFALSE;
481 : }
482 : else {
483 0 : HLTError("unknown argument for highLevelOutput: %s", toCompareTo.Data());
484 : iResult=-EINVAL;
485 0 : break;
486 : }
487 0 : continue;
488 0 : }
489 0 : else if (argument.CompareTo("-emulateHLToutput")==0) {
490 0 : if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
491 0 : TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
492 0 : if (toCompareTo.CompareTo("yes")==0){
493 0 : HLTWarning("Setting emulateHLToutput to: %s", toCompareTo.Data());
494 0 : fEmulateHLTClusters=kTRUE;
495 0 : }
496 0 : else if (toCompareTo.CompareTo("no")==0){
497 0 : HLTInfo("Setting emulateHLToutput to: %s", toCompareTo.Data());
498 0 : fEmulateHLTClusters=kFALSE;
499 : }
500 : else {
501 0 : HLTError("unknown argument for emulateHLToutput: %s", toCompareTo.Data());
502 : iResult=-EINVAL;
503 0 : break;
504 : }
505 0 : continue;
506 0 : }
507 0 : else if (argument.CompareTo("-yPosMethod")==0) {
508 0 : if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
509 0 : TString toCompareTo=((TObjString*)pTokens->At(i))->GetString();
510 0 : if (toCompareTo.CompareTo("COG")==0){
511 0 : HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
512 0 : fyPosMethod=0;
513 0 : }
514 0 : else if (toCompareTo.CompareTo("LUT")==0){
515 0 : HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
516 0 : fyPosMethod=1;
517 0 : }
518 0 : else if (toCompareTo.CompareTo("Gauss")==0){
519 0 : HLTInfo("Setting yPosMethod method to: %s", toCompareTo.Data());
520 0 : fyPosMethod=2;
521 : }
522 : else {
523 0 : HLTError("unknown argument for yPosMethod: %s", toCompareTo.Data());
524 : iResult=-EINVAL;
525 0 : break;
526 : }
527 0 : continue;
528 0 : }
529 :
530 : else {
531 0 : HLTError("unknown argument: %s", argument.Data());
532 : iResult=-EINVAL;
533 0 : break;
534 : }
535 : }
536 0 : delete pTokens;
537 : }
538 0 : if (bMissingParam) {
539 0 : HLTError("missing parameter for argument %s", argument.Data());
540 : iResult=-EINVAL;
541 0 : }
542 0 : if(iResult>=0){
543 0 : iResult=SetParams();
544 0 : }
545 : return iResult;
546 0 : }
547 :
548 : int AliHLTTRDClusterizerComponent::SetParams()
549 : {
550 : Int_t iResult=0;
551 0 : if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
552 0 : HLTError("DefaultStorage is not set in CDBManager");
553 0 : return -EINVAL;
554 : }
555 0 : if(AliCDBManager::Instance()->GetRun()<0){
556 0 : HLTError("Run Number is not set in CDBManager");
557 0 : return -EINVAL;
558 : }
559 0 : HLTInfo("CDB default storage: %s; RunNo: %i", (AliCDBManager::Instance()->GetDefaultStorage()->GetBaseFolder()).Data(), AliCDBManager::Instance()->GetRun());
560 :
561 0 : if(!AliGeomManager::GetGeometry()){
562 0 : if(fgeometryFileName.CompareTo("")==0 || !TFile::Open(fgeometryFileName.Data())){
563 0 : HLTInfo("Loading standard geometry file");
564 0 : AliGeomManager::LoadGeometry();
565 0 : }else{
566 0 : HLTWarning("Loading NON-standard geometry file");
567 0 : AliGeomManager::LoadGeometry(fgeometryFileName.Data());
568 : }
569 0 : if(!AliGeomManager::GetGeometry()){
570 0 : HLTError("Could not load geometry");
571 0 : return -EINVAL;
572 : }
573 0 : HLTInfo("Applying Alignment from CDB object");
574 0 : AliGeomManager::ApplyAlignObjsFromCDB("TRD");
575 0 : }
576 : else{
577 0 : HLTInfo("Geometry Already Loaded!");
578 : }
579 :
580 0 : if(fReconstructor->GetRecoParam()){
581 0 : fRecoParam = new AliTRDrecoParam(*fReconstructor->GetRecoParam());
582 0 : HLTInfo("RecoParam already set!");
583 : }else{
584 0 : if(fRecoParamType == 0){
585 : #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
586 0 : if(fHLTflag){
587 0 : HLTInfo("Low flux HLT params init.");
588 0 : fRecoParam = AliTRDrecoParam::GetLowFluxHLTParam();
589 0 : }else
590 : #endif
591 : {
592 0 : HLTInfo("Low flux params init.");
593 0 : fRecoParam = AliTRDrecoParam::GetLowFluxParam();
594 : }
595 : }
596 0 : if(fRecoParamType == 1){
597 : #ifndef HAVE_NOT_ALITRD_RECOPARAM_r41621
598 0 : if(fHLTflag){
599 0 : HLTInfo("High flux HLT params init.");
600 0 : fRecoParam = AliTRDrecoParam::GetHighFluxHLTParam();
601 0 : }else
602 : #endif
603 : {
604 0 : HLTInfo("High flux params init.");
605 0 : fRecoParam = AliTRDrecoParam::GetHighFluxParam();
606 : }
607 : }
608 0 : if(fRecoParamType == 2){
609 0 : HLTInfo("Cosmic Test params init.");
610 0 : fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
611 0 : }
612 : }
613 :
614 0 : if (!fRecoParam)
615 : {
616 0 : HLTError("No reco params initialized. Sniffing big trouble!");
617 0 : return -EINVAL;
618 : }
619 :
620 0 : if(fTC){fRecoParam->SetTailCancelation(kTRUE); HLTDebug("Enableing Tail Cancelation"); }
621 0 : else{fRecoParam->SetTailCancelation(kFALSE); HLTDebug("Disableing Tail Cancelation"); }
622 :
623 0 : switch(fyPosMethod){
624 0 : case 0: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kFALSE); break;
625 0 : case 1: fRecoParam->SetGAUS(kFALSE); fRecoParam->SetLUT(kTRUE); break;
626 0 : case 2: fRecoParam->SetGAUS(kTRUE); fRecoParam->SetLUT(kFALSE); break;
627 : }
628 :
629 0 : fRecoParam->SetStreamLevel(AliTRDrecoParam::kClusterizer, 0);
630 0 : fReconstructor->SetRecoParam(fRecoParam);
631 :
632 0 : if(!fClusterizer){
633 0 : fClusterizer = new AliHLTTRDClusterizer("TRDCclusterizer", "TRDCclusterizer");
634 : HLTDebug("TRDClusterizer at 0x%x", fClusterizer);
635 0 : }
636 :
637 0 : TString recoOptions="!cw";
638 0 : if(fHLTflag){
639 0 : recoOptions += ",hlt";
640 :
641 : // we transfer clusters that do no contain the XYZ coodrinates (AliHLTTRDCluster),
642 : // thus this coordinate transformation ist useless
643 : #ifndef HAVE_NOT_ALITRD_CLUSTERIZER_r42837
644 0 : fClusterizer->SetSkipTransform();
645 : #endif
646 : }
647 0 : if(fProcessTracklets) recoOptions += ",tp";
648 0 : else recoOptions += ",!tp";
649 :
650 : HLTDebug("Reconstructor options are: %s",recoOptions.Data());
651 0 : fReconstructor->SetOption(recoOptions.Data());
652 :
653 0 : if (fRecoDataType < 0 || fRecoDataType > 1)
654 : {
655 0 : HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
656 0 : fRecoDataType = 0;
657 0 : }
658 :
659 0 : if (fRecoDataType == 0)
660 : {
661 : HLTDebug("Data type expected is SIMULATION!");
662 : }
663 :
664 : if (fRecoDataType == 1)
665 : {
666 : HLTDebug("Data type expected is EXPERIMENT!");
667 : }
668 :
669 0 : fClusterizer->SetRawVersion(fRawDataVersion);
670 :
671 : return iResult;
672 0 : }
673 :
674 : int AliHLTTRDClusterizerComponent::Reconfigure(const char* cdbEntry, const char* chainId)
675 : {
676 : // see header file for class documentation
677 :
678 : int iResult=0;
679 : const char* path="HLT/ConfigTRD/ClusterizerComponent";
680 : const char* defaultNotify="";
681 0 : if (cdbEntry) {
682 : path=cdbEntry;
683 : defaultNotify=" (default)";
684 0 : }
685 0 : if (path) {
686 0 : HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
687 0 : AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
688 0 : if (pEntry) {
689 0 : TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
690 0 : if (pString) {
691 0 : HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
692 0 : iResult=Configure(pString->GetString().Data());
693 0 : } else {
694 0 : HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
695 : }
696 0 : } else {
697 0 : HLTError("cannot fetch object \"%s\" from CDB", path);
698 : }
699 0 : }
700 :
701 0 : return iResult;
702 0 : }
703 :
704 : void AliHLTTRDClusterizerComponent::GetOCDBObjectDescription(TMap* const targetMap){
705 : // Get a list of OCDB object description needed for the particular component
706 0 : if (!targetMap) return;
707 0 : targetMap->Add(new TObjString("HLT/ConfigTRD/ClusterizerComponent"), new TObjString("component arguments"));
708 0 : targetMap->Add(new TObjString("TRD/Calib/ChamberGainFactor"), new TObjString("gain factor of chambers"));
709 0 : targetMap->Add(new TObjString("TRD/Calib/ChamberT0"), new TObjString("T0 of chambers"));
710 0 : targetMap->Add(new TObjString("TRD/Calib/ChamberVdrift"), new TObjString("drift velocity of chambers"));
711 0 : targetMap->Add(new TObjString("TRD/Calib/DetNoise"), new TObjString("noise of chambers"));
712 0 : targetMap->Add(new TObjString("TRD/Calib/LocalGainFactor"), new TObjString("per pad gain factor"));
713 0 : targetMap->Add(new TObjString("TRD/Calib/LocalT0"), new TObjString("per pad T0"));
714 0 : targetMap->Add(new TObjString("TRD/Calib/LocalVdrift"), new TObjString("per pad drift velocity"));
715 0 : targetMap->Add(new TObjString("TRD/Calib/PadNoise"), new TObjString("per pad noise"));
716 0 : targetMap->Add(new TObjString("TRD/Calib/PadStatus"), new TObjString("pad status"));
717 0 : targetMap->Add(new TObjString("TRD/Calib/PRFWidth"), new TObjString("pad response function"));
718 0 : }
|