]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalEsdConverterComponent.cxx
Adding component to generate histograms for clusters.
[u/mrichter/AliRoot.git] / HLT / global / AliHLTGlobalEsdConverterComponent.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   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 "AliESDEvent.h"
38 #include "AliESDtrack.h"
39 #include "AliESDMuonTrack.h"
40 #include "AliCDBEntry.h"
41 #include "AliCDBManager.h"
42 #include "AliPID.h"
43 #include "TTree.h"
44 #include "TList.h"
45 #include "TClonesArray.h"
46 #include "AliHLTESDCaloClusterMaker.h"
47 #include "AliHLTCaloClusterDataStruct.h"
48 #include "AliHLTCaloClusterReader.h"
49 #include "AliESDCaloCluster.h"
50
51 /** ROOT macro for the implementation of ROOT specific class methods */
52 ClassImp(AliHLTGlobalEsdConverterComponent)
53
54 AliHLTGlobalEsdConverterComponent::AliHLTGlobalEsdConverterComponent()
55   : AliHLTProcessor()
56   , fWriteTree(0)
57   , fVerbosity(0)
58   , fESD(NULL)
59   , fSolenoidBz(-5.00668)
60 {
61   // see header file for class documentation
62   // or
63   // refer to README to build package
64   // or
65   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
66 }
67
68 AliHLTGlobalEsdConverterComponent::~AliHLTGlobalEsdConverterComponent()
69 {
70   // see header file for class documentation
71   if (fESD) delete fESD;
72   fESD=NULL;
73 }
74
75 int AliHLTGlobalEsdConverterComponent::Configure(const char* arguments)
76 {
77   // see header file for class documentation
78   int iResult=0;
79   if (!arguments) return iResult;
80
81   TString allArgs=arguments;
82   TString argument;
83   int bMissingParam=0;
84
85   TObjArray* pTokens=allArgs.Tokenize(" ");
86   if (pTokens) {
87     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
88       argument=((TObjString*)pTokens->At(i))->GetString();      
89       if (argument.IsNull()) continue;
90       
91       if (argument.CompareTo("-solenoidBz")==0) {
92         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
93         HLTWarning("argument -solenoidBz is deprecated, magnetic field set up globally (%f)", GetBz());
94         continue;
95       } else {
96         HLTError("unknown argument %s", argument.Data());
97         iResult=-EINVAL;
98         break;
99       }
100     }
101     delete pTokens;
102   }
103   if (bMissingParam) {
104     HLTError("missing parameter for argument %s", argument.Data());
105     iResult=-EINVAL;
106   }
107
108   return iResult;
109 }
110
111 int AliHLTGlobalEsdConverterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
112 {
113   // see header file for class documentation
114   int iResult=0;
115   const char* path=NULL;
116   const char* defaultNotify="";
117   if (cdbEntry) {
118     path=cdbEntry;
119     defaultNotify=" (default)";
120   }
121   if (path) {
122     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
123     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
124     if (pEntry) {
125       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
126       if (pString) {
127         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
128         iResult=Configure(pString->GetString().Data());
129       } else {
130         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
131       }
132     } else {
133       HLTError("can not fetch object \"%s\" from CDB", path);
134     }
135   }
136   
137   return iResult;
138 }
139
140 void AliHLTGlobalEsdConverterComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
141 {
142   // see header file for class documentation
143   list.push_back(kAliHLTDataTypeTrack);
144   list.push_back(kAliHLTDataTypeTrackMC);
145   list.push_back(kAliHLTDataTypeCaloCluster);
146   list.push_back(kAliHLTDataTypedEdx );
147   list.push_back(kAliHLTDataTypeESDVertex );
148   list.push_back(kAliHLTDataTypeESDObject);
149   list.push_back(kAliHLTDataTypeTObject);
150 }
151
152 AliHLTComponentDataType AliHLTGlobalEsdConverterComponent::GetOutputDataType()
153 {
154   // see header file for class documentation
155   return kAliHLTDataTypeESDObject|kAliHLTDataOriginOut;
156 }
157
158 void AliHLTGlobalEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
159 {
160   // see header file for class documentation
161   constBase=2000000;
162   inputMultiplier=10.0;
163 }
164
165 int AliHLTGlobalEsdConverterComponent::DoInit(int argc, const char** argv)
166 {
167   // see header file for class documentation
168   int iResult=0;
169   TString argument="";
170   int bMissingParam=0;
171
172   iResult=Reconfigure(NULL, NULL);
173   TString allArgs = "";
174   for ( int i = 0; i < argc; i++ ) {
175     if ( !allArgs.IsNull() ) allArgs += " ";
176     allArgs += argv[i];
177   }
178
179   TObjArray* pTokens=allArgs.Tokenize(" ");
180   if (pTokens) {
181     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
182       argument=((TObjString*)pTokens->At(i))->GetString();      
183       if (argument.IsNull()) continue;
184
185       // -notree
186       if (argument.CompareTo("-notree")==0) {
187         fWriteTree=0;
188         
189         // -tree
190       } else if (argument.CompareTo("-tree")==0) {
191         fWriteTree=1;
192       } else if (argument.CompareTo("-solenoidBz")==0) {
193         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
194         HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
195         fSolenoidBz=((TObjString*)pTokens->At(i))->GetString().Atof();
196         continue;
197       } else {
198         HLTError("unknown argument %s", argument.Data());
199         break;
200       }
201     }
202   }
203   if (bMissingParam) {
204     HLTError("missing parameter for argument %s", argument.Data());
205     iResult=-EINVAL;
206   }
207
208   fSolenoidBz=GetBz();
209
210   if (iResult>=0) {
211     fESD = new AliESDEvent;
212     if (fESD) {
213       fESD->CreateStdContent();
214     } else {
215       iResult=-ENOMEM;
216     }
217
218     SetupCTPData();
219   }
220
221   return iResult;
222 }
223
224 int AliHLTGlobalEsdConverterComponent::DoDeinit()
225 {
226   // see header file for class documentation
227   if (fESD) delete fESD;
228   fESD=NULL;
229
230   return 0;
231 }
232
233 int AliHLTGlobalEsdConverterComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, 
234                                                AliHLTComponentTriggerData& trigData)
235 {
236   // see header file for class documentation
237   int iResult=0;
238   if (!fESD) return -ENODEV;
239
240   AliESDEvent* pESD = fESD;
241
242   pESD->Reset(); 
243   pESD->SetMagneticField(fSolenoidBz);
244   pESD->SetRunNumber(GetRunNo());
245   pESD->SetPeriodNumber(GetPeriodNumber());
246   pESD->SetOrbitNumber(GetOrbitNumber());
247   pESD->SetBunchCrossNumber(GetBunchCrossNumber());
248   pESD->SetTimeStamp(GetTimeStamp());
249
250   const AliHLTCTPData* pCTPData=CTPData();
251   if (pCTPData) {
252     AliHLTUInt64_t mask=pCTPData->ActiveTriggers(trigData);
253     for (int index=0; index<gkNCTPTriggerClasses; index++) {
254       if ((mask&((AliHLTUInt64_t)0x1<<index)) == 0) continue;
255       pESD->SetTriggerClass(pCTPData->Name(index), index);
256     }
257     pESD->SetTriggerMask(mask);
258   }
259
260   TTree* pTree = NULL;
261   if (fWriteTree)
262     pTree = new TTree("esdTree", "Tree with HLT ESD objects");
263  
264   if (pTree) {
265     pTree->SetDirectory(0);
266   }
267
268   if ((iResult=ProcessBlocks(pTree, pESD))>=0) {
269     // TODO: set the specification correctly
270     if (pTree) {
271       // the esd structure is written to the user info and is
272       // needed in te ReadFromTree method to read all objects correctly
273       pTree->GetUserInfo()->Add(pESD);
274       pESD->WriteToTree(pTree);
275       iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginOut, 0);
276     } else {
277       iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
278     }
279   }
280   if (pTree) {
281     // clear user info list to prevent objects from being deleted
282     pTree->GetUserInfo()->Clear();
283     delete pTree;
284   }
285   return iResult;
286 }
287
288 int AliHLTGlobalEsdConverterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD)
289 {
290   // see header file for class documentation
291
292   int iResult=0;
293   int iAddedDataBlocks=0;
294
295   // Barrel tracking
296
297   // in the first attempt this component reads the TPC tracks and updates in the
298   // second step from the ITS tracks
299
300
301   // first read MC information (if present)
302   std::map<int,int> mcLabels;
303
304   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
305        pBlock!=NULL; pBlock=GetNextInputBlock()) {
306     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
307     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
308       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
309         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
310         mcLabels[lab.fTrackID] = lab.fMCLabel;
311       }
312     } else {
313       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
314                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
315                  dataPtr->fCount, pBlock->fSize);
316     }
317   }
318
319   // read dEdx information (if present)
320
321   AliHLTFloat32_t *dEdxTPC = 0; 
322   Int_t ndEdxTPC = 0;
323   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypedEdx|kAliHLTDataOriginTPC);
324        pBlock!=NULL; pBlock=GetNextInputBlock()) {
325     dEdxTPC = reinterpret_cast<AliHLTFloat32_t*>( pBlock->fPtr );
326     ndEdxTPC = pBlock->fSize / sizeof(AliHLTFloat32_t);
327     break;
328   }
329
330   // convert the TPC tracks to ESD tracks
331   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
332        pBlock!=NULL; pBlock=GetNextInputBlock()) {
333     vector<AliHLTGlobalBarrelTrack> tracks;
334     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
335       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
336            element!=tracks.end(); element++) {
337         Float_t points[4] = {
338           element->GetX(),
339           element->GetY(),
340           element->GetLastPointX(),
341           element->GetLastPointY()
342         };
343
344         Int_t mcLabel = -1;
345         if( mcLabels.find(element->TrackID())!=mcLabels.end() )
346           mcLabel = mcLabels[element->TrackID()];
347         element->SetLabel( mcLabel );
348
349         AliESDtrack iotrack;
350
351         // for the moment, the number of clusters are not set when processing the
352         // kTPCin update, only at kTPCout
353         // there ar emainly three parameters updated for kTPCout
354         //   number of clusters
355         //   chi2
356         //   pid signal
357         // The first one can be updated already at that stage here, while the two others
358         // eventually require to update from the ITS tracks before. The exact scheme
359         // needs to be checked 
360         iotrack.SetID( element->TrackID() );
361         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
362         {
363           AliHLTGlobalBarrelTrack outPar(*element);       
364           outPar.AliExternalTrackParam::PropagateTo( element->GetLastPointX(), fSolenoidBz );
365           iotrack.UpdateTrackParams(&outPar,AliESDtrack::kTPCout);
366         }
367         iotrack.SetTPCPoints(points);
368         if( element->TrackID()<ndEdxTPC ){
369           iotrack.SetTPCsignal( dEdxTPC[element->TrackID()], 0, 0 ); 
370           //AliTPCseed s;
371           //s.Set( element->GetX(), element->GetAlpha(),
372           //element->GetParameter(), element->GetCovariance() );
373           //s.SetdEdx( dEdxTPC[element->TrackID()] );
374           //s.CookPID();
375           //iotrack.SetTPCpid(s.TPCrPIDs() );
376         } else {
377           if( dEdxTPC ) HLTWarning("Wrong number of dEdx TPC labels");
378         }
379         pESD->AddTrack(&iotrack);
380         if (fVerbosity>0) element->Print();
381       }
382       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
383       iAddedDataBlocks++;
384     } else if (iResult<0) {
385       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
386                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
387     }
388   }
389
390
391   // Get ITS SPD vertex
392
393   for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); iter != NULL; iter = GetNextInputObject() ) {
394     AliESDVertex *vtx = dynamic_cast<AliESDVertex*>(const_cast<TObject*>( iter ) );
395     pESD->SetPrimaryVertexSPD( vtx );
396   }
397
398   // now update ESD tracks with the ITSOut info
399   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITSOut);
400        pBlock!=NULL; pBlock=GetNextInputBlock()) {
401     vector<AliHLTGlobalBarrelTrack> tracks;
402     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
403       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
404            element!=tracks.end(); element++) {
405         int tpcID=element->TrackID();
406         // the ITS tracker assigns the TPC track used as seed for a certain track to
407         // the trackID
408         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
409         AliESDtrack *tESD = pESD->GetTrack( tpcID );
410         if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSout );
411       }
412     }
413   }
414
415   // now update ESD tracks with the ITS info
416   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
417        pBlock!=NULL; pBlock=GetNextInputBlock()) {
418     vector<AliHLTGlobalBarrelTrack> tracks;
419     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
420       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
421            element!=tracks.end(); element++) {
422         int tpcID=element->TrackID();
423         // the ITS tracker assigns the TPC track used as seed for a certain track to
424         // the trackID
425         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
426         AliESDtrack *tESD = pESD->GetTrack( tpcID );
427         if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
428       }
429     }
430   }
431
432
433   // convert the HLT TRD tracks to ESD tracks                        
434   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
435        pBlock!=NULL; pBlock=GetNextInputBlock()) {
436     vector<AliHLTGlobalBarrelTrack> tracks;
437     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
438       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
439            element!=tracks.end(); element++) {
440         
441         Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
442         for(Int_t i=0; i<AliPID::kSPECIES; i++){
443           switch(i){
444           case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
445           default: TRDpid[i]=restProb; break;
446           }
447         }
448         
449         AliESDtrack iotrack;
450         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDin);
451         iotrack.SetTRDpid(TRDpid);
452         
453         pESD->AddTrack(&iotrack);
454         if (fVerbosity>0) element->Print();
455       }
456       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
457       iAddedDataBlocks++;
458     } else if (iResult<0) {
459       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
460                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
461     }
462   }
463   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeCaloCluster | kAliHLTDataOriginPHOS); pBlock!=NULL; pBlock=GetNextInputBlock()) 
464     {
465       AliHLTCaloClusterHeaderStruct *caloClusterHeaderPtr = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
466
467       HLTDebug("%d HLT clusters from spec: 0x%X", caloClusterHeaderPtr->fNClusters, pBlock->fSpecification);
468
469       AliHLTCaloClusterReader reader;
470       reader.SetMemory(caloClusterHeaderPtr);
471
472       AliHLTESDCaloClusterMaker clusterMaker;
473
474       int nClusters = clusterMaker.FillESD(pESD, caloClusterHeaderPtr);
475
476       HLTInfo("converted %d cluster(s) to AliESDCaloCluster and added to ESD", nClusters);
477       iAddedDataBlocks++;
478     }
479   
480   // Add tracks from MUON.
481   for (const TObject* obj = GetFirstInputObject(kAliHLTAnyDataType | kAliHLTDataOriginMUON);
482        obj != NULL;
483        obj = GetNextInputObject()
484       )
485   {
486     const TClonesArray* tracklist = NULL;
487     if (obj->IsA() == AliESDEvent::Class())
488     {
489       const AliESDEvent* event = static_cast<const AliESDEvent*>(obj);
490       HLTDebug("Received a MUON ESD with specification: 0x%X", GetSpecification(obj));
491       if (event->GetList() == NULL) continue;
492       tracklist = dynamic_cast<const TClonesArray*>(event->GetList()->FindObject("MuonTracks"));
493       if (tracklist == NULL) continue;
494     }
495     else if (obj->IsA() == TClonesArray::Class())
496     {
497       tracklist = static_cast<const TClonesArray*>(obj);
498       HLTDebug("Received a MUON TClonesArray of tracks with specification: 0x%X", GetSpecification(obj));
499     }
500     else
501     {
502       // Cannot handle this object type.
503       continue;
504     }
505     HLTDebug("Received %d MUON tracks.", tracklist->GetEntriesFast());
506     if (tracklist->GetEntriesFast() > 0)
507     {
508       const AliESDMuonTrack* track = dynamic_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(0));
509       if (track == NULL)
510       {
511         HLTError(Form("%s from MUON does not contain AliESDMuonTrack objects.", obj->ClassName()));
512         continue;
513       }
514     }
515     for (Int_t i = 0; i < tracklist->GetEntriesFast(); ++i)
516     {
517       const AliESDMuonTrack* track = static_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(i));
518       pESD->AddMuonTrack(track);
519     }
520   }
521   
522   if (iAddedDataBlocks>0 && pTree) {
523     pTree->Fill();
524   }
525   
526   if (iResult>=0) iResult=iAddedDataBlocks;
527   return iResult;
528 }