]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/BASE/AliHLTTrackGeometry.cxx
bugfix: suppress points without associated track point in residual calculation
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTTrackGeometry.cxx
CommitLineData
e1e03704 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 AliHLTTrackGeometry.cxx
20/// @author Matthias Richter
21/// @date 2011-05-20
22/// @brief Desciption of a track by a sequence of track points
23///
24
25#include "AliHLTTrackGeometry.h"
26#include "AliHLTSpacePointContainer.h"
27#include "TObjArray.h"
28#include "TMarker.h"
29#include "TMath.h"
54ff4c01 30#include "TH2.h"
e1e03704 31#include <memory>
32#include <iostream>
33#include <algorithm>
34
35/** ROOT macro for the implementation of ROOT specific class methods */
36ClassImp(AliHLTTrackGeometry)
37
38AliHLTTrackGeometry::AliHLTTrackGeometry()
39 : TObject(), AliHLTLogging()
40 , fTrackPoints()
54ff4c01 41 , fSelectionMasks()
e1e03704 42 , fTrackId(-1)
54ff4c01 43 , fVerbosity(0)
e1e03704 44{
45 /// standard constructor
46}
47
48AliHLTTrackGeometry::AliHLTTrackGeometry(const AliHLTTrackGeometry& src)
49 : TObject(src), AliHLTLogging()
50 , fTrackPoints(src.fTrackPoints)
54ff4c01 51 , fSelectionMasks(src.fSelectionMasks)
e1e03704 52 , fTrackId(src.fTrackId)
54ff4c01 53 , fVerbosity(src.fVerbosity)
e1e03704 54{
55 /// copy constructor
56}
57
58AliHLTTrackGeometry& AliHLTTrackGeometry::operator=(const AliHLTTrackGeometry& src)
59{
60 /// assignment operator
61 if (this!=&src) {
62 fTrackPoints.assign(src.fTrackPoints.begin(), src.fTrackPoints.end());
54ff4c01 63 fSelectionMasks.assign(src.fSelectionMasks.begin(), src.fSelectionMasks.end());
e1e03704 64 fTrackId=src.fTrackId;
54ff4c01 65 fVerbosity=src.fVerbosity;
e1e03704 66 }
67 return *this;
68}
69
70AliHLTTrackGeometry::~AliHLTTrackGeometry()
71{
72 /// destructor
73}
74
54ff4c01 75int AliHLTTrackGeometry::AddTrackPoint(const AliHLTTrackPoint& point, AliHLTUInt32_t selectionMask)
e1e03704 76{
77 /// add a track point to the list
78 vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), point);
79 if (element==fTrackPoints.end()) {
80 fTrackPoints.push_back(point);
54ff4c01 81 if (std::find(fSelectionMasks.begin(), fSelectionMasks.end(), selectionMask)==fSelectionMasks.end()) {
82 fSelectionMasks.push_back(selectionMask);
83 }
e1e03704 84 } else {
85 HLTError("track point of id %08x already existing", point.GetId());
86 return -EEXIST;
87 }
88 return 0;
89}
90
91void AliHLTTrackGeometry::Clear(Option_t * /*option*/)
92{
93 // internal cleanup
94}
95
96void AliHLTTrackGeometry::Print(Option_t *option) const
97{
98 // print info
99 Print(cout, option);
100}
101
102void AliHLTTrackGeometry::Print(ostream& out, Option_t */*option*/) const
103{
104 // print to stream
105 out << "AliHLTTrackGeometry::Print" << endl;
106}
107
108void AliHLTTrackGeometry::Draw(Option_t *option)
109{
110 /// Inherited from TObject, draw the track
111 float scale=250;
112 float center[2]={0.5,0.5};
113 int markerColor=1;
114 int markerSize=1;
3488f70b 115 int verbosity=0;
e1e03704 116
117 TString strOption(option);
118 std::auto_ptr<TObjArray> tokens(strOption.Tokenize(" "));
119 if (!tokens.get()) return;
120 for (int i=0; i<tokens->GetEntriesFast(); i++) {
121 if (!tokens->At(i)) continue;
122 const char* key="";
123 TString arg=tokens->At(i)->GetName();
124
125 key="scale=";
126 if (arg.BeginsWith(key)) {
127 arg.ReplaceAll(key, "");
128 scale=arg.Atof();
129 continue;
130 }
131 key="centerx=";
132 if (arg.BeginsWith(key)) {
133 arg.ReplaceAll(key, "");
134 center[0]=arg.Atof();
135 continue;
136 }
137 key="centery=";
138 if (arg.BeginsWith(key)) {
139 arg.ReplaceAll(key, "");
140 center[1]=arg.Atof();
141 continue;
142 }
143
144 key="markercolor=";
145 if (arg.BeginsWith(key)) {
146 arg.ReplaceAll(key, "");
147 markerColor=arg.Atoi();
148 continue;
149 }
150
151 key="markersize=";
152 if (arg.BeginsWith(key)) {
153 arg.ReplaceAll(key, "");
154 markerSize=arg.Atoi();
155 continue;
156 }
3488f70b 157
158 key="verbosity=";
159 if (arg.BeginsWith(key)) {
160 arg.ReplaceAll(key, "");
161 verbosity=arg.Atoi();
162 continue;
163 }
e1e03704 164 }
165
166 bool bFirstPoint=true;
167 float firstalpha=0.0;
168 for (vector<AliHLTTrackPoint>::const_iterator point=fTrackPoints.begin();
169 point!=fTrackPoints.end();
170 point++) {
171 float alpha=GetPlaneAlpha(point->GetId());
172 float r=GetPlaneR(point->GetId());
173 float cosa=TMath::Cos(alpha);
174 float sina=TMath::Sin(alpha);
175 float x = r*sina + point->GetU()*cosa;
176 float y =-r*cosa + point->GetU()*sina;
3488f70b 177 if (verbosity>0) {
178 HLTInfo("ID 0x%08x: x=% .4f y=% .4f alpha=% .4f", point->GetId(), r, point->GetU(), alpha);
179 }
e1e03704 180 int color=markerColor;
181 if (bFirstPoint) {
182 bFirstPoint=false;
183 TMarker* m=new TMarker(x/(2*scale)+center[0], y/(2*scale)+center[1], 29);
184 m->SetMarkerSize(2);
185 m->SetMarkerColor(2);
186 m->Draw("same");
187 firstalpha=alpha;
188 } else {
409faf4b 189 color+=int(9*TMath::Abs(alpha-firstalpha)/TMath::Pi());
e1e03704 190 }
191 TMarker* m=new TMarker(x/(2*scale)+center[0], y/(2*scale)+center[1], point->GetV()>0?2:5);
192 m->SetMarkerColor(color);
193 m->SetMarkerSize(markerSize);
194 m->Draw("same");
195 }
196}
197
54ff4c01 198int AliHLTTrackGeometry::SetAssociatedSpacePoint(UInt_t planeId, UInt_t spacepointId, int status, float fdU, float fdV)
e1e03704 199{
200 /// set the spacepoint associated with a track point
201 vector<AliHLTTrackPoint>::iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
202 if (element==fTrackPoints.end()) return -ENOENT;
203 element->SetAssociatedSpacePoint(spacepointId, status);
54ff4c01 204 element->SetResidual(0, fdU);
205 element->SetResidual(1, fdV);
e1e03704 206 return 0;
207}
208
209int AliHLTTrackGeometry::GetAssociatedSpacePoint(UInt_t planeId, UInt_t& spacepointId) const
210{
211 /// get the spacepoint associated with a track point
212 /// return status flag if found, -ENOENT if no associated spacepoint found
213 vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
214 if (element==fTrackPoints.end()) return -ENOENT;
215 if (!element->HaveAssociatedSpacePoint()) return -ENODATA;
216 return element->GetAssociatedSpacePoint(spacepointId);
217}
218
219const AliHLTTrackGeometry::AliHLTTrackPoint* AliHLTTrackGeometry::GetTrackPoint(AliHLTUInt32_t id) const
220{
221 /// get const pointer to track point
222 vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), id);
223 if (element==fTrackPoints.end()) return NULL;
224 return &(*element);
225}
226
227AliHLTTrackGeometry::AliHLTTrackPoint* AliHLTTrackGeometry::GetTrackPoint(AliHLTUInt32_t id)
228{
229 /// get const pointer to track point
230 vector<AliHLTTrackPoint>::iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), id);
231 if (element==fTrackPoints.end()) return NULL;
232 return &(*element);
233}
234
54ff4c01 235AliHLTSpacePointContainer* AliHLTTrackGeometry::ConvertToSpacePoints(bool /*bAssociated*/) const
3488f70b 236{
237 /// create a collection of all points
238 HLTError("implementation of child method missing");
239 return NULL;
240}
241
e1e03704 242int AliHLTTrackGeometry::AssociateSpacePoints(AliHLTSpacePointContainer& points)
243{
244 /// associate the track space points to the calculated track points
245 vector<AliHLTUInt32_t> ids;
246 points.GetClusterIDs(ids);
247 if (ids.size()>0) return 0;
248 int result=AssociateSpacePoints(&ids[0], ids.size(), points);
249 if (result>0) {
54ff4c01 250 HLTInfo("associated %d of %d space point(s) to track points", result, ids.size());
e1e03704 251 }
252 return result;
253}
254
255int AliHLTTrackGeometry::AssociateSpacePoints(const AliHLTUInt32_t* trackpoints, AliHLTUInt32_t nofPoints, AliHLTSpacePointContainer& points)
256{
257 /// associate the track space points to the calculated track points
258 if (nofPoints==0) return 0;
259 if (trackpoints==NULL) return -EINVAL;
260 int count=0;
261 for (int i=nofPoints-1; i>=0; i--) {
262 if (!points.Check(trackpoints[i])) {
263 HLTWarning("can not find point id %08x", trackpoints[i]);
264 continue;
265 }
266 float xyz[3]={points.GetX(trackpoints[i]), points.GetY(trackpoints[i]), points.GetZ(trackpoints[i])};
267 AliHLTUInt32_t planeId=0;
268 int result=FindMatchingTrackPoint(trackpoints[i], xyz, planeId);
269 if (result<0) {
0ca33cc6 270 if (GetVerbosity()>0) HLTWarning("no associated track point found for space point id %08x x=%f y=%f z=%f", trackpoints[i], xyz[0], xyz[1], xyz[2]);
e1e03704 271 continue;
272 } else if (result==0) {
273 HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", trackpoints[i], xyz[0], xyz[1], xyz[2]);
274 continue;
275 }
54ff4c01 276 vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
277 SetAssociatedSpacePoint(planeId, trackpoints[i], 1, xyz[1]-element->GetU(), xyz[2]-element->GetV());
e1e03704 278 if (points.GetTrackID(trackpoints[i])<0 && GetTrackId()>=0) {
279 points.SetTrackID(GetTrackId(), trackpoints[i]);
280 HLTDebug("associating unused cluster %08x with track %d", trackpoints[i], GetTrackId());
281 }
282 count++;
283 }
284 return count;
285}
286
287int AliHLTTrackGeometry::AssociateUnusedSpacePoints(AliHLTSpacePointContainer& points)
288{
289 /// associate the track space points to the calculated track points
290 int count=0;
54ff4c01 291 for (vector<AliHLTUInt32_t>::iterator mask=fSelectionMasks.begin();
292 mask!=fSelectionMasks.end(); mask++) {
293 int subcount=0;
294 const vector<AliHLTUInt32_t>* selectedPoints=points.GetClusterIDs(*mask);
295 if (!selectedPoints) {
296 HLTWarning("space point collection does not contain data for mask 0x%08x", *mask);
e1e03704 297 continue;
298 }
54ff4c01 299 for (vector<AliHLTUInt32_t>::const_iterator id=selectedPoints->begin();
300 id!=selectedPoints->end(); id++) {
301 if (points.GetTrackID(*id)>=0) continue;
302 float xyz[3]={points.GetX(*id), points.GetY(*id), points.GetZ(*id)};
303 AliHLTUInt32_t planeId=0;
304 int result=FindMatchingTrackPoint(*id, xyz, planeId);
305 if (result<0) {
306 //HLTWarning("no associated track point found for space point id %08x x=%f y=%f z=%f", *id, xyz[0], xyz[1], xyz[2]);
307 continue;
308 } else if (result==0) {
309 //HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", *id, xyz[0], xyz[1], xyz[2]);
310 continue;
311 }
312 SetAssociatedSpacePoint(planeId, *id, 1);
313 if (points.GetTrackID(*id)<0 && GetTrackId()>=0) {
314 points.SetTrackID(GetTrackId(), *id);
315 HLTDebug("associating unused cluster %08x with track %d", *id, GetTrackId());
316 }
317 subcount++;
e1e03704 318 }
54ff4c01 319 if (fVerbosity>0) {
320 HLTInfo("associated %d of %d spacepoint(s) from selection 0x%08x to track %d",
321 subcount, selectedPoints->size(), *mask, GetTrackId());
322 }
323 count+=subcount;
e1e03704 324 }
325 return count;
326}
327
54ff4c01 328int AliHLTTrackGeometry::FillResidual(int coordinate, TH2* histo) const
329{
330 // fill residual histogram
331 const vector<AliHLTTrackPoint>& trackPoints=TrackPoints();
332 for (vector<AliHLTTrackPoint>::const_iterator trackpoint=trackPoints.begin();
333 trackpoint!=trackPoints.end(); trackpoint++) {
0ca33cc6 334 if (!trackpoint->HaveAssociatedSpacePoint()) continue;
54ff4c01 335 histo->Fill(GetPlaneR(trackpoint->GetId()), trackpoint->GetResidual(coordinate));
336 }
337 return 0;
338}
339
e1e03704 340ostream& operator<<(ostream &out, const AliHLTTrackGeometry& p)
341{
342 p.Print(out);
343 return out;
344}