3 /**************************************************************************
4 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
6 * Authors: Matthias Richter <Matthias.Richter@ift.uib.no> *
7 * Timm Steinbeck <timm@kip.uni-heidelberg.de> *
8 * for The ALICE Off-line Project. *
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 **************************************************************************/
19 /** @file AliHLTTRDClusterizerComponent.cxx
20 @author Timm Steinbeck, Matthias Richter
22 @brief A TRDClusterizer processing component for the HLT. */
24 // see header file for class documentation //
26 // refer to README to build package //
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt //
38 #include "AliHLTTRDClusterizerComponent.h"
39 #include "AliHLTTRDDefinitions.h"
41 #include "AliGeomManager.h"
42 #include "AliTRDReconstructor.h"
43 #include "AliCDBManager.h"
44 #include "AliTRDclusterizerHLT.h"
45 #include "AliTRDrecoParam.h"
46 #include "AliTRDrawStreamBase.h"
48 #include "AliRawReaderMemory.h"
54 // this is a global object used for automatic component registration, do not use this
55 AliHLTTRDClusterizerComponent gAliHLTTRDClusterizerComponent;
57 ClassImp(AliHLTTRDClusterizerComponent);
59 AliHLTTRDClusterizerComponent::AliHLTTRDClusterizerComponent()
61 , fOutputPercentage(100) // By default we copy to the output exactly what we got as input
62 , fStrorageDBpath("local://$ALICE_ROOT")
67 , fReconstructor(NULL)
68 , fGeometryFileName("")
71 // Default constructor
73 fGeometryFileName = getenv("ALICE_ROOT");
74 fGeometryFileName += "/HLT/TRD/geometry.root";
77 AliHLTTRDClusterizerComponent::~AliHLTTRDClusterizerComponent()
80 // Work is Done in DoDeInit()
82 const char* AliHLTTRDClusterizerComponent::GetComponentID()
84 // Return the component ID const char *
85 return "TRDClusterizer"; // The ID of this component
88 void AliHLTTRDClusterizerComponent::GetInputDataTypes( vector<AliHLTComponent_DataType>& list)
90 // Get the list of input data
91 list.clear(); // We do not have any requirements for our input data type(s).
92 list.push_back( AliHLTTRDDefinitions::fgkDDLRawDataType );
95 AliHLTComponent_DataType AliHLTTRDClusterizerComponent::GetOutputDataType()
97 // Get the output data type
98 return AliHLTTRDDefinitions::fgkClusterDataType;
101 void AliHLTTRDClusterizerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
103 // Get the output data size
105 inputMultiplier = ((double)fOutputPercentage)/100.0;
108 AliHLTComponent* AliHLTTRDClusterizerComponent::Spawn()
110 // Spawn function, return new instance of this class
111 return new AliHLTTRDClusterizerComponent;
114 int AliHLTTRDClusterizerComponent::DoInit( int argc, const char** argv )
116 // perform initialization. We check whether our relative output size is specified in the arguments.
117 fOutputPercentage = 100;
118 Int_t iRawDataVersion = 2;
121 Bool_t bWriteClusters = kFALSE;
124 Int_t iRecoParamType = -1; // default will be the low flux
126 // the data type will become obsolete as soon as the formats are established
127 Int_t iRecoDataType = -1; // default will be simulation
131 HLTDebug("argv[%d] == %s", i, argv[i] );
132 if ( !strcmp( argv[i], "output_percentage" ) )
136 HLTError("Missing output_percentage parameter");
139 HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
140 fOutputPercentage = strtoul( argv[i+1], &cpErr, 0 );
143 HLTError("Cannot convert output_percentage parameter '%s'", argv[i+1] );
146 HLTInfo("Output percentage set to %lu %%", fOutputPercentage );
149 else if ( strcmp( argv[i], "-cdb" ) == 0)
153 HLTError("Missing -cdb argument");
156 fStrorageDBpath = argv[i+1];
157 HLTInfo("DB storage is %s", fStrorageDBpath.c_str() );
162 else if ( strcmp( argv[i], "-lowflux" ) == 0)
165 HLTDebug("Low flux reco selected.");
170 else if ( strcmp( argv[i], "-highflux" ) == 0)
173 HLTDebug("High flux reco selected.");
178 else if ( strcmp( argv[i], "-cosmics" ) == 0)
181 HLTDebug("Cosmic test reco selected.");
186 // raw data type - sim or experiment
187 else if ( strcmp( argv[i], "-simulation" ) == 0)
194 else if ( strcmp( argv[i], "-experiment" ) == 0)
201 else if ( strcmp( argv[i], "-rawver" ) == 0)
205 HLTError("Missing -rawver argument");
208 iRawDataVersion = atoi( argv[i+1] );
209 HLTInfo("Raw data version is %d", iRawDataVersion );
214 else if ( strcmp( argv[i], "-geometry" ) == 0)
218 HLTError("Missing -geometry argument");
221 fGeometryFileName = argv[i+1];
222 HLTInfo("GeomFile storage is %s", fGeometryFileName.c_str() );
226 else if ( strcmp( argv[i], "-writeClusters" ) == 0)
228 bWriteClusters = kTRUE;
232 HLTError("Unknown option '%s'", argv[i] );
238 // THE "REAL" INIT COMES HERE
240 if (iRecoParamType < 0 || iRecoParamType > 2)
242 HLTWarning("No reco param selected. Use -lowflux or -highflux flag. Defaulting to low flux.");
246 if (iRecoParamType == 0)
248 fRecoParam = AliTRDrecoParam::GetLowFluxParam();
249 HLTDebug("Low flux params init.");
252 if (iRecoParamType == 1)
254 fRecoParam = AliTRDrecoParam::GetHighFluxParam();
255 HLTDebug("High flux params init.");
258 if (iRecoParamType == 2)
260 fRecoParam = AliTRDrecoParam::GetCosmicTestParam();
261 HLTDebug("Cosmic Test params init.");
266 HLTError("No reco params initialized. Sniffing big trouble!");
270 fReconstructor = new AliTRDReconstructor();
271 fReconstructor->SetRecoParam(fRecoParam);
272 fReconstructor->SetStreamLevel(0, AliTRDReconstructor::kClusterizer); // default value
275 HLTInfo("Writing clusters. I.e. output is a TTree with clusters");
276 fReconstructor->SetOption("cw,sl_cf_0");
280 HLTInfo("Not witing clusters. I.e. output is a TClonesArray of clusters");
281 fReconstructor->SetOption("!cw,sl_cf_0");
285 // init the raw data type to be used...
286 // the switch here will become obsolete as soon as the data structures is fixed
287 // both: in sim and reality
288 if (iRecoDataType < 0 || iRecoDataType > 1)
290 HLTWarning("No data type selected. Use -simulation or -experiment flag. Defaulting to simulation.");
294 if (iRecoDataType == 0)
296 AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDsimStream);
297 HLTDebug("Data type expected is SIMULATION!");
300 if (iRecoDataType == 1)
302 AliTRDrawStreamBase::SetRawStreamVersion(AliTRDrawStreamBase::kTRDrealStream);
303 HLTDebug("Data type expected is EXPERIMENT!");
306 // the DATA BASE STUFF
307 fCDB = AliCDBManager::Instance();
310 HLTError("Could not get CDB instance", "fCDB 0x%x", fCDB);
314 fCDB->SetRun(0); // THIS HAS TO BE RETRIEVED !!!
315 fCDB->SetDefaultStorage(fStrorageDBpath.c_str());
316 HLTDebug("CDB instance; fCDB 0x%x", fCDB);
319 fGeometryFile = TFile::Open(fGeometryFileName.c_str());
322 AliGeomManager::LoadGeometry(fGeometryFileName.c_str());
326 HLTError("Unable to open file. FATAL!");
330 fMemReader = new AliRawReaderMemory;
332 fClusterizer = new AliTRDclusterizerHLT("TRDCclusterizer", "TRDCclusterizer");
333 fClusterizer->SetReconstructor(fReconstructor);
335 fClusterizer->SetRawVersion(iRawDataVersion);
336 fClusterizer->InitClusterTree();
340 int AliHLTTRDClusterizerComponent::DoDeinit()
342 // Deinitialization of the component
348 fReconstructor->SetClusters(0x0);
349 delete fReconstructor;
350 fReconstructor = 0x0;
355 fGeometryFile->Close();
356 delete fGeometryFile;
362 HLTDebug("destroy fCDB");
369 HLTDebug("Deleting fRecoParam");
375 int AliHLTTRDClusterizerComponent::DoEvent( const AliHLTComponentEventData& evtData,
376 const AliHLTComponentBlockData* blocks,
377 AliHLTComponent_TriggerData& /*trigData*/,
378 AliHLTUInt8_t* /*outputPtr*/,
379 AliHLTUInt32_t& size,
380 vector<AliHLTComponent_BlockData>& /*outputBlocks*/ )
383 HLTDebug( "NofBlocks %lu", evtData.fBlockCnt );
385 Bool_t bWriteClusters = fReconstructor->IsWritingClusters();
387 //implement a usage of the following
388 // AliHLTUInt32_t triggerDataStructSize = trigData.fStructSize;
389 // AliHLTUInt32_t triggerDataSize = trigData.fDataSize;
390 // void *triggerData = trigData.fData;
391 //HLTDebug( "Trigger data received. Struct size %d Data size %d Data location 0x%x", trigData.fStructSize, trigData.fDataSize, (UInt_t*)trigData.fData);
393 // Loop over all input blocks in the event
394 AliHLTComponentDataType expectedDataType = (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
395 for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
397 // lets not use the internal TRD data types here : AliHLTTRDDefinitions::fgkDDLRawDataType
398 // which is depreciated - we use HLT global defs instead
399 // if ( blocks[i].fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD) )
400 AliHLTComponentDataType inputDataType = blocks[i].fDataType;
401 if ( inputDataType != expectedDataType)
403 HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s; Skipping",
404 i, evtData.fBlockCnt,
405 evtData.fEventID, evtData.fEventID,
406 DataType2Text(inputDataType).c_str(),
407 DataType2Text(expectedDataType).c_str());
412 HLTDebug("We get the right data type!");
415 // fMemReader->Reset();
416 fMemReader->SetMemory((UChar_t*) blocks[i].fPtr, blocks[i].fSize);
418 AliHLTUInt32_t spec = blocks[i].fSpecification;
422 for ( Int_t ii = 0; ii < 18 ; ii++ ) {
423 if ( spec & 0x00000001 ) {
430 fMemReader->SetEquipmentID( id );
432 fClusterizer->ResetTree();
433 Bool_t iclustered = fClusterizer->Raw2ClustersChamber(fMemReader);
434 if (iclustered == kTRUE)
436 HLTDebug( "Clustered successfully");
440 HLTError("Clustering ERROR");
444 // put the tree into output blocks
448 TTree *fcTree = fClusterizer->GetClusterTree();
451 HLTDebug("fcTree: Entries - %i; Size - %i",fcTree->GetEntriesFast(),sizeof(fcTree));
452 PushBack(fcTree, AliHLTTRDDefinitions::fgkClusterDataType, blocks[i].fSpecification);
454 fClusterizer->RecPoints()->Delete();
459 TClonesArray *clustersArray = fClusterizer->RecPoints();
460 fClusterizer->SetClustersOwner(kFALSE);
463 clustersArray->BypassStreamer(kFALSE);
464 HLTDebug("clustersArray: Entries - %i; Size - %i",clustersArray->GetEntriesFast(),sizeof(clustersArray));
465 PushBack(clustersArray, AliHLTTRDDefinitions::fgkClusterDataType, blocks[i].fSpecification);
466 clustersArray->Delete();
467 delete clustersArray;
470 HLTWarning("Array of clusters is empty!");
476 }// for ( unsigned long i = 0; i < evtData.fBlockCnt; i++ )
481 size=0; // this function did not write data to the buffer directly