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