]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
adding functionality to indicate regions of spacepoints, residual calculation; more...
authorrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 1 Aug 2011 10:37:25 +0000 (10:37 +0000)
committerrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 1 Aug 2011 10:37:25 +0000 (10:37 +0000)
HLT/BASE/AliHLTSpacePointContainer.cxx
HLT/BASE/AliHLTSpacePointContainer.h
HLT/BASE/AliHLTTrackGeometry.cxx
HLT/BASE/AliHLTTrackGeometry.h

index 0c6080fbcac2a393ff5541f11fc0ace295575012..9bc103b6226aacfbf481ce04ce1e7eeff5640ba3 100644 (file)
@@ -172,6 +172,13 @@ AliHLTUInt8_t* AliHLTSpacePointContainer::Alloc(int size)
   return reinterpret_cast<AliHLTUInt8_t*>(buffer->GetArray());
 }
 
+AliHLTSpacePointContainer* AliHLTSpacePointContainer::SelectByMask(AliHLTUInt32_t /*mask*/, bool /*bAlloc*/) const
+{
+  /// create a collection of clusters for a space point mask
+  /// default implementation, nothing to do
+  return NULL;
+}
+
 AliHLTSpacePointContainer* AliHLTSpacePointContainer::SelectByTrack(int /*trackId*/, bool /*bAlloc*/) const
 {
   /// create a collection of clusters for a specific track
@@ -323,6 +330,8 @@ void AliHLTSpacePointContainer::Draw(Option_t *option)
   /// Inherited from TObject, draw clusters
   float scale=250;
   float center[2]={0.5,0.5};
+  int markersize=1;
+  int markercolor=5;
 
   TString strOption(option);
   std::auto_ptr<TObjArray> tokens(strOption.Tokenize(" "));
@@ -347,6 +356,16 @@ void AliHLTSpacePointContainer::Draw(Option_t *option)
       arg.ReplaceAll(key, "");
       center[1]=arg.Atof();
     }
+    key="markersize=";
+    if (arg.BeginsWith(key)) {
+      arg.ReplaceAll(key, "");
+      markersize=arg.Atoi();
+    }
+    key="markercolor=";
+    if (arg.BeginsWith(key)) {
+      arg.ReplaceAll(key, "");
+      markercolor=arg.Atoi();
+    }
   }
 
   vector<AliHLTUInt32_t> clusterIDs;
@@ -359,13 +378,14 @@ void AliHLTSpacePointContainer::Draw(Option_t *option)
     float sinphi=TMath::Sin(clusterphi);
     float clusterx=GetX(*clusterID);
     float clustery=GetY(*clusterID);
-    cout << " " << *clusterID << " " << clusterx << " " << clustery << " " << clusterphi << endl;
+    //cout << " " << *clusterID << " " << clusterx << " " << clustery << " " << clusterphi << endl;
     // rotate
     float pointx= clusterx*sinphi + clustery*cosphi;
     float pointy=-clusterx*cosphi + clustery*sinphi;
 
     TMarker* m=new TMarker(pointx/(2*scale)+center[0], pointy/(2*scale)+center[1], 3);
-    m->SetMarkerColor(5);
+    m->SetMarkerSize(markersize);
+    m->SetMarkerColor(markercolor);
     m->Draw("same");    
   }  
 }
index 7276c74db64db3e230e49102d2f1c07dace11994..e39672dde55394b9cf58fbbb506b111cd70de32c 100644 (file)
@@ -48,6 +48,7 @@ class AliHLTSpacePointContainer : public TObject, public AliHLTLogging
   virtual int GetNumberOfSpacePoints() const;
   virtual bool Check(AliHLTUInt32_t clusterID) const;
   virtual int GetClusterIDs(vector<AliHLTUInt32_t>& tgt) const = 0;
+  virtual const vector<AliHLTUInt32_t>* GetClusterIDs(AliHLTUInt32_t /*mask*/) {return NULL;}
   virtual float GetX(AliHLTUInt32_t clusterID) const = 0;
   virtual float GetXWidth(AliHLTUInt32_t clusterID) const = 0;
   virtual float GetY(AliHLTUInt32_t clusterID) const = 0;
@@ -57,6 +58,9 @@ class AliHLTSpacePointContainer : public TObject, public AliHLTLogging
   virtual float GetCharge(AliHLTUInt32_t clusterID) const = 0;
   virtual float GetPhi(AliHLTUInt32_t /*clusterID*/) const {return 0.0;}
 
+  /// create a collection of clusters for a space point mask
+  virtual AliHLTSpacePointContainer* SelectByMask(AliHLTUInt32_t mask, bool bAlloc=false) const;
+
   /// create a collection of clusters for a specific track
   virtual AliHLTSpacePointContainer* SelectByTrack(int trackId, bool bAlloc=false) const;
 
