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 The TPC conformal mapping tracker component.
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()
63 , fMeanMultiplicity(NULL)
68 //, fNClusterVsXY(NULL)
70 //, fClustersArray(NULL)
71 //, fNSpacePoints(NULL)
73 // see header file for class documentation
75 // refer to README to build package
77 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
80 const char* AliHLTTPCTrackHistoComponent::fgkOCDBEntry="HLT/ConfigTPC/TPCTrackHisto";
82 AliHLTTPCTrackHistoComponent::~AliHLTTPCTrackHistoComponent(){
83 // see header file for class documentation
86 // Public functions to implement AliHLTComponent's interface.
87 // These functions are required for the registration process
89 const char* AliHLTTPCTrackHistoComponent::GetComponentID(){
90 // see header file for class documentation
92 return "TPCTrackHisto";
95 void AliHLTTPCTrackHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list){
96 // see header file for class documentation
99 list.push_back(AliHLTTPCDefinitions::fgkClustersDataType|kAliHLTDataOriginTPC);
100 list.push_back(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
101 //list.push_back(kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC);
104 AliHLTComponentDataType AliHLTTPCTrackHistoComponent::GetOutputDataType(){
105 // see header file for class documentation
106 return kAliHLTMultipleDataType;
109 int AliHLTTPCTrackHistoComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList){
110 // see header file for class documentation
112 tgtList.push_back(kAliHLTDataTypeTNtuple|kAliHLTDataOriginTPC);
113 tgtList.push_back(kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC);
114 return tgtList.size();
117 void AliHLTTPCTrackHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ){
118 // see header file for class documentation
121 inputMultiplier = 1;// XXX TODO: Find more realistic value
124 AliHLTComponent* AliHLTTPCTrackHistoComponent::Spawn(){
125 // see header file for class documentation
126 return new AliHLTTPCTrackHistoComponent;
129 int AliHLTTPCTrackHistoComponent::DoInit( int argc, const char** argv ){
130 // see header file for class documentation
132 fClusters = new TNtuple("fClusters", "fClusters", "charge:qmax:residualY:residualZ");
133 fTracks = new TNtuple("fTracks", "fTracks", "pt:eta:psi:nclusters");
135 fClusters->SetCircular(5000);
136 fTracks->SetCircular(5000);
138 //fTracksArray = new AliHLTTPCTrackArray();
140 fMultiplicity = new TH1F("fMultiplicity", "Track multiplicity per event", 1000, 0, 1000);
141 fMeanMultiplicity = new TH1F("fMeanMultiplicity", "Mean track multiplicity vs. #evt", 10000/fEvtMod, 0, 10000);
142 fDeDxVsP = new TProfile("fDeDxVsP", "E deposition per unit length vs. p",100, 0, 100);
143 fDeDxVsP->SetXTitle("p (GeV/c)");
145 // first configure the default
147 if (iResult>=0) iResult=ConfigureFromCDBTObjString(fgkOCDBEntry);
149 // configure from the command line parameters if specified
150 if (iResult>=0 && argc>0)
151 iResult=ConfigureFromArgumentString(argc, argv);
156 int AliHLTTPCTrackHistoComponent::DoDeinit(){
157 // see header file for class documentation
161 //delete fTracksArray;
163 delete fMultiplicity;
164 delete fMeanMultiplicity;
170 int AliHLTTPCTrackHistoComponent::Reconfigure(const char* cdbEntry, const char* /*chainId*/){
171 // see header file for class documentation
173 // configure from the specified antry or the default one
174 const char* entry=cdbEntry;
175 if (!entry || entry[0]==0) {
179 return ConfigureFromCDBTObjString(entry);
183 int AliHLTTPCTrackHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/){
184 // see header file for class documentation
186 if(GetFirstInputBlock(kAliHLTDataTypeSOR) || GetFirstInputBlock(kAliHLTDataTypeEOR)) return 0;
190 //if(!fTracksArray) fTracksArray = new AliHLTTPCTrackArray();
192 const AliHLTComponentBlockData *iter = NULL;
195 //----------------- loop over cluster blocks ---------------------//
197 Int_t totalSpacePoints = 0;
199 for(iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); iter != NULL; iter = GetNextInputBlock()){
201 if(iter->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType) continue;
203 AliHLTUInt8_t minSlice = AliHLTTPCDefinitions::GetMinSliceNr(*iter);
204 AliHLTUInt8_t minPartition = AliHLTTPCDefinitions::GetMinPatchNr(*iter);
205 //HLTDebug("Input Data - TPC cluster - slice/partition: %d/%d.", minSlice, minPartition);
207 const AliHLTTPCClusterData* clusterData = (const AliHLTTPCClusterData*)iter->fPtr;
208 Int_t nSpacepoint = (Int_t)clusterData->fSpacePointCnt;
209 totalSpacePoints += nSpacepoint;
210 HLTDebug("TrackHisto component found %d spacepoints in slice %d partition %d", nSpacepoint, minSlice, minPartition);
212 AliHLTTPCSpacePointData *clusters = (AliHLTTPCSpacePointData*)clusterData->fSpacePoints;
214 if(fClustersArray[minSlice][minPartition] != NULL){
215 //delete(fClustersArray[minSlice][minPartition]);
216 fClustersArray[minSlice][minPartition] = NULL;
219 // fill the array with AliHLTTPCSpacePointData pointers
220 // it will be used in the track loop to access information
221 // for the used clusters only
222 fClustersArray[minSlice][minPartition] = clusters;
223 fNSpacePoints[minSlice][minPartition] = nSpacepoint;
225 if(nSpacepoint==0) fClustersArray[minSlice][minPartition] = NULL;
227 } // end of loop over cluster data blocks
229 HLTDebug("TrackHisto found %d spacepoints",totalSpacePoints);
234 //----------------- loop over merged tracks -------------------//
236 Int_t totalTracks = 0;
238 for(iter = GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC); iter != NULL; iter = GetNextInputBlock()){
239 if(iter->fDataType != (kAliHLTDataTypeTrack|kAliHLTDataOriginTPC)) continue;
240 ReadTracks(iter,totalTracks);
243 HLTDebug("TrackHisto found %d tracks", totalTracks);
245 fMultiplicity->Fill(totalTracks);
247 fNtotTracks += totalTracks;
249 if(fNEvents%fEvtMod==0){
250 fMeanMultiplicity->Fill(fNEvents, Float_t(fNtotTracks)/(fEvtMod));
251 //HLTInfo("-------------- Event number: %d, total tracks accummulated %d", fNEvents, fNtotTracks);
260 void AliHLTTPCTrackHistoComponent::ReadTracks(const AliHLTComponentBlockData* iter,Int_t &tt){
261 // see header file for class documentation
263 AliHLTUInt8_t slice = AliHLTTPCDefinitions::GetMinSliceNr(*iter);
264 AliHLTUInt8_t partition = AliHLTTPCDefinitions::GetMinPatchNr(*iter);
266 if( slice < fMinSlice ) fMinSlice = slice;
267 if( slice > fMaxSlice ) fMaxSlice = slice;
268 if( partition < fMinPartition ) fMinPartition = partition;
269 if( partition > fMaxPartition ) fMaxPartition = partition;
271 Int_t usedSpacePoints = 0;
273 vector<AliHLTGlobalBarrelTrack> tracksVector;
274 AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(iter->fPtr), iter->fSize, tracksVector);
276 tt = tracksVector.size();
278 for(vector<AliHLTGlobalBarrelTrack>::iterator element=tracksVector.begin(); element!=tracksVector.end(); element++){
280 Double_t trackLength = 0.;
281 if(fdEdx==kTRUE) trackLength = element->GetPathLengthTo( element->GetLastPointX(), 5.0);
283 UInt_t nHits = element->GetNumberOfPoints();
284 fTracks->Fill( 1./element->OneOverPt(), element->GetSnp(), element->GetTgl(), nHits );
286 Double_t totCharge = 0;
287 const UInt_t *hitnum = element->GetPoints();
288 for(UInt_t i=0; i<element->GetNumberOfPoints(); i++){
290 UInt_t idTrack = hitnum[i];
291 Int_t sliceTrack = (idTrack>>25) & 0x7f;
292 Int_t patchTrack = (idTrack>>22) & 0x7;
293 UInt_t pos = idTrack&0x3fffff;
295 if( !fClustersArray[sliceTrack][patchTrack] ) continue;
296 if( fNSpacePoints[sliceTrack][patchTrack]<pos ) HLTError("Space point array out of boundaries!");
298 totCharge += (fClustersArray[sliceTrack][patchTrack])[pos].fCharge;
300 Float_t xyz[3]; xyz[0] = xyz[1] = xyz[2] = 0.;
302 xyz[0] = (fClustersArray[sliceTrack][patchTrack])[pos].fX;
303 xyz[1] = (fClustersArray[sliceTrack][patchTrack])[pos].fY;
304 xyz[2] = (fClustersArray[sliceTrack][patchTrack])[pos].fZ;
306 AliHLTTPCTransform::Local2Global(xyz,slice);
308 Double_t p[2] = { xyz[1], xyz[2]};
309 Double_t cov[3] = { (fClustersArray[sliceTrack][patchTrack])[pos].fSigmaY2, 0., (fClustersArray[sliceTrack][patchTrack])[pos].fSigmaZ2};
310 Double_t *res = element->GetResiduals(p,cov,kFALSE);
312 //HLTInfo("resy: %f, resz: %f", res[0], res[1]);
314 if(!res) res[0] = res[1] = -1000.;
315 else fClusters->Fill( (fClustersArray[sliceTrack][patchTrack])[pos].fCharge, (fClustersArray[sliceTrack][patchTrack])[pos].fQMax, res[0], res[1]);
319 if(fdEdx==kTRUE && trackLength > 0) fDeDxVsP->Fill(element->OneOverPt()*TMath::Sqrt(1.+element->GetTgl()*element->GetTgl()), totCharge/trackLength);
323 void AliHLTTPCTrackHistoComponent::PushHisto(){
324 // see header file for class documentation
326 AliHLTUInt32_t fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(fMinSlice,fMaxSlice,fMinPartition,fMaxPartition);
328 PushBack( (TObject*)fTracks, kAliHLTDataTypeTNtuple |kAliHLTDataOriginTPC, fSpecification);
329 PushBack( (TObject*)fClusters, kAliHLTDataTypeTNtuple |kAliHLTDataOriginTPC, fSpecification);
330 PushBack( (TObject*)fMultiplicity, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, fSpecification);
331 PushBack( (TObject*)fMeanMultiplicity, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, fSpecification);
332 PushBack( (TObject*)fDeDxVsP, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, fSpecification);
335 int AliHLTTPCTrackHistoComponent::ScanConfigurationArgument(int argc, const char** argv){
336 // see header file for class documentation
338 if (argc<=0) return 0;
340 TString argument=argv[i];
343 if (argument.CompareTo("-event-modulo")==0) {
344 if (++i>=argc) return -EPROTO;
346 fEvtMod=argument.Atof();
351 if (argument.CompareTo("-buffer-size")==0) {
352 if (++i>=argc) return -EPROTO;
354 fBufferSize=argument.Atof();
359 if (argument.CompareTo("-dEdx")==0) {