Bug fix: AliHLTTPCRawSpacePointContainer
[u/mrichter/AliRoot.git] / HLT / global / AliHLTGlobalEsdConverterComponent.cxx
... / ...
CommitLineData
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//* 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 AliHLTGlobalEsdConverterComponent.cxx
20// @author Matthias Richter
21// @date
22// @brief Global ESD converter component.
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 <cassert>
32#include "AliHLTGlobalEsdConverterComponent.h"
33#include "AliHLTGlobalBarrelTrack.h"
34#include "AliHLTExternalTrackParam.h"
35#include "AliHLTTrackMCLabel.h"
36#include "AliHLTCTPData.h"
37#include "AliHLTErrorGuard.h"
38#include "AliESDEvent.h"
39#include "AliESDtrack.h"
40#include "AliESDMuonTrack.h"
41#include "AliESDMuonCluster.h"
42#include "AliCDBEntry.h"
43#include "AliCDBManager.h"
44#include "AliPID.h"
45#include "TTree.h"
46#include "TList.h"
47#include "TClonesArray.h"
48#include "AliHLTESDCaloClusterMaker.h"
49#include "AliHLTCaloClusterDataStruct.h"
50#include "AliHLTCaloClusterReader.h"
51#include "AliESDCaloCluster.h"
52#include "AliESDVZERO.h"
53#include "AliHLTGlobalVertexerComponent.h"
54#include "AliHLTVertexFinderBase.h"
55
56/** ROOT macro for the implementation of ROOT specific class methods */
57ClassImp(AliHLTGlobalEsdConverterComponent)
58
59AliHLTGlobalEsdConverterComponent::AliHLTGlobalEsdConverterComponent()
60 : AliHLTProcessor()
61 , fWriteTree(0)
62 , fVerbosity(0)
63 , fESD(NULL)
64 , fSolenoidBz(-5.00668)
65 , fBenchmark("EsdConverter")
66{
67 // see header file for class documentation
68 // or
69 // refer to README to build package
70 // or
71 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
72}
73
74AliHLTGlobalEsdConverterComponent::~AliHLTGlobalEsdConverterComponent()
75{
76 // see header file for class documentation
77 if (fESD) delete fESD;
78 fESD=NULL;
79}
80
81int AliHLTGlobalEsdConverterComponent::Configure(const char* arguments)
82{
83 // see header file for class documentation
84 int iResult=0;
85 if (!arguments) return iResult;
86
87 TString allArgs=arguments;
88 TString argument;
89 int bMissingParam=0;
90
91 TObjArray* pTokens=allArgs.Tokenize(" ");
92 if (pTokens) {
93 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
94 argument=((TObjString*)pTokens->At(i))->String();
95 if (argument.IsNull()) continue;
96
97 if (argument.CompareTo("-solenoidBz")==0) {
98 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
99 HLTWarning("argument -solenoidBz is deprecated, magnetic field set up globally (%f)", GetBz());
100 continue;
101 } else {
102 HLTError("unknown argument %s", argument.Data());
103 iResult=-EINVAL;
104 break;
105 }
106 }
107 delete pTokens;
108 }
109 if (bMissingParam) {
110 HLTError("missing parameter for argument %s", argument.Data());
111 iResult=-EINVAL;
112 }
113
114 return iResult;
115}
116
117int AliHLTGlobalEsdConverterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
118{
119 // see header file for class documentation
120 int iResult=0;
121 const char* path=NULL;
122 const char* defaultNotify="";
123 if (cdbEntry) {
124 path=cdbEntry;
125 defaultNotify=" (default)";
126 }
127 if (path) {
128 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
129 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
130 if (pEntry) {
131 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
132 if (pString) {
133 HLTInfo("received configuration object string: \'%s\'", pString->String().Data());
134 iResult=Configure(pString->String().Data());
135 } else {
136 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
137 }
138 } else {
139 HLTError("can not fetch object \"%s\" from CDB", path);
140 }
141 }
142
143 return iResult;
144}
145
146void AliHLTGlobalEsdConverterComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
147{
148 // see header file for class documentation
149 list.push_back(kAliHLTDataTypeTrack);
150 list.push_back(kAliHLTDataTypeTrackMC);
151 list.push_back(kAliHLTDataTypeCaloCluster);
152 list.push_back(kAliHLTDataTypedEdx );
153 list.push_back(kAliHLTDataTypeESDVertex );
154 list.push_back(kAliHLTDataTypeESDObject);
155 list.push_back(kAliHLTDataTypeTObject);
156 list.push_back(kAliHLTDataTypeGlobalVertexer);
157 list.push_back(kAliHLTDataTypeV0Finder); // array of track ids for V0s
158 list.push_back(kAliHLTDataTypeKFVertex); // KFVertex object from vertexer
159 list.push_back(kAliHLTDataTypePrimaryFinder); // array of track ids for prim vertex
160 list.push_back(kAliHLTDataTypeESDContent);
161}
162
163AliHLTComponentDataType AliHLTGlobalEsdConverterComponent::GetOutputDataType()
164{
165 // see header file for class documentation
166 return kAliHLTDataTypeESDObject|kAliHLTDataOriginOut;
167}
168
169void AliHLTGlobalEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
170{
171 // see header file for class documentation
172 constBase=2000000;
173 inputMultiplier=10.0;
174}
175
176int AliHLTGlobalEsdConverterComponent::DoInit(int argc, const char** argv)
177{
178 // see header file for class documentation
179 int iResult=0;
180 TString argument="";
181 int bMissingParam=0;
182
183 // default list of skiped ESD objects
184 TString skipObjects=
185 // "AliESDRun,"
186 // "AliESDHeader,"
187 // "AliESDZDC,"
188 "AliESDFMD,"
189 // "AliESDVZERO,"
190 // "AliESDTZERO,"
191 // "TPCVertex,"
192 // "SPDVertex,"
193 // "PrimaryVertex,"
194 // "AliMultiplicity,"
195 // "PHOSTrigger,"
196 // "EMCALTrigger,"
197 // "SPDPileupVertices,"
198 // "TrkPileupVertices,"
199 "Cascades,"
200 "Kinks,"
201 "AliRawDataErrorLogs,"
202 "AliESDACORDE";
203
204 iResult=Reconfigure(NULL, NULL);
205 TString allArgs = "";
206 for ( int i = 0; i < argc; i++ ) {
207 if ( !allArgs.IsNull() ) allArgs += " ";
208 allArgs += argv[i];
209 }
210
211 TObjArray* pTokens=allArgs.Tokenize(" ");
212 if (pTokens) {
213 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
214 argument=((TObjString*)pTokens->At(i))->String();
215 if (argument.IsNull()) continue;
216
217 // -notree
218 if (argument.CompareTo("-notree")==0) {
219 fWriteTree=0;
220
221 // -tree
222 } else if (argument.CompareTo("-tree")==0) {
223 fWriteTree=1;
224 } else if (argument.CompareTo("-solenoidBz")==0) {
225 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
226 HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->String().Data());
227 HLTWarning("argument '-solenoidBz' is deprecated, solenoid field initiaized from CDB settings");
228 continue;
229 } else if (argument.Contains("-skipobject=")) {
230 argument.ReplaceAll("-skipobject=", "");
231 skipObjects=argument;
232 } else {
233 HLTError("unknown argument %s", argument.Data());
234 iResult=-EINVAL;
235 break;
236 }
237 }
238 }
239 if (bMissingParam) {
240 HLTError("missing parameter for argument %s", argument.Data());
241 iResult=-EINVAL;
242 }
243
244 fSolenoidBz=GetBz();
245
246 if (iResult>=0) {
247 fESD = new AliESDEvent;
248 if (fESD) {
249 fESD->CreateStdContent();
250
251 // remove some of the objects which are not needed
252 if (fESD->GetList() && !skipObjects.IsNull()) {
253 pTokens=skipObjects.Tokenize(",");
254 if (pTokens) {
255 const char* id=NULL;
256 TIter next(pTokens);
257 TObject* pObject=NULL;
258 while ((pObject=next())!=NULL) {
259 id=((TObjString*)pObject)->String().Data();
260 TObject* pObj=fESD->GetList()->FindObject(id);
261 if (pObj) {
262 HLTDebug("removing object %s", id);
263 fESD->GetList()->Remove(pObj);
264 delete pObj;
265 } else {
266 HLTWarning("failed to remove object '%s' from ESD", id);
267 }
268 }
269 fESD->GetStdContent();
270 delete pTokens;
271 }
272 }
273 } else {
274 iResult=-ENOMEM;
275 }
276
277 SetupCTPData();
278 }
279
280 fBenchmark.SetTimer(0,"total");
281
282 return iResult;
283}
284
285int AliHLTGlobalEsdConverterComponent::DoDeinit()
286{
287 // see header file for class documentation
288 if (fESD) delete fESD;
289 fESD=NULL;
290
291 return 0;
292}
293
294int AliHLTGlobalEsdConverterComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
295 AliHLTComponentTriggerData& trigData)
296{
297 // see header file for class documentation
298 int iResult=0;
299 if (!fESD) return -ENODEV;
300
301 if (IsDataEvent()) fBenchmark.StartNewEvent();
302 fBenchmark.Start(0);
303
304 AliESDEvent* pESD = fESD;
305
306 pESD->Reset();
307 pESD->SetMagneticField(fSolenoidBz);
308 pESD->SetRunNumber(GetRunNo());
309 pESD->SetPeriodNumber(GetPeriodNumber());
310 pESD->SetOrbitNumber(GetOrbitNumber());
311 pESD->SetBunchCrossNumber(GetBunchCrossNumber());
312 pESD->SetTimeStamp(GetTimeStamp());
313
314 const AliHLTCTPData* pCTPData=CTPData();
315 if (pCTPData) {
316 AliHLTUInt64_t mask=pCTPData->ActiveTriggers(trigData);
317 for (int index=0; index<gkNCTPTriggerClasses; index++) {
318 if ((mask&((AliHLTUInt64_t)0x1<<index)) == 0) continue;
319 pESD->SetTriggerClass(pCTPData->Name(index), index);
320 }
321 pESD->SetTriggerMask(mask);
322 }
323
324 TTree* pTree = NULL;
325 if (fWriteTree)
326 pTree = new TTree("esdTree", "Tree with HLT ESD objects");
327
328 if (pTree) {
329 pTree->SetDirectory(0);
330 }
331
332 if ((iResult=ProcessBlocks(pTree, pESD))>=0) {
333 // TODO: set the specification correctly
334 if (pTree) {
335 // the esd structure is written to the user info and is
336 // needed in te ReadFromTree method to read all objects correctly
337 pTree->GetUserInfo()->Add(pESD);
338 pESD->WriteToTree(pTree);
339 iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginOut, 0);
340 } else {
341 iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
342 }
343 fBenchmark.AddOutput(GetLastObjectSize());
344 }
345 if (pTree) {
346 // clear user info list to prevent objects from being deleted
347 pTree->GetUserInfo()->Clear();
348 delete pTree;
349 }
350
351 fBenchmark.Stop(0);
352 HLTInfo( fBenchmark.GetStatistics() );
353
354 return iResult;
355}
356
357int AliHLTGlobalEsdConverterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD)
358{
359 // see header file for class documentation
360
361 int iResult=0;
362 int iAddedDataBlocks=0;
363
364 // check if there is an ESD block in the input and copy its content first to the
365 // ESD
366 const TObject* pEsdObj = GetFirstInputObject(kAliHLTDataTypeESDObject, "AliESDEvent");
367 AliESDEvent* pInputESD=NULL;
368 if (pEsdObj) pInputESD=dynamic_cast<AliESDEvent*>(const_cast<TObject*>(pEsdObj));
369 if (pInputESD) {
370 pInputESD->GetStdContent();
371 *pESD=*pInputESD;
372 }
373 if (GetNextInputObject()!=NULL) {
374 HLTWarning("handling of multiple ESD input objects not implemented, skipping input");
375 }
376
377 // Barrel tracking
378 // tracks are based on the TPC tracks, and only updated from the ITS information
379 // Sequence:
380 // 1) extract MC information for TPC and ITS from specific data blocks and store in
381 // intermediate vector arrays
382 // 2) extract TPC tracks, update with MC labels if available, the track parameters
383 // are estimated at the first cluster position
384 // 2.1) propagate to last cluster position and update kTPCout, sets also outer param (fOp)
385 // 2.2) update kTPCin, sets also inner param (fIp) and TPC inner param (fTPCInner)
386 // 2.3) update kTPCrefit using the same parameters at the first cluster position
387 // HLT has strictly spoking no refit, but we want the flag to be set
388 // can be changed to be done after all the individual barrel detector parameters
389 // have been updated by looping over the tracks again
390 // 3) extract ITS tracks, the tracks are actually TPC tracks updated from the ITS
391 // tracking information
392 // 3.1) TODO 2010-07-12: handle ITS standalone tracks by updating kITSout before kITSin
393 // 3.2) update with kITSin
394 // TODO 2010-07-12 find out if the kITSrefit has to be set as well
395 // 4) extract TRD tracks and add to ESD
396 // TODO 2010-07-12 at the moment there is no matching or merging of TPC and TRD tracks
397 // 5) Add Trigger Detectors
398 // VZERO, ZDC
399
400 // 1) first read MC information (if present)
401 std::map<int,int> mcLabelsTPC;
402 std::map<int,int> mcLabelsITS;
403
404 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
405 pBlock!=NULL; pBlock=GetNextInputBlock()) {
406
407 fBenchmark.AddInput(pBlock->fSize);
408
409 AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
410 if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
411 for( unsigned int il=0; il<dataPtr->fCount; il++ ){
412 AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
413 mcLabelsTPC[lab.fTrackID] = lab.fMCLabel;
414 }
415 } else {
416 HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information",
417 DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification,
418 dataPtr->fCount, pBlock->fSize);
419 }
420 }
421
422 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginITS);
423 pBlock!=NULL; pBlock=GetNextInputBlock()) {
424
425 fBenchmark.AddInput(pBlock->fSize);
426
427 AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
428 if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
429 for( unsigned int il=0; il<dataPtr->fCount; il++ ){
430 AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
431 mcLabelsITS[lab.fTrackID] = lab.fMCLabel;
432 }
433 } else {
434 HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information",
435 DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification,
436 dataPtr->fCount, pBlock->fSize);
437 }
438 }
439
440
441 // read dEdx information (if present)
442 // TODO 2010-07-12 this needs to be verified
443
444 AliHLTFloat32_t *dEdxTPC = 0;
445 Int_t ndEdxTPC = 0;
446 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypedEdx|kAliHLTDataOriginTPC);
447 pBlock!=NULL; pBlock=NULL/*GetNextInputBlock() there is only one block*/) {
448 fBenchmark.AddInput(pBlock->fSize);
449 dEdxTPC = reinterpret_cast<AliHLTFloat32_t*>( pBlock->fPtr );
450 ndEdxTPC = pBlock->fSize / (3*sizeof(AliHLTFloat32_t));
451 }
452
453 // 2) convert the TPC tracks to ESD tracks
454 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
455 pBlock!=NULL; pBlock=GetNextInputBlock()) {
456 if (pInputESD && pInputESD->GetNumberOfTracks()>0) {
457 HLTWarning("Tracks array already filled from the input esd block, additional filling from TPC tracks block might cause inconsistent content");
458 }
459 fBenchmark.AddInput(pBlock->fSize);
460 vector<AliHLTGlobalBarrelTrack> tracks;
461 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
462 for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
463 element!=tracks.end(); element++) {
464 Float_t points[4] = {
465 static_cast<Float_t>(element->GetX()),
466 static_cast<Float_t>(element->GetY()),
467 static_cast<Float_t>(element->GetLastPointX()),
468 static_cast<Float_t>(element->GetLastPointY())
469 };
470
471 Int_t mcLabel = -1;
472 if( mcLabelsTPC.find(element->TrackID())!=mcLabelsTPC.end() )
473 mcLabel = mcLabelsTPC[element->TrackID()];
474 element->SetLabel( mcLabel );
475
476 AliESDtrack iotrack;
477
478 // for the moment, the number of clusters are not set when processing the
479 // kTPCin update, only at kTPCout
480 // there ar emainly three parameters updated for kTPCout
481 // number of clusters
482 // chi2
483 // pid signal
484 // The first one can be updated already at that stage here, while the two others
485 // eventually require to update from the ITS tracks before. The exact scheme
486 // needs to be checked
487 iotrack.SetID( element->TrackID() );
488
489 // 2.1 set kTPCout
490 // TODO 2010-07-12 disabled for the moment because of bug
491 // https://savannah.cern.ch/bugs/index.php?69875
492 // the propagation does not work, there is an offset in z
493 // propagate to last cluster and update parameters with flag kTPCout
494 // Note: updating with flag kITSout further down will overwrite the outer
495 // parameter again so this can be done only for ITS standalone tracks, meaning
496 // tracks in the ITS not associated with any TPC track
497 // HLT does not provide such standalone tracking
498 {
499 AliHLTGlobalBarrelTrack outPar(*element);
500 //outPar.AliExternalTrackParam::PropagateTo( element->GetLastPointX(), fSolenoidBz );
501 const Int_t N=10; // number of steps.
502 const Float_t xRange = element->GetLastPointX() - element->GetX();
503 const Float_t xStep = xRange / N ;
504 for(int i = 1; i <= N; ++i) {
505 if(!outPar.AliExternalTrackParam::PropagateTo(element->GetX() + xStep * i, fSolenoidBz)) break;
506 }
507 outPar.SetLabel(element->GetLabel());
508 iotrack.UpdateTrackParams(&outPar,AliESDtrack::kTPCout);
509 }
510 // 2.2 TPC tracking estimates parameters at first cluster
511 iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
512
513 // 2.3 use the same parameters also for kTPCrefit, there is no proper refit in HLT
514 // maybe this can be done later, however also the inner parameter is changed which
515 // is used as a reference point in the display. That point should be at the first
516 // TPC cluster
517 iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCrefit);
518 iotrack.SetTPCPoints(points);
519 if( element->TrackID()<ndEdxTPC ){
520 AliHLTFloat32_t *val = &(dEdxTPC[3*element->TrackID()]);
521 iotrack.SetTPCsignal( val[0], val[1], (UChar_t) val[2] );
522 //AliTPCseed s;
523 //s.Set( element->GetX(), element->GetAlpha(),
524 //element->GetParameter(), element->GetCovariance() );
525 //s.SetdEdx( val[0] );
526 //s.CookPID();
527 //iotrack.SetTPCpid(s.TPCrPIDs() );
528 } else {
529 if( dEdxTPC ) HLTWarning("Wrong number of dEdx TPC labels");
530 }
531 iotrack.SetLabel(mcLabel);
532 pESD->AddTrack(&iotrack);
533 if (fVerbosity>0) element->Print();
534 }
535 HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
536 iAddedDataBlocks++;
537 } else if (iResult<0) {
538 HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d",
539 DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
540 }
541 }
542
543
544 // Get ITS SPD vertex
545 for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); i!=NULL; i=GetNextInputBlock() ){
546 fBenchmark.AddInput(i->fSize);
547 }
548
549 for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); iter != NULL; iter = GetNextInputObject() ) {
550 AliESDVertex *vtx = dynamic_cast<AliESDVertex*>(const_cast<TObject*>( iter ) );
551 pESD->SetPrimaryVertexSPD( vtx );
552 }
553
554 // 3.1. now update ESD tracks with the ITSOut info
555 // updating track parameters with flag kITSout will overwrite parameter set above with flag kTPCout
556 // TODO 2010-07-12 there are some issues with this updating sequence, for the moment update with
557 // flag kITSout is disabled, would require
558 // a) to update kTPCout after kITSout
559 // b) update only for ITS standalone tracks. The HLT ITS tracker does not provide standalone
560 // tracking, every track is related to a TPC track
561 // Furthermore there seems to be a bug as the data blocks describe track parameters around the
562 // vertex, not at the outer ITS
563 // bug https://savannah.cern.ch/bugs/index.php?69872
564 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITSOut);
565 pBlock!=NULL; pBlock=GetNextInputBlock()) {
566 fBenchmark.AddInput(pBlock->fSize);
567 vector<AliHLTGlobalBarrelTrack> tracks;
568 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
569 for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
570 element!=tracks.end(); element++) {
571 int tpcID=element->TrackID();
572 // the ITS tracker assigns the TPC track used as seed for a certain track to
573 // the trackID
574 if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
575 AliESDtrack *tESD = pESD->GetTrack( tpcID );
576 element->SetLabel(tESD->GetLabel());
577 // 2010-07-12 disabled, see above, bugfix https://savannah.cern.ch/bugs/index.php?69872
578 //if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSout );
579 }
580 }
581 }
582
583 // 3.2. now update ESD tracks with the ITS info
584 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
585 pBlock!=NULL; pBlock=GetNextInputBlock()) {
586 fBenchmark.AddInput(pBlock->fSize);
587 vector<AliHLTGlobalBarrelTrack> tracks;
588 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
589 for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
590 element!=tracks.end(); element++) {
591 int tpcID=element->TrackID();
592 // the ITS tracker assigns the TPC track used as seed for a certain track to
593 // the trackID
594 if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
595 Int_t mcLabel = -1;
596 if( mcLabelsITS.find(element->TrackID())!=mcLabelsITS.end() )
597 mcLabel = mcLabelsITS[element->TrackID()];
598 AliESDtrack *tESD = pESD->GetTrack( tpcID );
599 if (!tESD) continue;
600 // the labels for the TPC and ITS tracking params can be different, e.g.
601 // there can be a decay. The ITS label should then be the better one, the
602 // TPC label is saved in a member of AliESDtrack
603 if (mcLabel>=0) {
604 // upadte only if the ITS label is available, otherwise keep TPC label
605 element->SetLabel( mcLabel );
606 } else {
607 // bugfix https://savannah.cern.ch/bugs/?69713
608 element->SetLabel( tESD->GetLabel() );
609 }
610 tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
611
612 // TODO: add a proper refit
613 //tESD->UpdateTrackParams( &(*element), AliESDtrack::kTPCrefit );
614 }
615 }
616 }
617
618 // update with vertices and vertex-fitted tracks
619 // output of the GlobalVertexerComponent
620 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeGlobalVertexer);
621 pBlock!=NULL; pBlock=GetNextInputBlock()) {
622 fBenchmark.AddInput(pBlock->fSize);
623 AliHLTGlobalVertexerComponent::FillESD( pESD, reinterpret_cast<AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData* >(pBlock->fPtr) );
624 }
625
626 // update with vertices and vertex-fitted tracks
627 // output of PrimaryVertexer and V0Finder components
628 TObject* pBase = (TObject*)GetFirstInputObject(kAliHLTDataTypeKFVertex | kAliHLTDataOriginOut);
629 if (pBase) {
630 AliKFVertex* kfVertex = dynamic_cast<AliKFVertex *>(pBase);
631 if (kfVertex) {
632 const AliHLTComponentBlockData* pP = GetFirstInputBlock(kAliHLTDataTypePrimaryFinder | kAliHLTDataOriginOut);
633 if (pP && pP->fSize && pP->fPtr) {
634 const AliHLTComponentBlockData* pV0 = GetFirstInputBlock(kAliHLTDataTypeV0Finder | kAliHLTDataOriginOut);
635 if (pV0 && pV0->fPtr && pInputESD && pInputESD->GetNumberOfV0s()>0) {
636 const int* v0s = static_cast<const int*>(pV0->fPtr);
637 HLTWarning("V0 array already filled from the input esd block, additional filling from V0 block of %d entries might cause inconsistent content", v0s[0]);
638 }
639 AliHLTVertexFinderBase::FillESD(pESD, kfVertex, pP->fPtr, pV0?pV0->fPtr:NULL);
640 } else
641 HLTWarning("Problem with primary finder's data block");
642 } else {
643 HLTWarning("primary vertex block of wrong type, expecting AliKFVertex instead of %s", pBase->GetName());
644 }
645 } else {
646 // throw an error if there is a V0 data block which can not be handled without
647 // the AliKFVertex object
648 if (GetFirstInputBlock(kAliHLTDataTypeV0Finder | kAliHLTDataOriginOut)!=NULL) {
649 ALIHLTERRORGUARD(1, "missing AliKFVertex object ignoring V0 data block of type %s",
650 DataType2Text(kAliHLTDataTypeV0Finder|kAliHLTDataOriginOut).c_str());
651 }
652 }
653
654 // Fill DCA parameters for TPC tracks
655 // TODO 2010-07-12 this propagates also the TPC inner param to beamline
656 // sounds not very reasonable
657 // https://savannah.cern.ch/bugs/index.php?69873
658 for (int i=0; i<pESD->GetNumberOfTracks(); i++) {
659 if (!pESD->GetTrack(i) ||
660 !pESD->GetTrack(i)->GetTPCInnerParam() ) continue;
661 pESD->GetTrack(i)->RelateToVertexTPC(pESD->GetPrimaryVertexTracks(), fSolenoidBz, 1000 );
662 }
663
664 // loop over all tracks and set the TPC refit flag by updating with the
665 // original TPC inner parameter if not yet set
666 // TODO: replace this by a proper refit
667 // code is comented for the moment as it does not fully solve the problems with
668 // the display
669 // - would set the main parameters to the TPC inner wall again, or
670 // - changes the inner param if the parameters are propagated, so we loose the track
671 // reference point for the display
672 // with the current sequence we have the latter case as the DCA operations above
673 // change the TPC inner parameters
674 /*
675 for (int i=0; i<pESD->GetNumberOfTracks(); i++) {
676 if (!pESD->GetTrack(i) ||
677 !pESD->GetTrack(i)->GetTPCInnerParam() ||
678 pESD->GetTrack(i)->IsOn(AliESDtrack::kTPCrefit)) continue;
679 AliESDtrack* tESD=pESD->GetTrack(i);
680 AliHLTGlobalBarrelTrack inner(*tESD->GetTPCInnerParam());
681 inner.SetLabel(tESD->GetLabel());
682 tESD->UpdateTrackParams(&inner, AliESDtrack::kTPCrefit);
683 }
684 */
685
686 // 4. convert the HLT TRD tracks to ESD tracks
687 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
688 pBlock!=NULL; pBlock=GetNextInputBlock()) {
689 fBenchmark.AddInput(pBlock->fSize);
690 vector<AliHLTGlobalBarrelTrack> tracks;
691 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
692 for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
693 element!=tracks.end(); element++) {
694
695 Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
696 for(Int_t i=0; i<AliPID::kSPECIES; i++){
697 switch(i){
698 case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
699 default: TRDpid[i]=restProb; break;
700 }
701 }
702
703 AliESDtrack iotrack;
704 iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDout);
705 iotrack.SetStatus(AliESDtrack::kTRDin);
706 iotrack.SetTRDpid(TRDpid);
707
708 pESD->AddTrack(&iotrack);
709 if (fVerbosity>0) element->Print();
710 }
711 HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
712 iAddedDataBlocks++;
713 } else if (iResult<0) {
714 HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d",
715 DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
716 }
717 }
718 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeCaloCluster | kAliHLTDataOriginAny); pBlock!=NULL; pBlock=GetNextInputBlock())
719 {
720 fBenchmark.AddInput(pBlock->fSize);
721 AliHLTCaloClusterHeaderStruct *caloClusterHeaderPtr = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
722
723 HLTDebug("%d HLT clusters from spec: 0x%X", caloClusterHeaderPtr->fNClusters, pBlock->fSpecification);
724
725 //AliHLTCaloClusterReader reader;
726 //reader.SetMemory(caloClusterHeaderPtr);
727
728 AliHLTESDCaloClusterMaker clusterMaker;
729
730 int nClusters = clusterMaker.FillESD(pESD, caloClusterHeaderPtr);
731
732 HLTInfo("converted %d cluster(s) to AliESDCaloCluster and added to ESD", nClusters);
733 iAddedDataBlocks++;
734 }
735
736 // 5) Add Trigger Detectors
737 // VZERO, ZDC
738
739 // FIXME: the size of all input blocks can be added in one loop
740 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeESDContent|kAliHLTDataOriginVZERO);
741 pBlock!=NULL; pBlock=GetNextInputBlock()) {
742 fBenchmark.AddInput(pBlock->fSize);
743 }
744
745 for ( const TObject *pObject = GetFirstInputObject(kAliHLTDataTypeESDContent|kAliHLTDataOriginVZERO);
746 pObject != NULL; pObject = GetNextInputObject() ) {
747 AliESDVZERO *esdVZERO = dynamic_cast<AliESDVZERO*>(const_cast<TObject*>( pObject ) );
748 if (esdVZERO) {
749 pESD->SetVZEROData( esdVZERO );
750 break;
751 } else {
752 ALIHLTERRORGUARD(1, "input object of data type %s is not of class AliESDVZERO",
753 DataType2Text(kAliHLTDataTypeESDContent|kAliHLTDataOriginVZERO).c_str());
754 }
755 }
756
757 // FIXME: the size of all input blocks can be added in one loop
758 for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeESDContent|kAliHLTDataOriginZDC);
759 pBlock!=NULL; pBlock=GetNextInputBlock()) {
760 fBenchmark.AddInput(pBlock->fSize);
761 }
762 for ( const TObject *pObject = GetFirstInputObject(kAliHLTDataTypeESDContent|kAliHLTDataOriginZDC);
763 pObject != NULL; pObject = GetNextInputObject() ) {
764 AliESDZDC *esdZDC = dynamic_cast<AliESDZDC*>(const_cast<TObject*>( pObject ) );
765 if (esdZDC) {
766#ifndef HAVE_NOT_ALIZDCRECONSTRUCTOR_r43770
767 pESD->SetZDCData( esdZDC );
768#else
769 ALIHLTERRORGUARD(1, "Processing of ZDC data requires AliRoot r43770m skipping data block of type %s",
770 DataType2Text(kAliHLTDataTypeESDContent|kAliHLTDataOriginZDC).c_str());
771#endif
772 break;
773 } else {
774 ALIHLTERRORGUARD(1, "input object of data type %s is not of class AliESDZDC",
775 DataType2Text(kAliHLTDataTypeESDContent|kAliHLTDataOriginZDC).c_str());
776 }
777 }
778
779 // Add tracks and clusters from MUON.
780 for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTAnyDataType | kAliHLTDataOriginMUON); i!=NULL; i=GetNextInputBlock() ){
781 fBenchmark.AddInput(i->fSize);
782 }
783
784 for (const TObject* obj = GetFirstInputObject(kAliHLTAnyDataType | kAliHLTDataOriginMUON);
785 obj != NULL;
786 obj = GetNextInputObject()
787 )
788 {
789 const TClonesArray* tracklist = NULL;
790 const TClonesArray* clusterlist = NULL;
791 if (obj->IsA() == AliESDEvent::Class())
792 {
793 const AliESDEvent* event = static_cast<const AliESDEvent*>(obj);
794 HLTDebug("Received a MUON ESD with specification: 0x%X", GetSpecification(obj));
795 if (event->GetList() == NULL) continue;
796 tracklist = dynamic_cast<const TClonesArray*>(event->GetList()->FindObject("MuonTracks"));
797 if (tracklist == NULL) continue;
798 clusterlist = dynamic_cast<const TClonesArray*>(event->GetList()->FindObject("MuonClusters"));
799 if (clusterlist == NULL) continue;
800 }
801 else if (obj->IsA() == TClonesArray::Class())
802 {
803 if (!strcmp(obj->GetName(), "MuonTracks")) {
804 tracklist = static_cast<const TClonesArray*>(obj);
805 HLTDebug("Received a MUON TClonesArray of tracks with specification: 0x%X", GetSpecification(obj));
806 } else {
807 clusterlist = static_cast<const TClonesArray*>(obj);
808 HLTDebug("Received a MUON TClonesArray of clusters with specification: 0x%X", GetSpecification(obj));
809 }
810 }
811 else
812 {
813 // Cannot handle this object type.
814 continue;
815 }
816 // copy tracks
817 if (tracklist) {
818 HLTDebug("Received %d MUON tracks.", tracklist->GetEntriesFast());
819 if (tracklist->GetEntriesFast() > 0)
820 {
821 const AliESDMuonTrack* track = dynamic_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(0));
822 if (track == NULL)
823 {
824 HLTError(Form("%s from MUON does not contain AliESDMuonTrack objects.", obj->ClassName()));
825 continue;
826 }
827 }
828 for (Int_t i = 0; i < tracklist->GetEntriesFast(); ++i)
829 {
830 AliESDMuonTrack* track = pESD->NewMuonTrack();
831 *track = *(static_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(i)));
832 }
833 }
834 // copy clusters
835 if (clusterlist) {
836 HLTDebug("Received %d MUON clusters.", clusterlist->GetEntriesFast());
837 if (clusterlist->GetEntriesFast() > 0)
838 {
839 const AliESDMuonCluster* cluster = dynamic_cast<const AliESDMuonCluster*>(clusterlist->UncheckedAt(0));
840 if (cluster == NULL)
841 {
842 HLTError(Form("%s from MUON does not contain AliESDMuonCluster objects.", obj->ClassName()));
843 continue;
844 }
845 }
846 for (Int_t i = 0; i < clusterlist->GetEntriesFast(); ++i)
847 {
848 AliESDMuonCluster* cluster = pESD->NewMuonCluster();
849 *cluster = *(static_cast<const AliESDMuonCluster*>(clusterlist->UncheckedAt(i)));
850 }
851 }
852 }
853
854 if (iAddedDataBlocks>0 && pTree) {
855 pTree->Fill();
856 }
857
858 if (iResult>=0) iResult=iAddedDataBlocks;
859 return iResult;
860}