index e5b5c2c78ae99b07c321da13d31cc02fa6e69cc9..459e8d8d4845b0761d6fbb7f47fd3ab0c090bef6 100644 (file)
@@ -27,6 +27,7 @@
 #include "TObjArray.h"
 #include "TMarker.h"
 #include "TMath.h"
+#include "TH2.h"
 #include <memory>
 #include <iostream>
 #include <algorithm>
@@ -37,7 +38,9 @@ ClassImp(AliHLTTrackGeometry)
 AliHLTTrackGeometry::AliHLTTrackGeometry()
   : TObject(), AliHLTLogging()
   , fTrackPoints()
+  , fSelectionMasks()
   , fTrackId(-1)
+  , fVerbosity(0)
 {
   /// standard constructor
 }
@@ -45,7 +48,9 @@ AliHLTTrackGeometry::AliHLTTrackGeometry()
 AliHLTTrackGeometry::AliHLTTrackGeometry(const AliHLTTrackGeometry& src)
   : TObject(src), AliHLTLogging()
   , fTrackPoints(src.fTrackPoints)
+  , fSelectionMasks(src.fSelectionMasks)
   , fTrackId(src.fTrackId)
+  , fVerbosity(src.fVerbosity)
 {
   /// copy constructor
 }
@@ -55,7 +60,9 @@ AliHLTTrackGeometry& AliHLTTrackGeometry::operator=(const AliHLTTrackGeometry& s
   /// assignment operator
   if (this!=&src) {
     fTrackPoints.assign(src.fTrackPoints.begin(), src.fTrackPoints.end());
+    fSelectionMasks.assign(src.fSelectionMasks.begin(), src.fSelectionMasks.end());
     fTrackId=src.fTrackId;
+    fVerbosity=src.fVerbosity;
   }
   return *this;
 }
@@ -65,12 +72,15 @@ AliHLTTrackGeometry::~AliHLTTrackGeometry()
   /// destructor
 }
 
-int AliHLTTrackGeometry::AddTrackPoint(const AliHLTTrackPoint& point)
+int AliHLTTrackGeometry::AddTrackPoint(const AliHLTTrackPoint& point, AliHLTUInt32_t selectionMask)
 {
   /// add a track point to the list
   vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), point);
   if (element==fTrackPoints.end()) {
     fTrackPoints.push_back(point);
+    if (std::find(fSelectionMasks.begin(), fSelectionMasks.end(), selectionMask)==fSelectionMasks.end()) {
+      fSelectionMasks.push_back(selectionMask);
+    }
   } else {
     HLTError("track point of id %08x already existing", point.GetId());
     return -EEXIST;
@@ -185,12 +195,14 @@ void AliHLTTrackGeometry::Draw(Option_t *option)
   }
 }
 
-int AliHLTTrackGeometry::SetAssociatedSpacePoint(UInt_t planeId, UInt_t spacepointId, int status)
+int AliHLTTrackGeometry::SetAssociatedSpacePoint(UInt_t planeId, UInt_t spacepointId, int status, float fdU, float fdV)
 {
   /// set the spacepoint associated with a track point
   vector<AliHLTTrackPoint>::iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
   if (element==fTrackPoints.end()) return -ENOENT;
   element->SetAssociatedSpacePoint(spacepointId, status);
+  element->SetResidual(0, fdU);
+  element->SetResidual(1, fdV);
   return 0;
 }
 
@@ -220,7 +232,7 @@ AliHLTTrackGeometry::AliHLTTrackPoint* AliHLTTrackGeometry::GetTrackPoint(AliHLT
   return &(*element);
 }
 
-AliHLTSpacePointContainer* AliHLTTrackGeometry::ConvertToSpacePoints() const
+AliHLTSpacePointContainer* AliHLTTrackGeometry::ConvertToSpacePoints(bool /*bAssociated*/) const
 {
   /// create a collection of all points
   HLTError("implementation of child method missing");
@@ -235,7 +247,7 @@ int AliHLTTrackGeometry::AssociateSpacePoints(AliHLTSpacePointContainer& points)
   if (ids.size()>0) return 0;
   int result=AssociateSpacePoints(&ids[0], ids.size(), points);
   if (result>0) {
-    HLTInfo("associated %d space point(s) to track points", result);
+    HLTInfo("associated %d of %d space point(s) to track points", result, ids.size());
   }
   return result;
 }
@@ -261,7 +273,8 @@ int AliHLTTrackGeometry::AssociateSpacePoints(const AliHLTUInt32_t* trackpoints,
       HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", trackpoints[i], xyz[0], xyz[1], xyz[2]);
       continue;
     }
