]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCEsdWriterComponent.cxx
- adding TPC-mcTrackMarker configuration (propagation of MC information
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCEsdWriterComponent.cxx
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   AliHLTTPCEsdWriterComponent.cxx
20     @author Matthias Richter
21     @date   
22     @brief  Writer component to store tracks of the HLT TPC conformal
23             mapping tracker in the AliESD format
24 */
25
26 // see header file for class documentation
27 // or
28 // refer to README to build package
29 // or
30 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
31
32 #include <cassert>
33 #include "AliHLTTPCEsdWriterComponent.h"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "AliCDBEntry.h"
37 #include "AliCDBManager.h"
38 #include "TTree.h"
39 #include "TList.h"
40 #include "AliHLTTPCTrack.h"
41 #include "AliHLTTPCTrackArray.h"
42 #include "AliHLTTPCTrackletDataFormat.h"
43 #include "AliHLTTPCDefinitions.h"
44 #include "AliHLTTPCTransform.h"
45 #include "AliHLTExternalTrackParam.h"
46 #include "AliHLTKalmanTrack.h"
47 #include "AliHLTTrackMCLabel.h"
48
49 #include <vector>
50
51 /** ROOT macro for the implementation of ROOT specific class methods */
52 ClassImp(AliHLTTPCEsdWriterComponent)
53
54 AliHLTTPCEsdWriterComponent::AliHLTTPCEsdWriterComponent()
55   :
56   fSolenoidBz(0)
57 {
58   // see header file for class documentation
59   // or
60   // refer to README to build package
61   // or
62   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
63 }
64
65 AliHLTTPCEsdWriterComponent::~AliHLTTPCEsdWriterComponent()
66 {
67   // see header file for class documentation
68 }
69
70 AliHLTTPCEsdWriterComponent::AliWriter::AliWriter()
71   :
72   fTree(NULL),
73   fESD(NULL),
74   fBase(new AliHLTTPCEsdWriterComponent)
75 {
76   // see header file for class documentation
77   // or
78   // refer to README to build package
79   // or
80   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
81 }
82
83 AliHLTTPCEsdWriterComponent::AliWriter::~AliWriter()
84 {
85   // see header file for class documentation
86   if (fBase) delete fBase;
87   fBase=NULL;
88 }
89
90 void AliHLTTPCEsdWriterComponent::AliWriter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
91 {
92   // see header file for class documentation
93   list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
94   list.push_back(AliHLTTPCDefinitions::fgkTracksDataType);
95   list.push_back( kAliHLTDataTypeTrack|kAliHLTDataOriginTPC );
96 }
97
98 int AliHLTTPCEsdWriterComponent::AliWriter::InitWriter()
99 {
100   // see header file for class documentation
101   int iResult=0;
102   fESD = new AliESDEvent;
103   if (fESD) {
104     fESD->CreateStdContent();
105     fTree = new TTree("esdTree", "Tree with HLT ESD objects");
106     if (fTree) {
107       fTree->SetDirectory(0);
108       fESD->WriteToTree(fTree);
109     }
110   }
111   if (fTree==NULL) {
112     iResult=-ENOMEM;
113   }
114
115   if (iResult>=0) {
116     iResult=fBase->Reconfigure(NULL, NULL);
117   }
118
119   return iResult;
120 }
121
122 int AliHLTTPCEsdWriterComponent::AliWriter::CloseWriter()
123 {
124   // see header file for class documentation
125   int iResult=0;
126   if (fTree) {
127     // the esd structure is written to the user info and is
128     // needed in te ReadFromTree method to read all objects correctly
129     if (fESD) fTree->GetUserInfo()->Add(fESD);
130     WriteObject(kAliHLTVoidEventID, fTree);
131     fTree->GetUserInfo()->Clear();
132     TTree* pTree=fTree;
133     fTree=NULL;
134     delete pTree;
135   } else {
136     HLTWarning("not initialized");
137   }
138
139   if (fESD) {
140     delete fESD;
141   }
142   iResult=AliHLTRootFileWriterComponent::CloseWriter();
143   return iResult;
144 }
145
146 int AliHLTTPCEsdWriterComponent::AliWriter::DumpEvent( const AliHLTComponentEventData& evtData,
147                                             const AliHLTComponentBlockData* blocks, 
148                                             AliHLTComponentTriggerData& /*trigData*/ )
149 {
150   // see header file for class documentation
151   int iResult=0;
152   TTree* pTree=fTree;
153   assert(fBase);
154   if (pTree && fBase) {
155     if (fESD) {
156       AliESDEvent* pESD=fESD;
157
158       iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt);
159
160     } else {
161       iResult=-ENOMEM;
162     }
163   }
164   return iResult;
165 }
166
167 int AliHLTTPCEsdWriterComponent::AliWriter::ScanArgument(int argc, const char** argv)
168 {
169   // see header file for class documentation
170   int iResult=AliHLTRootFileWriterComponent::ScanArgument(argc, argv);
171   return iResult;
172 }
173
174 int AliHLTTPCEsdWriterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD,
175                                                const AliHLTComponentBlockData* blocks,
176                                                int nBlocks, int* pMinSlice,
177                                                int* pMaxSlice)
178 {
179   // see header file for class documentation
180
181   int iResult=0;
182   int iAddedDataBlocks=0;
183   
184   if (pESD && blocks) {
185       pESD->Reset(); 
186       pESD->SetMagneticField(fSolenoidBz);
187       const AliHLTComponentBlockData* iter = NULL;
188       int bIsTrackSegs=0;
189
190       // first read MC information (if present)
191
192       std::map<int,int> mcLabels;
193
194       for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
195         iter = blocks+ndx;
196         if(iter->fDataType == (kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC) ) {
197           AliHLTTrackMCData* dataPtr = ( AliHLTTrackMCData* )( iter->fPtr );
198           for( unsigned int il=0; il<dataPtr->fCount; il++ ){
199             AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
200             mcLabels[lab.fTrackID] = lab.fMCLabel;
201           }
202         }
203       }
204
205
206       // do the conversion of tracks
207
208       for (int ndx=0; ndx<nBlocks && iResult>=0; ndx++) {
209         iter = blocks+ndx;
210         
211         if( iter->fDataType == ( kAliHLTDataTypeTrack|kAliHLTDataOriginTPC ) ){   
212           AliHLTTracksData* dataPtr = ( AliHLTTracksData* ) iter->fPtr;
213           int nTracks = dataPtr->fCount;
214           AliHLTExternalTrackParam* currOutTrack = dataPtr->fTracklets;
215
216           for( int itr=0; itr<nTracks; itr++ ){
217             AliHLTKalmanTrack t(*currOutTrack);
218             Float_t points[4] = {currOutTrack->fX, currOutTrack->fY, currOutTrack->fLastX, currOutTrack->fLastY };
219             
220             Int_t mcLabel = -1;
221             if( mcLabels.find(currOutTrack->fTrackID)!=mcLabels.end() )
222               mcLabel = mcLabels[currOutTrack->fTrackID];
223             
224             t.SetLabel( mcLabel );
225             
226             AliESDtrack iotrack;
227             iotrack.UpdateTrackParams( &t,AliESDtrack::kTPCin);
228             iotrack.SetTPCPoints(points);
229             pESD->AddTrack(&iotrack);
230             unsigned int dSize = sizeof( AliHLTExternalTrackParam ) + currOutTrack->fNPoints * sizeof( unsigned int );
231             currOutTrack = ( AliHLTExternalTrackParam* )( (( Byte_t * )currOutTrack) + dSize );
232           }
233           iAddedDataBlocks++;
234         }
235
236         
237         if ( (bIsTrackSegs=(iter->fDataType == AliHLTTPCDefinitions::fgkTrackSegmentsDataType))==1 ||
238              iter->fDataType == AliHLTTPCDefinitions::fgkTracksDataType ) {
239           Int_t minslice=AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
240           Int_t maxslice=AliHLTTPCDefinitions::GetMaxSliceNr(iter->fSpecification);
241           if (bIsTrackSegs==0) {
242             // slice parameter and data specification ignored, tracks already in global coordinates
243             minslice=-1;
244             maxslice=-1;
245             if (pMinSlice) *pMinSlice=0;
246             if (pMaxSlice) *pMaxSlice=AliHLTTPCTransform::GetNSlice()-1;
247           } else {
248             if (pMinSlice && (*pMinSlice==-1 || *pMinSlice>minslice)) *pMinSlice=minslice;
249             if (pMaxSlice && (*pMaxSlice==-1 || *pMaxSlice<maxslice)) *pMaxSlice=maxslice;
250           }
251           //HLTDebug("dataspec %#x minslice %d", iter->fSpecification, minslice);
252           if (minslice >=-1 && minslice<AliHLTTPCTransform::GetNSlice()) {
253             if (minslice!=maxslice) {
254               HLTWarning("data from multiple sectors in one block: "
255                          "possible mismatch in treatment of local coordinate system");
256             }
257             AliHLTTPCTrackArray tracks;
258             AliHLTTPCTrackletData* inPtr = (AliHLTTPCTrackletData*) iter->fPtr;     
259             HLTDebug("reading block %d (slice %d): %d tracklets", ndx, minslice, inPtr->fTrackletCnt);
260             if ((iResult=tracks.FillTracksChecked(inPtr->fTracklets, inPtr->fTrackletCnt, iter->fSize, minslice, 0/*don't rotate*/))>=0) {
261               if ((iResult=Tracks2ESD(&tracks, pESD ))>=0) {
262                 iAddedDataBlocks++;
263               }
264             }
265           } else {
266             HLTError("invalid sector number");
267             iResult=-EBADF;
268           }
269         }      
270       }
271
272       if (iAddedDataBlocks>0 && pTree) {
273         pTree->Fill();
274       }
275   
276   } else {
277     iResult=-EINVAL;
278   }
279   if (iResult>=0) iResult=iAddedDataBlocks;
280   return iResult;
281 }
282
283
284 int AliHLTTPCEsdWriterComponent::Tracks2ESD(AliHLTTPCTrackArray* pTracks, AliESDEvent* pESD )
285 {
286   // see header file for class documentation
287   int iResult=0;
288
289   if (pTracks && pESD) {    
290  
291     for (int i=0; i<pTracks->GetNTracks() && iResult>=0; i++) {
292       AliHLTTPCTrack* pTrack=(*pTracks)[i];
293       if (pTrack) {
294         
295         if( pTrack->Convert2AliKalmanTrack() ){   
296           HLTError("conversion to AliKalmanTrack failed for track %d of %d", i, pTracks->GetNTracks()); 
297           continue;
298         }
299
300         Float_t points[4] = {pTrack->GetFirstPointX(), pTrack->GetFirstPointY(), pTrack->GetLastPointX(), pTrack->GetLastPointY() };
301
302         if(pTrack->GetSector() == -1){ // Set first and last points for global tracks
303           Double_t s = TMath::Sin( pTrack->GetAlpha() );
304           Double_t c = TMath::Cos( pTrack->GetAlpha() );
305           points[0] =  pTrack->GetFirstPointX()*c + pTrack->GetFirstPointY()*s;
306           points[1] = -pTrack->GetFirstPointX()*s + pTrack->GetFirstPointY()*c;   
307           points[2] =  pTrack->GetLastPointX() *c + pTrack->GetLastPointY() *s;
308           points[3] = -pTrack->GetLastPointX() *s + pTrack->GetLastPointY() *c;   
309         }
310
311         AliESDtrack iotrack;
312         iotrack.UpdateTrackParams(pTrack,AliESDtrack::kTPCin);
313         iotrack.SetTPCPoints(points);
314
315         pESD->AddTrack(&iotrack);
316       } else {
317         HLTError("internal mismatch in array");
318         iResult=-EFAULT;
319       }
320     }
321     
322   } else {
323     iResult=-EINVAL;
324   }
325   return iResult;
326 }
327
328 int AliHLTTPCEsdWriterComponent::Configure(const char* arguments)
329 {
330   // see header file for class documentation
331   int iResult=0;
332   if (!arguments) return iResult;
333
334   TString allArgs=arguments;
335   TString argument;
336   int bMissingParam=0;
337
338   TObjArray* pTokens=allArgs.Tokenize(" ");
339   if (pTokens) {
340     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
341       argument=((TObjString*)pTokens->At(i))->GetString();
342       if (argument.IsNull()) continue;
343       
344       if (argument.CompareTo("-solenoidBz")==0) {
345         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
346         HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
347         fSolenoidBz=((TObjString*)pTokens->At(i))->GetString().Atof();
348         continue;
349       } else {
350         HLTError("unknown argument %s", argument.Data());
351         iResult=-EINVAL;
352         break;
353       }
354     }
355     delete pTokens;
356   }
357   if (bMissingParam) {
358     HLTError("missing parameter for argument %s", argument.Data());
359     iResult=-EINVAL;
360   }
361
362   return iResult;
363 }
364
365 int AliHLTTPCEsdWriterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
366 {
367   // see header file for class documentation
368   int iResult=0;
369   const char* path=kAliHLTCDBSolenoidBz;
370   const char* defaultNotify="";
371   if (cdbEntry) {
372     path=cdbEntry;
373     defaultNotify=" (default)";
374   }
375   if (path) {
376     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
377     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
378     if (pEntry) {
379       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
380       if (pString) {
381         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
382         iResult=Configure(pString->GetString().Data());
383       } else {
384         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
385       }
386     } else {
387       HLTError("can not fetch object \"%s\" from CDB", path);
388     }
389   }
390   
391   return iResult;
392 }
393
394 AliHLTTPCEsdWriterComponent::AliConverter::AliConverter()
395   :
396   fESD(NULL),
397   fBase(new AliHLTTPCEsdWriterComponent),
398   fWriteTree(0)
399 {
400   // see header file for class documentation
401   // or
402   // refer to README to build package
403   // or
404   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
405 }
406
407 AliHLTTPCEsdWriterComponent::AliConverter::~AliConverter()
408 {
409   // see header file for class documentation
410   if (fBase) delete fBase;
411   fBase=NULL;
412
413   if (fESD) delete fESD;
414   fESD=NULL;
415 }
416
417 void AliHLTTPCEsdWriterComponent::AliConverter::GetInputDataTypes(AliHLTComponentDataTypeList& list)
418 {
419   // see header file for class documentation
420   list.push_back(AliHLTTPCDefinitions::fgkTrackSegmentsDataType);
421   list.push_back(AliHLTTPCDefinitions::fgkTracksDataType);
422   list.push_back( kAliHLTDataTypeTrack|kAliHLTDataOriginTPC );
423 }
424
425 AliHLTComponentDataType AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataType()
426 {
427   // see header file for class documentation
428   return kAliHLTDataTypeESDTree;
429 }
430
431 void AliHLTTPCEsdWriterComponent::AliConverter::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
432 {
433   // see header file for class documentation
434   constBase=2000000;
435   inputMultiplier=10.0;
436 }
437
438 int AliHLTTPCEsdWriterComponent::AliConverter::DoInit(int argc, const char** argv)
439 {
440   // see header file for class documentation
441   int iResult=0;
442   TString argument="";
443   int bMissingParam=0;
444   for (int i=0; i<argc && iResult>=0; i++) {
445     argument=argv[i];
446     if (argument.IsNull()) continue;
447
448     // -notree
449     if (argument.CompareTo("-notree")==0) {
450       fWriteTree=0;
451
452       // -tree
453     } else if (argument.CompareTo("-tree")==0) {
454       fWriteTree=1;
455
456       // -solenoidBz
457     } else if (argument.CompareTo("-solenoidBz")==0) {
458       TString tmp="-solenoidBz ";
459       tmp+=argv[++i];
460       fBase->Configure(tmp.Data());
461     } else {
462       HLTError("unknown argument %s", argument.Data());
463       break;
464     }
465   }
466   if (bMissingParam) {
467     HLTError("missing parameter for argument %s", argument.Data());
468     iResult=-EINVAL;
469   }
470
471   if (iResult>=0) {
472     iResult=fBase->Reconfigure(NULL, NULL);
473   }
474
475   return iResult;
476 }
477
478 int AliHLTTPCEsdWriterComponent::AliConverter::DoDeinit()
479 {
480   // see header file for class documentation
481   return 0;
482 }
483
484 int AliHLTTPCEsdWriterComponent::AliConverter::DoEvent(const AliHLTComponentEventData& evtData, 
485                                                        const AliHLTComponentBlockData* blocks, 
486                                                        AliHLTComponentTriggerData& /*trigData*/,
487                                                        AliHLTUInt8_t* /*outputPtr*/, 
488                                                        AliHLTUInt32_t& size,
489                                                        AliHLTComponentBlockDataList& /*outputBlocks*/ )
490 {
491   // see header file for class documentation
492   int iResult=0;
493   // no direct writing to the output buffer
494   size=0;
495
496   assert(fBase);
497   if (!fESD) {
498     fESD = new AliESDEvent;
499     if (fESD) {
500       fESD->CreateStdContent();
501     } else {
502       iResult=-ENOMEM;
503     }
504   }
505
506   AliESDEvent* pESD = fESD;
507
508   if (pESD && fBase) {
509   
510     TTree* pTree = NULL;
511     // TODO: Matthias 06.12.2007
512     // Tried to write the ESD directly instead to a tree, but this did not work
513     // out. Information in the ESD is different, needs investigation.
514     
515     if (fWriteTree)
516       pTree = new TTree("esdTree", "Tree with HLT ESD objects");
517  
518     if (pTree) {
519       pTree->SetDirectory(0);
520     }
521
522     if ((iResult=fBase->ProcessBlocks(pTree, pESD, blocks, (int)evtData.fBlockCnt))>0) {
523       // TODO: set the specification correctly
524       if (pTree) {
525         // the esd structure is written to the user info and is
526         // needed in te ReadFromTree method to read all objects correctly
527         pTree->GetUserInfo()->Add(pESD);
528         pESD->WriteToTree(pTree);
529         iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginTPC, 0);
530       } else {
531         iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginTPC, 0);
532       }
533     }
534     if (pTree) {
535       // clear user info list to prevent objects from being deleted
536       pTree->GetUserInfo()->Clear();
537       delete pTree;
538     }
539   }
540   return iResult;
541 }
542