CMS 3D CMS Logo

List of all members | Static Public Member Functions | Private Member Functions | Static Private Member Functions
SiStripFineDelayTOF Class Reference

#include <SiStripFineDelayTOF.h>

Static Public Member Functions

static double timeOfFlight (bool cosmics, bool field, double *trackParameters, double *hit, double *phit, bool onDisk)
 
static void trackParameters (const reco::Track &tk, double *trackParameters)
 

Private Member Functions

 SiStripFineDelayTOF ()=delete
 
virtual ~SiStripFineDelayTOF ()=delete
 

Static Private Member Functions

static double getPhi (double *trackParameters, double *hit, bool onDisk)
 
static double timeOfFlightBeam (double *hit, double *phit)
 
static double timeOfFlightBeamB (double *trackParameters, double *hit, double *phit, bool onDisk)
 
static double timeOfFlightCosmic (double *hit, double *phit)
 
static double timeOfFlightCosmicB (double *trackParameters, double *hit, double *phit, bool onDisk)
 
static double x (double *trackParameters, double phi)
 
static double y (double *trackParameters, double phi)
 
static double z (double *trackParameters, double phi)
 

Detailed Description

Definition at line 10 of file SiStripFineDelayTOF.h.

Constructor & Destructor Documentation

◆ SiStripFineDelayTOF()

SiStripFineDelayTOF::SiStripFineDelayTOF ( )
privatedelete

◆ ~SiStripFineDelayTOF()

virtual SiStripFineDelayTOF::~SiStripFineDelayTOF ( )
privatevirtualdelete

Member Function Documentation

◆ getPhi()

double SiStripFineDelayTOF::getPhi ( double *  trackParameters,
double *  hit,
bool  onDisk 
)
staticprivate

Definition at line 114 of file SiStripFineDelayTOF.cc.

114  {
115  if (onDisk) //use z coordinate to find phi
116  return trackParameters[2] + trackParameters[0] * tan(trackParameters[1]) * (hit[2] - trackParameters[4]);
117  else // use x,y coordinate to find phi
118  {
119  double phi = 0;
120  if (trackParameters[2] > -2.35 && trackParameters[2] < -0.78) {
121  // first use y
122  phi = acos(-trackParameters[0] *
123  (hit[1] - (trackParameters[3] - 1 / trackParameters[0]) * (-cos(trackParameters[2]))));
124  // use x to resolve the ambiguity
125  if (fabs(x(trackParameters, phi) - hit[0]) > fabs(x(trackParameters, -phi) - hit[0]))
126  phi *= -1.;
127  } else {
128  // first use x
129  phi =
130  asin(trackParameters[0] * (hit[0] - (trackParameters[3] - 1 / trackParameters[0]) * sin(trackParameters[2])));
131  // use y to resolve the ambiguity
132  if ((fabs(y(trackParameters, phi) - hit[1]) > fabs(y(trackParameters, TMath::Pi() - phi) - hit[1])))
133  phi = phi > 0 ? TMath::Pi() - phi : -TMath::Pi() - phi;
134  }
135  return phi;
136  }
137  return 0.;
138 }

References funct::cos(), phi, Pi, funct::sin(), funct::tan(), trackParameters(), x(), and y().

Referenced by timeOfFlightCosmicB().

◆ timeOfFlight()

double SiStripFineDelayTOF::timeOfFlight ( bool  cosmics,
bool  field,
double *  trackParameters,
double *  hit,
double *  phit,
bool  onDisk 
)
static

Definition at line 10 of file SiStripFineDelayTOF.cc.

11  {
12  // case 1: cosmics with no field.
13  if (cosmics && !field) {
14  return timeOfFlightCosmic(hit, phit);
15  }
16  // case 2: cosmics with field.
17  else if (cosmics && field) {
18  return timeOfFlightCosmicB(trackParameters, hit, phit, onDisk);
19  }
20  // case 3: beam with no field.
21  else if (!cosmics && !field) {
22  return timeOfFlightBeam(hit, phit);
23  }
24  // case 4: beam with field.
25  else {
26  return timeOfFlightBeamB(trackParameters, hit, phit, onDisk);
27  }
28 }

References GlobalTrackerMuonAlignment_cfi::cosmics, timeOfFlightBeam(), timeOfFlightBeamB(), timeOfFlightCosmic(), timeOfFlightCosmicB(), and trackParameters().

Referenced by SiStripFineDelayHit::detId().

◆ timeOfFlightBeam()

double SiStripFineDelayTOF::timeOfFlightBeam ( double *  hit,
double *  phit 
)
staticprivate

Definition at line 85 of file SiStripFineDelayTOF.cc.

85  {
86  // constants
87  const double c = 30; // cm/ns
88  TVector3 r0(hit[0], hit[1], hit[2]);
89  return r0.Mag() / c;
90 }

References HltBtagPostValidation_cff::c.

Referenced by timeOfFlight().

◆ timeOfFlightBeamB()

double SiStripFineDelayTOF::timeOfFlightBeamB ( double *  trackParameters,
double *  hit,
double *  phit,
bool  onDisk 
)
staticprivate

Definition at line 92 of file SiStripFineDelayTOF.cc.

92  {
93  // constants
94  const double c = 30; // cm/ns
95  // track parameters
96  double& theta = trackParameters[1];
97  // returns the time of flight from the origin
98  return onDisk ? fabs(hit[2] / (c * cos(theta))) : fabs(sqrt(hit[0] * hit[0] + hit[1] * hit[1]) / (c * sin(theta)));
99 }

References HltBtagPostValidation_cff::c, funct::cos(), funct::sin(), mathSSE::sqrt(), theta(), and trackParameters().

Referenced by timeOfFlight().

◆ timeOfFlightCosmic()

double SiStripFineDelayTOF::timeOfFlightCosmic ( double *  hit,
double *  phit 
)
staticprivate

Definition at line 30 of file SiStripFineDelayTOF.cc.

30  {
31  // constants
32  const double c = 30; // cm/ns
33 #ifndef TIF_COSMIC_SETUP
34  const double Rmu = 385; // cm
35  const double zmu = 560; // cm
36  // estimate the time for crossing the barrel
37  TVector3 r0(hit[0], hit[1], hit[2]);
38  TVector3 pr(phit[0], phit[1], phit[2]);
39  TVector2 r0_xy = r0.XYvector();
40  TVector2 pr_xy = pr.Unit().XYvector();
41  double t_barrel =
42  ((r0_xy * pr_xy) + sqrt((r0_xy * pr_xy) * (r0_xy * pr_xy) - r0_xy.Mod2() + Rmu * Rmu)) / c * pr.Mag() / pr.Perp();
43  // estimate the time for crossing endcaps
44  double t_endcap = fabs(((phit[2] / fabs(phit[2]) * zmu) + hit[2]) / c * pr.Mag() / pr.Pz());
45  // take the minimum
46  return t_barrel < t_endcap ? t_barrel : t_endcap;
47 #else
48  const double y_trigger = 100; //cm
49  double p = sqrt(phit[0] * phit[0] + phit[1] * phit[1] + phit[2] * phit[2]);
50  // LogDebug("TOF") << "momentum:" << phit[0] << " " << phit[1] << " " << phit[2];
51  // LogDebug("TOF") << "p/py=" << p/phit[1] << " Y0,Y,dY = " << y_trigger << " " << hit[1] << " " << (y_trigger-hit[1]);
52  // LogDebug("TOF") << "d, t=d/c : " << ((y_trigger-hit[1])*(p/phit[1])) << " " << ((y_trigger-hit[1])*(p/phit[1])/c);
53  return fabs((y_trigger - hit[1]) * (p / phit[1]) / c);
54 #endif
55 }

References HltBtagPostValidation_cff::c, AlCaHLTBitMon_ParallelJobs::p, mathSSE::sqrt(), and SiPixelPI::t_barrel.

Referenced by timeOfFlight().

◆ timeOfFlightCosmicB()

double SiStripFineDelayTOF::timeOfFlightCosmicB ( double *  trackParameters,
double *  hit,
double *  phit,
bool  onDisk 
)
staticprivate

Definition at line 57 of file SiStripFineDelayTOF.cc.

57  {
58  // constants
59  const double Rmu = 385; // cm
60  const double zmu = 560; // cm
61  const double c = 30; // cm/ns
62  // track parameters
63  double& kappa = trackParameters[0];
64  double& theta = trackParameters[1];
65  double& phi_0 = trackParameters[2];
66  double& d_0 = trackParameters[3];
67  double& z_0 = trackParameters[4];
68  double invkappa = 1 / kappa;
69  // computes the value of the track parameter that correspond to the hit, relative to phi_0
70  //double phi = kappa*tan(theta)*(hit[2]-z_0);
71  double phi = getPhi(trackParameters, hit, onDisk) - phi_0;
72  // computes the value of the track parameter that correspond to the muon system, relative to phi_0
73  // phi_mu = phi_mu0 - phi_0
74  double phi_mu_b = (kappa > 0 ? -1 : 1) * acos((Rmu * Rmu - d_0 * d_0 - 2 * invkappa * invkappa + 2 * invkappa * d_0) /
75  (2 * invkappa * d_0 - 2 * invkappa * invkappa));
76  double phi_mu_e = kappa * tan(theta) * ((phit[2] < 0 ? 1 : -1) * zmu - z_0);
77  // estimate the time for crossing the barrel
78  double t_barrel = fabs(invkappa / (c * sin(theta)) * (phi - phi_mu_b));
79  // estimate the time for crossing endcaps
80  double t_endcap = fabs(invkappa / (c * sin(theta)) * (phi - phi_mu_e));
81  // take the minimum
82  return t_barrel < t_endcap ? t_barrel : t_endcap;
83 }

References HltBtagPostValidation_cff::c, getPhi(), kappa, phi, funct::sin(), SiPixelPI::t_barrel, funct::tan(), theta(), and trackParameters().

Referenced by timeOfFlight().

◆ trackParameters()

void SiStripFineDelayTOF::trackParameters ( const reco::Track tk,
double *  trackParameters 
)
static

