CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
CSCCondSegFit Class Reference

#include <CSCCondSegFit.h>

Inheritance diagram for CSCCondSegFit:
CSCSegFit

Public Member Functions

 CSCCondSegFit (const edm::ParameterSet &ps, const CSCChamber *csc, const CSCSetOfHits &hits)
 
void fit (bool condpass1=false, bool condpass2=false)
 
int worstHit (void)
 
 ~CSCCondSegFit () override
 
- Public Member Functions inherited from CSCSegFit
const CSCChamberchamber () const
 
double chi2 (void) const
 
AlgebraicSymMatrix covarianceMatrix (void)
 
 CSCSegFit (const CSCChamber *csc, CSCSetOfHits hits)
 
void fit (void)
 
bool fitdone () const
 
CSCSetOfHits hits (void) const
 
LocalPoint intercept () const
 
LocalVector localdir () const
 
int ndof (void) const
 
size_t nhits (void) const
 
float Rdev (float x, float y, float z) const
 
double scaleXError (void) const
 
void setScaleXError (double factor)
 
float xdev (float x, float z) const
 
float xfit (float z) const
 
float ydev (float y, float z) const
 
float yfit (float z) const
 
virtual ~CSCSegFit ()
 

Private Member Functions

void correctTheCovMatrix (CSCSegFit::SMatrixSym2 &IC)
 
void correctTheCovX (void)
 
void setChi2 (bool condpass1, bool condpass2)
 

Private Attributes

double chi2Norm_
 
double condSeed1_
 
double condSeed2_
 
double covAnyNumber_
 Allow to use any number for covariance for all RecHits. More...
 
bool covToAnyNumber_
 The correction parameters. More...
 
bool covToAnyNumberAll_
 Allow to use any number for covariance (by hand) More...
 
std::vector< double > lex_
 
int worstHit_
 

Additional Inherited Members

- Public Types inherited from CSCSegFit
typedef std::vector< const CSCRecHit2D * > CSCSetOfHits
 
typedef ROOT::Math::SMatrix< double, 12, 4 > SMatrix12by4
 
typedef ROOT::Math::SMatrix< double, 4 > SMatrix4
 
typedef ROOT::Math::SMatrix< double, 12, 12, ROOT::Math::MatRepSym< double, 12 > > SMatrixSym12
 
typedef ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepSym< double, 2 > > SMatrixSym2
 
typedef ROOT::Math::SMatrix< double, 4, 4, ROOT::Math::MatRepSym< double, 4 > > SMatrixSym4
 
typedef ROOT::Math::SVector< double, 4 > SVector4
 
- Protected Member Functions inherited from CSCSegFit
SMatrix12by4 derivativeMatrix (void)
 
AlgebraicSymMatrix flipErrors (const SMatrixSym4 &)
 
void setOutFromIP (void)
 
SMatrixSym12 weightMatrix (void)
 
- Protected Attributes inherited from CSCSegFit
const CSCChamberchamber_
 
double chi2_
 
bool fitdone_
 
CSCSetOfHits hits_
 
LocalPoint intercept_
 
LocalVector localdir_
 
int ndof_
 
double scaleXError_
 
float uslope_
 
float vslope_
 

Detailed Description

Definition at line 18 of file CSCCondSegFit.h.

Constructor & Destructor Documentation

◆ CSCCondSegFit()

CSCCondSegFit::CSCCondSegFit ( const edm::ParameterSet ps,
const CSCChamber csc,
const CSCSetOfHits hits 
)
inline

Definition at line 20 of file CSCCondSegFit.h.

