#include <EdbSigma.h>
|
double | DAL (double t, double sxy, double sz, double dz) |
|
double | DALbt (double t) |
|
double | DALmt (double t) |
|
double | DAT (double t, double sxy, double dz) |
|
double | DATbt (double t) |
|
double | DATmt (double t) |
|
double | DP (double t, double sxy, double da, double dz) |
|
double | DPLbt (double t) |
|
double | DPLmt (double t) |
|
double | DPTbt (double t) |
|
double | DPTmt (double t) |
|
void | Draw () |
|
| EdbSigma () |
| RS agreement: 0 - in-view (default); 1 - in-zone; 2-in-plate; 3-in-brick. More...
|
|
double | FDAL (double *x, double *par) |
|
double | FDALbt (double *x, double *par) |
|
double | FDAT (double *x, double *par) |
|
double | FDPLmt (double *x, double *par) |
|
void | Set0 () |
|
double | SqSum (double a, double b) |
|
virtual | ~EdbSigma () |
|
|
double | eDZbase |
| = 210. emulsion base More...
|
|
double | eDZcell |
| =1300. plate-to plate distance More...
|
|
double | eDZem |
| geometry defines the "in-view" mt/bt accuracy: More...
|
|
double | eSaPlate |
| "in-plate" basetrack angular tolrance (nonplanarity) More...
|
|
double | eSaZone |
| "in-zone" microtrack angular tolrance (incorrect shrinkage+offsets) More...
|
|
double | eSxy |
| grains measurement accuracy defines the "in view" bt and mt errors More...
|
|
double | eSxyPlate |
| additional terms defines "in-plate" resolution More...
|
|
double | eSxyZone |
| additional terms defines "in-zone" resolution More...
|
|
double | eSz |
| = 2. grains z-uncertainty "field depth" More...
|
|
◆ EdbSigma()
RS agreement: 0 - in-view (default); 1 - in-zone; 2-in-plate; 3-in-brick.
◆ ~EdbSigma()
virtual EdbSigma::~EdbSigma |
( |
| ) |
|
|
inlinevirtual |
◆ DAL()
double EdbSigma::DAL |
( |
double |
t, |
|
|
double |
sxy, |
|
|
double |
sz, |
|
|
double |
dz |
|
) |
| |
return the estimation of the longitudinal angular error of segment (theta error) due to points uncertanty
38{
40 return Sqrt( sxy*Cos(
t)*sxy*Cos(
t) + sz*Sin(
t)*sz*Sin(
t)) / (
dz/Cos(
t));
41}
brick dz
Definition: RecDispMC.C:107
TTree * t
Definition: check_shower.C:4
◆ DALbt()
double EdbSigma::DALbt |
( |
double |
t | ) |
|
return the estimation of the longitudinal angular error of basetrack (theta error)
89{
92}
double eSz
= 2. grains z-uncertainty "field depth"
Definition: EdbSigma.h:13
double DPLmt(double t)
Definition: EdbSigma.cxx:73
double DAL(double t, double sxy, double sz, double dz)
Definition: EdbSigma.cxx:37
double eDZbase
= 210. emulsion base
Definition: EdbSigma.h:17
◆ DALmt()
double EdbSigma::DALmt |
( |
double |
t | ) |
|
return the estimation of the longitudinal angular error of microtrack (theta error)
60{
63}
double eDZem
geometry defines the "in-view" mt/bt accuracy:
Definition: EdbSigma.h:16
double eSxy
grains measurement accuracy defines the "in view" bt and mt errors
Definition: EdbSigma.h:12
◆ DAT()
double EdbSigma::DAT |
( |
double |
t, |
|
|
double |
sxy, |
|
|
double |
dz |
|
) |
| |
return the estimation of the transverse angular error of segment due to points uncertanty
45{
47 return sxy/(
dz/Cos(
t));
48}
◆ DATbt()
double EdbSigma::DATbt |
( |
double |
t | ) |
|
return the estimation of the transverse angular error of basetrack
96{
99}
double DATmt(double t)
Definition: EdbSigma.cxx:66
double DAT(double t, double sxy, double dz)
Definition: EdbSigma.cxx:44
◆ DATmt()
double EdbSigma::DATmt |
( |
double |
t | ) |
|
return the estimation of the transverse angular error of microtrack
◆ DP()
double EdbSigma::DP |
( |
double |
t, |
|
|
double |
sxy, |
|
|
double |
da, |
|
|
double |
dz |
|
) |
| |
return the estimation of the position error of segment propagated to dz
52{
54 double dpt =
dz/Cos(
t)*da /Cos(
t);
55 return Sqrt(sxy*sxy + dpt*dpt);
56}
◆ DPLbt()
double EdbSigma::DPLbt |
( |
double |
t | ) |
|
return the estimation of the longitudinal position error of basetrack
◆ DPLmt()
double EdbSigma::DPLmt |
( |
double |
t | ) |
|
return the estimation of the longitudinal position error of microtrack
74{
77}
double DP(double t, double sxy, double da, double dz)
Definition: EdbSigma.cxx:51
double DALmt(double t)
Definition: EdbSigma.cxx:59
◆ DPTbt()
double EdbSigma::DPTbt |
( |
double |
t | ) |
|
return the estimation of the transverse position error of basetrack
111{
114}
double DPTmt(double t)
Definition: EdbSigma.cxx:80
◆ DPTmt()
double EdbSigma::DPTmt |
( |
double |
t | ) |
|
return the estimation of the transverse position error of microtrack
◆ Draw()
118{
119 double fmin=0, fmax=1.;
120 TCanvas *cmt =
new TCanvas(
"mt",
"Intrinsic emulsion resolutions",800,600);
121 cmt->Divide(2,2);
122
124
125 cmt->cd(1);
126 TF1 *dal =
new TF1(
"dal",pts,&
EdbSigma::FDAL,fmin,fmax,3,
"EdbSigma",
"FDAL");
127 dal->SetParameter(0,
eSxy);
128 dal->SetParameter(1,
eSz);
129 dal->SetParameter(2,
eDZem);
130 dal->SetTitle("angle longitudinal microtrack resolution");
131 dal->Draw();
132
133 cmt->cd(2);
134 TF1 *dat =
new TF1(
"dat",pts,&
EdbSigma::FDAT,fmin,fmax,2,
"EdbSigma",
"FDAT");
135 dat->SetParameter(0,
eSxy);
136 dat->SetParameter(1,
eDZem);
137 dat->SetTitle("angle transverse microtrack resolution");
138 dat->Draw();
139
140 cmt->cd(4);
141 TF1 *dplmt =
new TF1(
"dplmt",pts,&
EdbSigma::FDPLmt,fmin,fmax,0,
"EdbSigma",
"FDPLmt");
142 dplmt->SetTitle("position longitudinal microtrack resolution");
143 dplmt->Draw();
144
145
146 cmt->cd(3);
147 TF1 *dalbt =
new TF1(
"dalbt",pts,&
EdbSigma::FDALbt,fmin,fmax,0,
"EdbSigma",
"FDALbt");
148 dalbt->SetTitle("angle longitudinal basetrack resolution");
149 dalbt->Draw();
150
151}
EdbSigma()
RS agreement: 0 - in-view (default); 1 - in-zone; 2-in-plate; 3-in-brick.
Definition: EdbSigma.h:31
double FDALbt(double *x, double *par)
Definition: EdbSigma.h:55
double FDPLmt(double *x, double *par)
Definition: EdbSigma.h:54
double FDAL(double *x, double *par)
Definition: EdbSigma.h:52
double FDAT(double *x, double *par)
Definition: EdbSigma.h:53
◆ FDAL()
double EdbSigma::FDAL |
( |
double * |
x, |
|
|
double * |
par |
|
) |
| |
|
inline |
52{
return DAL(*x, par[0], par[1], par[2]);}
◆ FDALbt()
double EdbSigma::FDALbt |
( |
double * |
x, |
|
|
double * |
par |
|
) |
| |
|
inline |
double DALbt(double t)
Definition: EdbSigma.cxx:88
◆ FDAT()
double EdbSigma::FDAT |
( |
double * |
x, |
|
|
double * |
par |
|
) |
| |
|
inline |
53{
return DAT(*x, par[0], par[1]);}
◆ FDPLmt()
double EdbSigma::FDPLmt |
( |
double * |
x, |
|
|
double * |
par |
|
) |
| |
|
inline |
◆ Set0()
◆ SqSum()
double EdbSigma::SqSum |
( |
double |
a, |
|
|
double |
b |
|
) |
| |
|
inline |
36{
return Sqrt(
a*
a + b*b);}
void a()
Definition: check_aligned.C:59
◆ eDZbase
◆ eDZcell
=1300. plate-to plate distance
◆ eDZem
geometry defines the "in-view" mt/bt accuracy:
= 45.; emulsion thickness
◆ eSaPlate
double EdbSigma::eSaPlate |
"in-plate" basetrack angular tolrance (nonplanarity)
◆ eSaZone
"in-zone" microtrack angular tolrance (incorrect shrinkage+offsets)
◆ eSxy
grains measurement accuracy defines the "in view" bt and mt errors
= 0.35 mean in-view grain accuracy xy
◆ eSxyPlate
double EdbSigma::eSxyPlate |
additional terms defines "in-plate" resolution
"in-plate" basetrack position tolerance (XY plate deformations)
◆ eSxyZone
double EdbSigma::eSxyZone |
additional terms defines "in-zone" resolution
"in-zone" microtrack position tolerance (local deformations)
◆ eSz
= 2. grains z-uncertainty "field depth"
The documentation for this class was generated from the following files: