]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTTrackGeometry.cxx
use explicite cast from double to int
[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 <memory>
31 #include <iostream>
32 #include <algorithm>
33
34 /** ROOT macro for the implementation of ROOT specific class methods */
35 ClassImp(AliHLTTrackGeometry)
36
37 AliHLTTrackGeometry::AliHLTTrackGeometry()
38   : TObject(), AliHLTLogging()
39   , fTrackPoints()
40   , fTrackId(-1)
41 {
42   /// standard constructor
43 }
44
45 AliHLTTrackGeometry::AliHLTTrackGeometry(const AliHLTTrackGeometry& src)
46   : TObject(src), AliHLTLogging()
47   , fTrackPoints(src.fTrackPoints)
48   , fTrackId(src.fTrackId)
49 {
50   /// copy constructor
51 }
52
53 AliHLTTrackGeometry& AliHLTTrackGeometry::operator=(const AliHLTTrackGeometry& src)
54 {
55   /// assignment operator
56   if (this!=&src) {
57     fTrackPoints.assign(src.fTrackPoints.begin(), src.fTrackPoints.end());
58     fTrackId=src.fTrackId;
59   }
60   return *this;
61 }
62
63 AliHLTTrackGeometry::~AliHLTTrackGeometry()
64 {
65   /// destructor
66 }
67
68 int AliHLTTrackGeometry::AddTrackPoint(const AliHLTTrackPoint& point)
69 {
70   /// add a track point to the list
71   vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), point);
72   if (element==fTrackPoints.end()) {
73     fTrackPoints.push_back(point);
74   } else {
75     HLTError("track point of id %08x already existing", point.GetId());
76     return -EEXIST;
77   }
78   return 0;
79 }
80
81 void AliHLTTrackGeometry::Clear(Option_t * /*option*/)
82 {
83   // internal cleanup
84 }
85
86 void AliHLTTrackGeometry::Print(Option_t *option) const
87 {
88   // print info
89   Print(cout, option);
90 }
91
92 void AliHLTTrackGeometry::Print(ostream& out, Option_t */*option*/) const
93 {
94   // print to stream
95   out << "AliHLTTrackGeometry::Print" << endl;
96 }
97
98 void AliHLTTrackGeometry::Draw(Option_t *option)
99 {
100   /// Inherited from TObject, draw the track
101   float scale=250;
102   float center[2]={0.5,0.5};
103   int markerColor=1;
104   int markerSize=1;
105
106   TString strOption(option);
107   std::auto_ptr<TObjArray> tokens(strOption.Tokenize(" "));
108   if (!tokens.get()) return;
109   for (int i=0; i<tokens->GetEntriesFast(); i++) {
110     if (!tokens->At(i)) continue;
111     const char* key="";
112     TString arg=tokens->At(i)->GetName();
113
114     key="scale=";
115     if (arg.BeginsWith(key)) {
116       arg.ReplaceAll(key, "");
117       scale=arg.Atof();
118       continue;
119     }
120     key="centerx=";
121     if (arg.BeginsWith(key)) {
122       arg.ReplaceAll(key, "");
123       center[0]=arg.Atof();
124       continue;
125     }
126     key="centery=";
127     if (arg.BeginsWith(key)) {
128       arg.ReplaceAll(key, "");
129       center[1]=arg.Atof();
130       continue;
131     }
132
133     key="markercolor=";
134     if (arg.BeginsWith(key)) {
135       arg.ReplaceAll(key, "");
136       markerColor=arg.Atoi();
137       continue;
138     }
139
140     key="markersize=";
141     if (arg.BeginsWith(key)) {
142       arg.ReplaceAll(key, "");
143       markerSize=arg.Atoi();
144       continue;
145     }
146   }
147
148   bool bFirstPoint=true;
149   float firstalpha=0.0;
150   for (vector<AliHLTTrackPoint>::const_iterator point=fTrackPoints.begin();
151        point!=fTrackPoints.end(); 
152        point++) {
153     float alpha=GetPlaneAlpha(point->GetId());
154     float r=GetPlaneR(point->GetId());
155     float cosa=TMath::Cos(alpha);
156     float sina=TMath::Sin(alpha);
157     float x = r*sina + point->GetU()*cosa;
158     float y =-r*cosa + point->GetU()*sina;
159     int color=markerColor;
160     if (bFirstPoint) {
161       bFirstPoint=false;
162       TMarker* m=new TMarker(x/(2*scale)+center[0], y/(2*scale)+center[1], 29);
163       m->SetMarkerSize(2);
164       m->SetMarkerColor(2);
165       m->Draw("same");
166       firstalpha=alpha;
167     } else {
168       color+=int(9*TMath::Abs(alpha-firstalpha)/TMath::Pi());
169     }
170     TMarker* m=new TMarker(x/(2*scale)+center[0], y/(2*scale)+center[1], point->GetV()>0?2:5);
171     m->SetMarkerColor(color);
172     m->SetMarkerSize(markerSize);
173     m->Draw("same");
174   }
175 }
176
177 int AliHLTTrackGeometry::SetAssociatedSpacePoint(UInt_t planeId, UInt_t spacepointId, int status)
178 {
179   /// set the spacepoint associated with a track point
180   vector<AliHLTTrackPoint>::iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
181   if (element==fTrackPoints.end()) return -ENOENT;
182   element->SetAssociatedSpacePoint(spacepointId, status);
183   return 0;
184 }
185
186 int AliHLTTrackGeometry::GetAssociatedSpacePoint(UInt_t planeId, UInt_t& spacepointId) const
187 {
188   /// get the spacepoint associated with a track point
189   /// return status flag if found, -ENOENT if no associated spacepoint found
190   vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
191   if (element==fTrackPoints.end()) return -ENOENT;
192   if (!element->HaveAssociatedSpacePoint()) return -ENODATA;
193   return element->GetAssociatedSpacePoint(spacepointId);
194 }
195
196 const AliHLTTrackGeometry::AliHLTTrackPoint* AliHLTTrackGeometry::GetTrackPoint(AliHLTUInt32_t id) const
197 {
198   /// get const pointer to track point
199   vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), id);
200   if (element==fTrackPoints.end()) return NULL;
201   return &(*element);
202 }
203
204 AliHLTTrackGeometry::AliHLTTrackPoint* AliHLTTrackGeometry::GetTrackPoint(AliHLTUInt32_t id)
205 {
206   /// get const pointer to track point
207   vector<AliHLTTrackPoint>::iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), id);
208   if (element==fTrackPoints.end()) return NULL;
209   return &(*element);
210 }
211
212 int AliHLTTrackGeometry::AssociateSpacePoints(AliHLTSpacePointContainer& points)
213 {
214   /// associate the track space points to the calculated track points
215   vector<AliHLTUInt32_t> ids;
216   points.GetClusterIDs(ids);
217   if (ids.size()>0) return 0;
218   int result=AssociateSpacePoints(&ids[0], ids.size(), points);
219   if (result>0) {
220     HLTInfo("associated %d space point(s) to track points", result);
221   }
222   return result;
223 }
224
225 int AliHLTTrackGeometry::AssociateSpacePoints(const AliHLTUInt32_t* trackpoints, AliHLTUInt32_t nofPoints, AliHLTSpacePointContainer& points)
226 {
227   /// associate the track space points to the calculated track points
228   if (nofPoints==0) return 0;
229   if (trackpoints==NULL) return -EINVAL;
230   int count=0;
231   for (int i=nofPoints-1; i>=0; i--) {
232     if (!points.Check(trackpoints[i])) {
233       HLTWarning("can not find point id %08x", trackpoints[i]);
234       continue;
235     }
236     float xyz[3]={points.GetX(trackpoints[i]), points.GetY(trackpoints[i]), points.GetZ(trackpoints[i])};
237     AliHLTUInt32_t planeId=0;
238     int result=FindMatchingTrackPoint(trackpoints[i], xyz, planeId);
239     if (result<0) {
240       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]);
241       continue;
242     } else if (result==0) {
243       HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", trackpoints[i], xyz[0], xyz[1], xyz[2]);
244       continue;
245     }
246     SetAssociatedSpacePoint(planeId, trackpoints[i], 1);
247     if (points.GetTrackID(trackpoints[i])<0 && GetTrackId()>=0) {
248       points.SetTrackID(GetTrackId(), trackpoints[i]);
249       HLTDebug("associating unused cluster %08x with track %d", trackpoints[i], GetTrackId());
250     }
251     count++;
252   }
253   return count;
254 }
255
256 int AliHLTTrackGeometry::AssociateUnusedSpacePoints(AliHLTSpacePointContainer& points)
257 {
258   /// associate the track space points to the calculated track points
259   int count=0;
260   vector<AliHLTUInt32_t> ids;
261   points.GetClusterIDs(ids);
262   for (vector<AliHLTUInt32_t>::iterator id=ids.begin();
263        id!=ids.end(); id++) {
264     float xyz[3]={points.GetX(*id), points.GetY(*id), points.GetZ(*id)};
265     AliHLTUInt32_t planeId=0;
266     int result=FindMatchingTrackPoint(*id, xyz, planeId);
267     if (result<0) {
268       //HLTWarning("no associated track point found for space point id %08x x=%f y=%f z=%f", *id, xyz[0], xyz[1], xyz[2]);
269       continue;
270     } else if (result==0) {
271       //HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", *id, xyz[0], xyz[1], xyz[2]);
272       continue;
273     }
274     SetAssociatedSpacePoint(planeId, *id, 1);
275     if (points.GetTrackID(*id)<0 && GetTrackId()>=0) {
276       points.SetTrackID(GetTrackId(), *id);
277       HLTDebug("associating unused cluster %08x with track %d", *id, GetTrackId());
278     }
279     count++;
280   }
281   return count;
282 }
283
284 ostream& operator<<(ostream &out, const AliHLTTrackGeometry& p)
285 {
286   p.Print(out);
287   return out;
288 }