3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project *
5 //* ALICE Experiment at CERN, All rights reserved. *
7 //* Primary Authors: Gaute Ovrebekk <ovrebekk@ift.uib.no> *
8 //* for The ALICE HLT Project. *
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 //**************************************************************************
19 /** @file AliHLTTPCTrackHistoComponent.cxx
20 @author Gaute Ovrebekk, Matthias Richter
22 @brief A histogram component with track and associated cluster properties
29 #include "AliHLTTPCTrackHistoComponent.h"
30 #include "AliHLTTPCTransform.h"
31 #include "AliHLTTPCClusterDataFormat.h"
32 #include "AliHLTTPCTrackletDataFormat.h"
33 #include "AliHLTTPCMemHandler.h"
34 #include "AliHLTTPCDefinitions.h"
35 #include "AliHLTGlobalBarrelTrack.h"
36 #include "AliHLTExternalTrackParam.h"
37 #include "AliHLTDataTypes.h"
44 #include "TObjString.h"
45 #include "TObjArray.h"
48 /** ROOT macro for the implementation of ROOT specific class methods */
49 ClassImp(AliHLTTPCTrackHistoComponent)
51 AliHLTTPCTrackHistoComponent::AliHLTTPCTrackHistoComponent()
59 , fMeanMultiplicity(NULL)
65 //, fNClusterVsXY(NULL)
67 //, fClustersArray(NULL)
68 //, fNSpacePoints(NULL)
70 // see header file for class documentation
72 // refer to README to build package
74 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
77 const char* AliHLTTPCTrackHistoComponent::fgkOCDBEntry="HLT/ConfigTPC/TPCTrackHisto";
79 AliHLTTPCTrackHistoComponent::~AliHLTTPCTrackHistoComponent(){
80 // see header file for class documentation
83 // Public functions to implement AliHLTComponent's interface.
84 // These functions are required for the registration process
86 const char* AliHLTTPCTrackHistoComponent::GetComponentID(){
87 // see header file for class documentation
89 return "TPCTrackHisto";
92 void AliHLTTPCTrackHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list){
93 // see header file for class documentation
96 list.push_back(AliHLTTPCDefinitions::fgkClustersDataType|kAliHLTDataOriginTPC);
97 list.push_back(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
98 //list.push_back(kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC);
101 AliHLTComponentDataType AliHLTTPCTrackHistoComponent::GetOutputDataType(){
102 // see header file for class documentation
103 return kAliHLTMultipleDataType;
106 int AliHLTTPCTrackHistoComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList){
107 // see header file for class documentation
109 tgtList.push_back(kAliHLTDataTypeTNtuple|kAliHLTDataOriginTPC);
110 tgtList.push_back(kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
111 return tgtList.size();
114 void AliHLTTPCTrackHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ){
115 // see header file for class documentation
118 inputMultiplier = 1;// XXX TODO: Find more realistic value
121 AliHLTComponent* AliHLTTPCTrackHistoComponent::Spawn(){
122 // see header file for class documentation
123 return new AliHLTTPCTrackHistoComponent;
126 int AliHLTTPCTrackHistoComponent::DoInit( int argc, const char** argv ){
127 // see header file for class documentation
129 fClusters = new TNtuple("fClusters", "fClusters", "charge:qmax:residualY:residualZ");
130 fTracks = new TNtuple("fTracks", "fTracks", "pt:eta:psi:nclusters");
132 fClusters->SetCircular(fBufferSize);
133 fTracks->SetCircular(fBufferSize);
135 fMultiplicity = new TH1F("fMultiplicity", "Track multiplicity per event", 1000, 0, 1000);
136 fMeanMultiplicity = new TH1F("fMeanMultiplicity", "Mean track multiplicity vs. #evt", 10000/fEvtMod, 0, 10000);
137 fDeDxVsP = new TProfile("fDeDxVsP", "E deposition per unit length vs. p",100, 0, 100);
138 fDeDxVsP->SetXTitle("p (GeV/c)");
139 //fdNdEta = new TH1F("fdNdEta", "dN/d#eta",100,-10,10);
141 // first configure the default
143 if (iResult>=0) iResult=ConfigureFromCDBTObjString(fgkOCDBEntry);
145 // configure from the command line parameters if specified
146 if (iResult>=0 && argc>0) iResult=ConfigureFromArgumentString(argc, argv);
151 int AliHLTTPCTrackHistoComponent::DoDeinit(){
152 // see header file for class documentation
157 delete fMultiplicity;
158 delete fMeanMultiplicity;
165 int AliHLTTPCTrackHistoComponent::Reconfigure(const char* cdbEntry, const char* /*chainId*/){
166 // see header file for class documentation
168 // configure from the specified antry or the default one
169 const char* entry=cdbEntry;
170 if (!entry || entry[0]==0) {
174 return ConfigureFromCDBTObjString(entry);
178 int AliHLTTPCTrackHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/){
179 // see header file for class documentation
181 if(GetFirstInputBlock(kAliHLTDataTypeSOR) || GetFirstInputBlock(kAliHLTDataTypeEOR)) return 0;
185 const AliHLTComponentBlockData *iter = NULL;
188 //----------------- loop over cluster blocks ---------------------//
190 Int_t totalSpacePoints = 0;
192 for(iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); iter != NULL; iter = GetNextInputBlock()){
194 if(iter->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType) continue;
196 AliHLTUInt8_t minSlice = AliHLTTPCDefinitions::GetMinSliceNr(*iter);
197 AliHLTUInt8_t minPartition = AliHLTTPCDefinitions::GetMinPatchNr(*iter);
198 //HLTDebug("Input Data - TPC cluster - slice/partition: %d/%d.", minSlice, minPartition);
200 const AliHLTTPCClusterData* clusterData = (const AliHLTTPCClusterData*)iter->fPtr;
201 Int_t nSpacepoint = (Int_t)clusterData->fSpacePointCnt;
202 totalSpacePoints += nSpacepoint;
203 HLTDebug("TrackHisto component found %d spacepoints in slice %d partition %d", nSpacepoint, minSlice, minPartition);
205 AliHLTTPCSpacePointData *clusters = (AliHLTTPCSpacePointData*)clusterData->fSpacePoints;
207 if(fClustersArray[minSlice][minPartition] != NULL){
208 //delete(fClustersArray[minSlice][minPartition]);
209 fClustersArray[minSlice][minPartition] = NULL;
212 // fill the array with AliHLTTPCSpacePointData pointers
213 // it will be used in the track loop to access information
214 // for the used clusters only
215 fClustersArray[minSlice][minPartition] = clusters;
216 fNSpacePoints[minSlice][minPartition] = nSpacepoint;
218 if(nSpacepoint==0) fClustersArray[minSlice][minPartition] = NULL;
220 } // end of loop over cluster data blocks
222 HLTDebug("TrackHisto found %d spacepoints",totalSpacePoints);
227 //----------------- loop over merged tracks -------------------//
229 Int_t totalTracks = 0;
231 for(iter = GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC); iter != NULL; iter = GetNextInputBlock()){
232 if(iter->fDataType != (kAliHLTDataTypeTrack|kAliHLTDataOriginTPC)) continue;
233 ReadTracks(iter,totalTracks);
236 HLTDebug("TrackHisto found %d tracks", totalTracks);
238 fMultiplicity->Fill(totalTracks);
240 fNtotTracks += totalTracks;
242 if(fNEvents%fEvtMod==0){
243 fMeanMultiplicity->Fill(fNEvents, Float_t(fNtotTracks)/(fEvtMod));
244 //HLTInfo("-------------- Event number: %d, total tracks accummulated %d", fNEvents, fNtotTracks);
253 void AliHLTTPCTrackHistoComponent::ReadTracks(const AliHLTComponentBlockData* iter,Int_t &tt){
254 // see header file for class documentation
256 AliHLTUInt8_t slice = AliHLTTPCDefinitions::GetMinSliceNr(*iter);
258 Int_t usedSpacePoints = 0;
260 vector<AliHLTGlobalBarrelTrack> tracksVector;
261 AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(iter->fPtr), iter->fSize, tracksVector);
263 tt = tracksVector.size();
265 for(vector<AliHLTGlobalBarrelTrack>::iterator element=tracksVector.begin(); element!=tracksVector.end(); element++){
267 Double_t trackLength = 0.;
268 if(fdEdx==kTRUE) trackLength = element->GetPathLengthTo( element->GetLastPointX(), 5.0);
270 UInt_t nHits = element->GetNumberOfPoints();
271 fTracks->Fill( 1./element->OneOverPt(), element->GetSnp(), element->GetTgl(), nHits );
272 //fdNdEta->Fill(element->GetSnp());
274 Double_t totCharge = 0;
275 const UInt_t *hitnum = element->GetPoints();
276 for(UInt_t i=0; i<element->GetNumberOfPoints(); i++){
278 UInt_t idTrack = hitnum[i];
279 Int_t sliceTrack = (idTrack>>25) & 0x7f;
280 Int_t patchTrack = (idTrack>>22) & 0x7;
281 UInt_t pos = idTrack&0x3fffff;
283 if( !fClustersArray[sliceTrack][patchTrack] ) continue;
285 if(sliceTrack<0 || sliceTrack>36 || patchTrack<0 || patchTrack>5 ){
286 HLTError("Corrupted TPC cluster Id: slice %d, patch %d, cluster %d", sliceTrack, patchTrack, idTrack);
290 if(fNSpacePoints[sliceTrack][patchTrack]<=pos ){
291 HLTError("Space point array out of boundaries!");
295 totCharge += (fClustersArray[sliceTrack][patchTrack])[pos].fCharge;
297 //Float_t xyz[3]; xyz[0] = xyz[1] = xyz[2] = 0.;
299 //xyz[0] = (fClustersArray[sliceTrack][patchTrack])[pos].fX;
300 //xyz[1] = (fClustersArray[sliceTrack][patchTrack])[pos].fY;
301 //xyz[2] = (fClustersArray[sliceTrack][patchTrack])[pos].fZ;
303 //AliHLTTPCTransform::Local2Global(xyz,slice);
305 //Double_t p[2] = { xyz[1], xyz[2] };
306 //Double_t cov[3] = { (fClustersArray[sliceTrack][patchTrack])[pos].fSigmaY2, 0., (fClustersArray[sliceTrack][patchTrack])[pos].fSigmaZ2};
307 //Double_t *res = element->GetResiduals(p,cov,kFALSE);
309 //HLTInfo("resy: %f, resz: %f", res[0], res[1]);
311 //if(!res) res[0] = res[1] = -1000.;
312 //else fClusters->Fill( (fClustersArray[sliceTrack][patchTrack])[pos].fCharge, (fClustersArray[sliceTrack][patchTrack])[pos].fQMax, res[0], res[1]);
314 fClusters->Fill( (fClustersArray[sliceTrack][patchTrack])[pos].fCharge, (fClustersArray[sliceTrack][patchTrack])[pos].fQMax, -1000., -1000.);
318 if(fdEdx==kTRUE && trackLength > 0) fDeDxVsP->Fill(element->OneOverPt()*TMath::Sqrt(1.+element->GetTgl()*element->GetTgl()), totCharge/trackLength);
322 void AliHLTTPCTrackHistoComponent::PushHisto(){
323 // see header file for class documentation
325 PushBack( (TObject*)fTracks, kAliHLTDataTypeTNtuple |kAliHLTDataOriginTPC, 0x0);
326 PushBack( (TObject*)fClusters, kAliHLTDataTypeTNtuple |kAliHLTDataOriginTPC, 0x0);
327 PushBack( (TObject*)fMultiplicity, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0);
328 PushBack( (TObject*)fMeanMultiplicity, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0);
329 PushBack( (TObject*)fDeDxVsP, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0);
330 //PushBack( (TObject*)fdNdEta, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, 0x0);
333 int AliHLTTPCTrackHistoComponent::ScanConfigurationArgument(int argc, const char** argv){
334 // see header file for class documentation
336 if (argc<=0) return 0;
338 TString argument=argv[i];
341 if (argument.CompareTo("-event-modulo")==0) {
342 if (++i>=argc) return -EPROTO;
344 fEvtMod=argument.Atof();
349 if (argument.CompareTo("-buffer-size")==0) {
350 if (++i>=argc) return -EPROTO;
352 fBufferSize=argument.Atof();
357 if (argument.CompareTo("-dEdx")==0) {