2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project *
4 //* ALICE Experiment at CERN, All rights reserved. *
6 //* Primary Authors: Gaute Ovrebekk <ovrebekk@ift.uib.no> *
7 //* for The ALICE HLT Project. *
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 //**************************************************************************
18 /** @file AliHLTITSClusterHistoComponent.cxx
19 @author Gaute Ovrebekk
20 @brief Component for ploting charge in clusters
27 #include "AliHLTITSClusterHistoComponent.h"
28 #include "AliHLTITSClusterDataFormat.h"
29 #include "AliCDBEntry.h"
30 #include "AliCDBManager.h"
31 #include "AliITSRecPoint.h"
35 #include "TObjString.h"
36 #include "TObjArray.h"
41 /** ROOT macro for the implementation of ROOT specific class methods */
42 ClassImp(AliHLTITSClusterHistoComponent)
44 AliHLTITSClusterHistoComponent::AliHLTITSClusterHistoComponent()
52 // see header file for class documentation
54 // refer to README to build package
56 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
60 AliHLTITSClusterHistoComponent::~AliHLTITSClusterHistoComponent()
62 // see header file for class documentation
65 // Public functions to implement AliHLTComponent's interface.
66 // These functions are required for the registration process
68 const char* AliHLTITSClusterHistoComponent::GetComponentID()
70 // see header file for class documentation
72 return "ITSClusterHisto";
75 void AliHLTITSClusterHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
77 // see header file for class documentation
79 list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD );
80 list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSDD );
81 list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD );
84 AliHLTComponentDataType AliHLTITSClusterHistoComponent::GetOutputDataType()
86 // see header file for class documentation
87 return kAliHLTDataTypeHistogram;
91 void AliHLTITSClusterHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
93 // see header file for class documentation
94 // XXX TODO: Find more realistic values.
99 AliHLTComponent* AliHLTITSClusterHistoComponent::Spawn()
101 // see header file for class documentation
102 return new AliHLTITSClusterHistoComponent;
105 int AliHLTITSClusterHistoComponent::DoInit( int argc, const char** argv )
110 if(fPlotCharge){fCharge = new TH1F("fCharge","Total Charge of clusters",2000,0,2000);}
112 fXY = new TH2F("fXY","Global XY of ITS clusters",1600,-80,80,1600,-80,80);
115 fPhieta = new TH2F*[6];
116 for (Int_t iLay=0;iLay<6;iLay++) {
117 sprintf(name,"Phi_vs_Eta_ITS_Layer%d",iLay+1);
118 sprintf(title,"Phi vs Eta - ITS Layer %d",iLay+1);
119 fPhieta[iLay]=new TH2F(name,title,60,-1.5,1.5,60,0.,2*TMath::Pi());
120 fPhieta[iLay]->GetXaxis()->SetTitle("Pseudorapidity");
121 fPhieta[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
127 TString configuration="";
129 for (int i=0; i<argc && iResult>=0; i++) {
131 if (!configuration.IsNull()) configuration+=" ";
132 configuration+=argument;
135 if (!configuration.IsNull()) {
136 iResult=Configure(configuration.Data());
142 int AliHLTITSClusterHistoComponent::DoDeinit()
144 // see header file for class documentation
145 if(fCharge!=NULL) delete fCharge;
148 for(Int_t i=0;i<6;i++) delete fPhieta[5-i];
153 int AliHLTITSClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
156 static Int_t event = 0;
158 int TotalSpacePoint = 0;
160 const AliHLTComponentBlockData* iter = NULL;
165 for ( iter = GetFirstInputBlock(kAliHLTDataTypeClusters); iter != NULL; iter = GetNextInputBlock() ) {
168 if(iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITSSPD) &&
169 iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITSSDD) &&
170 iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITSSSD))
173 const AliHLTITSClusterData* clusterData = (const AliHLTITSClusterData*) iter->fPtr;
174 Int_t nSpacepoint = (Int_t) clusterData->fSpacePointCnt;
175 TotalSpacePoint += nSpacepoint;
176 AliHLTITSSpacePointData *clusters = (AliHLTITSSpacePointData*) clusterData->fSpacePoints;
178 for(int i=0;i<nSpacepoint;i++){
179 Int_t lab[4]={0,0,0,0};
180 Float_t hit[6]={0,0,0,0,0,0};
181 Int_t info[3]={0,0,0};
183 lab[0]=clusters[i].fTracks[0];
184 lab[1]=clusters[i].fTracks[1];
185 lab[2]=clusters[i].fTracks[2];
186 lab[3]=clusters[i].fIndex;
187 hit[0]=clusters[i].fY;
188 hit[1]=clusters[i].fZ;
189 hit[2]=clusters[i].fSigmaY2;
190 hit[3]=clusters[i].fSigmaZ2;
191 hit[4]=clusters[i].fQ;
192 hit[5]=clusters[i].fSigmaYZ;
193 info[0]=clusters[i].fNy;
194 info[1]=clusters[i].fNz;
195 info[2]=clusters[i].fLayer;
198 AliITSRecPoint recpoint(lab,hit,info);
199 recpoint.GetGlobalXYZ(xyz);
200 Int_t layer = recpoint.GetLayer();
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);
211 fCharge->Fill(recpoint.GetQ());
217 AliHLTUInt32_t fSpecification = 0x0;
218 PushBack( (TObject*) fCharge,kAliHLTDataTypeHistogram,fSpecification);
221 AliHLTUInt32_t fSpecification = 0x0;
222 PushBack( (TObject*) fXY,kAliHLTDataTypeHistogram,fSpecification);
223 for(Int_t ii=0;ii<6;ii++)
224 PushBack( (TObject*) fPhieta[ii],kAliHLTDataTypeHistogram,fSpecification);
227 HLTInfo("ITSClusterHisto found %d Total Spacepoints", TotalSpacePoint);
232 int AliHLTITSClusterHistoComponent::Configure(const char* arguments)
237 if (!arguments) return iResult;
239 TString allArgs=arguments;
242 TObjArray* pTokens=allArgs.Tokenize(" ");
245 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
246 argument=((TObjString*)pTokens->At(i))->GetString();
247 if (argument.IsNull()) continue;
249 if (argument.CompareTo("-plot-all")==0) {
250 HLTInfo("Ploting all historgams");
251 fPlotXYPhiEta = kTRUE;
256 else if (argument.CompareTo("-plot-xy")==0) {
257 HLTInfo("Ploting Global XY");
258 fPlotXYPhiEta = kTRUE;
262 else if (argument.CompareTo("-plot-charge")==0) {
263 HLTInfo("Ploting charge of clusters");
269 HLTError("unknown argument %s", argument.Data());
277 if(!fCharge && fPlotCharge){fCharge = new TH1F("fCharge","Total Charge of clusters",2000,0,2000);}
278 if(!fXY && fPlotXYPhiEta){
279 fXY = new TH2F("fXY","Global XY of ITS clusters",1600,-80,80,1600,-80,80);
282 fPhieta = new TH2F*[6];
283 for (Int_t iLay=0;iLay<6;iLay++) {
284 sprintf(name,"Phi_vs_Eta_ITS_Layer%d",iLay+1);
285 sprintf(title,"Phi vs Eta - ITS Layer %d",iLay+1);
286 fPhieta[iLay]=new TH2F(name,title,30,-1.5,1.5,200,0.,2*TMath::Pi());
287 fPhieta[iLay]->GetXaxis()->SetTitle("Pseudorapidity");
288 fPhieta[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
295 int AliHLTITSClusterHistoComponent::Reconfigure(const char* cdbEntry, const char* chainId)
297 // see header file for class documentation
300 const char* path="HLT/ConfigITS/HistoComponent";
301 const char* defaultNotify="";
304 defaultNotify=" (default)";
307 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
308 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path);
310 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
312 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
313 iResult=Configure(pString->GetString().Data());
315 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
318 HLTError("can not fetch object \"%s\" from CDB", path);