Definition at line 140 of file SiStripFineDelayTOF.cc.

140  {
142  const math::XYZVector& momentum = tk.outerMomentum();
143  LogDebug("trackParameters") << "outer position: " << position.x() << " " << position.y() << " " << position.z();
144  LogDebug("trackParameters") << "outer momentum: " << momentum.x() << " " << momentum.y() << " " << momentum.z();
145  math::XYZVector fieldDirection(0, 0, 1);
146  // computes the center of curvature
147  math::XYZVector radius = momentum.Cross(fieldDirection).Unit();
149  LogDebug("trackParameters") << "center of curvature: " << position.x() << " " << position.y() << " " << position.z();
150  // the transverse IP
151  trackParameters[3] = position.rho() - fabs(1 / trackParameters[0]);
152  if (trackParameters[0] > 0)
153  trackParameters[3] *= -1.;
154  // phi_0
155  double phi_out = trackParameters[2];
156  double phi_0 = position.phi() - TMath::Pi() / 2;
157  if (trackParameters[0] < 0)
158  phi_0 -= TMath::Pi();
159  if (phi_0 < -2 * TMath::Pi())
160  phi_0 += 2 * TMath::Pi();
161  if (phi_0 > 2 * TMath::Pi())
162  phi_0 -= 2 * TMath::Pi();
163  if (phi_0 > 0)
164  phi_0 -= 2 * TMath::Pi();
165  LogDebug("trackParameters") << "phi_0: " << phi_0;
166  trackParameters[2] = phi_0;
167  // z_0
168  trackParameters[4] = tk.outerPosition().z() - (phi_out - phi_0) / TMath::Tan(trackParameters[1]) / trackParameters[0];
169  LogDebug("trackParameters") << "z_0: " << tk.outerPosition().z() << " " << trackParameters[4] << " "
170  << tk.innerPosition().z();
171 }

References reco::Track::innerPosition(), LogDebug, reco::Track::outerMomentum(), reco::Track::outerPosition(), Pi, position, and CosmicsPD_Skims::radius.

Referenced by SiStripFineDelayHit::detId(), getPhi(), timeOfFlight(), timeOfFlightBeamB(), timeOfFlightCosmicB(), x(), y(), and z().

◆ x()

double SiStripFineDelayTOF::x ( double *  trackParameters,
double  phi 
)
staticprivate

◆ y()

double SiStripFineDelayTOF::y ( double *  trackParameters,
double  phi 
)
staticprivate

◆ z()

double SiStripFineDelayTOF::z ( double *  trackParameters,
double  phi 
)
staticprivate

Definition at line 110 of file SiStripFineDelayTOF.cc.

110  {
111  return trackParameters[4] + (phi - trackParameters[2]) / (trackParameters[0] * tan(trackParameters[1]));
112 }

References phi, funct::tan(), and trackParameters().

Referenced by geometryXMLparser.Alignable::pos(), and ntupleDataFormat._HitObject::r3D().

reco::Track::outerPosition
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:62
reco::Track::outerMomentum
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:65
SiStripFineDelayTOF::trackParameters
static void trackParameters(const reco::Track &tk, double *trackParameters)
Definition: SiStripFineDelayTOF.cc:140
SiPixelPI::t_barrel
Definition: SiPixelPayloadInspectorHelper.h:571
SiStripFineDelayTOF::y
static double y(double *trackParameters, double phi)
Definition: SiStripFineDelayTOF.cc:105
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
SiStripFineDelayTOF::timeOfFlightCosmicB
static double timeOfFlightCosmicB(double *trackParameters, double *hit, double *phit, bool onDisk)
Definition: SiStripFineDelayTOF.cc:57
SiStripFineDelayTOF::timeOfFlightBeam
static double timeOfFlightBeam(double *hit, double *phit)
Definition: SiStripFineDelayTOF.cc:85
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
sipixeldigitoraw
Definition: SiPixelDigiToRaw.cc:39
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
funct::tan
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
reco::Track::innerPosition
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:56
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
DDAxes::phi
SiStripFineDelayTOF::getPhi
static double getPhi(double *trackParameters, double *hit, bool onDisk)
Definition: SiStripFineDelayTOF.cc:114
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
Pi
const double Pi
Definition: CosmicMuonParameters.h:18
kappa
static const G4double kappa
Definition: UrbanMscModel93.cc:35
SiStripFineDelayTOF::timeOfFlightBeamB
static double timeOfFlightBeamB(double *trackParameters, double *hit, double *phit, bool onDisk)
Definition: SiStripFineDelayTOF.cc:92
SiStripFineDelayTOF::timeOfFlightCosmic
static double timeOfFlightCosmic(double *hit, double *phit)
Definition: SiStripFineDelayTOF.cc:30
GlobalTrackerMuonAlignment_cfi.cosmics
cosmics
Definition: GlobalTrackerMuonAlignment_cfi.py:5
hit
Definition: SiStripHitEffFromCalibTree.cc:88
SiStripFineDelayTOF::x
static double x(double *trackParameters, double phi)
Definition: SiStripFineDelayTOF.cc:101