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