1 /**************************************************************************
2 * This file is property of and copyright by the ALICE HLT Project *
3 * All rights reserved. *
6 * Artur Szostak <artursz@iafrica.com> *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
20 /// @file AliHLTMUONESDMaker.cxx
21 /// @author Artur Szostak <artursz@iafrica.com>
22 /// @date 30 June 2008
23 /// @brief Implementation of the AliHLTMUONESDMaker component.
25 /// The ESD maker component converts dHLT raw internal reconstructed information
26 /// into AliESDEvent objects.
29 #include "AliHLTMUONESDMaker.h"
30 #include "AliHLTMUONEvent.h"
31 #include "AliHLTMUONConstants.h"
32 #include "AliHLTMUONUtils.h"
33 #include "AliHLTMUONRecHit.h"
34 #include "AliHLTMUONTriggerRecord.h"
35 #include "AliHLTMUONMansoTrack.h"
36 #include "AliHLTMUONDecision.h"
37 #include "AliMUONConstants.h"
38 #include "AliMUONVCluster.h"
39 #include "AliESDEvent.h"
40 #include "AliESDRun.h"
41 #include "AliESDMuonTrack.h"
42 #include "AliESDMuonCluster.h"
43 #include "TClonesArray.h"
49 ClassImp(AliHLTMUONESDMaker);
52 AliHLTMUONESDMaker::AliHLTMUONESDMaker() :
53 AliHLTMUONProcessor(),
54 fWarnForUnexpecedBlock(false),
55 fMakeMinimalESD(false),
56 fAddCustomData(false),
57 fMakeClonesArray(false),
58 fMakeESDDataBlock(true),
61 /// Default constructor.
65 AliHLTMUONESDMaker::~AliHLTMUONESDMaker()
67 /// Default destructor.
71 bool AliHLTMUONESDMaker::IgnoreArgument(const char* arg) const
73 /// Return true if the argument is one of -cdbpath -run or -delaysetup
74 /// to prevent the parent class from parsing these arguments in DoInit.
76 if (strcmp(arg, "-cdbpath") == 0 or strcmp(arg, "-run") == 0 or
77 strcmp(arg, "-delaysetup") == 0)
88 int AliHLTMUONESDMaker::DoInit(int argc, const char** argv)
90 /// Inherited from AliHLTComponent.
91 /// Parses the command line parameters and initialises the component.
93 HLTInfo("Initialising dHLT ESD maker component.");
95 // Inherit the parents functionality.
96 int result = AliHLTMUONProcessor::DoInit(argc, argv);
97 if (result != 0) return result;
99 fWarnForUnexpecedBlock = false;
100 fMakeMinimalESD = false;
101 fMakeClonesArray = false;
102 fMakeESDDataBlock = true;
104 for (int i = 0; i < argc; i++)
106 if (ArgumentAlreadyHandled(i, argv[i])) continue;
108 if (strcmp(argv[i], "-make_minimal_esd") == 0)
110 fMakeMinimalESD = true;
114 if (strcmp(argv[i], "-warn_on_unexpected_block") == 0)
116 fWarnForUnexpecedBlock = true;
120 if (strcmp(argv[i], "-add_rootified_objects") == 0)
122 fAddCustomData = true;
126 if (strcmp(argv[i], "-makeclonesarray") == 0)
128 fMakeClonesArray = true;
132 if (strcmp(argv[i], "-makeonlyclonesarray") == 0)
134 fMakeMinimalESD = true;
135 fMakeClonesArray = true;
136 fMakeESDDataBlock = false;
140 HLTError("Unknown option '%s'.", argv[i]);
148 int AliHLTMUONESDMaker::DoDeinit()
150 /// Inherited from AliHLTComponent. Performs a cleanup of the component.
152 HLTInfo("Deinitialising dHLT ESD maker component.");
157 const char* AliHLTMUONESDMaker::GetComponentID()
159 /// Inherited from AliHLTComponent. Returns the component ID.
161 return AliHLTMUONConstants::ESDMakerId();
165 AliHLTComponentDataType AliHLTMUONESDMaker::GetOutputDataType()
167 /// Inherited from AliHLTComponent. Returns kAliHLTMultipleDataType.
168 /// Refer to GetOutputDataTypes for all returned data types.
170 return kAliHLTMultipleDataType;
174 int AliHLTMUONESDMaker::GetOutputDataTypes(AliHLTComponentDataTypeList& list)
176 /// Inherited from AliHLTComponent.
177 /// Returns the ESD object and kAliHLTDataTypeTObject data type with MUON origin.
179 list.push_back(AliHLTMUONConstants::ESDDataType());
180 list.push_back(kAliHLTDataTypeTObject | kAliHLTDataOriginMUON);
185 void AliHLTMUONESDMaker::GetInputDataTypes(AliHLTComponentDataTypeList& list)
187 /// Inherited from AliHLTProcessor.
188 /// Returns the list of expected input data types.
190 list.push_back(AliHLTMUONConstants::TriggerRecordsBlockDataType());
191 list.push_back(AliHLTMUONConstants::MansoTracksBlockDataType());
192 list.push_back(AliHLTMUONConstants::TracksBlockDataType());
196 void AliHLTMUONESDMaker::GetOutputDataSize(
197 unsigned long& constBase, double& inputMultiplier
200 /// Inherited from AliHLTComponent.
201 /// Returns an estimate of the expected output data size.
203 constBase = sizeof(AliESDEvent) + 1024*1024; // The extra 1 MByte is for auxilary objects created in AliESDEvent.
204 inputMultiplier = 10;
208 AliHLTComponent* AliHLTMUONESDMaker::Spawn()
210 /// Inherited from AliHLTComponent. Creates a new object instance.
212 return new AliHLTMUONESDMaker();
216 int AliHLTMUONESDMaker::DoEvent(
217 const AliHLTComponentEventData& evtData,
218 AliHLTComponentTriggerData& trigData
221 /// Inherited from AliHLTProcessor. Processes the new event data.
224 fClusterIndex = 0; // for the cluster unique ID.
226 // Create and fill in the standard ESD objects or just create the muon
227 // tracks array if so requested.
230 TClonesArray* muonTracks = new TClonesArray("AliESDMuonTrack",2);
231 muonTracks->SetName("MuonTracks");
232 event.AddObject(muonTracks);
233 TClonesArray* muonClusters = new TClonesArray("AliESDMuonCluster",0);
234 muonClusters->SetName("MuonClusters");
235 event.AddObject(muonClusters);
236 event.GetStdContent();
240 event.CreateStdContent();
241 event.SetRunNumber(GetRunNo());
244 const AliHLTComponentBlockData* block = NULL;
245 AliHLTUInt32_t specification = 0; // Contains the output data block spec bits.
246 std::vector<const AliHLTMUONTriggerRecordStruct*> triggerRecords;
248 // First process the blocks of trigger records. We simply mark each trigger
249 // record in the triggerRecords array.
250 for (int i = 0; i < GetNumberOfInputBlocks(); i++)
252 block = GetInputBlock(i);
253 assert( block != NULL );
255 HLTDebug("Handling block: %u, with fDataType = '%s', fPtr = %p and fSize = %u bytes.",
256 i, DataType2Text(block->fDataType).c_str(), block->fPtr, block->fSize
259 if (block->fDataType == AliHLTMUONConstants::TriggerRecordsBlockDataType())
261 specification |= block->fSpecification;
262 AliHLTMUONTriggerRecordsBlockReader inblock(block->fPtr, block->fSize);
263 if (not BlockStructureOk(inblock))
265 if (DumpDataOnError()) DumpEvent(evtData, trigData);
269 for (AliHLTUInt32_t n = 0; n < inblock.Nentries(); n++)
271 triggerRecords.push_back(&inblock[n]);
274 else if (block->fDataType == AliHLTMUONConstants::RootifiedEventDataType() and fAddCustomData)
276 // Do nothing for now, will handle this later.
280 if (block->fDataType != AliHLTMUONConstants::MansoTracksBlockDataType() and
281 block->fDataType != AliHLTMUONConstants::TracksBlockDataType() )
283 // Log a message indicating that we got a data block that we
284 // do not know how to handle.
285 if (fWarnForUnexpecedBlock)
286 HLTWarning("Received a data block of a type we cannot handle: '%s', spec: 0x%X",
287 DataType2Text(block->fDataType).c_str(), block->fSpecification
291 HLTDebug("Received a data block of a type we cannot handle: '%s', spec: 0x%X",
292 DataType2Text(block->fDataType).c_str(), block->fSpecification
299 // If we were requested to add all dHLT rootified data objects then do so.
302 const AliHLTComponentDataType& type = AliHLTMUONConstants::RootifiedEventDataType();
303 const char* classname = AliHLTMUONEvent::Class_Name();
304 const TObject* obj = NULL;
305 for (obj = GetFirstInputObject(type, classname); obj != NULL; obj = GetNextInputObject())
307 // Clone the object since the ESD takes ownership of it.
308 event.AddObject(obj->Clone());
312 // Now we can look for tracks to add. We needed the ROOT trigger records
313 // and reco hits created before we can create track objects.
314 for (block = GetFirstInputBlock(AliHLTMUONConstants::MansoTracksBlockDataType());
316 block = GetNextInputBlock()
319 specification |= block->fSpecification;
320 AliHLTMUONMansoTracksBlockReader inblock(block->fPtr, block->fSize);
321 if (not BlockStructureOk(inblock))
323 if (DumpDataOnError()) DumpEvent(evtData, trigData);
327 for (AliHLTUInt32_t n = 0; n < inblock.Nentries(); n++)
329 const AliHLTMUONMansoTrackStruct& t = inblock[n];
330 AliESDMuonTrack *muTrack = event.NewMuonTrack();
332 AliHLTMUONParticleSign sign;
334 AliHLTMUONUtils::UnpackMansoTrackFlags(t.fFlags, sign, hitset);
339 case kSignMinus: signVal = +1.; break;
340 case kSignUnknown: signVal = 0.; break;
341 case kSignPlus: signVal = -1.; break;
343 HLTWarning("Got a track with an invalid sign value: %d", int(sign));
346 TVector3 mom(t.fPx, t.fPy, t.fPz);
348 muTrack->SetInverseBendingMomentum(signVal/mom.Mag());
350 muTrack->SetInverseBendingMomentum(0.);
351 muTrack->SetThetaX(atan2(t.fPx, t.fPz));
352 muTrack->SetThetaY(atan2(t.fPy, t.fPz));
354 muTrack->SetBendingCoor(0.);
355 muTrack->SetNonBendingCoor(0.);
357 // The Manso algorithm assumes the information at the
358 // Distance of Closest Approach and chamber 1 is the same
361 muTrack->SetInverseBendingMomentumAtDCA(1./mom.Mag());
363 muTrack->SetInverseBendingMomentumAtDCA(0.);
364 muTrack->SetThetaXAtDCA(atan2(t.fPx, t.fPz));
365 muTrack->SetThetaYAtDCA(atan2(t.fPy, t.fPz));
366 muTrack->SetBendingCoorAtDCA(0.);
367 muTrack->SetNonBendingCoorAtDCA(0.);
370 muTrack->SetInverseBendingMomentumUncorrected(1./mom.Mag());
372 muTrack->SetInverseBendingMomentumUncorrected(0.);
373 muTrack->SetThetaXUncorrected(atan2(t.fPx, t.fPz));
374 muTrack->SetThetaYUncorrected(atan2(t.fPy, t.fPz));
375 muTrack->SetZUncorrected(0.);
376 muTrack->SetBendingCoorUncorrected(0.);
377 muTrack->SetNonBendingCoorUncorrected(0.);
379 muTrack->SetChi2(t.fChi2);
381 // Fill in the track hit points.
382 for (int i = 0; i < 4; i++)
384 if (not hitset[i]) continue;
386 AliHLTUInt8_t chamber;
387 AliHLTUInt16_t detElemId;
388 AliHLTMUONUtils::UnpackRecHitFlags(t.fHit[i].fFlags, chamber, detElemId);
390 AliESDMuonCluster *cluster = event.NewMuonCluster();
391 cluster->SetUniqueID(AliMUONVCluster::BuildUniqueID(chamber, detElemId, fClusterIndex++));
392 cluster->SetXYZ(t.fHit[i].fX, t.fHit[i].fY, t.fHit[i].fZ);
393 cluster->SetErrXY( // Use nominal values.
394 AliHLTMUONConstants::DefaultNonBendingReso(),
395 AliHLTMUONConstants::DefaultBendingReso()
397 cluster->SetCharge(-1.); // Indicate no total charge calculated.
398 cluster->SetChi2(-1.); // Indicate no fit made.
399 muTrack->AddClusterId(cluster->GetUniqueID());
400 muTrack->AddInMuonClusterMap(i+6);
403 FillTriggerInfo(triggerRecords, t.fTrigRec, t.fId, *muTrack, event);
407 for (block = GetFirstInputBlock(AliHLTMUONConstants::TracksBlockDataType());
409 block = GetNextInputBlock()
412 specification |= block->fSpecification;
413 AliHLTMUONTracksBlockReader inblock(block->fPtr, block->fSize);
414 if (not BlockStructureOk(inblock))
416 if (DumpDataOnError()) DumpEvent(evtData, trigData);
420 for (AliHLTUInt32_t n = 0; n < inblock.Nentries(); n++)
422 const AliHLTMUONTrackStruct& t = inblock[n];
423 AliESDMuonTrack *muTrack = event.NewMuonTrack();
425 AliHLTMUONParticleSign sign;
427 AliHLTMUONUtils::UnpackTrackFlags(t.fFlags, sign, hitset);
429 muTrack->SetInverseBendingMomentum(t.fInverseBendingMomentum);
430 muTrack->SetThetaX(t.fThetaX);
431 muTrack->SetThetaY(t.fThetaY);
433 muTrack->SetBendingCoor(t.fY);
434 muTrack->SetNonBendingCoor(t.fX);
436 // The full tracker assumes the information at the
437 // Distance of Closest Approach and chamber 1 is the same
439 muTrack->SetInverseBendingMomentumAtDCA(t.fInverseBendingMomentum);
440 muTrack->SetThetaXAtDCA(t.fThetaX);
441 muTrack->SetThetaYAtDCA(t.fThetaY);
442 muTrack->SetBendingCoorAtDCA(t.fY);
443 muTrack->SetNonBendingCoorAtDCA(t.fX);
445 muTrack->SetInverseBendingMomentumUncorrected(t.fInverseBendingMomentum);
446 muTrack->SetThetaXUncorrected(t.fThetaX);
447 muTrack->SetThetaYUncorrected(t.fThetaY);
448 muTrack->SetZUncorrected(t.fZ);
449 muTrack->SetBendingCoorUncorrected(t.fY);
450 muTrack->SetNonBendingCoorUncorrected(t.fX);
452 muTrack->SetChi2(t.fChi2);
454 // Fill in the track hit points.
455 for (int i = 0; i < 16; i++)
457 if (not hitset[i]) continue;
459 AliHLTUInt8_t chamber;
460 AliHLTUInt16_t detElemId;
461 AliHLTMUONUtils::UnpackRecHitFlags(t.fHit[i].fFlags, chamber, detElemId);
463 AliESDMuonCluster *cluster = event.NewMuonCluster();
464 cluster->SetUniqueID(AliMUONVCluster::BuildUniqueID(chamber, detElemId, fClusterIndex++));
465 cluster->SetXYZ(t.fHit[i].fX, t.fHit[i].fY, t.fHit[i].fZ);
466 cluster->SetErrXY( // Use nominal values.
467 AliHLTMUONConstants::DefaultNonBendingReso(),
468 AliHLTMUONConstants::DefaultBendingReso()
470 cluster->SetCharge(-1.); // Indicate no total charge calculated.
471 cluster->SetChi2(-1.); // Indicate no fit made.
472 muTrack->AddClusterId(cluster->GetUniqueID());
473 muTrack->AddInMuonClusterMap(chamber);
476 FillTriggerInfo(triggerRecords, t.fTrigRec, t.fId, *muTrack, event);
480 if (fMakeClonesArray and event.GetList() != NULL)
483 event.GetList()->FindObject("MuonTracks"),
484 kAliHLTDataTypeTObject | kAliHLTDataOriginMUON,
488 event.GetList()->FindObject("MuonClusters"),
489 kAliHLTDataTypeTObject | kAliHLTDataOriginMUON,
493 if (fMakeESDDataBlock)
495 PushBack(&event, AliHLTMUONConstants::ESDDataType(), specification);
501 void AliHLTMUONESDMaker::FillTriggerInfo(
502 const AliTriggerRecordList& triggerRecords,
503 AliHLTInt32_t trigRecId, AliHLTInt32_t trackId,
504 AliESDMuonTrack& muTrack, AliESDEvent& event
507 /// Finds the trigger record with ID = 'id' and fills the output ESD track structure.
509 // Find the corresponding trigger record.
510 const AliHLTMUONTriggerRecordStruct* trigrec = NULL;
511 for (size_t k = 0; k < triggerRecords.size(); k++)
513 if (triggerRecords[k]->fId == trigRecId)
515 trigrec = triggerRecords[k];
519 // If the trigger record was found then fill its hit information also.
522 AliHLTMUONParticleSign trsign;
524 AliHLTMUONUtils::UnpackTriggerRecordFlags(
525 trigrec->fFlags, trsign, trhitset
528 for (int i = 0; i < 4; i++)
530 if (not trhitset[i]) continue;
532 AliHLTUInt8_t chamber;
533 AliHLTUInt16_t detElemId;
534 AliHLTMUONUtils::UnpackRecHitFlags(trigrec->fHit[i].fFlags, chamber, detElemId);
536 AliESDMuonCluster *cluster = event.NewMuonCluster();
537 cluster->SetUniqueID(AliMUONVCluster::BuildUniqueID(chamber, detElemId, fClusterIndex++));
543 cluster->SetErrXY( // Use nominal values.
544 AliMUONConstants::TriggerNonBendingReso(),
545 AliMUONConstants::TriggerBendingReso()
547 cluster->SetCharge(-1.); // Indicate no total charge calculated.
548 cluster->SetChi2(-1.); // Indicate no fit made.
549 muTrack.AddClusterId(cluster->GetUniqueID());
550 muTrack.AddInMuonClusterMap(i+10);
555 HLTWarning("Trigger record (ID = %d) not found for track ID = %d.",