21  : CSCSegFit(csc, hits),
22  worstHit_(0),
23  chi2Norm_(ps.getParameter<double>("NormChi2Cut2D")),
24  condSeed1_(ps.getParameter<double>("SeedSmall")),
25  condSeed2_(ps.getParameter<double>("SeedBig")),
26  covToAnyNumber_(ps.getParameter<bool>("ForceCovariance")),
27  covToAnyNumberAll_(ps.getParameter<bool>("ForceCovarianceAll")),
28  covAnyNumber_(ps.getParameter<double>("Covariance")) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
CSCSegFit(const CSCChamber *csc, CSCSetOfHits hits)
Definition: CSCSegFit.h:55
bool covToAnyNumber_
The correction parameters.
Definition: CSCCondSegFit.h:61
double covAnyNumber_
Allow to use any number for covariance for all RecHits.
Definition: CSCCondSegFit.h:63
bool covToAnyNumberAll_
Allow to use any number for covariance (by hand)
Definition: CSCCondSegFit.h:62
double chi2Norm_
Definition: CSCCondSegFit.h:55
double condSeed2_
Definition: CSCCondSegFit.h:60
Definition: L1Track.h:19
double condSeed1_
Definition: CSCCondSegFit.h:60
CSCSetOfHits hits(void) const
Definition: CSCSegFit.h:79

◆ ~CSCCondSegFit()

CSCCondSegFit::~CSCCondSegFit ( )
inlineoverride

Definition at line 30 of file CSCCondSegFit.h.

30 {}

Member Function Documentation

◆ correctTheCovMatrix()

void CSCCondSegFit::correctTheCovMatrix ( CSCSegFit::SMatrixSym2 IC)
private

Definition at line 239 of file CSCCondSegFit.cc.

References condSeed1_, condSeed2_, covAnyNumber_, covToAnyNumber_, covToAnyNumberAll_, and mathSSE::sqrt().

Referenced by fit(), and setChi2().

239  {
240  //double condNumberCorr1=0.0;
241  double condNumberCorr2 = 0.0;
242  double detCov = 0.0;
243  double diag1 = 0.0;
244  double diag2 = 0.0;
245  double IC_01_corr = 0.0;
246  double IC_00_corr = 0.0;
247  if (!covToAnyNumberAll_) {
248  //condNumberCorr1=condSeed1_*IC(1,1);
249  condNumberCorr2 = condSeed2_ * IC(1, 1);
250  diag1 = IC(0, 0) * IC(1, 1);
251  diag2 = IC(0, 1) * IC(0, 1);
252  detCov = fabs(diag1 - diag2);
253  if ((diag1 < condNumberCorr2) && (diag2 < condNumberCorr2)) {
254  if (covToAnyNumber_)
255  IC(0, 1) = covAnyNumber_;
256  else {
257  IC_00_corr = condSeed1_ + fabs(IC(0, 1)) / IC(1, 1);
258  IC(0, 0) = IC_00_corr;
259  }
260  }
261 
262  if (((detCov < condNumberCorr2) && (diag1 > condNumberCorr2)) ||
263  ((diag2 > condNumberCorr2) && (detCov < condNumberCorr2))) {
264  if (covToAnyNumber_)
265  IC(0, 1) = covAnyNumber_;
266  else {
267  IC_01_corr = sqrt(fabs(diag1 - condNumberCorr2));
268  if (IC(0, 1) < 0)
269  IC(0, 1) = (-1) * IC_01_corr;
270  else
271  IC(0, 1) = IC_01_corr;
272  }
273  }
274  } else {
275  IC(0, 1) = covAnyNumber_;
276  }
277  // With SMatrix formulation. the following might be unnecessary since it might
278  // enforce the symmetry. But can't hurt to leave it (it WAS a real bug fix for
279  // the original CLHEP formulation.)
280  IC(1, 0) = IC(0, 1); //@@ ADDED TO MAINTAIN SYMMETRY OF IC
281 }
bool covToAnyNumber_
The correction parameters.
Definition: CSCCondSegFit.h:61
double covAnyNumber_
Allow to use any number for covariance for all RecHits.
Definition: CSCCondSegFit.h:63
bool covToAnyNumberAll_
Allow to use any number for covariance (by hand)
Definition: CSCCondSegFit.h:62
T sqrt(T t)
Definition: SSEVec.h:19
double condSeed2_
Definition: CSCCondSegFit.h:60
double condSeed1_
Definition: CSCCondSegFit.h:60

◆ correctTheCovX()

void CSCCondSegFit::correctTheCovX ( void  )
private

