]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalEsdConverterComponent.cxx
Change of a convention of the sign of the Bz field, following the change in the off...
[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 }
146
147 AliHLTComponentDataType AliHLTGlobalEsdConverterComponent::GetOutputDataType()
148 {
149   // see header file for class documentation
150   return kAliHLTDataTypeESDObject|kAliHLTDataOriginOut;
151 }
152
153 void AliHLTGlobalEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
154 {
155   // see header file for class documentation
156   constBase=2000000;
157   inputMultiplier=10.0;
158 }
159
160 int AliHLTGlobalEsdConverterComponent::DoInit(int argc, const char** argv)
161 {
162   // see header file for class documentation
163   int iResult=0;
164   TString argument="";
165   int bMissingParam=0;
166   for (int i=0; i<argc && iResult>=0; i++) {
167     argument=argv[i];
168     if (argument.IsNull()) continue;
169
170     // -notree
171     if (argument.CompareTo("-notree")==0) {
172       fWriteTree=0;
173
174       // -tree
175     } else if (argument.CompareTo("-tree")==0) {
176       fWriteTree=1;
177
178     } else {
179       HLTError("unknown argument %s", argument.Data());
180       break;
181     }
182   }
183   if (bMissingParam) {
184     HLTError("missing parameter for argument %s", argument.Data());
185     iResult=-EINVAL;
186   }
187
188   if (iResult>=0) {
189     iResult=Reconfigure(NULL, NULL);
190   }
191
192   if (iResult>=0) {
193     fESD = new AliESDEvent;
194     if (fESD) {
195       fESD->CreateStdContent();
196     } else {
197       iResult=-ENOMEM;
198     }
199
200     SetupCTPData();
201   }
202
203   return iResult;
204 }
205
206 int AliHLTGlobalEsdConverterComponent::DoDeinit()
207 {
208   // see header file for class documentation
209   if (fESD) delete fESD;
210   fESD=NULL;
211
212   return 0;
213 }
214
215 int AliHLTGlobalEsdConverterComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, 
216                                                AliHLTComponentTriggerData& trigData)
217 {
218   // see header file for class documentation
219   int iResult=0;
220   if (!fESD) return -ENODEV;
221
222   AliESDEvent* pESD = fESD;
223
224   pESD->Reset(); 
225   pESD->SetMagneticField(fSolenoidBz);
226   pESD->SetRunNumber(GetRunNo());
227   pESD->SetPeriodNumber(GetPeriodNumber());
228   pESD->SetOrbitNumber(GetOrbitNumber());
229   pESD->SetBunchCrossNumber(GetBunchCrossNumber());
230   pESD->SetTimeStamp(GetTimeStamp());
231
232   const AliHLTCTPData* pCTPData=CTPData();
233   if (pCTPData) {
234     AliHLTUInt64_t mask=pCTPData->ActiveTriggers(trigData);
235     for (int index=0; index<gkNCTPTriggerClasses; index++) {
236       if ((mask&((AliHLTUInt64_t)0x1<<index)) == 0) continue;
237       pESD->SetTriggerClass(pCTPData->Name(index), index);
238     }
239     pESD->SetTriggerMask(mask);
240   }
241
242   TTree* pTree = NULL;
243   if (fWriteTree)
244     pTree = new TTree("esdTree", "Tree with HLT ESD objects");
245  
246   if (pTree) {
247     pTree->SetDirectory(0);
248   }
249
250   if ((iResult=ProcessBlocks(pTree, pESD))>=0) {
251     // TODO: set the specification correctly
252     if (pTree) {
253       // the esd structure is written to the user info and is
254       // needed in te ReadFromTree method to read all objects correctly
255       pTree->GetUserInfo()->Add(pESD);
256       pESD->WriteToTree(pTree);
257       iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginOut, 0);
258     } else {
259       iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
260     }
261   }
262   if (pTree) {
263     // clear user info list to prevent objects from being deleted
264     pTree->GetUserInfo()->Clear();
265     delete pTree;
266   }
267   return iResult;
268 }
269
270 int AliHLTGlobalEsdConverterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD)
271 {
272   // see header file for class documentation
273
274   int iResult=0;
275   int iAddedDataBlocks=0;
276
277   // Barrel tracking
278
279   // in the first attempt this component reads the TPC tracks and updates in the
280   // second step from the ITS tracks
281
282
283   // first read MC information (if present)
284   std::map<int,int> mcLabels;
285
286   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
287        pBlock!=NULL; pBlock=GetNextInputBlock()) {
288     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
289     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
290       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
291         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
292         mcLabels[lab.fTrackID] = lab.fMCLabel;
293       }
294     } else {
295       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
296                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
297                  dataPtr->fCount, pBlock->fSize);
298     }
299   }
300
301   // convert the TPC tracks to ESD tracks
302   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
303        pBlock!=NULL; pBlock=GetNextInputBlock()) {
304     vector<AliHLTGlobalBarrelTrack> tracks;
305     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
306       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
307            element!=tracks.end(); element++) {
308         Float_t points[4] = {
309           element->GetX(),
310           element->GetY(),
311           element->GetLastPointX(),
312           element->GetLastPointY()
313         };
314
315         Int_t mcLabel = -1;
316         if( mcLabels.find(element->TrackID())!=mcLabels.end() )
317           mcLabel = mcLabels[element->TrackID()];
318         element->SetLabel( mcLabel );
319
320         AliESDtrack iotrack;
321
322         // for the moment, the number of clusters are not set when processing the
323         // kTPCin update, only at kTPCout
324         // there ar emainly three parameters updated for kTPCout
325         //   number of clusters
326         //   chi2
327         //   pid signal
328         // The first one can be updated already at that stage here, while the two others
329         // eventually require to update from the ITS tracks before. The exact scheme
330         // needs to be checked 
331         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
332         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCout);
333         iotrack.SetTPCPoints(points);
334
335         pESD->AddTrack(&iotrack);
336         if (fVerbosity>0) element->Print();
337       }
338       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
339       iAddedDataBlocks++;
340     } else if (iResult<0) {
341       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
342                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
343     }
344   }
345
346   // now update ESD tracks with the ITS info
347   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
348        pBlock!=NULL; pBlock=GetNextInputBlock()) {
349     vector<AliHLTGlobalBarrelTrack> tracks;
350     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
351       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
352            element!=tracks.end(); element++) {
353         int ncl=0;
354         const UInt_t* pointsArray=element->GetPoints();
355         for( unsigned il=0; il<element->GetNumberOfPoints(); il++ ){
356           // TODO: check what needs to be done with the clusters
357           if( pointsArray[il]<~(UInt_t)0 ) {/*tITS.SetClusterIndex(ncl,  tr.fClusterIds[il]);*/}
358           ncl++;
359         }
360         //tITS.SetNumberOfClusters( ncl );
361         int tpcID=element->TrackID();
362         // the ITS tracker assigns the TPC track used as seed for a certain track to
363         // the trackID
364         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
365
366         AliESDtrack *tESD = pESD->GetTrack( tpcID );
367         if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
368       }
369     }
370   }
371
372   // convert the HLT TRD tracks to ESD tracks                        
373   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
374        pBlock!=NULL; pBlock=GetNextInputBlock()) {
375     vector<AliHLTGlobalBarrelTrack> tracks;
376     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
377       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
378            element!=tracks.end(); element++) {
379         
380         Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
381         for(Int_t i=0; i<AliPID::kSPECIES; i++){
382           switch(i){
383           case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
384           default: TRDpid[i]=restProb; break;
385           }
386         }
387         
388         AliESDtrack iotrack;
389         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDin);
390         iotrack.SetTRDpid(TRDpid);
391         
392         pESD->AddTrack(&iotrack);
393         if (fVerbosity>0) element->Print();
394       }
395       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
396       iAddedDataBlocks++;
397     } else if (iResult<0) {
398       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
399                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
400     }
401   }
402   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeCaloCluster | kAliHLTDataOriginPHOS); pBlock!=NULL; pBlock=GetNextInputBlock()) 
403     {
404       AliHLTCaloClusterHeaderStruct *caloClusterHeaderPtr = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
405
406       HLTDebug("%d HLT clusters from spec: 0x%X", caloClusterHeaderPtr->fNClusters, pBlock->fSpecification);
407
408       AliHLTCaloClusterReader reader;
409       reader.SetMemory(caloClusterHeaderPtr);
410
411       AliHLTESDCaloClusterMaker clusterMaker;
412
413       int nClusters = clusterMaker.FillESD(pESD, caloClusterHeaderPtr);
414
415       HLTInfo("converted %d cluster(s) to AliESDCaloCluster and added to ESD", nClusters);
416       iAddedDataBlocks++;
417     }
418       
419       // primary vertex & V0's 
420       /*
421       //AliHLTVertexer vertexer;
422       //vertexer.SetESD( pESD );
423       //vertexer.FindPrimaryVertex();
424       //vertexer.FindV0s();
425       */
426       if (iAddedDataBlocks>0 && pTree) {
427        pTree->Fill();
428       }
429   
430   if (iResult>=0) iResult=iAddedDataBlocks;
431   return iResult;
432 }