FairRoot/PandaRoot
Public Member Functions | Private Attributes | List of all members
PndTrkCombiLegendreTransform Class Reference

#include <PndTrkCombiLegendreTransform.h>

Inheritance diagram for PndTrkCombiLegendreTransform:

Public Member Functions

 PndTrkCombiLegendreTransform ()
 
 ~PndTrkCombiLegendreTransform ()
 
void SetUpLegendreHisto ()
 
void SetUpLegendreHisto (double thetaNofBin, double thetaMin, double thetaMax, double rNofBin, double rMin, double rMax)
 
void ResetLegendreHisto ()
 
TH2F * GetLegendreHisto ()
 
void FillLegendreHisto (double x1, double y1, double radius1, double x2, double y2, double radius2)
 
void FillHisto (TH2F *histo, double x1, double y1, double r1, double x2, double y2, double r2)
 
void ComputeThetaR (double x1, double y1, double r1, double x2, double y2, double r2, double &theta, double &r)
 
int ExtractLegendreMaximum (double &theta_max, double &r_max)
 
int ExtractMaximumFromHisto (TH2F *histo, double &theta_max, double &r_max)
 
void ExtractLegendreSingleLineParameters (double &slope, double &intercept)
 
void ExtractLineParameters (double theta, double r, double &slope, double &intercept)
 
void Draw ()
 

Private Attributes

Double_t fThetaNofBin
 
Double_t fThetaMin
 
Double_t fThetaMax
 
Double_t fRNofBin
 
Double_t fRMin
 
Double_t fRMax
 
TH2F * fhLegendre
 

Detailed Description

Definition at line 7 of file PndTrkCombiLegendreTransform.h.

Constructor & Destructor Documentation

PndTrkCombiLegendreTransform::PndTrkCombiLegendreTransform ( )
PndTrkCombiLegendreTransform::~PndTrkCombiLegendreTransform ( )

Definition at line 40 of file PndTrkCombiLegendreTransform.cxx.

40  {
41 
42 }

Member Function Documentation

void PndTrkCombiLegendreTransform::ComputeThetaR ( double  x1,
double  y1,
double  r1,
double  x2,
double  y2,
double  r2,
double &  theta,
double &  r 
)

Definition at line 74 of file PndTrkCombiLegendreTransform.cxx.

References alpha, CAMath::ATan2(), CAMath::Cos(), dx, dy, fabs(), ir, Pi, r1, r2, CAMath::Sin(), and CAMath::Sqrt().

Referenced by PndTrkCombiLegendreTask::ComputePlaneExtremities().

