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