PandaRoot
PndRiemannTrack.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 //-----------------------------------------------------------
14 // File and Version Information:
15 // $Id$
16 //
17 // Description:
18 // Track on Riemann Sphere
19 // Circle parameters can be calculated from plane parameters
20 // plane(c,nx,ny,nz);
21 //
22 //
23 // Environment:
24 // Software developed for the PANDA Detector at FAIR.
25 //
26 // Author List:
27 // Sebastian Neubert TUM (original author)
28 //
29 //
30 //-----------------------------------------------------------
31 
32 #ifndef PNDRIEMANNTRACK_HH
33 #define PNDRIEMANNTRACK_HH
34 
35 // Base Class Headers ----------------
36 #include "TObject.h"
37 
38 // Collaborating Class Headers -------
39 #include <vector>
40 #include "TVectorD.h"
41 #include "TVector3.h"
42 #include "TMatrixD.h"
43 #include "TMath.h"
44 
45 // Collaborating Class Declarations --
46 #include "PndRiemannHit.h"
47 #include "PndTrack.h"
48 #include "PndTrackCand.h"
49 
50 #include "PndSttHit.h"
51 #include "PndSttTube.h"
52 
53 #include <iostream>
54 #include <iomanip>
55 #include <algorithm>
56 
57 class PndRiemannTrack : public TObject {
58  public:
59  // Constructors/Destructors ---------
61  PndRiemannTrack(PndTrackCand *trackCand);
63 
64  // Accessors -----------------------
65  const TVectorD &n() const { return fn; }
66  double c() const { return fc; }
67  const TVectorD &av() const { return fav; }
68 
69  TVectorD orig() const;
70  double dX() { return TMath::Sqrt(fcovRXY[1][1]); }
71  double dY() { return TMath::Sqrt(fcovRXY[2][2]); }
72  double r() const;
73  double dR();
74  double dip();
75  double dipangle(); //< dipangle theta
76  double dDip();
77  double Pt(double B); //< transvers momentum Pt calculated by the magnetic field B [tesla]
78  double Pl(double B);
79  double P(double B);
80  double sign() const;
81  double getSZm() const { return fm; }
82  double getSZt() const { return ft; }
83  unsigned int getNumHits() { return fHits.size(); }
84  PndRiemannHit *getHit(unsigned int i)
85  {
86  PndRiemannHit *myHit = &(fHits[i]);
87  return myHit;
88  }
89  PndRiemannHit *getLastHit() { return getHit(getNumHits() - 1); }
90  std::vector<PndRiemannHit> getHits() const { return fHits; };
91  TVector3 getPforHit(int i, double B);
92  PndTrack getPndTrack(Double_t B);
93  FairTrackParP getTrackParPForHit(Int_t i, Double_t B);
94  Int_t getCharge(Double_t B);
95 
96  double calcZPosByS(double s);
97  TVector3 calcPosByS(double s);
98  void calcStartStopAlpha();
99  double calcAlpha(PndRiemannHit *myHit);
100 
101  void calcSForHits();
102 
103  int calcIntersection(PndRiemannTrack &track, TVector3 &p1, TVector3 &p2);
104 
105  void sortHits() { std::sort(fHits.begin(), fHits.end()); }
106 
107  double weight() const { return fweight; };
108  TMatrixD covPlane() const { return fcovPlane; };
109  TMatrixD jacRXY() const { return fjacRXY; };
110  TMatrixD covRXY() const { return fcovRXY; };
111  double m() const { return fm; }
112  double mError() const { return fmError; }
113  double t() const { return ft; }
114  double tError() const { return ftError; }
115 
116  // Modifiers -----------------------
117  void addHit(PndRiemannHit &hit);
118  void addPndTrackCand(PndTrackCand *trackCand);
119  void init(double x0, double y0, double R, double dip, double z0);
120 
121  // Operations ----------------------
122  void refit(bool withErrorCalc = true);
123  double dist(PndRiemannHit *hit);
124  double distError(PndRiemannHit *hit);
125  double distCircle(PndRiemannHit *hit);
126  double ChiSquareDistCircle();
127  void szFit(bool withErrorCalc = true);
128  double calcChi2Plane();
129  double calcSZChi2(PndRiemannHit *hit); // calculates the chi2 of the track plus the additional hit
130  double szDist(PndRiemannHit *hit);
131  double szError(PndRiemannHit *hit);
132  double szChi2() const { return fChi2; };
133 
134  void SetVerbose(int i) { fVerbose = i; }
135  void SetVertexCut(double cut) { fVertexCut = cut; }
136  TVector3 calcErrorPosByS(Double_t s, Double_t dS);
137 
138  void correctSttHits();
141  void PrintHits();
142 
143  virtual void Print(std::ostream &out = std::cout)
144  {
145 
146  out << std::setprecision(6) << "Riemann Track: Radius " << r() << " +/- " << dR() << " Origin: " << orig()[0] << " +/- " << dX() << " / " << orig()[1] << " +/- " << dY()
147  << std::endl;
148  out << "RiemannTrack: Normal: " << n()[0] << "/" << n()[1] << "/" << n()[2] << " c: " << c() << std::endl;
149  out << "Dip: " << dip() << " +/- " << dDip() << " StartAlpha: " << fStartAlpha << " StopAlpha: " << fStopAlpha << std::endl;
150  out << "Chi2Plane: " << calcChi2Plane() << std::endl;
151  PrintHits();
152  }
153 
154  friend std::ostream &operator<<(std::ostream &out, PndRiemannTrack &track)
155  {
156  track.Print(out);
157  return out;
158  }
159 
160  private:
161  // Private Data Members ------------
162  TVectorD fn;
163  TVectorD fav;
164  double fc;
165 
166  double fm;
167  double ft;
168  double fmError;
169  double ftError;
170  double fChi2;
171  TMatrixD fcovPlane;
172  TMatrixD fjacRXY;
173  TMatrixD fcovRXY;
174 
175  int fVerbose;
176 
177  bool fFitDone;
178  bool fSZFitDone;
179  bool fErrorCalcDone;
180  double fweight;
181  bool ftrefit;
182  double fVertexCut;
183 
184  std::vector<PndRiemannHit> fHits;
185  Double_t fStartAlpha;
186  Double_t fStopAlpha;
187 
188  std::map<TString, Int_t> fBranchNameMap;
191 
192  void calcJacRXY();
193 
194  TVector3 calcErrorLineNorm(PndRiemannTrack &track);
195  TVector3 calcErrorLineOffset(PndRiemannTrack &track);
196  TVectorD calcErrorXY1XY2(TVector3 &line, TVector3 &dLine, TVector3 &offset, TVector3 &dOffset);
197  Double_t calcErrorS(TVector2 &XY, TVector2 &dXY, PndRiemannTrack *track);
198 
199  Int_t GetBranchId(TString branchName);
200 
201  // Private Methods -----------------
202 
203  public:
204  ClassDef(PndRiemannTrack, 3)
205 };
206 
207 #endif
208 
209 //--------------------------------------------------------------
210 // $Log$
211 //--------------------------------------------------------------
PndRiemannHit * getLastHit()
unsigned int getNumHits()
double calcSZChi2(PndRiemannHit *hit)
double distError(PndRiemannHit *hit)
PndRiemannHit correctSttSkewedHit(PndSttHit *mySttHit, PndSttTube *myTube)
TMatrixD jacRXY() const
void init(double x0, double y0, double R, double dip, double z0)
double calcZPosByS(double s)
std::vector< PndRiemannHit > getHits() const
PndRiemannHit correctSttHit(PndSttHit *mySttHit)
double t() const
double dipangle()
static T Sqrt(const T &x)
Definition: PndCAMath.h:57
double distCircle(PndRiemannHit *hit)
double szChi2() const
double szDist(PndRiemannHit *hit)
friend std::ostream & operator<<(std::ostream &out, PndRiemannTrack &track)
double tError() const
double sign() const
const TVectorD & av() const
double P(double B)
double dist(PndRiemannHit *hit)
unsigned int i
Definition: P4_F32vec4.h:33
TMatrixD covPlane() const
TVector3 calcErrorPosByS(Double_t s, Double_t dS)
FairTrackParP getTrackParPForHit(Int_t i, Double_t B)
double c() const
Int_t getCharge(Double_t B)
const TVectorD & n() const
void addPndTrackCand(PndTrackCand *trackCand)
void SetVertexCut(double cut)
double getSZt() const
void refit(bool withErrorCalc=true)
double calcChi2Plane()
void SetVerbose(int i)
double Pt(double B)
void calcStartStopAlpha()
double ChiSquareDistCircle()
double mError() const
int calcIntersection(PndRiemannTrack &track, TVector3 &p1, TVector3 &p2)
TVectorD orig() const
virtual void Print(std::ostream &out=std::cout)
double m() const
double r() const
void szFit(bool withErrorCalc=true)
void addHit(PndRiemannHit &hit)
TVector3 getPforHit(int i, double B)
double getSZm() const
PndTrack getPndTrack(Double_t B)
void correctSttHits()
TMatrixD covRXY() const
double calcAlpha(PndRiemannHit *myHit)
double weight() const
TVector3 calcPosByS(double s)
PndRiemannHit * getHit(unsigned int i)
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:64
double szError(PndRiemannHit *hit)
double Pl(double B)