75 {
76 
77  double dx = x1 - x2;
78  double dy = y1 - y2;
79 
80  if(r1 == 0 && r2 == 0) { // two points
81  theta = TMath::ATan2(-dx, dy);
82  r = x1 * TMath::Cos(theta) + y1 * TMath::Sin(theta);
83  // cout << "2points " << theta * TMath::RadToDeg() << "/ r " << r << endl;
84  }
85  else if(r1 > 0 && r2 > 0) { // two circles
86  double dr[4] = {r1 - r2, r1 + r2, -r1 - r2, -r1 + r2};
87  for(int ir = 0; ir < 4; ir++) {
88 
89  double alpha = TMath::ATan2(dy, dx);
90  if(dy < 0) alpha += (2 * TMath::Pi());
91 
92  double distance = TMath::Sqrt(dx * dx + dy * dy);
93  theta = alpha + TMath::ACos(-dr[ir]/distance);
94 
95  while(theta < 0) theta += TMath::Pi();
96  while(theta > TMath::Pi()) theta -= TMath::Pi();
97 
98  // try all the combinations (++ +- -+ --) and...
99  double rA[2], rB[2];
100  rA[0] = x1 * TMath::Cos(theta) + y1 * TMath::Sin(theta) - r1;
101  rA[1] = x1 * TMath::Cos(theta) + y1 * TMath::Sin(theta) + r1;
102  rB[0] = x2 * TMath::Cos(theta) + y2 * TMath::Sin(theta) - r2;
103  rB[1] = x2 * TMath::Cos(theta) + y2 * TMath::Sin(theta) + r2;
104 
105  // ...get the one where curve A and curve B have the same r --> they cross
106  if(fabs(rA[0] - rB[0]) < 1e-10 || fabs(rA[0] - rB[1]) < 1e-10) r = rA[0];
107  else if(fabs(rA[1] - rB[0]) < 1e-10 || fabs(rA[1] - rB[1]) < 1e-10) r = rA[1];
108 
109  // cout << "2circles " << theta * TMath::RadToDeg() << "/ r " << r << endl;
110  // cout << "ir " << ir << " theta " << theta * TMath::RadToDeg() << "/ r " << r << endl;
111  }
112  }
113  else { // one point/one circle
114  double rj;
115  if(r1 > 0) rj = r1;
116  else rj = r2;
117 
118  double dr[2] = {-rj, rj};
119  for(int ir = 0; ir < 2; ir++) {
120 
121  double alpha = TMath::ATan2(dy, dx);
122  if(dy < 0) alpha += (2 * TMath::Pi());
123 
124  double distance = TMath::Sqrt(dx * dx + dy * dy);
125  theta = alpha + TMath::ACos(dr[ir]/distance);
126 
127  while(theta < 0) theta += TMath::Pi();
128  while(theta > TMath::Pi()) theta -= TMath::Pi();
129 
130  if(r1 > 0) r = x2 * TMath::Cos(theta) + y2 * TMath::Sin(theta);
131  else r = x1 * TMath::Cos(theta) + y1 * TMath::Sin(theta);
132  // cout << "1pt/1cir " << theta * TMath::RadToDeg() << "/ r " << r << endl;
133 
134  // cout << "ir " << ir << " theta " << theta * TMath::RadToDeg() << "/ r " << r << endl;
135  }
136  }
137 }
double dy
double r
Definition: RiemannTest.C:14
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
static T Sin(const T &x)
Definition: PndCAMath.h:42
double r1
static T Cos(const T &x)
Definition: PndCAMath.h:43
static T ATan2(const T &y, const T &x)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
double dx
double alpha
Definition: f_Init.h:9
Double_t Pi
double r2
static int ir
Definition: ranlxd.cxx:374
void PndTrkCombiLegendreTransform::Draw ( )
int PndTrkCombiLegendreTransform::ExtractLegendreMaximum ( double &  theta_max,
double &  r_max 
)

Definition at line 148 of file PndTrkCombiLegendreTransform.cxx.

Referenced by PndTrkCombiLegendreTask::LegendreFit().

148  {
149  ExtractMaximumFromHisto(fhLegendre, theta_max, r_max);
150 }
int ExtractMaximumFromHisto(TH2F *histo, double &theta_max, double &r_max)
void PndTrkCombiLegendreTransform::ExtractLegendreSingleLineParameters ( double &  slope,
double &  intercept 
)

Definition at line 162 of file PndTrkCombiLegendreTransform.cxx.

References r, and theta.

Referenced by PndTrkCombiLegendreTask::LegendreFit().

162  {
163  double theta, r;
164  ExtractLegendreMaximum(theta, r);
165  ExtractLineParameters(theta, r, slope, intercept);
166 }
double r
Definition: RiemannTest.C:14
void ExtractLineParameters(double theta, double r, double &slope, double &intercept)
int ExtractLegendreMaximum(double &theta_max, double &r_max)
void PndTrkCombiLegendreTransform::ExtractLineParameters ( double  theta,
double  r,
double &  slope,
double &  intercept 
)

Definition at line 168 of file PndTrkCombiLegendreTransform.cxx.

References CAMath::Sin(), and CAMath::Tan().