Make a one dimensional fit in U-Z plane (i.e. chamber local x-z) U=U0+UZ*Z fit parameters
Calculate the fit line trace Calculate one dimensional chi^2 and normalize the errors if needed

Definition at line 174 of file CSCCondSegFit.cc.

References CSCSegFit::chamber(), chi2Norm_, makePileupJSON::denom, runTauDisplay::gp, CSCSegFit::hits_, mps_fire::i, nano_mu_digi_cff::layer, CSCChamber::layer(), lex_, funct::pow(), CSCSegFit::setScaleXError(), GeomDet::toLocal(), findQualityFiles::v, worstHit_, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), z, PV3DBase< T, PVType, FrameType >::z(), and geometryCSVtoXML::zz.

Referenced by fit().

174  {
175  std::vector<double> uu, vv, zz; // Vectors of coordinates
176 
177  lex_.clear(); // x component of local error for each hit
178 
179  double sum_U_err = 0.0;
180  double sum_Z_U_err = 0.0;
181  double sum_Z2_U_err = 0.0;
182  double sum_U_U_err = 0.0;
183  double sum_UZ_U_err = 0.0;
184  std::vector<double> chiUZind;
185  std::vector<double>::iterator chiContribution;
186  double chiUZ = 0.0;
187  CSCSetOfHits::const_iterator ih = hits_.begin();
188  for (ih = hits_.begin(); ih != hits_.end(); ++ih) {
189  const CSCRecHit2D& hit = (**ih);
190  lex_.push_back(hit.localPositionError().xx());
191 
192  const CSCLayer* layer = chamber()->layer(hit.cscDetId().layer());
193  GlobalPoint gp = layer->toGlobal(hit.localPosition());
194  LocalPoint lp = chamber()->toLocal(gp);
195  // Local position of hit w.r.t. chamber
196  double u = lp.x();
197  double v = lp.y();
198  double z = lp.z();
199  uu.push_back(u);
200  vv.push_back(v);
201  zz.push_back(z);
202  // Prepare the sums for the standard linear fit
203  sum_U_err += 1. / lex_.back();
204  sum_Z_U_err += z / lex_.back();
205  sum_Z2_U_err += (z * z) / lex_.back();
206  sum_U_U_err += u / lex_.back();
207  sum_UZ_U_err += (u * z) / lex_.back();
208  }
209 
212 
213  double denom = sum_U_err * sum_Z2_U_err - pow(sum_Z_U_err, 2);
214  double U0 = (sum_Z2_U_err * sum_U_U_err - sum_Z_U_err * sum_UZ_U_err) / denom;
215  double UZ = (sum_U_err * sum_UZ_U_err - sum_Z_U_err * sum_U_U_err) / denom;
216 
219 
220  for (unsigned i = 0; i < uu.size(); ++i) {
221  double uMean = U0 + UZ * zz[i];
222  chiUZind.push_back((pow((uMean - uu[i]), 2)) / lex_[i]);
223  chiUZ += (pow((uMean - uu[i]), 2)) / lex_[i];
224  }
225  chiUZ = chiUZ / (uu.size() - 2);
226 
227  if (chiUZ >= chi2Norm_) {
228  double chi2uCorrection = chiUZ / chi2Norm_;
229  for (unsigned i = 0; i < uu.size(); ++i)
230  lex_[i] = lex_[i] * chi2uCorrection;
231  setScaleXError(chi2uCorrection);
232  }
233 
234  // Find most deviant hit
235  chiContribution = max_element(chiUZind.begin(), chiUZind.end());
236  worstHit_ = chiContribution - chiUZind.begin();
237 }
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:30
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
T z() const
Definition: PV3DBase.h:61
std::vector< double > lex_
Definition: CSCCondSegFit.h:54
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const CSCChamber * chamber() const
Definition: CSCSegFit.h:86
double chi2Norm_
Definition: CSCCondSegFit.h:55
CSCSetOfHits hits_
Definition: CSCSegFit.h:109
void setScaleXError(double factor)
Definition: CSCSegFit.h:67
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ fit()

