3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project *
5 //* ALICE Experiment at CERN, All rights reserved. *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no> *
8 //* for The ALICE HLT 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 AliHLTGlobalEsdConverterComponent.cxx
20 @author Matthias Richter
22 @brief Global ESD converter component.
25 // see header file for class documentation
27 // refer to README to build package
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
32 #include "AliHLTGlobalEsdConverterComponent.h"
33 #include "AliHLTGlobalBarrelTrack.h"
34 #include "AliHLTExternalTrackParam.h"
35 #include "AliHLTTrackMCLabel.h"
36 #include "AliESDEvent.h"
37 #include "AliESDtrack.h"
38 #include "AliCDBEntry.h"
39 #include "AliCDBManager.h"
44 /** ROOT macro for the implementation of ROOT specific class methods */
45 ClassImp(AliHLTGlobalEsdConverterComponent)
47 AliHLTGlobalEsdConverterComponent::AliHLTGlobalEsdConverterComponent()
54 // see header file for class documentation
56 // refer to README to build package
58 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
61 AliHLTGlobalEsdConverterComponent::~AliHLTGlobalEsdConverterComponent()
63 // see header file for class documentation
64 if (fESD) delete fESD;
68 int AliHLTGlobalEsdConverterComponent::Configure(const char* arguments)
70 // see header file for class documentation
72 if (!arguments) return iResult;
74 TString allArgs=arguments;
78 TObjArray* pTokens=allArgs.Tokenize(" ");
80 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
81 argument=((TObjString*)pTokens->At(i))->GetString();
82 if (argument.IsNull()) continue;
84 if (argument.CompareTo("-solenoidBz")==0) {
85 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
86 HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
87 fSolenoidBz=((TObjString*)pTokens->At(i))->GetString().Atof();
90 HLTError("unknown argument %s", argument.Data());
98 HLTError("missing parameter for argument %s", argument.Data());
105 int AliHLTGlobalEsdConverterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
107 // see header file for class documentation
109 const char* path=kAliHLTCDBSolenoidBz;
110 const char* defaultNotify="";
113 defaultNotify=" (default)";
116 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
117 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
119 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
121 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
122 iResult=Configure(pString->GetString().Data());
124 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
127 HLTError("can not fetch object \"%s\" from CDB", path);
134 void AliHLTGlobalEsdConverterComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
136 // see header file for class documentation
137 list.push_back(kAliHLTDataTypeTrack);
138 list.push_back(kAliHLTDataTypeTrackMC);
141 AliHLTComponentDataType AliHLTGlobalEsdConverterComponent::GetOutputDataType()
143 // see header file for class documentation
144 return kAliHLTDataTypeESDObject|kAliHLTDataOriginOut;
147 void AliHLTGlobalEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
149 // see header file for class documentation
151 inputMultiplier=10.0;
154 int AliHLTGlobalEsdConverterComponent::DoInit(int argc, const char** argv)
156 // see header file for class documentation
160 for (int i=0; i<argc && iResult>=0; i++) {
162 if (argument.IsNull()) continue;
165 if (argument.CompareTo("-notree")==0) {
169 } else if (argument.CompareTo("-tree")==0) {
173 HLTError("unknown argument %s", argument.Data());
178 HLTError("missing parameter for argument %s", argument.Data());
183 iResult=Reconfigure(NULL, NULL);
187 fESD = new AliESDEvent;
189 fESD->CreateStdContent();
198 int AliHLTGlobalEsdConverterComponent::DoDeinit()
200 // see header file for class documentation
201 if (fESD) delete fESD;
207 int AliHLTGlobalEsdConverterComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
208 AliHLTComponentTriggerData& /*trigData*/)
210 // see header file for class documentation
212 if (!fESD) return -ENODEV;
214 AliESDEvent* pESD = fESD;
217 pESD->SetMagneticField(fSolenoidBz);
218 pESD->SetPeriodNumber(GetPeriodNumber());
219 pESD->SetOrbitNumber(GetOrbitNumber());
220 pESD->SetBunchCrossNumber(GetBunchCrossNumber());
221 pESD->SetTimeStamp(GetTimeStamp());
225 pTree = new TTree("esdTree", "Tree with HLT ESD objects");
228 pTree->SetDirectory(0);
231 if ((iResult=ProcessBlocks(pTree, pESD))>=0) {
232 // TODO: set the specification correctly
234 // the esd structure is written to the user info and is
235 // needed in te ReadFromTree method to read all objects correctly
236 pTree->GetUserInfo()->Add(pESD);
237 pESD->WriteToTree(pTree);
238 iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginOut, 0);
240 iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
244 // clear user info list to prevent objects from being deleted
245 pTree->GetUserInfo()->Clear();
251 int AliHLTGlobalEsdConverterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD)
253 // see header file for class documentation
256 int iAddedDataBlocks=0;
260 // in the first attempt this component reads the TPC tracks and updates in the
261 // second step from the ITS tracks
264 // first read MC information (if present)
265 std::map<int,int> mcLabels;
267 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
268 pBlock!=NULL; pBlock=GetNextInputBlock()) {
269 AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
270 if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
271 for( unsigned int il=0; il<dataPtr->fCount; il++ ){
272 AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
273 mcLabels[lab.fTrackID] = lab.fMCLabel;
276 HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information",
277 DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification,
278 dataPtr->fCount, pBlock->fSize);
282 // convert the TPC tracks to ESD tracks
283 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
284 pBlock!=NULL; pBlock=GetNextInputBlock()) {
285 vector<AliHLTGlobalBarrelTrack> tracks;
286 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
287 for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
288 element!=tracks.end(); element++) {
289 Float_t points[4] = {
292 element->GetLastPointX(),
293 element->GetLastPointY()
297 if( mcLabels.find(element->TrackID())!=mcLabels.end() )
298 mcLabel = mcLabels[element->TrackID()];
299 element->SetLabel( mcLabel );
302 iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
303 iotrack.SetTPCPoints(points);
305 pESD->AddTrack(&iotrack);
306 if (fVerbosity>0) element->Print();
308 HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
310 } else if (iResult<0) {
311 HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d",
312 DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
316 // now update ESD tracks with the ITS info
317 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
318 pBlock!=NULL; pBlock=GetNextInputBlock()) {
319 vector<AliHLTGlobalBarrelTrack> tracks;
320 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
321 for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
322 element!=tracks.end(); element++) {
324 const UInt_t* pointsArray=element->GetPoints();
325 for( unsigned il=0; il<element->GetNumberOfPoints(); il++ ){
326 // TODO: check what needs to be done with the clusters
327 if( pointsArray[il]<~(UInt_t)0 ) {/*tITS.SetClusterIndex(ncl, tr.fClusterIds[il]);*/}
330 //tITS.SetNumberOfClusters( ncl );
331 int tpcID=element->TrackID();
332 // the ITS tracker assigns the TPC track used as seed for a certain track to
334 if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
336 AliESDtrack *tESD = pESD->GetTrack( tpcID );
337 if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
342 // convert the HLT TRD tracks to ESD tracks
343 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
344 pBlock!=NULL; pBlock=GetNextInputBlock()) {
345 vector<AliHLTGlobalBarrelTrack> tracks;
346 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
347 for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
348 element!=tracks.end(); element++) {
350 Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
351 for(Int_t i=0; i<AliPID::kSPECIES; i++){
353 case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
354 default: TRDpid[i]=restProb; break;
359 iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDin);
360 iotrack.SetTRDpid(TRDpid);
362 pESD->AddTrack(&iotrack);
363 if (fVerbosity>0) element->Print();
365 HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
367 } else if (iResult<0) {
368 HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d",
369 DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
373 // primary vertex & V0's
375 //AliHLTVertexer vertexer;
376 //vertexer.SetESD( pESD );
377 //vertexer.FindPrimaryVertex();
378 //vertexer.FindV0s();
380 if (iAddedDataBlocks>0 && pTree) {
384 if (iResult>=0) iResult=iAddedDataBlocks;