ATLAS Offline Software
TrigInDetTrackFitPar.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #ifndef TRIGINDETTRACKFITPAR_H
6 #define TRIGINDETTRACKFITPAR_H
9 #include <vector>
10 
68 
69  public:
70 
72 
75  m_a0(0.0), m_phi0(0.0), m_z0(0.0),m_eta(0.0), m_pT(0.0),
76  m_ea0(0.0), m_ephi0(0.0), m_ez0(0.0), m_eeta(0.0), m_epT(0.0),
78  {};
79 
81  TrigInDetTrackFitPar (const double a0,
82  const double phi0,
83  const double z0,
84  const double eta,
85  const double pT,
86  const double ea0,
87  const double ephi0,
88  const double ez0,
89  const double eeta,
90  const double epT,
91  const std::vector<double>* cov=0) :
92  m_a0(a0), m_phi0(phi0), m_z0(z0),
96  {};
97 
98 
100  TrigInDetTrackFitPar (const double a0,
101  const double phi0,
102  const double z0,
103  const double eta,
104  const double pT,
105  const double ea0,
106  const double ephi0,
107  const double ez0,
108  const double eeta,
109  const double epT,
110  const TrigSurfaceType t,
111  const double c,
112  const std::vector<double>* cov=0) :
113  m_a0(a0), m_phi0(phi0), m_z0(z0),
116  {};
117 
119  TrigInDetTrackFitPar (const double a0,
120  const double phi0,
121  const double z0,
122  const double eta,
123  const double pT,
124  const std::vector<double>* cov=0) :
125  m_a0(a0), m_phi0(phi0), m_z0(z0),
126  m_eta(eta), m_pT(pT),
127  m_ea0(0.0), m_ephi0(0.0), m_ez0(0.0), m_eeta(0.0),
128  m_epT(0.0), m_cov(cov), m_surfaceType(PERIGEE),
129  m_surfaceCoordinate(0.0){
130  if(cov!=NULL) {
131  if((*cov)[0]>=0.0) m_ea0 = sqrt((*cov)[0]);
132  if((*cov)[5]>=0.0) m_ephi0 = sqrt((*cov)[5]);
133  if((*cov)[9]>=0.0) m_ez0 = sqrt((*cov)[9]);
134  if((*cov)[12]>=0.0) m_eeta = sqrt((*cov)[12]);
135  if((*cov)[14]>=0.0) m_epT = sqrt((*cov)[14]);
136  }
137  }
138 
140  TrigInDetTrackFitPar (const double a0,
141  const double phi0,
142  const double z0,
143  const double eta,
144  const double pT,
145  const TrigSurfaceType t,
146  const double c,
147  const std::vector<double>* cov=0) :
148  m_a0(a0), m_phi0(phi0), m_z0(z0),
149  m_eta(eta), m_pT(pT),
150  m_ea0(0.0), m_ephi0(0.0), m_ez0(0.0), m_eeta(0.0), m_epT(0.0), m_cov(cov), m_surfaceType(t),
152  {
153  if(cov!=NULL) {
154  if((*cov)[0]>=0.0) m_ea0 = sqrt((*cov)[0]);
155  if((*cov)[5]>=0.0) m_ephi0 = sqrt((*cov)[5]);
156  if((*cov)[9]>=0.0) m_ez0 = sqrt((*cov)[9]);
157  if((*cov)[12]>=0.0) m_eeta = sqrt((*cov)[12]);
158  if((*cov)[14]>=0.0) m_epT = sqrt((*cov)[14]);
159  }
160  }
161 
164  : m_a0 (rhs.m_a0),
165  m_phi0 (rhs.m_phi0),
166  m_z0 (rhs.m_z0),
167  m_eta (rhs.m_eta),
168  m_pT (rhs.m_pT),
169  m_ea0 (rhs.m_ea0),
170  m_ephi0 (rhs.m_ephi0),
171  m_ez0 (rhs.m_ez0),
172  m_eeta (rhs.m_eeta),
173  m_epT (rhs.m_epT),
176  {
177  if (rhs.m_cov)
178  m_cov = new std::vector<double> (*rhs.m_cov);
179  else
180  m_cov = nullptr;
181  }
182 
185  {
186  if (this != &rhs) {
187  m_a0 = rhs.m_a0;
188  m_phi0 = rhs.m_phi0;
189  m_z0 = rhs.m_z0;
190  m_eta = rhs.m_eta;
191  m_pT = rhs.m_pT;
192  m_ea0 = rhs.m_ea0;
193  m_ephi0 = rhs.m_ephi0;
194  m_ez0 = rhs.m_ez0;
195  m_eeta = rhs.m_eeta;
196  m_epT = rhs.m_epT;
197  m_surfaceType = rhs.m_surfaceType;
198  m_surfaceCoordinate = rhs.m_surfaceCoordinate;
199 
200  delete m_cov;
201  m_cov = rhs.m_cov;
202  rhs.m_cov = nullptr;
203  }
204  return *this;
205  }
206 
209 
210 
211  // Methods to set data members
213  void a0 ( const double a0 ) { m_a0 = a0; }
215  void z0 ( const double z0 ) { m_z0 = z0; }
217  void phi0( const double phi0 ) { m_phi0 = phi0; }
219  void eta ( const double eta ) { m_eta = eta; }
221  void pT ( const double pT ) { m_pT = pT; }
223  void cov ( const std::vector<double>* cov ) { m_cov = cov; }
228 
229  // Methods to retrieve data members
231  double a0() const { return m_a0; }
233  double z0() const { return m_z0; }
235  double phi0() const { return m_phi0; }
237  double eta() const { return m_eta; }
239  double pT() const { return m_pT; }
241  double ea0() const { return m_ea0; }
243  double ez0() const { return m_ez0; }
245  double ephi0() const { return m_ephi0; }
247  double eeta() const { return m_eeta; }
249  double epT() const { return m_epT; }
251  const std::vector<double>* cov() const { return m_cov; }
255  double surfaceCoordinate() const { return m_surfaceCoordinate; }
256 
257  private:
258  double m_a0;
259  double m_phi0;
260  double m_z0;
261  double m_eta;
262  double m_pT;
263  double m_ea0;
264  double m_ephi0;
265  double m_ez0;
266  double m_eeta;
267  double m_epT;
268  const std::vector<double>* m_cov;
271 };
272 
273 CLASS_DEF( TrigInDetTrackFitPar , 159201573 , 1 )
274 CLASS_DEF( DataVector<TrigInDetTrackFitPar> , 99002908 , 1 )
275 
276 #endif
277 
278 
279 
280 
281 
282 
TrigInDetTrackFitPar::operator=
TrigInDetTrackFitPar & operator=(TrigInDetTrackFitPar &&rhs)
Move assignment.
Definition: TrigInDetTrackFitPar.h:184
TrigInDetTrackFitPar::pT
double pT() const
transverse momentum
Definition: TrigInDetTrackFitPar.h:239
TrigInDetTrackFitPar::ez0
double ez0() const
variance of longitudinal impact parameter
Definition: TrigInDetTrackFitPar.h:243
python.SystemOfUnits.s
int s
Definition: SystemOfUnits.py:131
TrigInDetTrackFitPar::m_epT
double m_epT
Definition: TrigInDetTrackFitPar.h:267
TrigInDetTrackFitPar::PERIGEE
@ PERIGEE
Definition: TrigInDetTrackFitPar.h:71
TrigInDetTrackFitPar::m_a0
double m_a0
see detailed description below
Definition: TrigInDetTrackFitPar.h:258
TrigInDetTrackFitPar::m_phi0
double m_phi0
see detailed description below
Definition: TrigInDetTrackFitPar.h:259
TrigInDetTrackFitPar::surfaceCoordinate
void surfaceCoordinate(double c)
Setter: surface reference coordinate for non-perigee surfaces.
Definition: TrigInDetTrackFitPar.h:227
TrigInDetTrackFitPar::TrigInDetTrackFitPar
TrigInDetTrackFitPar(const double a0, const double phi0, const double z0, const double eta, const double pT, const TrigSurfaceType t, const double c, const std::vector< double > *cov=0)
Constructor for non-PERIGEE parameters without errors or covariance.
Definition: TrigInDetTrackFitPar.h:140
TrigInDetTrackFitPar::m_eta
double m_eta
pseudorapidity
Definition: TrigInDetTrackFitPar.h:261
TrigInDetTrackFitPar
Definition: TrigInDetTrackFitPar.h:67
TrigInDetTrackFitPar::surfaceType
TrigSurfaceType surfaceType() const
surface type
Definition: TrigInDetTrackFitPar.h:253
TrigInDetTrackFitPar::TrigInDetTrackFitPar
TrigInDetTrackFitPar(const double a0, const double phi0, const double z0, const double eta, const double pT, const std::vector< double > *cov=0)
Constructor for PERIGEE parameters without errors or covariance.
Definition: TrigInDetTrackFitPar.h:119
TrigInDetTrackFitPar::z0
double z0() const
longitudinal impact parameter
Definition: TrigInDetTrackFitPar.h:233
TrigInDetTrackFitPar::eta
void eta(const double eta)
Setter: pseudorapidity.
Definition: TrigInDetTrackFitPar.h:219
TrigInDetTrackFitPar::epT
double epT() const
variance of transverse momentum
Definition: TrigInDetTrackFitPar.h:249
TrigInDetTrackFitPar::z0
void z0(const double z0)
Setter: longitudinal impact parameter.
Definition: TrigInDetTrackFitPar.h:215
TrigInDetTrackFitPar::TrigInDetTrackFitPar
TrigInDetTrackFitPar(const double a0, const double phi0, const double z0, const double eta, const double pT, const double ea0, const double ephi0, const double ez0, const double eeta, const double epT, const TrigSurfaceType t, const double c, const std::vector< double > *cov=0)
Constructor for parameters on non-PERIGEE surface.
Definition: TrigInDetTrackFitPar.h:100
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
TrigInDetTrackFitPar::cov
const std::vector< double > * cov() const
covariance (packed) of track parameters
Definition: TrigInDetTrackFitPar.h:251
TrigInDetTrackFitPar::m_pT
double m_pT
transverse momentum
Definition: TrigInDetTrackFitPar.h:262
TrigInDetTrackFitPar::surfaceCoordinate
double surfaceCoordinate() const
surface reference coordinate (radius or z-position) for non-perigee parameters
Definition: TrigInDetTrackFitPar.h:255
TrigInDetTrackFitPar::a0
double a0() const
transverse impact parameter
Definition: TrigInDetTrackFitPar.h:231
TrigInDetTrackFitPar::m_surfaceType
TrigSurfaceType m_surfaceType
type of track parameters - perigee, barrel, or endcap
Definition: TrigInDetTrackFitPar.h:269
TrigInDetTrackFitPar::m_surfaceCoordinate
double m_surfaceCoordinate
barrel radius or z of endcap disk
Definition: TrigInDetTrackFitPar.h:270
TrigInDetTrackFitPar::UNDEFINED
@ UNDEFINED
Definition: TrigInDetTrackFitPar.h:71
TrigInDetTrackFitPar::m_eeta
double m_eeta
Definition: TrigInDetTrackFitPar.h:266
TrigInDetTrackFitPar::TrigSurfaceType
TrigSurfaceType
Definition: TrigInDetTrackFitPar.h:71
TrigInDetTrackFitPar::surfaceType
void surfaceType(TrigSurfaceType s)
Setter: surface type PERIGEE=0, BARREL=1, ENDCAP=2.
Definition: TrigInDetTrackFitPar.h:225
TrigInDetTrackFitPar::m_ez0
double m_ez0
Definition: TrigInDetTrackFitPar.h:265
TrigInDetTrackFitPar::cov
void cov(const std::vector< double > *cov)
Setter: covariance matrix of track parameters.
Definition: TrigInDetTrackFitPar.h:223
TrigInDetTrackFitPar::~TrigInDetTrackFitPar
~TrigInDetTrackFitPar()
Destructor.
Definition: TrigInDetTrackFitPar.h:208
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
TrigInDetTrackFitPar::eeta
double eeta() const
variance of pseudorapidity
Definition: TrigInDetTrackFitPar.h:247
TrigInDetTrackFitPar::ea0
double ea0() const
variance of transverse impact parameter
Definition: TrigInDetTrackFitPar.h:241
TrigInDetTrackFitPar::TrigInDetTrackFitPar
TrigInDetTrackFitPar()
Constructor for POOL only.
Definition: TrigInDetTrackFitPar.h:74
TrigInDetTrackFitPar::ephi0
double ephi0() const
variance of azimuthal angle of the momentum
Definition: TrigInDetTrackFitPar.h:245
TrigInDetTrackFitPar::m_ephi0
double m_ephi0
Definition: TrigInDetTrackFitPar.h:264
TrigInDetTrackFitPar::m_cov
const std::vector< double > * m_cov
covariance matrix packed as described below
Definition: TrigInDetTrackFitPar.h:268
DataVector.h
An STL vector of pointers that by default owns its pointed-to elements.
TrigInDetTrackFitPar::a0
void a0(const double a0)
Setter: transverse impact parameter.
Definition: TrigInDetTrackFitPar.h:213
TrigInDetTrackFitPar::TrigInDetTrackFitPar
TrigInDetTrackFitPar(const double a0, const double phi0, const double z0, const double eta, const double pT, const double ea0, const double ephi0, const double ez0, const double eeta, const double epT, const std::vector< double > *cov=0)
Constructor for parameters on PERIGEE surface.
Definition: TrigInDetTrackFitPar.h:81
TrigInDetTrackFitPar::ENDCAP
@ ENDCAP
Definition: TrigInDetTrackFitPar.h:71
CLASS_DEF
#define CLASS_DEF(NAME, CID, VERSION)
associate a clid and a version to a type eg
Definition: Control/AthenaKernel/AthenaKernel/CLASS_DEF.h:64
TrigInDetTrackFitPar::m_z0
double m_z0
see detailed description below
Definition: TrigInDetTrackFitPar.h:260
TrigInDetTrackFitPar::phi0
double phi0() const
azimuthal angle of the momentum
Definition: TrigInDetTrackFitPar.h:235
TrigInDetTrackFitPar::m_ea0
double m_ea0
Definition: TrigInDetTrackFitPar.h:263
TrigInDetTrackFitPar::eta
double eta() const
pseudorapidity
Definition: TrigInDetTrackFitPar.h:237
TrigInDetTrackFitPar::BARREL
@ BARREL
Definition: TrigInDetTrackFitPar.h:71
TrigInDetTrackFitPar::phi0
void phi0(const double phi0)
Setter: azimuthal angle of the momentum.
Definition: TrigInDetTrackFitPar.h:217
python.compressB64.c
def c
Definition: compressB64.py:93
CLASS_DEF.h
macros to associate a CLID to a type
TrigInDetTrackFitPar::pT
void pT(const double pT)
Setter: transverse momentum.
Definition: TrigInDetTrackFitPar.h:221
TrigInDetTrackFitPar::TrigInDetTrackFitPar
TrigInDetTrackFitPar(const TrigInDetTrackFitPar &rhs)
Copy constructor.
Definition: TrigInDetTrackFitPar.h:163