void CSCCondSegFit::fit ( bool  condpass1 = false,
bool  condpass2 = false 
)

Definition at line 12 of file CSCCondSegFit.cc.

References B, CSCSegFit::chamber(), correctTheCovMatrix(), correctTheCovX(), runTauDisplay::gp, CSCSegFit::hits_, CSCSegFit::intercept_, nano_mu_digi_cff::layer, CSCChamber::layer(), lex_, LogTrace, convertSQLiteXML::ok, AlCaHLTBitMon_ParallelJobs::p, CSCSegFit::setChi2(), CSCSegFit::setOutFromIP(), GeomDet::toLocal(), CSCSegFit::uslope_, findQualityFiles::v, CSCSegFit::vslope_, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by CSCSegAlgoST::buildSegments(), trackingPlots.Iteration::modules(), and CSCSegAlgoST::prune_bad_hits().

12  {
13  // See base class CSCSegFit::fit for detailed explanation of algorithm.
14  // This is exactly the same, except it adds two conditioning passes
15  // which can improve the robustness of the calculations.
16 
17  lex_.clear(); // Vector of the error matrix (only xx)
18 
19  if (condpass1 && !condpass2) {
21  if (lex_.size() != hits_.size()) {
22  LogTrace("CSCSegment|CSC") << "[CSCSegFit::fit] lex_.size()!=hits_.size() ALARM! Please inform CSC DPG "
23  << std::endl;
24  }
25  }
26 
27  SMatrix4 M; // 4x4, init to 0
28  SVector4 B; // 4x1, init to 0;
29 
30  CSCSetOfHits::const_iterator ih = hits_.begin();
31 
32  for (ih = hits_.begin(); ih != hits_.end(); ++ih) {
33  const CSCRecHit2D& hit = (**ih);
34  const CSCLayer* layer = chamber()->layer(hit.cscDetId().layer());
35  GlobalPoint gp = layer->toGlobal(hit.localPosition());
36  LocalPoint lp = chamber()->toLocal(gp);
37 
38  // Local position of hit w.r.t. chamber
39  double u = lp.x();
40  double v = lp.y();
41  double z = lp.z();
42 
43  // Covariance matrix of local errors
44  SMatrixSym2 IC; // 2x2, init to 0
45 
46  // Adjust covariance matrix elements to improve numerical stability?
47  if (condpass1 && !condpass2) {
48  IC(0, 0) = lex_.at(ih - hits_.begin());
49  } else {
50  IC(0, 0) = hit.localPositionError().xx();
51  }
52  IC(1, 1) = hit.localPositionError().yy();
53  //@@ NOT SURE WHICH OFF-DIAGONAL ELEMENT MUST BE DEFINED BUT (1,0) WORKS
54  //@@ (and SMatrix enforces symmetry)
55  IC(1, 0) = hit.localPositionError().xy();
56  // IC(0,1) = IC(1,0);
57 
58  // Correct the covariance matrix to improve numerical stability?
59  if (condpass2) {
61  }
62 
63  // Invert covariance matrix (and trap if it fails!)
64  bool ok = IC.Invert();
65  if (!ok) {
66  LogTrace("CSCSegment|CSC") << "[CSCSegFit::fit] Failed to invert covariance matrix: \n" << IC;
67  // return ok; //@@ SHOULD PASS THIS BACK TO CALLER?
68  }
69 
70  M(0, 0) += IC(0, 0);
71  M(0, 1) += IC(0, 1);
72  M(0, 2) += IC(0, 0) * z;
73  M(0, 3) += IC(0, 1) * z;
74  B(0) += u * IC(0, 0) + v * IC(0, 1);
75 
76  M(1, 0) += IC(1, 0);
77  M(1, 1) += IC(1, 1);
78  M(1, 2) += IC(1, 0) * z;
79  M(1, 3) += IC(1, 1) * z;
80  B(1) += u * IC(1, 0) + v * IC(1, 1);
81 
82  M(2, 0) += IC(0, 0) * z;
83  M(2, 1) += IC(0, 1) * z;
84  M(2, 2) += IC(0, 0) * z * z;
85  M(2, 3) += IC(0, 1) * z * z;
86  B(2) += (u * IC(0, 0) + v * IC(0, 1)) * z;
87 
88  M(3, 0) += IC(1, 0) * z;
89  M(3, 1) += IC(1, 1) * z;
90  M(3, 2) += IC(1, 0) * z * z;
91  M(3, 3) += IC(1, 1) * z * z;
92  B(3) += (u * IC(1, 0) + v * IC(1, 1)) * z;
93  }
94 
95  SVector4 p;
96  bool ok = M.Invert();
97  if (!ok) {
98  LogTrace("CSCSegment|CSC") << "[CSCSegFit::fit] Failed to invert matrix: \n" << M;
99  // return ok; //@@ SHOULD PASS THIS BACK TO CALLER?
100  } else {
101  p = M * B;
102  }
103 
104  // LogTrace("CSCSegFit") << "[CSCSegFit::fit] p = "
105  // << p(0) << ", " << p(1) << ", " << p(2) << ", " << p(3);
106 
107  // fill member variables (note origin has local z = 0)
108  intercept_ = LocalPoint(p(0), p(1), 0.);
109 
110  // localdir_
111  uslope_ = p(2);
112  vslope_ = p(3);
113  setOutFromIP();
114 
115  // calculate chi2 of fit
116  setChi2(condpass1, condpass2);
117 }
LocalPoint intercept_
Definition: CSCSegFit.h:112
void correctTheCovX(void)
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:30
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
Definition: APVGainStruct.h:7
void correctTheCovMatrix(CSCSegFit::SMatrixSym2 &IC)
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
T z() const
Definition: PV3DBase.h:61
float vslope_
Definition: CSCSegFit.h:111
ROOT::Math::SVector< double, 4 > SVector4
Definition: CSCSegFit.h:50
std::vector< double > lex_
Definition: CSCCondSegFit.h:54
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepSym< double, 2 > > SMatrixSym2
Definition: CSCSegFit.h:47
#define LogTrace(id)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const CSCChamber * chamber() const
Definition: CSCSegFit.h:86
ROOT::Math::SMatrix< double, 4 > SMatrix4
Definition: CSCSegFit.h:43
CSCSetOfHits hits_
Definition: CSCSegFit.h:109
void setOutFromIP(void)
Definition: CSCSegFit.cc:331
float uslope_
Definition: CSCSegFit.h:110
void setChi2(void)
Definition: CSCSegFit.cc:236