168  {
169  // r = x cost + y sint --> y = r/sint - x/tant
170 
171  if(theta == 0.) theta = 1.e-6; // CHECK
172 
173  slope = - 1./TMath::Tan(TMath::DegToRad() * theta);
174  intercept = r/TMath::Sin(TMath::DegToRad() * theta);
175 
176 }
double r
Definition: RiemannTest.C:14
static T Sin(const T &x)
Definition: PndCAMath.h:42
float Tan(float x)
Definition: PndCAMath.h:165
int PndTrkCombiLegendreTransform::ExtractMaximumFromHisto ( TH2F *  histo,
double &  theta_max,
double &  r_max 
)

Definition at line 152 of file PndTrkCombiLegendreTransform.cxx.

152  {
153  int bin = histo->GetMaximumBin();
154  int binx, biny, binz;
155  histo->GetBinXYZ(bin, binx, biny, binz);
156  theta_max = histo->GetXaxis()->GetBinCenter(binx);
157  r_max = histo->GetYaxis()->GetBinCenter(biny);
158  return histo->GetMaximum();
159 
160 }
void PndTrkCombiLegendreTransform::FillHisto ( TH2F *  histo,
double  x1,
double  y1,
double  r1,
double  x2,
double  y2,
double  r2 
)

Definition at line 67 of file PndTrkCombiLegendreTransform.cxx.

References r, and theta.

68 {
69  double theta, r;
70  ComputeThetaR(x1, y1, r1, x2, y2, r2, theta, r);
71  histo->Fill(theta * TMath::RadToDeg(), r);
72 }
double r
Definition: RiemannTest.C:14
double r1
void ComputeThetaR(double x1, double y1, double r1, double x2, double y2, double r2, double &theta, double &r)
double r2
void PndTrkCombiLegendreTransform::FillLegendreHisto ( double  x1,
double  y1,
double  radius1,
double  x2,
double  y2,
double  radius2 
)

Definition at line 140 of file PndTrkCombiLegendreTransform.cxx.

Referenced by PndTrkCombiLegendreTask::FillPeakCouplesHisto(), and PndTrkCombiLegendreTask::FillPeakNeighCouplesHisto().

140  {
141  FillHisto(fhLegendre, x1, y1, radius1, x2, y2, radius2);
142 }
void FillHisto(TH2F *histo, double x1, double y1, double r1, double x2, double y2, double r2)
TH2F* PndTrkCombiLegendreTransform::GetLegendreHisto ( )
inline
void PndTrkCombiLegendreTransform::ResetLegendreHisto ( )

Definition at line 61 of file PndTrkCombiLegendreTransform.cxx.

Referenced by PndTrkCombiLegendreTask::LegendreFit().

61  {
62  fhLegendre->Reset();
63  fhLegendre->SetStats(kFALSE);
64 }
void PndTrkCombiLegendreTransform::SetUpLegendreHisto ( )
void PndTrkCombiLegendreTransform::SetUpLegendreHisto ( double  thetaNofBin,
double  thetaMin,
double  thetaMax,
double  rNofBin,
double  rMin,
double  rMax 
)

Member Data Documentation

TH2F* PndTrkCombiLegendreTransform::fhLegendre
private

Definition at line 31 of file PndTrkCombiLegendreTransform.h.

Referenced by GetLegendreHisto().

Double_t PndTrkCombiLegendreTransform::fRMax
private

Definition at line 30 of file PndTrkCombiLegendreTransform.h.

Double_t PndTrkCombiLegendreTransform::fRMin
private

Definition at line 30 of file PndTrkCombiLegendreTransform.h.

Double_t PndTrkCombiLegendreTransform::fRNofBin
private

Definition at line 30 of file PndTrkCombiLegendreTransform.h.

Double_t PndTrkCombiLegendreTransform::fThetaMax
private

Definition at line 30 of file PndTrkCombiLegendreTransform.h.

Double_t PndTrkCombiLegendreTransform::fThetaMin
private

Definition at line 30 of file PndTrkCombiLegendreTransform.h.

Double_t PndTrkCombiLegendreTransform::fThetaNofBin
private

Definition at line 30 of file PndTrkCombiLegendreTransform.h.


The documentation for this class was generated from the following files: