]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/EVE/AliHLTTPCEVE.cxx
minor compilation warning fixed
[u/mrichter/AliRoot.git] / HLT / TPCLib / EVE / AliHLTTPCEVE.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
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 //**************************************************************************
18
19 //  @file   AliHLTTPCEVE.cxx
20 //  @author Matthias Richter
21 //  @date   2008-11-22
22 //  @brief  AliEVE bindings for the HLT TPC.
23 //  @note   
24
25 #include <cerrno>
26 #include <cassert>
27 #include "AliHLTTPCEVE.h"
28 #include "TEveManager.h"
29 #include "TEvePointSet.h"
30 #include "TEveElement.h"
31 #include "AliRawReader.h"
32 #include "AliHLTOUT.h"
33 #include "AliHLTTPCDefinitions.h"
34 #include "AliHLTTPCClusterDataFormat.h"
35 #include "AliHLTTPCSpacePointData.h"
36 #include "AliHLTTPCSpacePointContainer.h"
37 #include "TSystem.h"
38 #include "TClass.h"
39 #include "TString.h"
40 #include "TMath.h"
41
42 /** ROOT macro for the implementation of ROOT specific class methods */
43 ClassImp(AliHLTTPCEVE)
44
45 AliHLTTPCEVE::AliHLTTPCEVE()
46   : AliHLTLogging()
47 {
48   // see header file for class documentation
49   // or
50   // refer to README to build package
51   // or
52   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
53 }
54
55 AliHLTTPCEVE::~AliHLTTPCEVE()
56 {
57   // see header file for class documentation
58 }
59
60 TEvePointSet* AliHLTTPCEVE::MakePointSetFromHLTDigits(const char* path, int eventNo, TEveElement* cont, Float_t maxR) const
61 {
62   // see header file for class documentation
63   AliHLTOUT* pHLTOUT=AliHLTOUT::New(path, eventNo);
64   TEvePointSet* pointSet=MakePointSetFromHLTOUT(pHLTOUT, cont, maxR);
65   AliHLTOUT::Delete(pHLTOUT);
66   return pointSet;
67 }
68
69 TEvePointSet* AliHLTTPCEVE::MakePointSetFromHLTOUT(AliRawReader* pRawReader, TEveElement* cont, Float_t maxR) const
70 {
71   // see header file for class documentation
72   if (!pRawReader) return NULL;
73   AliHLTOUT* pHLTOUT=AliHLTOUT::New(pRawReader);
74   TEvePointSet* pointSet=MakePointSetFromHLTOUT(pHLTOUT, cont, maxR);
75   AliHLTOUT::Delete(pHLTOUT);
76   return pointSet;
77 }
78
79 TEvePointSet* AliHLTTPCEVE::MakePointSetFromHLTOUT(AliHLTOUT* pHLTOUT, TEveElement* /*cont*/, Float_t maxR) const
80 {
81   // see header file for class documentation
82   if (!pHLTOUT) return NULL;
83
84   const Int_t kMaxCl=100*160;
85   TEvePointSet* clusters = new TEvePointSet(kMaxCl);
86   if (clusters) {
87     clusters->SetOwnIds(kTRUE);
88     clusters->SetMarkerColor(2);
89     clusters->SetMarkerStyle(5);
90     
91     if (pHLTOUT->Init()>=0) {
92       for (int idx=pHLTOUT->SelectFirstDataBlock(AliHLTTPCDefinitions::fgkClustersDataType);
93            idx>=0;
94            idx=pHLTOUT->SelectNextDataBlock()) {
95         const AliHLTUInt8_t* pBuffer=NULL;
96         AliHLTUInt32_t size=0;
97         AliHLTComponentDataType dt=kAliHLTVoidDataType;
98         AliHLTUInt32_t spec=kAliHLTVoidDataSpec;
99         pHLTOUT->GetDataBlockDescription(dt, spec);
100
101         if (pHLTOUT->GetDataBuffer(pBuffer, size)>=0) {
102           int slice=AliHLTTPCDefinitions::GetMinSliceNr(spec);
103           if (slice!=AliHLTTPCDefinitions::GetMaxSliceNr(spec)) {
104             HLTWarning("cluster data array of multiple TPC slices, can not unambiguously determine phi; skipping data block of specification 0x%08x", spec);
105             continue;
106           }
107           if (size>0 && AddClusters(clusters, reinterpret_cast<const AliHLTTPCClusterData*>(pBuffer), size, slice, maxR)<0) {
108             // action if failed
109           }
110           pHLTOUT->ReleaseDataBuffer(pBuffer);
111         }
112       }
113     } else {
114       HLTError("initialization of HLTOUT handler failed");
115     }
116     TString name="HLT TPC Clusters";
117     clusters->SetName(name);
118     name.Form("N=%d", clusters->Size());
119     clusters->SetTitle(name);
120   }
121
122   return clusters;
123 }
124
125 int AliHLTTPCEVE::AddClusters(TEvePointSet* clusters, const AliHLTTPCClusterData* data, unsigned int sizeInByte, int slice, Float_t maxR) const
126 {
127   // see header file for class documentation
128   int iResult=0;
129   if (!clusters || !data) return -EINVAL;
130   Float_t phi     = ( slice + 0.5 ) * TMath::Pi() / 9.0;
131   Float_t cos     = TMath::Cos( phi );
132   Float_t sin     = TMath::Sin( phi );
133   Float_t maxRsqr = maxR*maxR; // maxR squared for a simple geometrical cut below
134
135   for (iResult=0; iResult<(int)data->fSpacePointCnt && iResult>=0; iResult++) {
136     if (reinterpret_cast<const AliHLTUInt8_t*>(data->fSpacePoints)+(iResult+1)*sizeof(AliHLTTPCSpacePointData)>reinterpret_cast<const AliHLTUInt8_t*>(data)+sizeInByte) {
137       HLTError("data missmatch: buffer of size %d does not match size of AliHLTTPCClusterData (%d) + %d tracks x AliHLTTPCSpacePointData (%d)", 
138                sizeInByte, sizeof(AliHLTTPCClusterData), data->fSpacePointCnt, sizeof(AliHLTTPCSpacePointData));
139       iResult=-ENOMSG;
140       break;
141     }
142     const AliHLTTPCSpacePointData* sp=&data->fSpacePoints[iResult];
143     if (sp->fX*sp->fX+sp->fY*sp->fY<=maxRsqr) {
144       clusters->SetNextPoint(cos*sp->fX - sin*sp->fY, sin*sp->fX + cos*sp->fY, sp->fZ);
145     }
146   }
147
148   if (iResult<0) {
149     // reset clusters
150   }
151
152   return iResult;
153 }
154
155 int AliHLTTPCEVE::AddClusters(TEvePointSet* clusters, const AliHLTSpacePointContainer* points, Float_t maxR) const
156 {
157   // add clusters from a space point collection
158   int iResult=0;
159   if (!clusters || !points) return -EINVAL;
160   Float_t maxRsqr = maxR*maxR; // maxR squared for a simple geometrical cut below
161   vector<AliHLTUInt32_t> clusterIDs;
162   points->GetClusterIDs(clusterIDs);
163   for (vector<AliHLTUInt32_t>::const_iterator element=clusterIDs.begin();
164        element!=clusterIDs.end(); element++) {
165     AliHLTUInt32_t clusterID=*element;
166     UInt_t clusterSlice =AliHLTTPCSpacePointData::GetSlice(clusterID);
167     //UInt_t clusterPart  =AliHLTTPCSpacePointData::GetPatch(clusterID);
168     //UInt_t clusterNumber=AliHLTTPCSpacePointData::GetNumber(clusterID);
169     Float_t phi     = ( clusterSlice + 0.5 ) * TMath::Pi() / 9.0;
170     Float_t cos     = TMath::Cos( phi );
171     Float_t sin     = TMath::Sin( phi );
172
173     Float_t X=points->GetX(clusterID);
174     Float_t Y=points->GetY(clusterID);
175     Float_t Z=points->GetZ(clusterID);
176     if (X*X+Y*Y<=maxRsqr) {
177       clusters->SetNextPoint(cos*X - sin*Y, sin*X + cos*Y, Z);
178     }
179   }
180
181   return iResult;
182 }
183
184 int AliHLTTPCEVE::AddPointSet(TEveManager* pEve, const AliHLTSpacePointContainer* points, Float_t maxR, const char* title) const
185 {
186   // create new point set and add to eve
187   if (!pEve) return -EINVAL;
188   int iResult=0;
189   TEvePointSet* clusters=new TEvePointSet(200000);
190   clusters->SetOwnIds(kTRUE);
191   clusters->SetMarkerColor(2);
192   clusters->SetMarkerStyle(5);
193   iResult=AddClusters(clusters, points, maxR);
194   TString name=title==NULL?"HLT TPC Clusters":title;
195   clusters->SetName(name);
196   name.Form("N=%d", clusters->Size());
197   clusters->SetTitle(name);
198   pEve->AddElement(clusters,NULL);
199   pEve->Redraw3D();
200
201   return iResult;
202 }