◆ setChi2()

void CSCCondSegFit::setChi2 ( bool  condpass1,
bool  condpass2 
)
private

Correct the cov matrix

Definition at line 119 of file CSCCondSegFit.cc.

References CSCSegFit::chamber(), CSCSegFit::chi2_, correctTheCovMatrix(), runTauDisplay::gp, CSCSegFit::hits_, CSCSegFit::intercept_, nano_mu_digi_cff::layer, CSCChamber::layer(), lex_, LogTrace, CSCSegFit::ndof_, convertSQLiteXML::ok, GeomDet::toLocal(), CSCSegFit::uslope_, findQualityFiles::v, CSCSegFit::vslope_, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

119  {
120  double chsq = 0.;
121 
122  CSCSetOfHits::const_iterator ih;
123  for (ih = hits_.begin(); ih != hits_.end(); ++ih) {
124  const CSCRecHit2D& hit = (**ih);
125  const CSCLayer* layer = chamber()->layer(hit.cscDetId().layer());
126  GlobalPoint gp = layer->toGlobal(hit.localPosition());
127  LocalPoint lp = chamber()->toLocal(gp);
128 
129  double u = lp.x();
130  double v = lp.y();
131  double z = lp.z();
132 
133  double du = intercept_.x() + uslope_ * z - u;
134  double dv = intercept_.y() + vslope_ * z - v;
135 
136  // LogTrace("CSCSegFit") << "[CSCSegFit::setChi2] u, v, z = " << u << ", " << v << ", " << z;
137 
138  SMatrixSym2 IC; // 2x2, init to 0
139 
140  if (condpass1 && !condpass2) {
141  IC(0, 0) = lex_.at(ih - hits_.begin());
142  } else {
143  IC(0, 0) = hit.localPositionError().xx();
144  }
145  // IC(0,1) = hit.localPositionError().xy();
146  IC(1, 0) = hit.localPositionError().xy();
147  IC(1, 1) = hit.localPositionError().yy();
148  // IC(1,0) = IC(0,1);
149 
150  // LogTrace("CSCSegFit") << "[CSCSegFit::setChi2] IC before = \n" << IC;
151 
153  if (condpass2) {
155  }
156 
157  // Invert covariance matrix
158  bool ok = IC.Invert();
159  if (!ok) {
160  LogTrace("CSCSegment|CSC") << "[CSCSegFit::setChi2] Failed to invert covariance matrix: \n" << IC;
161  // return ok;
162  }
163  // LogTrace("CSCSegFit") << "[CSCSegFit::setChi2] IC after = \n" << IC;
164  chsq += du * du * IC(0, 0) + 2. * du * dv * IC(0, 1) + dv * dv * IC(1, 1);
165  }
166 
167  // fill member variables
168  chi2_ = chsq;
169  ndof_ = 2. * hits_.size() - 4;
170 
171  // LogTrace("CSCSegFit") << "[CSCSegFit::setChi2] chi2 = " << chi2_ << "/" << ndof_ << " dof";
172 }
LocalPoint intercept_
Definition: CSCSegFit.h:112
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:30
void correctTheCovMatrix(CSCSegFit::SMatrixSym2 &IC)
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
T z() const
Definition: PV3DBase.h:61
double chi2_
Definition: CSCSegFit.h:114
float vslope_
Definition: CSCSegFit.h:111
std::vector< double > lex_
Definition: CSCCondSegFit.h:54
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepSym< double, 2 > > SMatrixSym2
Definition: CSCSegFit.h:47
#define LogTrace(id)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const CSCChamber * chamber() const
Definition: CSCSegFit.h:86
CSCSetOfHits hits_
Definition: CSCSegFit.h:109
float uslope_
Definition: CSCSegFit.h:110
int ndof_
Definition: CSCSegFit.h:115

