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 AliHLTTPCEsdWriterComponent.cxx
20 @author Matthias Richter
22 @brief Writer component to store tracks of the HLT TPC conformal
23 mapping tracker in the AliESD format
26 // see header file for class documentation
28 // refer to README to build package
30 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
33 #include "AliHLTTPCEsdWriterComponent.h"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "AliCDBEntry.h"
37 #include "AliCDBManager.h"
40 #include "AliHLTTPCTrack.h"
41 #include "AliHLTTPCTrackArray.h"
42 #include "AliHLTTPCTrackletDataFormat.h"
43 #include "AliHLTTPCDefinitions.h"
44 #include "AliHLTTPCTransform.h"
46 /** ROOT macro for the implementation of ROOT specific class methods */
47 ClassImp(AliHLTTPCEsdWriterComponent)
49 AliHLTTPCEsdWriterComponent::AliHLTTPCEsdWriterComponent()
53 // see header file for class documentation
55 // refer to README to build package
57 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
60 AliHLTTPCEsdWriterComponent::~AliHLTTPCEsdWriterComponent()
62 // see header file for class documentation
65 AliHLTTPCEsdWriterComponent::AliWriter::AliWriter()
69 fBase(new AliHLTTPCEsdWriterComponent)
71 // see header file for class documentation
73 // refer to README to build package
75 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
78 AliHLTTPCEsdWriterComponent::AliWriter::~AliWriter()
80 // see header file for class documentation
81 if (fBase) delete fBase;
85 void AliHLTTPCEsdWriterComponent::AliWriter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
87 // see header file for class documentation
88 list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
89 list.push_back(AliHLTTPCDefinitions::fgkTracksDataType);
92 int AliHLTTPCEsdWriterComponent::AliWriter::InitWriter()
94 // see header file for class documentation
96 fESD = new AliESDEvent;
98 fESD->CreateStdContent();
99 fTree = new TTree("esdTree", "Tree with HLT ESD objects");
101 fTree->SetDirectory(0);
102 fESD->WriteToTree(fTree);
110 iResult=fBase->Reconfigure(NULL, NULL);
116 int AliHLTTPCEsdWriterComponent::AliWriter::CloseWriter()
118 // see header file for class documentation
121 // the esd structure is written to the user info and is
122 // needed in te ReadFromTree method to read all objects correctly
123 if (fESD) fTree->GetUserInfo()->Add(fESD);
124 WriteObject(kAliHLTVoidEventID, fTree);
125 fTree->GetUserInfo()->Clear();
130 HLTWarning("not initialized");
136 iResult=AliHLTRootFileWriterComponent::CloseWriter();
140 int AliHLTTPCEsdWriterComponent::AliWriter::DumpEvent( const AliHLTComponentEventData& evtData,
141 const AliHLTComponentBlockData* blocks,
142 AliHLTComponentTriggerData& /*trigData*/ )
144 // see header file for class documentation
148 if (pTree && fBase) {
150 AliESDEvent* pESD=fESD;
152 iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt);
161 int AliHLTTPCEsdWriterComponent::AliWriter::ScanArgument(int argc, const char** argv)
163 // see header file for class documentation
164 int iResult=AliHLTRootFileWriterComponent::ScanArgument(argc, argv);
168 int AliHLTTPCEsdWriterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD,
169 const AliHLTComponentBlockData* blocks,
170 int nBlocks, int* pMinSlice,
173 // see header file for class documentation
175 if (pESD && blocks) {
176 pESD->SetMagneticField(fSolenoidBz);
177 const AliHLTComponentBlockData* iter = NULL;
178 AliHLTTPCTrackletData* inPtr=NULL;
181 for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
183 if ( (bIsTrackSegs=(iter->fDataType == AliHLTTPCDefinitions::fgkTrackSegmentsDataType))==1 ||
184 iter->fDataType == AliHLTTPCDefinitions::fgkTracksDataType ) {
185 Int_t minslice=AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
186 Int_t maxslice=AliHLTTPCDefinitions::GetMaxSliceNr(iter->fSpecification);
187 if (bIsTrackSegs==0) {
188 // slice parameter and data specification ignored, tracks already in global coordinates
191 if (pMinSlice) *pMinSlice=0;
192 if (pMaxSlice) *pMaxSlice=AliHLTTPCTransform::GetNSlice()-1;
194 if (pMinSlice && (*pMinSlice==-1 || *pMinSlice>minslice)) *pMinSlice=minslice;
195 if (pMaxSlice && (*pMaxSlice==-1 || *pMaxSlice<maxslice)) *pMaxSlice=maxslice;
197 //HLTDebug("dataspec %#x minslice %d", iter->fSpecification, minslice);
198 if (minslice >=-1 && minslice<AliHLTTPCTransform::GetNSlice()) {
199 if (minslice!=maxslice) {
200 HLTWarning("data from multiple sectors in one block: "
201 "possible mismatch in treatment of local coordinate system");
203 AliHLTTPCTrackArray tracks;
204 inPtr=(AliHLTTPCTrackletData*)iter->fPtr;
205 HLTDebug("reading block %d (slice %d): %d tracklets", ndx, minslice, inPtr->fTrackletCnt);
206 tracks.FillTracks(inPtr->fTrackletCnt, inPtr->fTracklets, minslice, 0/*don't rotate*/);
207 if ((iResult=Tracks2ESD(&tracks, pESD))>=0) {
210 HLTError("invalid sector number");
215 if (iResult>=0 && pTree) {
227 int AliHLTTPCEsdWriterComponent::Tracks2ESD(AliHLTTPCTrackArray* pTracks, AliESDEvent* pESD)
229 // see header file for class documentation
231 if (pTracks && pESD) {
232 HLTDebug("converting %d tracks from track array", pTracks->GetNTracks());
233 for (int i=0; i<pTracks->GetNTracks() && iResult>=0; i++) {
234 AliHLTTPCTrack* pTrack=(*pTracks)[i];
236 //HLTDebug("convert track %d", i);
238 int iLocal=pTrack->Convert2AliKalmanTrack();
241 iotrack.UpdateTrackParams(pTrack,AliESDtrack::kTPCin);
242 iotrack.SetTPCPoints(pTrack->GetPoints());
243 pESD->AddTrack(&iotrack);
245 HLTError("conversion to AliKalmanTrack failed for track %d of %d", i, pTracks->GetNTracks());
248 HLTError("internal mismatch in array");
259 int AliHLTTPCEsdWriterComponent::Configure(const char* arguments)
261 // see header file for class documentation
263 if (!arguments) return iResult;
265 TString allArgs=arguments;
269 TObjArray* pTokens=allArgs.Tokenize(" ");
271 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
272 argument=((TObjString*)pTokens->At(i))->GetString();
273 if (argument.IsNull()) continue;
275 if (argument.CompareTo("-solenoidBz")==0) {
276 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
277 HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
278 fSolenoidBz=((TObjString*)pTokens->At(i))->GetString().Atof();
281 HLTError("unknown argument %s", argument.Data());
289 HLTError("missing parameter for argument %s", argument.Data());
296 int AliHLTTPCEsdWriterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
298 // see header file for class documentation
300 const char* path="HLT/ConfigHLT/SolenoidBz";
301 const char* defaultNotify="";
304 defaultNotify=" (default)";
307 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
308 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
310 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
312 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
313 iResult=Configure(pString->GetString().Data());
315 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
318 HLTError("can not fetch object \"%s\" from CDB", path);
325 AliHLTTPCEsdWriterComponent::AliConverter::AliConverter()
327 fBase(new AliHLTTPCEsdWriterComponent),
330 // see header file for class documentation
332 // refer to README to build package
334 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
337 AliHLTTPCEsdWriterComponent::AliConverter::~AliConverter()
339 // see header file for class documentation
340 if (fBase) delete fBase;
344 void AliHLTTPCEsdWriterComponent::AliConverter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
346 // see header file for class documentation
347 list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
348 list.push_back(AliHLTTPCDefinitions::fgkTracksDataType);
351 AliHLTComponentDataType AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataType()
353 // see header file for class documentation
354 return kAliHLTDataTypeESDTree;
357 void AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
359 // see header file for class documentation
361 inputMultiplier=10.0;
364 int AliHLTTPCEsdWriterComponent::AliConverter::DoInit(int argc, const char** argv)
366 // see header file for class documentation
370 for (int i=0; i<argc && iResult>=0; i++) {
372 if (argument.IsNull()) continue;
375 if (argument.CompareTo("-notree")==0) {
379 } else if (argument.CompareTo("-tree")==0) {
383 HLTError("unknown argument %s", argument.Data());
388 HLTError("missing parameter for argument %s", argument.Data());
393 iResult=fBase->Reconfigure(NULL, NULL);
399 int AliHLTTPCEsdWriterComponent::AliConverter::DoDeinit()
401 // see header file for class documentation
405 int AliHLTTPCEsdWriterComponent::AliConverter::DoEvent(const AliHLTComponentEventData& evtData,
406 const AliHLTComponentBlockData* blocks,
407 AliHLTComponentTriggerData& /*trigData*/,
408 AliHLTUInt8_t* /*outputPtr*/,
409 AliHLTUInt32_t& /*size*/,
410 AliHLTComponentBlockDataList& /*outputBlocks*/ )
412 // see header file for class documentation
415 AliESDEvent* pESD = new AliESDEvent;
417 pESD->CreateStdContent();
419 // TODO: Matthias 06.12.2007
420 // Tried to write the ESD directly instead to a tree, but this did not work
421 // out. Information in the ESD is different, needs investigation.
423 pTree = new TTree("esdTree", "Tree with HLT ESD objects");
425 pTree->SetDirectory(0);
426 pESD->WriteToTree(pTree);
429 if ((iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt))>=0) {
430 // TODO: set the specification correctly
432 // the esd structure is written to the user info and is
433 // needed in te ReadFromTree method to read all objects correctly
434 pTree->GetUserInfo()->Add(pESD);
435 iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginTPC, 0);
437 iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginTPC, 0);
441 // clear user info list to prevent objects from being deleted
442 pTree->GetUserInfo()->Clear();