-    SetAssociatedSpacePoint(planeId, trackpoints[i], 1);
+    vector<AliHLTTrackPoint>::const_iterator element = find(fTrackPoints.begin(), fTrackPoints.end(), planeId);
+    SetAssociatedSpacePoint(planeId, trackpoints[i], 1, xyz[1]-element->GetU(), xyz[2]-element->GetV());
     if (points.GetTrackID(trackpoints[i])<0 && GetTrackId()>=0) {
       points.SetTrackID(GetTrackId(), trackpoints[i]);
       HLTDebug("associating unused cluster %08x with track %d", trackpoints[i], GetTrackId());
@@ -275,30 +288,54 @@ int AliHLTTrackGeometry::AssociateUnusedSpacePoints(AliHLTSpacePointContainer& p
 {
   /// associate the track space points to the calculated track points
   int count=0;
-  vector<AliHLTUInt32_t> ids;
-  points.GetClusterIDs(ids);
-  for (vector<AliHLTUInt32_t>::iterator id=ids.begin();
-       id!=ids.end(); id++) {
-    float xyz[3]={points.GetX(*id), points.GetY(*id), points.GetZ(*id)};
-    AliHLTUInt32_t planeId=0;
-    int result=FindMatchingTrackPoint(*id, xyz, planeId);
-    if (result<0) {
-      //HLTWarning("no associated track point found for space point id %08x x=%f y=%f z=%f", *id, xyz[0], xyz[1], xyz[2]);
-      continue;
-    } else if (result==0) {
-      //HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", *id, xyz[0], xyz[1], xyz[2]);
+  for (vector<AliHLTUInt32_t>::iterator mask=fSelectionMasks.begin();
+       mask!=fSelectionMasks.end(); mask++) {
+    int subcount=0;
+    const vector<AliHLTUInt32_t>* selectedPoints=points.GetClusterIDs(*mask);
+    if (!selectedPoints) {
+      HLTWarning("space point collection does not contain data for mask 0x%08x", *mask);
       continue;
     }
-    SetAssociatedSpacePoint(planeId, *id, 1);
-    if (points.GetTrackID(*id)<0 && GetTrackId()>=0) {
-      points.SetTrackID(GetTrackId(), *id);
-      HLTDebug("associating unused cluster %08x with track %d", *id, GetTrackId());
+    for (vector<AliHLTUInt32_t>::const_iterator id=selectedPoints->begin();
+        id!=selectedPoints->end(); id++) {
+      if (points.GetTrackID(*id)>=0) continue;
+      float xyz[3]={points.GetX(*id), points.GetY(*id), points.GetZ(*id)};
+      AliHLTUInt32_t planeId=0;
+      int result=FindMatchingTrackPoint(*id, xyz, planeId);
+      if (result<0) {
+       //HLTWarning("no associated track point found for space point id %08x x=%f y=%f z=%f", *id, xyz[0], xyz[1], xyz[2]);
+       continue;
+      } else if (result==0) {
+       //HLTWarning("associated track point for space pointid %08x x=%f y=%f z=%f occupied", *id, xyz[0], xyz[1], xyz[2]);
+       continue;
+      }
+      SetAssociatedSpacePoint(planeId, *id, 1);
+      if (points.GetTrackID(*id)<0 && GetTrackId()>=0) {
+       points.SetTrackID(GetTrackId(), *id);
+       HLTDebug("associating unused cluster %08x with track %d", *id, GetTrackId());
+      }
+      subcount++;
     }
-    count++;
+    if (fVerbosity>0) {
+      HLTInfo("associated %d of %d spacepoint(s) from selection 0x%08x to track %d",
+             subcount, selectedPoints->size(), *mask, GetTrackId());
+    }
+    count+=subcount;
   }
   return count;
 }
 
+int AliHLTTrackGeometry::FillResidual(int coordinate, TH2* histo) const
+{
+  // fill residual histogram
+  const vector<AliHLTTrackPoint>& trackPoints=TrackPoints();
+  for (vector<AliHLTTrackPoint>::const_iterator trackpoint=trackPoints.begin();
+       trackpoint!=trackPoints.end(); trackpoint++) {
+    histo->Fill(GetPlaneR(trackpoint->GetId()), trackpoint->GetResidual(coordinate));
+  }
+  return 0;
+}
+
 ostream& operator<<(ostream &out, const AliHLTTrackGeometry& p)
 {
   p.Print(out);
index 96354102bf44be6fb9394d1b85c10ebad2b63806..5171afb43dbfb16c1d37dc12ddae8c85a12371c7 100644 (file)
@@ -21,6 +21,7 @@
 
 class AliHLTGlobalBarrelTrack;
 class AliHLTSpacePointContainer;
+class TH2;
 
 /**
  * @class AliHLTTrackGeometry
@@ -127,13 +128,13 @@ class AliHLTTrackGeometry : public TObject, public AliHLTLogging
   public:
     // constructor
     AliHLTTrackPoint(AliHLTUInt32_t id, float u, float v)
-      : fId(id), fU(u), fV(v), fSpacePoint(0), fSpacePointStatus(-1) {}
+      : fId(id), fU(u), fV(v), fSpacePoint(0), fdU(0.), fdV(0.), fSpacePointStatus(-1) {}
     // copy constructor
     AliHLTTrackPoint(const AliHLTTrackPoint& src)
-      : fId(src.fId), fU(src.fU), fV(src.fV), fSpacePoint(src.fSpacePoint), fSpacePointStatus(src.fSpacePointStatus) {}
+      : fId(src.fId), fU(src.fU), fV(src.fV), fSpacePoint(src.fSpacePoint), fdU(src.fdU), fdV(src.fdV), fSpacePointStatus(src.fSpacePointStatus) {}
     // assignment operator
     AliHLTTrackPoint& operator=(const AliHLTTrackPoint& src) {
-      if (this!=&src) {fId=src.fId; fU=src.fU; fV=src.fV; fSpacePoint=src.fSpacePoint; fSpacePointStatus=src.fSpacePointStatus;}
+      if (this!=&src) {fId=src.fId; fU=src.fU; fV=src.fV; fSpacePoint=src.fSpacePoint; fdU=src.fdU; fdV=src.fdV; fSpacePointStatus=src.fSpacePointStatus;}
       return *this;
     }
     ~AliHLTTrackPoint() {}
@@ -166,6 +167,18 @@ class AliHLTTrackGeometry : public TObject, public AliHLTLogging
       spacepointId=fSpacePoint; return fSpacePointStatus<0?-ENOENT:fSpacePointStatus;
     }
 
+    int SetResidual(int coordinate, float value) {
+      if (coordinate==0) fdU=value;
+      else if (coordinate==1) fdV=value;
+      return 0;
+    }
+
+    float GetResidual(int coordinate) const {
+      if (coordinate==0) return fdU;
+      else if (coordinate==1) return fdV;
+      return 0.;
+    }
+
   private:
     // standard constructor prohibited
     AliHLTTrackPoint();
@@ -174,6 +187,8 @@ class AliHLTTrackGeometry : public TObject, public AliHLTLogging
     float fU;           // u coordinate
     float fV;           // v coordinate
     AliHLTUInt32_t fSpacePoint; // associated space point id
+    float fdU;          //  of the spacepoint
+    float fdV;          // residual of the spacepoint
     int fSpacePointStatus; // space point status
   };
 
@@ -211,14 +226,15 @@ class AliHLTTrackGeometry : public TObject, public AliHLTLogging
   AliHLTTrackPoint* GetTrackPoint(AliHLTUInt32_t id);
 
   /// create a collection of all points
-  virtual AliHLTSpacePointContainer* ConvertToSpacePoints() const;
+  virtual AliHLTSpacePointContainer* ConvertToSpacePoints() const {return ConvertToSpacePoints(false);}
+  virtual AliHLTSpacePointContainer* ConvertToSpacePoints(bool bAssociated) const;
 
   /// set the spacepoint associated with a track point
   /// @param  planeId       track point
   /// @param  spacepointId  space point id to be associated with track point
   /// @param  status        status flag
   /// @return 0 on success, -ENOENT planeId not found
-  int SetAssociatedSpacePoint(AliHLTUInt32_t planeId, AliHLTUInt32_t spacepointId, int status);
+  int SetAssociatedSpacePoint(AliHLTUInt32_t planeId, AliHLTUInt32_t spacepointId, int status, float fdU=0.0, float fdV=0.0);
 
   /// get the spacepoint associated with a track point
   /// @param  planeId       id of the track point
@@ -228,6 +244,11 @@ class AliHLTTrackGeometry : public TObject, public AliHLTLogging
 
   // services
 
+  int FillResidual(int coordinate, TH2* histo) const;
+
+  void SetVerbosity(int verbosity) {fVerbosity=verbosity;}
+  int GetVerbosity() const {return fVerbosity;}
+
   /// inherited from TObject: clear the object and reset pointer references
   virtual void Clear(Option_t * /*option*/ ="");
 
@@ -240,14 +261,16 @@ class AliHLTTrackGeometry : public TObject, public AliHLTLogging
   virtual void Draw(Option_t *option="");
 
  protected:
-  int AddTrackPoint(const AliHLTTrackPoint& point);
+  int AddTrackPoint(const AliHLTTrackPoint& point, AliHLTUInt32_t selectionMask);
 
   const vector<AliHLTTrackPoint>& TrackPoints() const {return fTrackPoints;}
 
  private:
   vector<AliHLTTrackPoint> fTrackPoints; // list of points
+  vector<AliHLTUInt32_t> fSelectionMasks; // selection masks
 
   int fTrackId; // track id
+  int fVerbosity; // verbosity
 
   ClassDef(AliHLTTrackGeometry, 0)
 };