]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTTrackGeometry.cxx
automatically added data sink components are now added directly to the internal insta...
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTTrackGeometry.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   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"
30 #include "TH2.h"
31 #include <memory>
32 #include <iostream>
33 #include <algorithm>
34
35 /** ROOT macro for the implementation of ROOT specific class methods */
36 ClassImp(AliHLTTrackGeometry)
37
38 AliHLTTrackGeometry::AliHLTTrackGeometry()
39   : TObject(), AliHLTLogging()
40   , fTrackPoints()
41   , fSelectionMasks()
42   , fTrackId(-1)
43   , fVerbosity(0)
44 {
45   /// standard constructor
46 }
47
48 AliHLTTrackGeometry::AliHLTTrackGeometry(const AliHLTTrackGeometry& src)
49   : TObject(src), AliHLTLogging()
50   , fTrackPoints(src.fTrackPoints)
51   , fSelectionMasks(src.fSelectionMasks)
52   , fTrackId(src.fTrackId)
53   , fVerbosity(src.fVerbosity)
54 {
55   /// copy constructor
56 }
57
58 AliHLTTrackGeometry& AliHLTTrackGeometry::operator=(const AliHLTTrackGeometry& src)
59 {
60   /// assignment operator
61   if (this!=&src) {
62     fTrackPoints.assign(src.fTrackPoints.begin(), src.fTrackPoints.end());
63     fSelectionMasks.assign(src.fSelectionMasks.begin(), src.fSelectionMasks.end());
64     fTrackId=src.fTrackId;
65     fVerbosity=src.fVerbosity;
66   }
67   return *this;
68 }
69
70 AliHLTTrackGeometry::~AliHLTTrackGeometry()
71 {
72   /// destructor
73 }
74
75 int AliHLTTrackGeometry::AddTrackPoint(const AliHLTTrackPoint& point, AliHLTUInt32_t selectionMask)
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);
81     if (std::find(fSelectionMasks.begin(), fSelectionMasks.end(), selectionMask)==fSelectionMasks.end()) {
82       fSelectionMasks.push_back(selectionMask);
83     }
84   } else {
85     HLTError("track point of id %08x already existing", point.GetId());
86     return -EEXIST;
87   }
88   return 0;
89 }
90
91 void AliHLTTrackGeometry::Clear(Option_t * /*option*/)
92 {
93   // internal cleanup
94 }
95
96 void AliHLTTrackGeometry::Print(Option_t *option) const
97 {
98   // print info
99   Print(cout, option);
100 }
101
102 void AliHLTTrackGeometry::Print(ostream& out, Option_t */*option*/) const
103 {
104   // print to stream
105   out << "AliHLTTrackGeometry::Print" << endl;
106 }
107
108 void 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;
115   int verbosity=0;
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     }
157
158     key="verbosity=";
159     if (arg.BeginsWith(key)) {
160       arg.ReplaceAll(key, "");
161       verbosity=arg.Atoi();
162       continue;
163     }
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;
177     if (verbosity>0) {
178       HLTInfo("ID 0x%08x: x=% .4f y=% .4f alpha=% .4f", point->GetId(), r, point->GetU(), alpha);
179     }
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 {
189       color+=int(9*TMath::Abs(alpha-firstalpha)/TMath::Pi());
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
198 int AliHLTTrackGeometry::SetAssociatedSpacePoint(UInt_t planeId, UInt_t spacepointId, int /*status*/, float dU, float dV)
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->AddAssociatedSpacePoint(spacepointId, dU, dV);
204   return 0;
205 }
206
207 int AliHLTTrackGeometry::GetAssociatedSpacePoint(UInt_t planeId, UInt_t& spacepointId) const
208 {
209   /// get the spacepoint associated with a track point
210   /// return status flag if found, -ENOENT if no associated spacepoint found
211   vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
212   if (element==fTrackPoints.end()) return -ENOENT;
213   if (!element->HaveAssociatedSpacePoint()) return -ENODATA;
214   spacepointId=(element->GetSpacepoints())[0].fId;
215   return 0;
216 }
217
218 int AliHLTTrackGeometry::RegisterTrackPoints(AliHLTTrackGrid* /*pGrid*/) const
219 {
220   /// default implementation, nothing to do
221   return -ENOSYS;
222 }
223
224 int AliHLTTrackGeometry::FillTrackPoints(AliHLTTrackGrid* /*pGrid*/) const
225 {
226   /// default implementation, nothing to do
227   return -ENOSYS;
228 }
229
230 const AliHLTTrackGeometry::AliHLTTrackPoint* AliHLTTrackGeometry::GetTrackPoint(AliHLTUInt32_t id) const
231 {
232   /// get const pointer to track point
233   vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), id);
234   if (element==fTrackPoints.end()) return NULL;
235   return &(*element);
236 }
237
238 AliHLTTrackGeometry::AliHLTTrackPoint* AliHLTTrackGeometry::GetTrackPoint(AliHLTUInt32_t id)
239 {
240   /// get const pointer to track point
241   vector<AliHLTTrackPoint>::iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), id);
242   if (element==fTrackPoints.end()) return NULL;
243   return &(*element);
244 }
245
246 AliHLTSpacePointContainer* AliHLTTrackGeometry::ConvertToSpacePoints(bool /*bAssociated*/) const
247 {
248   /// create a collection of all points
249   HLTError("implementation of child method missing");
250   return NULL;
251 }
252
253 int AliHLTTrackGeometry::AssociateSpacePoints(AliHLTSpacePointContainer& points)
254 {
255   /// associate the track space points to the calculated track points
256   vector<AliHLTUInt32_t> ids;
257   points.GetClusterIDs(ids);
258   if (ids.size()>0) return 0;
259   int result=AssociateSpacePoints(&ids[0], ids.size(), points);
260   if (result>0) {
261     HLTInfo("associated %d of %d space point(s) to track points", result, ids.size());
262   }
263   return result;
264 }
265
266 int AliHLTTrackGeometry::AssociateSpacePoints(const AliHLTUInt32_t* trackpoints, AliHLTUInt32_t nofPoints, AliHLTSpacePointContainer& points)
267 {
268   /// associate the track space points to the calculated track points
269   if (nofPoints==0) return 0;
270   if (trackpoints==NULL) return -EINVAL;
271   int count=0;
272   for (int i=nofPoints-1; i>=0; i--) {
273     if (!points.Check(trackpoints[i])) {
274       HLTWarning("can not find point id %08x", trackpoints[i]);
275       continue;
276     }
277     float xyz[3]={points.GetX(trackpoints[i]), points.GetY(trackpoints[i]), points.GetZ(trackpoints[i])};
278     AliHLTUInt32_t planeId=0;
279     int result=FindMatchingTrackPoint(trackpoints[i], xyz, planeId);
280     if (result<0) {
281       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]);
282       continue;
283     } else if (result==0) {
284       HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", trackpoints[i], xyz[0], xyz[1], xyz[2]);
285       continue;
286     }
287     vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
288     SetAssociatedSpacePoint(planeId, trackpoints[i], 1, xyz[1]-element->GetU(), xyz[2]-element->GetV());
289     if (points.GetTrackID(trackpoints[i])<0 && GetTrackId()>=0) {
290       points.SetTrackID(GetTrackId(), trackpoints[i]);
291       HLTDebug("associating unused cluster %08x with track %d", trackpoints[i], GetTrackId());
292     }
293     count++;
294   }
295   return count;
296 }
297
298 int AliHLTTrackGeometry::AssociateUnusedSpacePoints(AliHLTSpacePointContainer& points)
299 {
300   /// associate the track space points to the calculated track points
301   int count=0;
302   for (vector<AliHLTUInt32_t>::iterator mask=fSelectionMasks.begin();
303        mask!=fSelectionMasks.end(); mask++) {
304     int subcount=0;
305     const vector<AliHLTUInt32_t>* selectedPoints=points.GetClusterIDs(*mask);
306     if (!selectedPoints) {
307       HLTWarning("space point collection does not contain data for mask 0x%08x", *mask);
308       continue;
309     }
310     for (vector<AliHLTUInt32_t>::const_iterator id=selectedPoints->begin();
311          id!=selectedPoints->end(); id++) {
312       if (points.GetTrackID(*id)>=0) continue;
313       float xyz[3]={points.GetX(*id), points.GetY(*id), points.GetZ(*id)};
314       AliHLTUInt32_t planeId=0;
315       int result=FindMatchingTrackPoint(*id, xyz, planeId);
316       if (result<0) {
317         //HLTWarning("no associated track point found for space point id %08x x=%f y=%f z=%f", *id, xyz[0], xyz[1], xyz[2]);
318         continue;
319       } else if (result==0) {
320         //HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", *id, xyz[0], xyz[1], xyz[2]);
321         continue;
322       }
323       SetAssociatedSpacePoint(planeId, *id, 1);
324       if (points.GetTrackID(*id)<0 && GetTrackId()>=0) {
325         points.SetTrackID(GetTrackId(), *id);
326         HLTDebug("associating unused cluster %08x with track %d", *id, GetTrackId());
327       }
328       subcount++;
329     }
330     if (fVerbosity>0) {
331       HLTInfo("associated %d of %d spacepoint(s) from selection 0x%08x to track %d",
332               subcount, selectedPoints->size(), *mask, GetTrackId());
333     }
334     count+=subcount;
335   }
336   return count;
337 }
338
339 int AliHLTTrackGeometry::FillResidual(int coordinate, TH2* histo) const
340 {
341   // fill residual histogram
342   const vector<AliHLTTrackPoint>& trackPoints=TrackPoints();
343   for (vector<AliHLTTrackPoint>::const_iterator trackpoint=trackPoints.begin();
344        trackpoint!=trackPoints.end(); trackpoint++) {
345     if (!trackpoint->HaveAssociatedSpacePoint()) continue;
346     for (vector<AliHLTTrackSpacepoint>::const_iterator sp=(trackpoint->GetSpacepoints()).begin();
347          sp!=(trackpoint->GetSpacepoints()).end(); sp++) {
348       histo->Fill(GetPlaneR(trackpoint->GetId()), sp->GetResidual(coordinate));
349     }
350   }
351   return 0;
352 }
353
354 ostream& operator<<(ostream &out, const AliHLTTrackGeometry& p)
355 {
356   p.Print(out);
357   return out;
358 }