]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalEsdConverterComponent.cxx
- Adding merging of calorimeter clusters to the global ESD.
[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 "AliESDEvent.h"
37 #include "AliESDtrack.h"
38 #include "AliCDBEntry.h"
39 #include "AliCDBManager.h"
40 #include "AliPID.h"
41 #include "TTree.h"
42 #include "TList.h"
43 #include "AliHLTESDCaloClusterMaker.h"
44 #include "AliHLTCaloClusterDataStruct.h"
45 #include "AliHLTCaloClusterReader.h"
46 #include "AliESDCaloCluster.h"
47
48 /** ROOT macro for the implementation of ROOT specific class methods */
49 ClassImp(AliHLTGlobalEsdConverterComponent)
50
51 AliHLTGlobalEsdConverterComponent::AliHLTGlobalEsdConverterComponent()
52   : AliHLTProcessor()
53   , fWriteTree(0)
54   , fVerbosity(0)
55   , fESD(NULL)
56   , fSolenoidBz(5)
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 AliHLTGlobalEsdConverterComponent::~AliHLTGlobalEsdConverterComponent()
66 {
67   // see header file for class documentation
68   if (fESD) delete fESD;
69   fESD=NULL;
70 }
71
72 int AliHLTGlobalEsdConverterComponent::Configure(const char* arguments)
73 {
74   // see header file for class documentation
75   int iResult=0;
76   if (!arguments) return iResult;
77
78   TString allArgs=arguments;
79   TString argument;
80   int bMissingParam=0;
81
82   TObjArray* pTokens=allArgs.Tokenize(" ");
83   if (pTokens) {
84     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
85       argument=((TObjString*)pTokens->At(i))->GetString();
86       if (argument.IsNull()) continue;
87       
88       if (argument.CompareTo("-solenoidBz")==0) {
89         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
90         HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
91         fSolenoidBz=((TObjString*)pTokens->At(i))->GetString().Atof();
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=kAliHLTCDBSolenoidBz;
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 }
145
146 AliHLTComponentDataType AliHLTGlobalEsdConverterComponent::GetOutputDataType()
147 {
148   // see header file for class documentation
149   return kAliHLTDataTypeESDObject|kAliHLTDataOriginOut;
150 }
151
152 void AliHLTGlobalEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
153 {
154   // see header file for class documentation
155   constBase=2000000;
156   inputMultiplier=10.0;
157 }
158
159 int AliHLTGlobalEsdConverterComponent::DoInit(int argc, const char** argv)
160 {
161   // see header file for class documentation
162   int iResult=0;
163   TString argument="";
164   int bMissingParam=0;
165   for (int i=0; i<argc && iResult>=0; i++) {
166     argument=argv[i];
167     if (argument.IsNull()) continue;
168
169     // -notree
170     if (argument.CompareTo("-notree")==0) {
171       fWriteTree=0;
172
173       // -tree
174     } else if (argument.CompareTo("-tree")==0) {
175       fWriteTree=1;
176
177     } else {
178       HLTError("unknown argument %s", argument.Data());
179       break;
180     }
181   }
182   if (bMissingParam) {
183     HLTError("missing parameter for argument %s", argument.Data());
184     iResult=-EINVAL;
185   }
186
187   if (iResult>=0) {
188     iResult=Reconfigure(NULL, NULL);
189   }
190
191   if (iResult>=0) {
192     fESD = new AliESDEvent;
193     if (fESD) {
194       fESD->CreateStdContent();
195     } else {
196       iResult=-ENOMEM;
197     }
198   }
199
200   return iResult;
201 }
202
203 int AliHLTGlobalEsdConverterComponent::DoDeinit()
204 {
205   // see header file for class documentation
206   if (fESD) delete fESD;
207   fESD=NULL;
208
209   return 0;
210 }
211
212 int AliHLTGlobalEsdConverterComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, 
213                                                AliHLTComponentTriggerData& /*trigData*/)
214 {
215   // see header file for class documentation
216   int iResult=0;
217   if (!fESD) return -ENODEV;
218
219   AliESDEvent* pESD = fESD;
220
221   pESD->Reset(); 
222   pESD->SetMagneticField(fSolenoidBz);
223   pESD->SetPeriodNumber(GetPeriodNumber());
224   pESD->SetOrbitNumber(GetOrbitNumber());
225   pESD->SetBunchCrossNumber(GetBunchCrossNumber());
226   pESD->SetTimeStamp(GetTimeStamp());
227
228   TTree* pTree = NULL;
229   if (fWriteTree)
230     pTree = new TTree("esdTree", "Tree with HLT ESD objects");
231  
232   if (pTree) {
233     pTree->SetDirectory(0);
234   }
235
236   if ((iResult=ProcessBlocks(pTree, pESD))>=0) {
237     // TODO: set the specification correctly
238     if (pTree) {
239       // the esd structure is written to the user info and is
240       // needed in te ReadFromTree method to read all objects correctly
241       pTree->GetUserInfo()->Add(pESD);
242       pESD->WriteToTree(pTree);
243       iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginOut, 0);
244     } else {
245       iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
246     }
247   }
248   if (pTree) {
249     // clear user info list to prevent objects from being deleted
250     pTree->GetUserInfo()->Clear();
251     delete pTree;
252   }
253   return iResult;
254 }
255
256 int AliHLTGlobalEsdConverterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD)
257 {
258   // see header file for class documentation
259
260   int iResult=0;
261   int iAddedDataBlocks=0;
262
263   // Barrel tracking
264
265   // in the first attempt this component reads the TPC tracks and updates in the
266   // second step from the ITS tracks
267
268
269   // first read MC information (if present)
270   std::map<int,int> mcLabels;
271
272   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
273        pBlock!=NULL; pBlock=GetNextInputBlock()) {
274     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
275     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
276       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
277         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
278         mcLabels[lab.fTrackID] = lab.fMCLabel;
279       }
280     } else {
281       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
282                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
283                  dataPtr->fCount, pBlock->fSize);
284     }
285   }
286
287   // convert the TPC tracks to ESD tracks
288   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
289        pBlock!=NULL; pBlock=GetNextInputBlock()) {
290     vector<AliHLTGlobalBarrelTrack> tracks;
291     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
292       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
293            element!=tracks.end(); element++) {
294         Float_t points[4] = {
295           element->GetX(),
296           element->GetY(),
297           element->GetLastPointX(),
298           element->GetLastPointY()
299         };
300
301         Int_t mcLabel = -1;
302         if( mcLabels.find(element->TrackID())!=mcLabels.end() )
303           mcLabel = mcLabels[element->TrackID()];
304         element->SetLabel( mcLabel );
305
306         AliESDtrack iotrack;
307         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
308         iotrack.SetTPCPoints(points);
309
310         pESD->AddTrack(&iotrack);
311         if (fVerbosity>0) element->Print();
312       }
313       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
314       iAddedDataBlocks++;
315     } else if (iResult<0) {
316       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
317                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
318     }
319   }
320
321   // now update ESD tracks with the ITS info
322   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
323        pBlock!=NULL; pBlock=GetNextInputBlock()) {
324     vector<AliHLTGlobalBarrelTrack> tracks;
325     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
326       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
327            element!=tracks.end(); element++) {
328         int ncl=0;
329         const UInt_t* pointsArray=element->GetPoints();
330         for( unsigned il=0; il<element->GetNumberOfPoints(); il++ ){
331           // TODO: check what needs to be done with the clusters
332           if( pointsArray[il]<~(UInt_t)0 ) {/*tITS.SetClusterIndex(ncl,  tr.fClusterIds[il]);*/}
333           ncl++;
334         }
335         //tITS.SetNumberOfClusters( ncl );
336         int tpcID=element->TrackID();
337         // the ITS tracker assigns the TPC track used as seed for a certain track to
338         // the trackID
339         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
340
341         AliESDtrack *tESD = pESD->GetTrack( tpcID );
342         if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
343       }
344     }
345   }
346
347   // convert the HLT TRD tracks to ESD tracks                        
348   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
349        pBlock!=NULL; pBlock=GetNextInputBlock()) {
350     vector<AliHLTGlobalBarrelTrack> tracks;
351     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
352       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
353            element!=tracks.end(); element++) {
354         
355         Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
356         for(Int_t i=0; i<AliPID::kSPECIES; i++){
357           switch(i){
358           case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
359           default: TRDpid[i]=restProb; break;
360           }
361         }
362         
363         AliESDtrack iotrack;
364         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDin);
365         iotrack.SetTRDpid(TRDpid);
366         
367         pESD->AddTrack(&iotrack);
368         if (fVerbosity>0) element->Print();
369       }
370       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
371       iAddedDataBlocks++;
372     } else if (iResult<0) {
373       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
374                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
375     }
376   }
377   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeCaloCluster | kAliHLTDataOriginPHOS); pBlock!=NULL; pBlock=GetNextInputBlock()) 
378     {
379       AliHLTCaloClusterHeaderStruct *caloClusterHeaderPtr = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
380
381       HLTDebug("%d HLT clusters from spec: 0x%X", caloClusterHeaderPtr->fNClusters, pBlock->fSpecification);
382
383       AliHLTCaloClusterReader reader;
384       reader.SetMemory(caloClusterHeaderPtr);
385
386       AliHLTESDCaloClusterMaker clusterMaker;
387
388       int nClusters = clusterMaker.FillESD(pESD, caloClusterHeaderPtr);
389
390       HLTInfo("converted %d cluster(s) to AliESDCaloCluster and added to ESD", nClusters);
391       iAddedDataBlocks++;
392     }
393       
394       // primary vertex & V0's 
395       /*
396       //AliHLTVertexer vertexer;
397       //vertexer.SetESD( pESD );
398       //vertexer.FindPrimaryVertex();
399       //vertexer.FindV0s();
400       */
401       if (iAddedDataBlocks>0 && pTree) {
402        pTree->Fill();
403       }
404   
405   if (iResult>=0) iResult=iAddedDataBlocks;
406   return iResult;
407 }