]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSClusterHistoComponent.cxx
coverity warnings 15388 10083 10082 fixed
[u/mrichter/AliRoot.git] / HLT / ITS / AliHLTITSClusterHistoComponent.cxx
1 // $Id$
2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project        * 
4 //* ALICE Experiment at CERN, All rights reserved.                         *
5 //*                                                                        *
6 //* Primary Authors: Gaute Ovrebekk <ovrebekk@ift.uib.no>                  *
7 //*                  for The ALICE HLT Project.                            *
8 //*                                                                        *
9 //* Permission to use, copy, modify and distribute this software and its   *
10 //* documentation strictly for non-commercial purposes is hereby granted   *
11 //* without fee, provided that the above copyright notice appears in all   *
12 //* copies and that both the copyright notice and this permission notice   *
13 //* appear in the supporting documentation. The authors make no claims     *
14 //* about the suitability of this software for any purpose. It is          *
15 //* provided "as is" without express or implied warranty.                  *
16 //**************************************************************************
17
18 /** @file   AliHLTITSClusterHistoComponent.cxx
19     @author Gaute Ovrebekk
20     @brief  Component for plotting charge in clusters
21 */
22
23 #if __GNUC__>= 3
24 using namespace std;
25 #endif
26
27 #include "AliHLTITSClusterHistoComponent.h"
28 #include "AliHLTITSClusterDataFormat.h"
29 #include "AliCDBEntry.h"
30 #include "AliCDBManager.h"
31 #include "AliGeomManager.h"
32 #include "AliITSRecPoint.h"
33 #include "TMath.h"
34 #include <TFile.h>
35 #include <TString.h>
36 #include "TObjString.h"
37 #include "TObjArray.h"
38
39 /** ROOT macro for the implementation of ROOT specific class methods */
40 ClassImp(AliHLTITSClusterHistoComponent)
41
42 AliHLTITSClusterHistoComponent::AliHLTITSClusterHistoComponent()
43   :
44   fXY(NULL),                     
45   fPhieta(NULL),                   
46   fCharge(NULL),   
47   fPlotCharge(kFALSE),   
48   fPlotXYPhiEta(kTRUE),
49   fOutputSize(20000)
50 {
51   // see header file for class documentation
52   // or
53   // refer to README to build package
54   // or
55   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
56 }
57
58 AliHLTITSClusterHistoComponent::~AliHLTITSClusterHistoComponent(){
59   // see header file for class documentation
60 }
61
62 // Public functions to implement AliHLTComponent's interface.
63 // These functions are required for the registration process
64
65 const char* AliHLTITSClusterHistoComponent::GetComponentID(){
66   // see header file for class documentation
67   return "ITSClusterHisto";
68 }
69
70 void AliHLTITSClusterHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list){
71   // see header file for class documentation
72   list.clear();
73   list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD );
74   list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD );
75   list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD );
76   list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITS );
77 }
78
79 AliHLTComponentDataType AliHLTITSClusterHistoComponent::GetOutputDataType(){
80   // see header file for class documentation
81   return kAliHLTDataTypeHistogram;
82 }
83
84 void AliHLTITSClusterHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ){
85   // see header file for class documentation
86   // XXX TODO: Find more realistic values.
87   constBase = fOutputSize;
88   inputMultiplier = 20;
89 }
90
91 AliHLTComponent* AliHLTITSClusterHistoComponent::Spawn(){
92   // see header file for class documentation
93   return new AliHLTITSClusterHistoComponent;
94 }
95
96 int AliHLTITSClusterHistoComponent::DoInit( int argc, const char** argv ){  
97   // see header file for class documentation
98      
99   if(AliGeomManager::GetGeometry()==NULL){
100     AliGeomManager::LoadGeometry();
101   }
102   
103   fPlotCharge   = kFALSE;   
104   fPlotXYPhiEta = kTRUE;
105      
106   if(fPlotCharge) fCharge = new TH1F("fCharge","Total Charge of clusters",2000,0,2000);
107   if(fPlotXYPhiEta){
108      fXY = new TH2F("fXY","Global XY of ITS clusters",1600,-80,80,1600,-80,80);
109      TString name;
110      TString title;
111      fPhieta = new TH2F*[6];
112      for (Int_t iLay=0;iLay<6;iLay++){
113              name.Form("Phi_vs_Eta_ITS_Layer%d",iLay+1);
114              title.Form("Phi vs Eta - ITS Layer %d",iLay+1);
115              fPhieta[iLay]=new TH2F(name,title,60,-1.5,1.5,60,0.,2*TMath::Pi());
116              fPhieta[iLay]->GetXaxis()->SetTitle("Pseudorapidity");
117              fPhieta[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
118      }
119   }
120   
121   int iResult=0;
122   TString configuration="";
123   TString argument="";
124
125   for(int i=0; i<argc && iResult>=0; i++){
126     argument=argv[i];
127     if (!configuration.IsNull()) configuration+=" ";
128     configuration+=argument;
129   }
130   
131   if(!configuration.IsNull()){
132     iResult=Configure(configuration.Data());
133   }  
134
135   return iResult; 
136 }
137   
138 int AliHLTITSClusterHistoComponent::DoDeinit(){
139   // see header file for class documentation
140   if(fCharge!=NULL) delete fCharge;
141   if(fXY!=NULL) {
142      delete fXY;
143      for(Int_t i=0;i<6;i++) delete fPhieta[5-i];
144   }
145   return 0;
146 }
147
148 int AliHLTITSClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/){
149   // see header file for class documentation
150    
151   static Int_t event = 0;
152   event++;
153   int TotalSpacePoint = 0;
154   
155   const AliHLTComponentBlockData* iter = NULL;
156   
157   AliHLTUInt32_t eventType=gkAliEventTypeUnknown;
158   if(!IsDataEvent(&eventType) && eventType!=gkAliEventTypeEndOfRun) {
159     // publish the histograms at the end of run
160     return 0;
161   }
162   
163   for(iter = GetFirstInputBlock(kAliHLTDataTypeClusters); eventType!=gkAliEventTypeEndOfRun && iter != NULL; iter = GetNextInputBlock()){
164   
165       if(iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITSSPD) && 
166          iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITSSDD) && 
167          iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITSSSD) &&
168          iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITS))    
169          continue;
170
171       const AliHLTITSClusterData* clusterData = (const AliHLTITSClusterData*) iter->fPtr;
172       Int_t nSpacepoint = (Int_t) clusterData->fSpacePointCnt;
173       TotalSpacePoint += nSpacepoint;
174       AliHLTITSSpacePointData *clusters = (AliHLTITSSpacePointData*) clusterData->fSpacePoints;
175
176       for(int i=0;i<nSpacepoint;i++){
177          
178           Int_t lab[4]   = {0,0,0,0};
179           Float_t hit[6] = {0,0,0,0,0,0};
180           Int_t info[3]  = {0,0,0};
181           
182           lab[0]  = clusters[i].fTracks[0];
183           lab[1]  = clusters[i].fTracks[1];
184           lab[2]  = clusters[i].fTracks[2];
185           lab[3]  = clusters[i].fIndex;
186           hit[0]  = clusters[i].fY;
187           hit[1]  = clusters[i].fZ;
188           hit[2]  = clusters[i].fSigmaY2;
189           hit[3]  = clusters[i].fSigmaZ2;
190           hit[4]  = clusters[i].fQ;
191           hit[5]  = clusters[i].fSigmaYZ;
192           info[0] = clusters[i].fNy;
193           info[1] = clusters[i].fNz;
194           info[2] = clusters[i].fLayer;
195  
196           Float_t xyz[3];
197           AliITSRecPoint recpoint(lab,hit,info);
198           recpoint.GetGlobalXYZ(xyz);
199           Int_t layer = recpoint.GetLayer();
200           
201           if(fPlotXYPhiEta){
202              fXY->Fill(xyz[0],xyz[1]);
203              Float_t rad=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]); 
204              Float_t theta=TMath::ATan2(rad,xyz[2]);
205              Float_t eta=-1*TMath::Log(TMath::Tan(theta/2.0));
206              Float_t phi=TMath::ATan2(xyz[1],xyz[0]);
207              if(phi<0.0){phi=2 * TMath::Pi() - TMath::Abs(phi);} 
208              fPhieta[layer]->Fill(eta,phi);
209           }
210           if(fPlotCharge) fCharge->Fill(recpoint.GetQ());
211     }
212   }
213   
214   if(fPlotCharge){
215      AliHLTUInt32_t fSpecification = 0x0;
216      if (PushBack( (TObject*) fCharge,kAliHLTDataTypeHistogram|kAliHLTDataOriginITS,fSpecification)==-ENOSPC) {
217        fOutputSize+=GetLastObjectSize();
218      }
219   }
220   if(fPlotXYPhiEta){
221      AliHLTUInt32_t fSpecification = 0x0;
222      if (PushBack( (TObject*) fXY,kAliHLTDataTypeHistogram|kAliHLTDataOriginITS,fSpecification)==-ENOSPC) {
223        fOutputSize+=GetLastObjectSize();
224      }
225      for(Int_t ii=0;ii<6;ii++) 
226        if (PushBack( (TObject*) fPhieta[ii],kAliHLTDataTypeHistogram|kAliHLTDataOriginITS,fSpecification)==-ENOSPC) {
227          fOutputSize+=GetLastObjectSize();
228        }
229   }
230   
231   HLTDebug("ITSClusterHisto found %d Total Spacepoints", TotalSpacePoint);
232   
233   return 0;
234 }
235
236 int AliHLTITSClusterHistoComponent::Configure(const char* arguments){
237   // see header file for class documentation
238    
239   int iResult=0;
240   
241   if (!arguments) return iResult;
242   
243   TString allArgs=arguments;
244   TString argument;
245   
246   TObjArray* pTokens=allArgs.Tokenize(" ");
247   
248   if (pTokens) {
249     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
250       argument=((TObjString*)pTokens->At(i))->GetString();
251       if (argument.IsNull()) continue;
252       
253       if (argument.CompareTo("-plot-all")==0) {
254         HLTInfo("Ploting all historgams");
255         fPlotXYPhiEta = kTRUE;
256         fPlotCharge = kTRUE;
257         continue;
258       }
259       
260       else if (argument.CompareTo("-plot-xy")==0) {
261         HLTInfo("Ploting Global XY");
262         fPlotXYPhiEta = kTRUE;
263         continue;
264       }
265
266       else if (argument.CompareTo("-plot-charge")==0) {
267         HLTInfo("Ploting charge of clusters");
268         fPlotCharge = kTRUE;
269         continue;
270       }
271
272       else {
273         HLTError("unknown argument %s", argument.Data());
274         iResult=-EINVAL;
275         break;
276       }
277     }
278     delete pTokens;
279   }
280   
281   if(!fCharge && fPlotCharge){fCharge = new TH1F("fCharge","Total Charge of clusters",2000,0,2000);}
282   if(!fXY && fPlotXYPhiEta){
283         fXY = new TH2F("fXY","Global XY of ITS clusters",1600,-80,80,1600,-80,80);
284         TString name;
285         TString title;
286         fPhieta = new TH2F*[6];
287         for (Int_t iLay=0;iLay<6;iLay++) {
288                 name.Form("Phi_vs_Eta_ITS_Layer%d",iLay+1);
289                 title.Form("Phi vs Eta - ITS Layer %d",iLay+1);
290                 fPhieta[iLay]=new TH2F(name,title,30,-1.5,1.5,200,0.,2*TMath::Pi());
291                 fPhieta[iLay]->GetXaxis()->SetTitle("Pseudorapidity");
292                 fPhieta[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
293         }
294   }
295   
296   return iResult;
297 }
298
299 int AliHLTITSClusterHistoComponent::Reconfigure(const char* cdbEntry, const char* chainId){
300   // see header file for class documentation
301   int iResult=0;
302   
303   const char* path="HLT/ConfigITS/HistoComponent";
304   const char* defaultNotify="";
305   if (cdbEntry) {
306     path=cdbEntry;
307     defaultNotify=" (default)";
308   }
309   if (path) {
310     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
311     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
312     if (pEntry) {
313       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
314       if (pString) {
315         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
316         iResult=Configure(pString->GetString().Data());
317       } else {
318         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
319       }
320     } else {
321       HLTError("can not fetch object \"%s\" from CDB", path);
322     }
323   }
324   
325   return iResult;
326 }