◆ worstHit()

int CSCCondSegFit::worstHit ( void  )
inline

Definition at line 36 of file CSCCondSegFit.h.

References worstHit_.

Referenced by CSCSegAlgoST::buildSegments().

36 { return worstHit_; }

Member Data Documentation

◆ chi2Norm_

double CSCCondSegFit::chi2Norm_
private

Definition at line 55 of file CSCCondSegFit.h.

Referenced by correctTheCovX().

◆ condSeed1_

double CSCCondSegFit::condSeed1_
private

Definition at line 60 of file CSCCondSegFit.h.

Referenced by correctTheCovMatrix().

◆ condSeed2_

double CSCCondSegFit::condSeed2_
private

Definition at line 60 of file CSCCondSegFit.h.

Referenced by correctTheCovMatrix().

◆ covAnyNumber_

double CSCCondSegFit::covAnyNumber_
private

Allow to use any number for covariance for all RecHits.

Definition at line 63 of file CSCCondSegFit.h.

Referenced by correctTheCovMatrix().

◆ covToAnyNumber_

bool CSCCondSegFit::covToAnyNumber_
private

The correction parameters.

Definition at line 61 of file CSCCondSegFit.h.

Referenced by correctTheCovMatrix().

◆ covToAnyNumberAll_

bool CSCCondSegFit::covToAnyNumberAll_
private

Allow to use any number for covariance (by hand)

Definition at line 62 of file CSCCondSegFit.h.

Referenced by correctTheCovMatrix().

◆ lex_

std::vector<double> CSCCondSegFit::lex_
private

Definition at line 54 of file CSCCondSegFit.h.

Referenced by correctTheCovX(), fit(), and setChi2().

◆ worstHit_

int CSCCondSegFit::worstHit_
private

Definition at line 51 of file CSCCondSegFit.h.

Referenced by correctTheCovX(), and worstHit().