]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalHistoComponent.cxx
- coverity fix 10062
[u/mrichter/AliRoot.git] / HLT / global / AliHLTGlobalHistoComponent.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   AliHLTGlobalHistoComponent.cxx
20 /// @author Matthias Richter
21 /// @date   2010-09-16
22 /// @brief  A histogramming component for global ESD properties based
23 ///         on the AliHLTTTreeProcessor
24
25 #include "AliHLTGlobalHistoComponent.h"
26 #include "AliESDEvent.h"
27 #include "AliESDv0.h"
28 #include "AliKFParticle.h"
29 #include "AliKFVertex.h"
30
31 #include "TTree.h"
32 #include "TString.h"
33 #include <cassert>
34
35 /** ROOT macro for the implementation of ROOT specific class methods */
36 ClassImp(AliHLTGlobalHistoComponent)
37
38 AliHLTGlobalHistoComponent::AliHLTGlobalHistoComponent()
39   : AliHLTTTreeProcessor()
40   , fEvent(0)
41   , fNofTracks(0)
42   , fNofV0s(0)
43   , fNofContributors(0)
44   , fVertexX(-99)
45   , fVertexY(-99)
46   , fVertexZ(-99)
47   , fVertexStatus(kFALSE)
48   , fMaxTrackCount(20000)
49   , fMaxV0Count(1000)
50   , fFillV0(kFALSE)
51   , fTrackVariables()
52   , fTrackVariablesInt()
53   , fV0Variables()
54 {
55   // see header file for class documentation
56   // or
57   // refer to README to build package
58   // or
59   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
60 }
61
62 AliHLTGlobalHistoComponent::~AliHLTGlobalHistoComponent(){
63   // see header file for class documentation
64   fTrackVariables.Reset();
65   fTrackVariablesInt.Reset();
66   fV0Variables.Reset();
67 }
68
69 void AliHLTGlobalHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list){
70   // see header file for class documentation
71   list.push_back(kAliHLTAllDataTypes);
72 }
73
74 AliHLTComponentDataType AliHLTGlobalHistoComponent::GetOutputDataType(){
75   // see header file for class documentation
76   return kAliHLTDataTypeHistogram|kAliHLTDataOriginOut;
77 }
78
79 TTree* AliHLTGlobalHistoComponent::CreateTree(int /*argc*/, const char** /*argv*/){
80 // create the tree and branches
81  
82   int iResult=0;
83   TTree* pTree = new TTree("ESDproperties", "HLT ESD properties");
84   if (!pTree) return NULL;
85
86   const char* trackVariableNames = {
87     // Note the black at the end of each name!
88     "Track_pt "
89     "Track_phi "
90     "Track_eta "
91     "Track_p "
92     "Track_theta "
93     "Track_TPCclus "
94     "Track_ITSclus "
95     "Track_status "
96     "Track_charge "
97     "Track_DCAr "
98     "Track_DCAz "
99     "Track_dEdx "
100   };
101   
102   const char* trackIntVariableNames = {
103     // Note the black at the end of each name!
104     "Track_status "
105   };
106
107   const char* V0VariableNames = {
108     // Note the black at the end of each name!
109     "V0_AP "
110     "V0_pt "  
111     "clust1 "
112     "clust2 "
113     "dev1 "
114     "dev2 "
115     "devPrim "
116     "length "
117     "sigmaLength "
118     "r "
119   };
120      
121   if((iResult=fTrackVariables.Init(fMaxTrackCount, trackVariableNames))<0){
122     HLTError("failed to initialize internal structure for track properties (float)");
123   }
124   if((iResult=fTrackVariablesInt.Init(fMaxTrackCount, trackIntVariableNames))<0){
125     HLTError("failed to initialize internal structure for track properties (int)");
126   }
127   if((iResult=fV0Variables.Init(fMaxV0Count, V0VariableNames))<0){
128       HLTError("failed to initialize internal structure for V0 properties (float)");
129   }
130   
131   if(iResult>=0){     
132      pTree->Branch("event",        &fEvent,           "event/I");
133      pTree->Branch("trackcount",   &fNofTracks,       "trackcount/I");
134      pTree->Branch("vertexX",      &fVertexX,         "vertexX/F");
135      pTree->Branch("vertexY",      &fVertexY,         "vertexY/F");
136      pTree->Branch("vertexZ",      &fVertexZ,         "vertexZ/F");
137      pTree->Branch("nContributors",&fNofContributors, "nContributors/I");
138      pTree->Branch("vertexStatus", &fVertexStatus,    "vertexStatus/I");
139      if(fFillV0==kTRUE) pTree->Branch("V0", &fNofV0s, "V0/I");
140
141      int i=0;
142      // FIXME: this is a bit ugly since type 'f' and 'i' are specified
143      // explicitely. Would be better to use a function like
144      // AliHLTGlobalHistoVariables::GetType but could not get this working
145     
146      for(i=0; i<fTrackVariables.Variables(); i++){
147          TString specifier=fTrackVariables.GetKey(i);
148          float* pArray=fTrackVariables.GetArray(specifier);
149          specifier+="[trackcount]/f";
150          pTree->Branch(fTrackVariables.GetKey(i), pArray, specifier.Data());
151      }
152      for(i=0; i<fTrackVariablesInt.Variables(); i++){
153          TString specifier=fTrackVariablesInt.GetKey(i);
154          int* pArray=fTrackVariablesInt.GetArray(specifier);
155          specifier+="[trackcount]/i";
156          pTree->Branch(fTrackVariablesInt.GetKey(i), pArray, specifier.Data());
157      }    
158      if(fFillV0==kTRUE){
159         for(i=0; i<fV0Variables.Variables(); i++){
160             TString specifier=fV0Variables.GetKey(i);
161             float* pArray=fV0Variables.GetArray(specifier);
162             specifier+="[V0]/f";
163             pTree->Branch(fV0Variables.GetKey(i), pArray, specifier.Data());
164         }
165      }
166   } else {
167     delete pTree;
168     pTree = NULL;
169   }
170   
171   return pTree;
172 }
173
174 void AliHLTGlobalHistoComponent::FillHistogramDefinitions(){
175   /// default histogram definitions
176 }
177
178 int AliHLTGlobalHistoComponent::FillTree(TTree* pTree, const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/ ){
179
180   /// fill the tree from the ESD
181   int iResult=0;
182   if (!IsDataEvent()) return 0;
183
184   ResetVariables();
185
186   // fetch ESD from input stream
187   const TObject *obj = GetFirstInputObject(kAliHLTAllDataTypes, "AliESDEvent");
188   AliESDEvent *esd = dynamic_cast<AliESDEvent*>(const_cast<TObject*>(obj));
189   if(!esd) return 0;
190   esd->GetStdContent();
191
192   // fill track variables
193   fVertexX         = esd->GetPrimaryVertexTracks()->GetX();
194   fVertexY         = esd->GetPrimaryVertexTracks()->GetY();
195   fVertexZ         = esd->GetPrimaryVertexTracks()->GetZ();
196   fNofContributors = esd->GetPrimaryVertexTracks()->GetNContributors();
197   fVertexStatus    = esd->GetPrimaryVertexTracks()->GetStatus();
198   fNofV0s          = esd->GetNumberOfV0s();
199   
200   if(fFillV0==kTRUE && (fNofV0s > fMaxV0Count)){
201      HLTWarning("Found V0s are %d, while component argument is %d, the respective TTree branch is not filled properly. Need to reconfigure.\n", fNofV0s, fMaxV0Count);
202   }
203
204   Int_t nTracks = 0;
205
206   for(int i=0; i<esd->GetNumberOfTracks(); i++){    
207       AliESDtrack *esdTrack = esd->GetTrack(i);
208       if (!esdTrack) continue;
209       
210       Float_t DCAr, DCAz = -99;
211       esdTrack->GetImpactParametersTPC(DCAr, DCAz);
212       
213       fTrackVariables.Fill("Track_pt"        , esdTrack->Pt()                      );
214       fTrackVariables.Fill("Track_phi"       , esdTrack->Phi()*TMath::RadToDeg()   );
215       fTrackVariables.Fill("Track_eta"       , esdTrack->Eta()                     );
216       fTrackVariables.Fill("Track_p"         , esdTrack->P()                       );
217       fTrackVariables.Fill("Track_theta"     , esdTrack->Theta()*TMath::RadToDeg() );
218       fTrackVariables.Fill("Track_TPCclus"   , esdTrack->GetTPCNcls()              );
219       fTrackVariables.Fill("Track_ITSclus"   , esdTrack->GetNcls(0)                );
220       fTrackVariables.Fill("Track_status"    , esdTrack->GetStatus()               );
221       fTrackVariables.Fill("Track_charge"    , esdTrack->Charge()                  );
222       fTrackVariables.Fill("Track_DCAr"      , DCAr                                );
223       fTrackVariables.Fill("Track_DCAz"      , DCAz                                );   
224       fTrackVariables.Fill("Track_dEdx"      , esdTrack->GetTPCsignal()            );   
225       fTrackVariablesInt.Fill("Track_status" , esdTrack->GetStatus()               );   
226       
227       // selection of TPC tracks, the condition rejects e.g. TRD tracklets which would appear with 0 TPC clusters
228       if( esdTrack->GetTPCNcls()>0 ) nTracks++;
229   }
230   
231   fNofTracks = nTracks; 
232     
233   if(fFillV0==kTRUE){
234      for(int i=0; i<fNofV0s; i++){     
235          AliESDv0 *esdV0 = esd->GetV0(i);
236          if(!esdV0) continue;
237             
238         //fV0Variables.Fill("V0_AP", ap);
239         //fV0Variables.Fill("V0_pt", pt); 
240         //fV0Variables.Fill("clust1", t1->GetTPCNcls()); 
241         //fV0Variables.Fill("clust2", t2->GetTPCNcls()); 
242         //fV0Variables.Fill("dev1", dev1); 
243         //fV0Variables.Fill("dev2", dev2); 
244         //fV0Variables.Fill("devPrim", devPrim); 
245         //fV0Variables.Fill("length", length); 
246         //fV0Variables.Fill("sigmaLength", sigmaLength); 
247         //fV0Variables.Fill("r", r); 
248           
249      } // end of loop over V0s
250   } 
251   /* 
252   if(iResult<0){
253     // fill an empty event
254     ResetVariables();
255   }
256   */  
257   fEvent++;
258   pTree->Fill();
259   return iResult;
260 }
261
262 int AliHLTGlobalHistoComponent::ResetVariables(){
263 /// reset all filling variables
264   fNofTracks=0;
265   fNofV0s=0;
266   fTrackVariables.ResetCount();
267   fTrackVariablesInt.ResetCount();
268   fV0Variables.ResetCount();
269   return 0;
270 }
271
272 AliHLTComponentDataType AliHLTGlobalHistoComponent::GetOriginDataType() const{
273 // get the origin of the output data
274   return kAliHLTDataTypeHistogram|kAliHLTDataOriginHLT;
275 }
276
277 int AliHLTGlobalHistoComponent::ScanConfigurationArgument(int argc, const char** argv){
278 /// inherited from AliHLTComponent, scan argument
279
280   if (argv==NULL || argc<1) return 0;
281
282   int i=0;
283   TString argument = argv[i];
284
285   // -max-track-count
286   if(argument.CompareTo("-max-track-count")==0){    
287      if (++i>=argc) return -EPROTO; // missing parameter
288      argument = argv[i];
289      fMaxTrackCount = argument.Atoi();    
290      HLTInfo("got %s with parameter %s", argument.Data(), argv[i]);
291      return ++i; 
292   }
293  
294   // -max-V0-count 
295   if(argument.CompareTo("-max-V0-count")==0){    
296      if (++i>=argc) return -EPROTO; // missing parameter
297      argument = argv[i];
298      fMaxV0Count = argument.Atoi();    
299      HLTInfo("got %s with parameter %s", argument.Data(), argv[i]);
300      return ++i;
301   }
302   
303   // -fill-V0
304   if(argument.CompareTo("-fill-V0")==0){
305      fFillV0 = kTRUE;
306      HLTInfo("got %s", argument.Data());
307      return ++i;
308   }
309   // no recognized argument, forward to base class
310   return AliHLTTTreeProcessor::ScanConfigurationArgument(argc, argv);
311 }
312
313 int AliHLTGlobalHistoComponent::Reconfigure(const char* cdbEntry, const char* /*chainId*/){  
314 // see header file for class documentation
315
316   TString cdbPath;
317   if (cdbEntry) {
318     cdbPath=cdbEntry;
319   } else {
320     cdbPath="HLT/ConfigHLT/";
321     cdbPath+=GetComponentID();
322   }
323
324   return ConfigureFromCDBTObjString(cdbPath.Data());
325 }
326
327