CMS 3D CMS Logo

MuonResiduals6DOFrphiFitter.cc
Go to the documentation of this file.
1 // $Id:$
2 
3 #ifdef STANDALONE_FITTER
5 #else
7 //#include "DataFormats/MuonDetId/interface/CSCDetId.h"
8 #endif
9 
10 #include "TH2F.h"
11 #include "TMath.h"
12 #include "TTree.h"
13 #include "TFile.h"
14 
15 namespace {
16  TMinuit *minuit;
17 
18  double sum_of_weights;
19  double number_of_hits;
20  bool weight_alignment;
21 
22  //const CSCGeometry *cscGeometry;
23  //static CSCDetId cscDetId;
24  double csc_R;
25 
26  double getResidual(double delta_x,
27  double delta_y,
28  double delta_z,
29  double delta_phix,
30  double delta_phiy,
31  double delta_phiz,
32  double track_x,
33  double track_y,
34  double track_dxdz,
35  double track_dydz,
36  double R,
37  double alpha,
38  double resslope) {
39  return delta_x - (track_x / R - 3. * pow(track_x / R, 3)) * delta_y - track_dxdz * delta_z -
40  track_y * track_dxdz * delta_phix + track_x * track_dxdz * delta_phiy - track_y * delta_phiz +
41  resslope * alpha;
42  }
43 
44  double getResSlope(double delta_x,
45  double delta_y,
46  double delta_z,
47  double delta_phix,
48  double delta_phiy,
49  double delta_phiz,
50  double track_x,
51  double track_y,
52  double track_dxdz,
53  double track_dydz,
54  double R) {
55  return -track_dxdz / 2. / R * delta_y + (track_x / R - track_dxdz * track_dydz) * delta_phix +
56  (1. + track_dxdz * track_dxdz) * delta_phiy - track_dydz * delta_phiz;
57  }
58 
59  double effectiveR(double x, double y) {
60  // CSCDetId id = cscDetId;
61  // CSCDetId layerId(id.endcap(), id.station(), id.ring(), id.chamber(), 3);
62  // int strip = cscGeometry->layer(layerId)->geometry()->nearestStrip(align::LocalPoint(x, y, 0.));
63  // double angle = cscGeometry->layer(layerId)->geometry()->stripAngle(strip) - M_PI/2.;
64  // if (fabs(angle) < 1e-10) return csc_R;
65  // else return x/tan(angle) - y;
66  return csc_R;
67  }
68 
69  Double_t residual_trackx_TF1(Double_t *xx, Double_t *p) {
70  return 10. *
71  getResidual(
72  p[0], p[1], p[2], p[3], p[4], p[5], xx[0], p[7], p[8], p[9], effectiveR(xx[0], p[7]), p[10], p[11]);
73  }
74  Double_t residual_tracky_TF1(Double_t *xx, Double_t *p) {
75  return 10. *
76  getResidual(
77  p[0], p[1], p[2], p[3], p[4], p[5], p[6], xx[0], p[8], p[9], effectiveR(p[6], xx[0]), p[10], p[11]);
78  }
79  Double_t residual_trackdxdz_TF1(Double_t *xx, Double_t *p) {
80  return 10. * getResidual(
81  p[0], p[1], p[2], p[3], p[4], p[5], p[6], p[7], xx[0], p[9], effectiveR(p[6], p[7]), p[10], p[11]);
82  }
83  Double_t residual_trackdydz_TF1(Double_t *xx, Double_t *p) {
84  return 10. * getResidual(
85  p[0], p[1], p[2], p[3], p[4], p[5], p[6], p[7], p[8], xx[0], effectiveR(p[6], p[7]), p[10], p[11]);
86  }
87 
88  Double_t resslope_trackx_TF1(Double_t *xx, Double_t *p) {
89  return 1000. * getResSlope(p[0], p[1], p[2], p[3], p[4], p[5], xx[0], p[7], p[8], p[9], effectiveR(xx[0], p[7]));
90  }
91  Double_t resslope_tracky_TF1(Double_t *xx, Double_t *p) {
92  return 1000. * getResSlope(p[0], p[1], p[2], p[3], p[4], p[5], p[6], xx[0], p[8], p[9], effectiveR(p[6], xx[0]));
93  }
94  Double_t resslope_trackdxdz_TF1(Double_t *xx, Double_t *p) {
95  return 1000. * getResSlope(p[0], p[1], p[2], p[3], p[4], p[5], p[6], p[7], xx[0], p[9], effectiveR(p[6], p[7]));
96  }
97  Double_t resslope_trackdydz_TF1(Double_t *xx, Double_t *p) {
98  return 1000. * getResSlope(p[0], p[1], p[2], p[3], p[4], p[5], p[6], p[7], p[8], xx[0], effectiveR(p[6], p[7]));
99  }
100 } // namespace
101 
102 void MuonResiduals6DOFrphiFitter_FCN(int &npar, double *gin, double &fval, double *par, int iflag) {
103  MuonResidualsFitterFitInfo *fitinfo = (MuonResidualsFitterFitInfo *)(minuit->GetObjectFit());
104  MuonResidualsFitter *fitter = fitinfo->fitter();
105 
106  fval = 0.;
107  for (std::vector<double *>::const_iterator resiter = fitter->residuals_begin(); resiter != fitter->residuals_end();
108  ++resiter) {
109  const double residual = (*resiter)[MuonResiduals6DOFrphiFitter::kResid];
110  const double resslope = (*resiter)[MuonResiduals6DOFrphiFitter::kResSlope];
111  const double positionX = (*resiter)[MuonResiduals6DOFrphiFitter::kPositionX];
112  const double positionY = (*resiter)[MuonResiduals6DOFrphiFitter::kPositionY];
113  const double angleX = (*resiter)[MuonResiduals6DOFrphiFitter::kAngleX];
114  const double angleY = (*resiter)[MuonResiduals6DOFrphiFitter::kAngleY];
115  const double redchi2 = (*resiter)[MuonResiduals6DOFrphiFitter::kRedChi2];
116 
117  const double alignx = par[MuonResiduals6DOFrphiFitter::kAlignX];
118  const double aligny = par[MuonResiduals6DOFrphiFitter::kAlignY];
119  const double alignz = par[MuonResiduals6DOFrphiFitter::kAlignZ];
120  const double alignphix = par[MuonResiduals6DOFrphiFitter::kAlignPhiX];
121  const double alignphiy = par[MuonResiduals6DOFrphiFitter::kAlignPhiY];
122  const double alignphiz = par[MuonResiduals6DOFrphiFitter::kAlignPhiZ];
123  const double residsigma = par[MuonResiduals6DOFrphiFitter::kResidSigma];
124  const double resslopesigma = par[MuonResiduals6DOFrphiFitter::kResSlopeSigma];
125  const double alpha = par[MuonResiduals6DOFrphiFitter::kAlpha];
126  const double residgamma = par[MuonResiduals6DOFrphiFitter::kResidGamma];
127  const double resslopegamma = par[MuonResiduals6DOFrphiFitter::kResSlopeGamma];
128 
129  double coeff = alpha;
132  coeff = 0.;
133  double effr = effectiveR(positionX, positionY);
134  double residpeak = getResidual(alignx,
135  aligny,
136  alignz,
137  alignphix,
138  alignphiy,
139  alignphiz,
140  positionX,
141  positionY,
142  angleX,
143  angleY,
144  effr,
145  coeff,
146  resslope);
147  double resslopepeak = getResSlope(
148  alignx, aligny, alignz, alignphix, alignphiy, alignphiz, positionX, positionY, angleX, angleY, effr);
149 
150  double weight = (1. / redchi2) * number_of_hits / sum_of_weights;
151  if (!weight_alignment)
152  weight = 1.;
153 
154  if (!weight_alignment || TMath::Prob(redchi2 * 6, 6) < 0.99) // no spikes allowed
155  {
157  if (fitter->useRes() == MuonResidualsFitter::k1111 || fitter->useRes() == MuonResidualsFitter::k1110 ||
158  fitter->useRes() == MuonResidualsFitter::k1010) {
159  fval += -weight * MuonResidualsFitter_logPureGaussian(residual, residpeak, residsigma);
160  fval += -weight * MuonResidualsFitter_logPureGaussian(resslope, resslopepeak, resslopesigma);
161  } else if (fitter->useRes() == MuonResidualsFitter::k1100) {
162  fval += -weight * MuonResidualsFitter_logPureGaussian(residual, residpeak, residsigma);
163  } else if (fitter->useRes() == MuonResidualsFitter::k0010) {
164  fval += -weight * MuonResidualsFitter_logPureGaussian(resslope, resslopepeak, resslopesigma);
165  }
166  } else if (fitter->residualsModel() == MuonResidualsFitter::kPureGaussian2D) {
167  if (fitter->useRes() == MuonResidualsFitter::k1111 || fitter->useRes() == MuonResidualsFitter::k1110 ||
168  fitter->useRes() == MuonResidualsFitter::k1010) {
170  residual, resslope, residpeak, resslopepeak, residsigma, resslopesigma, alpha);
171  } else if (fitter->useRes() == MuonResidualsFitter::k1100) {
172  fval += -weight * MuonResidualsFitter_logPureGaussian(residual, residpeak, residsigma);
173  } else if (fitter->useRes() == MuonResidualsFitter::k0010) {
174  fval += -weight * MuonResidualsFitter_logPureGaussian(resslope, resslopepeak, resslopesigma);
175  }
176  } else if (fitter->residualsModel() == MuonResidualsFitter::kPowerLawTails) {
177  fval += -weight * MuonResidualsFitter_logPowerLawTails(residual, residpeak, residsigma, residgamma);
178  fval += -weight * MuonResidualsFitter_logPowerLawTails(resslope, resslopepeak, resslopesigma, resslopegamma);
179  } else if (fitter->residualsModel() == MuonResidualsFitter::kROOTVoigt) {
180  fval += -weight * MuonResidualsFitter_logROOTVoigt(residual, residpeak, residsigma, residgamma);
181  fval += -weight * MuonResidualsFitter_logROOTVoigt(resslope, resslopepeak, resslopesigma, resslopegamma);
182  } else if (fitter->residualsModel() == MuonResidualsFitter::kGaussPowerTails) {
183  fval += -weight * MuonResidualsFitter_logGaussPowerTails(residual, residpeak, residsigma);
184  fval += -weight * MuonResidualsFitter_logGaussPowerTails(resslope, resslopepeak, resslopesigma);
185  } else {
186  assert(false);
187  }
188  }
189  }
190 }
191 
193 
194 void MuonResiduals6DOFrphiFitter::inform(TMinuit *tMinuit) { minuit = tMinuit; }
195 
197  sum_of_weights = 0.;
198  number_of_hits = 0.;
199  weight_alignment = m_weightAlignment;
200  for (std::vector<double *>::const_iterator resiter = residuals_begin(); resiter != residuals_end(); ++resiter) {
201  if (m_weightAlignment) {
202  double redchi2 = (*resiter)[kRedChi2];
203  if (TMath::Prob(redchi2 * 6, 6) < 0.99) { // no spikes allowed
204  sum_of_weights += 1. / redchi2;
205  number_of_hits += 1.;
206  }
207  } else {
208  sum_of_weights += 1.;
209  number_of_hits += 1.;
210  }
211  }
212  return sum_of_weights;
213 }
214 
216  //cscGeometry = m_cscGeometry;
217  //cscDetId = CSCDetId(ali->geomDetId().rawId());
218 
219 #ifndef STANDALONE_FITTER
220  csc_R = sqrt(pow(ali->globalPosition().x(), 2) + pow(ali->globalPosition().y(), 2));
221 #else
222  csc_R = 200; // not important what number it is, as we usually not use the DOF where it matters
223 #endif
224 
225  initialize_table(); // if not already initialized
226  sumofweights();
227 
228  double res_std = 0.5;
229  double resslope_std = 0.002;
230 
231  int nums[10] = {kAlignX,
232  kAlignZ,
233  kAlignPhiX,
234  kAlignPhiY,
235  kAlignPhiZ,
236  kResidSigma,
238  kAlpha,
239  kResidGamma,
241  std::string names[10] = {"AlignX",
242  "AlignZ",
243  "AlignPhiX",
244  "AlignPhiY",
245  "AlignPhiZ",
246  "ResidSigma",
247  "ResSlopeSigma",
248  "Alpha",
249  "ResidGamma",
250  "ResSlopeGamma"};
251  double starts[10] = {0., 0., 0., 0., 0., res_std, resslope_std, 0., 0.1 * res_std, 0.1 * resslope_std};
252  double steps[10] = {
253  0.1, 0.1, 0.001, 0.001, 0.001, 0.001 * res_std, 0.001 * resslope_std, 0.001, 0.01 * res_std, 0.01 * resslope_std};
254  double lows[10] = {0., 0., 0., 0., 0., 0., 0., -1., 0., 0.};
255  double highs[10] = {0., 0., 0., 0., 0., 10., 0.1, 1., 0., 0.};
256 
257  std::vector<int> num(nums, nums + 5);
258  std::vector<std::string> name(names, names + 5);
259  std::vector<double> start(starts, starts + 5);
260  std::vector<double> step(steps, steps + 5);
261  std::vector<double> low(lows, lows + 5);
262  std::vector<double> high(highs, highs + 5);
263 
264  bool add_alpha = (residualsModel() == kPureGaussian2D);
265  bool add_gamma = (residualsModel() == kROOTVoigt || residualsModel() == kPowerLawTails);
266 
267  int idx[4], ni = 0;
268  if (useRes() == k1111 || useRes() == k1110 || useRes() == k1010) {
269  for (ni = 0; ni < 2; ni++)
270  idx[ni] = ni + 5;
271  if (add_alpha)
272  idx[ni++] = 7;
273  else if (add_gamma)
274  for (; ni < 4; ni++)
275  idx[ni] = ni + 6;
276  if (!add_alpha)
277  fix(kAlpha);
278  } else if (useRes() == k1100) {
279  idx[ni++] = 5;
280  if (add_gamma)
281  idx[ni++] = 8;
283  fix(kAlpha);
284  } else if (useRes() == k0010) {
285  idx[ni++] = 6;
286  if (add_gamma)
287  idx[ni++] = 9;
288  fix(kResidSigma);
289  fix(kAlpha);
290  }
291  for (int i = 0; i < ni; i++) {
292  num.push_back(nums[idx[i]]);
293  name.push_back(names[idx[i]]);
294  start.push_back(starts[idx[i]]);
295  step.push_back(steps[idx[i]]);
296  low.push_back(lows[idx[i]]);
297  high.push_back(highs[idx[i]]);
298  }
299 
301 }
302 
304  sumofweights();
305 
306  double mean_residual = 0., mean_resslope = 0.;
307  double mean_trackx = 0., mean_tracky = 0., mean_trackdxdz = 0., mean_trackdydz = 0.;
308  double sum_w = 0.;
309 
310  for (std::vector<double *>::const_iterator rit = residuals_begin(); rit != residuals_end(); ++rit) {
311  const double redchi2 = (*rit)[kRedChi2];
312  double weight = 1. / redchi2;
313  if (!m_weightAlignment)
314  weight = 1.;
315 
316  if (!m_weightAlignment || TMath::Prob(redchi2 * 6, 6) < 0.99) // no spikes allowed
317  {
318  double factor_w = 1. / (sum_w + weight);
319  mean_residual = factor_w * (sum_w * mean_residual + weight * (*rit)[kResid]);
320  mean_resslope = factor_w * (sum_w * mean_resslope + weight * (*rit)[kResSlope]);
321  mean_trackx = factor_w * (sum_w * mean_trackx + weight * (*rit)[kPositionX]);
322  mean_tracky = factor_w * (sum_w * mean_tracky + weight * (*rit)[kPositionY]);
323  mean_trackdxdz = factor_w * (sum_w * mean_trackdxdz + weight * (*rit)[kAngleX]);
324  mean_trackdydz = factor_w * (sum_w * mean_trackdydz + weight * (*rit)[kAngleY]);
325  sum_w += weight;
326  }
327  }
328 
329  std::string name_residual, name_resslope, name_residual_raw, name_resslope_raw, name_residual_cut, name_alpha;
330  std::string name_residual_trackx, name_resslope_trackx;
331  std::string name_residual_tracky, name_resslope_tracky;
332  std::string name_residual_trackdxdz, name_resslope_trackdxdz;
333  std::string name_residual_trackdydz, name_resslope_trackdydz;
334 
335  name_residual = name + "_residual";
336  name_resslope = name + "_resslope";
337  name_residual_raw = name + "_residual_raw";
338  name_resslope_raw = name + "_resslope_raw";
339  name_residual_cut = name + "_residual_cut";
340  name_alpha = name + "_alpha";
341  name_residual_trackx = name + "_residual_trackx";
342  name_resslope_trackx = name + "_resslope_trackx";
343  name_residual_tracky = name + "_residual_tracky";
344  name_resslope_tracky = name + "_resslope_tracky";
345  name_residual_trackdxdz = name + "_residual_trackdxdz";
346  name_resslope_trackdxdz = name + "_resslope_trackdxdz";
347  name_residual_trackdydz = name + "_residual_trackdydz";
348  name_resslope_trackdydz = name + "_resslope_trackdydz";
349 
350  double width = ali->surface().width();
351  double length = ali->surface().length();
352  int bins_residual = 100, bins_resslope = 100;
353  double min_residual = -50., max_residual = 50.;
354  double min_resslope = -50., max_resslope = 50.;
355  double min_trackx = -width / 2., max_trackx = width / 2.;
356  double min_tracky = -length / 2., max_tracky = length / 2.;
357  double min_trackdxdz = -1.5, max_trackdxdz = 1.5;
358  double min_trackdydz = -1.5, max_trackdydz = 1.5;
359 
360  TH1F *hist_residual = dir->make<TH1F>(name_residual.c_str(), "", bins_residual, min_residual, max_residual);
361  TH1F *hist_resslope = dir->make<TH1F>(name_resslope.c_str(), "", bins_resslope, min_resslope, max_resslope);
362  TH1F *hist_residual_raw = dir->make<TH1F>(name_residual_raw.c_str(), "", bins_residual, min_residual, max_residual);
363  TH1F *hist_resslope_raw = dir->make<TH1F>(name_resslope_raw.c_str(), "", bins_resslope, min_resslope, max_resslope);
364  TH1F *hist_residual_cut = dir->make<TH1F>(name_residual_cut.c_str(), "", bins_residual, min_residual, max_residual);
365  TH2F *hist_alpha = dir->make<TH2F>(name_alpha.c_str(), "", 50, min_resslope, max_resslope, 50, -50., 50.);
366 
367  TProfile *hist_residual_trackx =
368  dir->make<TProfile>(name_residual_trackx.c_str(), "", 50, min_trackx, max_trackx, min_residual, max_residual);
369  TProfile *hist_resslope_trackx =
370  dir->make<TProfile>(name_resslope_trackx.c_str(), "", 50, min_trackx, max_trackx, min_resslope, max_resslope);
371  TProfile *hist_residual_tracky =
372  dir->make<TProfile>(name_residual_tracky.c_str(), "", 50, min_tracky, max_tracky, min_residual, max_residual);
373  TProfile *hist_resslope_tracky =
374  dir->make<TProfile>(name_resslope_tracky.c_str(), "", 50, min_tracky, max_tracky, min_resslope, max_resslope);
375  TProfile *hist_residual_trackdxdz = dir->make<TProfile>(
376  name_residual_trackdxdz.c_str(), "", 500, min_trackdxdz, max_trackdxdz, min_residual, max_residual);
377  TProfile *hist_resslope_trackdxdz = dir->make<TProfile>(
378  name_resslope_trackdxdz.c_str(), "", 500, min_trackdxdz, max_trackdxdz, min_resslope, max_resslope);
379  TProfile *hist_residual_trackdydz = dir->make<TProfile>(
380  name_residual_trackdydz.c_str(), "", 500, min_trackdydz, max_trackdydz, min_residual, max_residual);
381  TProfile *hist_resslope_trackdydz = dir->make<TProfile>(
382  name_resslope_trackdydz.c_str(), "", 500, min_trackdydz, max_trackdydz, min_resslope, max_resslope);
383 
384  hist_residual_trackx->SetAxisRange(-10., 10., "Y");
385  hist_resslope_trackx->SetAxisRange(-10., 10., "Y");
386  hist_residual_tracky->SetAxisRange(-10., 10., "Y");
387  hist_resslope_tracky->SetAxisRange(-10., 10., "Y");
388  hist_residual_trackdxdz->SetAxisRange(-10., 10., "Y");
389  hist_resslope_trackdxdz->SetAxisRange(-10., 10., "Y");
390  hist_residual_trackdydz->SetAxisRange(-10., 10., "Y");
391  hist_resslope_trackdydz->SetAxisRange(-10., 10., "Y");
392 
393  name_residual += "_fit";
394  name_resslope += "_fit";
395  name_alpha += "_fit";
396  name_residual_trackx += "_fit";
397  name_resslope_trackx += "_fit";
398  name_residual_tracky += "_fit";
399  name_resslope_tracky += "_fit";
400  name_residual_trackdxdz += "_fit";
401  name_resslope_trackdxdz += "_fit";
402  name_residual_trackdydz += "_fit";
403  name_resslope_trackdydz += "_fit";
404 
405  TF1 *fit_residual = nullptr;
406  TF1 *fit_resslope = nullptr;
408  fit_residual = new TF1(name_residual.c_str(), MuonResidualsFitter_pureGaussian_TF1, min_residual, max_residual, 3);
409  fit_residual->SetParameters(
410  sum_of_weights * (max_residual - min_residual) / bins_residual, 10. * value(kAlignX), 10. * value(kResidSigma));
411  const double er_res[3] = {0., 10. * errorerror(kAlignX), 10. * errorerror(kResidSigma)};
412  fit_residual->SetParErrors(er_res);
413  fit_resslope = new TF1(name_resslope.c_str(), MuonResidualsFitter_pureGaussian_TF1, min_resslope, max_resslope, 3);
414  fit_resslope->SetParameters(sum_of_weights * (max_resslope - min_resslope) / bins_resslope,
415  1000. * value(kAlignPhiY),
416  1000. * value(kResSlopeSigma));
417  const double er_resslope[3] = {0., 1000. * errorerror(kAlignPhiY), 1000. * errorerror(kResSlopeSigma)};
418  fit_resslope->SetParErrors(er_resslope);
419  } else if (residualsModel() == kPowerLawTails) {
420  fit_residual = new TF1(name_residual.c_str(), MuonResidualsFitter_powerLawTails_TF1, min_residual, max_residual, 4);
421  fit_residual->SetParameters(sum_of_weights * (max_residual - min_residual) / bins_residual,
422  10. * value(kAlignX),
423  10. * value(kResidSigma),
424  10. * value(kResidGamma));
425  fit_resslope = new TF1(name_resslope.c_str(), MuonResidualsFitter_powerLawTails_TF1, min_resslope, max_resslope, 4);
426  fit_resslope->SetParameters(sum_of_weights * (max_resslope - min_resslope) / bins_resslope,
427  1000. * value(kAlignPhiY),
428  1000. * value(kResSlopeSigma),
429  1000. * value(kResSlopeGamma));
430  } else if (residualsModel() == kROOTVoigt) {
431  fit_residual = new TF1(name_residual.c_str(), MuonResidualsFitter_ROOTVoigt_TF1, min_residual, max_residual, 4);
432  fit_residual->SetParameters(sum_of_weights * (max_residual - min_residual) / bins_residual,
433  10. * value(kAlignX),
434  10. * value(kResidSigma),
435  10. * value(kResidGamma));
436  fit_resslope = new TF1(name_resslope.c_str(), MuonResidualsFitter_ROOTVoigt_TF1, min_resslope, max_resslope, 4);
437  fit_resslope->SetParameters(sum_of_weights * (max_resslope - min_resslope) / bins_resslope,
438  1000. * value(kAlignPhiY),
439  1000. * value(kResSlopeSigma),
440  1000. * value(kResSlopeGamma));
441  } else if (residualsModel() == kGaussPowerTails) {
442  fit_residual =
443  new TF1(name_residual.c_str(), MuonResidualsFitter_GaussPowerTails_TF1, min_residual, max_residual, 3);
444  fit_residual->SetParameters(
445  sum_of_weights * (max_residual - min_residual) / bins_residual, 10. * value(kAlignX), 10. * value(kResidSigma));
446  fit_resslope =
447  new TF1(name_resslope.c_str(), MuonResidualsFitter_GaussPowerTails_TF1, min_resslope, max_resslope, 3);
448  fit_resslope->SetParameters(sum_of_weights * (max_resslope - min_resslope) / bins_resslope,
449  1000. * value(kAlignPhiY),
450  1000. * value(kResSlopeSigma));
451  } else {
452  assert(false);
453  }
454 
455  fit_residual->SetLineColor(2);
456  fit_residual->SetLineWidth(2);
457  fit_residual->Write();
458  fit_resslope->SetLineColor(2);
459  fit_resslope->SetLineWidth(2);
460  fit_resslope->Write();
461 
462  TF1 *fit_alpha = new TF1(name_alpha.c_str(), "[0] + x*[1]", min_resslope, max_resslope);
463  double a = 10. * value(kAlignX), b = 10. * value(kAlpha) / 1000.;
464  if (residualsModel() == kPureGaussian2D) {
465  double sx = 10. * value(kResidSigma), sy = 1000. * value(kResSlopeSigma), r = value(kAlpha);
466  a = mean_residual;
467  b = 0.;
468  if (sx != 0.) {
469  b = 1. / (sy / sx * r);
470  a = -b * mean_resslope;
471  }
472  }
473  fit_alpha->SetParameters(a, b);
474  fit_alpha->SetLineColor(2);
475  fit_alpha->SetLineWidth(2);
476  fit_alpha->Write();
477 
478  TProfile *fit_residual_trackx = dir->make<TProfile>(name_residual_trackx.c_str(), "", 100, min_trackx, max_trackx);
479  TProfile *fit_resslope_trackx = dir->make<TProfile>(name_resslope_trackx.c_str(), "", 100, min_trackx, max_trackx);
480  TProfile *fit_residual_tracky = dir->make<TProfile>(name_residual_tracky.c_str(), "", 100, min_tracky, max_tracky);
481  TProfile *fit_resslope_tracky = dir->make<TProfile>(name_resslope_tracky.c_str(), "", 100, min_tracky, max_tracky);
482  TProfile *fit_residual_trackdxdz =
483  dir->make<TProfile>(name_residual_trackdxdz.c_str(), "", 500, min_trackdxdz, max_trackdxdz);
484  TProfile *fit_resslope_trackdxdz =
485  dir->make<TProfile>(name_resslope_trackdxdz.c_str(), "", 500, min_trackdxdz, max_trackdxdz);
486  TProfile *fit_residual_trackdydz =
487  dir->make<TProfile>(name_residual_trackdydz.c_str(), "", 500, min_trackdydz, max_trackdydz);
488  TProfile *fit_resslope_trackdydz =
489  dir->make<TProfile>(name_resslope_trackdydz.c_str(), "", 500, min_trackdydz, max_trackdydz);
490 
491  fit_residual_trackx->SetLineColor(2);
492  fit_residual_trackx->SetLineWidth(2);
493  fit_resslope_trackx->SetLineColor(2);
494  fit_resslope_trackx->SetLineWidth(2);
495  fit_residual_tracky->SetLineColor(2);
496  fit_residual_tracky->SetLineWidth(2);
497  fit_resslope_tracky->SetLineColor(2);
498  fit_resslope_tracky->SetLineWidth(2);
499  fit_residual_trackdxdz->SetLineColor(2);
500  fit_residual_trackdxdz->SetLineWidth(2);
501  fit_resslope_trackdxdz->SetLineColor(2);
502  fit_resslope_trackdxdz->SetLineWidth(2);
503  fit_residual_trackdydz->SetLineColor(2);
504  fit_residual_trackdydz->SetLineWidth(2);
505  fit_resslope_trackdydz->SetLineColor(2);
506  fit_resslope_trackdydz->SetLineWidth(2);
507 
508  name_residual_trackx += "line";
509  name_resslope_trackx += "line";
510  name_residual_tracky += "line";
511  name_resslope_tracky += "line";
512  name_residual_trackdxdz += "line";
513  name_resslope_trackdxdz += "line";
514  name_residual_trackdydz += "line";
515  name_resslope_trackdydz += "line";
516 
517  TF1 *fitline_residual_trackx = new TF1(name_residual_trackx.c_str(), residual_trackx_TF1, min_trackx, max_trackx, 12);
518  TF1 *fitline_resslope_trackx = new TF1(name_resslope_trackx.c_str(), resslope_trackx_TF1, min_trackx, max_trackx, 12);
519  TF1 *fitline_residual_tracky = new TF1(name_residual_tracky.c_str(), residual_tracky_TF1, min_tracky, max_tracky, 12);
520  TF1 *fitline_resslope_tracky = new TF1(name_resslope_tracky.c_str(), resslope_tracky_TF1, min_tracky, max_tracky, 12);
521  TF1 *fitline_residual_trackdxdz =
522  new TF1(name_residual_trackdxdz.c_str(), residual_trackdxdz_TF1, min_trackdxdz, max_trackdxdz, 12);
523  TF1 *fitline_resslope_trackdxdz =
524  new TF1(name_resslope_trackdxdz.c_str(), resslope_trackdxdz_TF1, min_trackdxdz, max_trackdxdz, 12);
525  TF1 *fitline_residual_trackdydz =
526  new TF1(name_residual_trackdydz.c_str(), residual_trackdydz_TF1, min_trackdydz, max_trackdydz, 12);
527  TF1 *fitline_resslope_trackdydz =
528  new TF1(name_resslope_trackdydz.c_str(), resslope_trackdydz_TF1, min_trackdydz, max_trackdydz, 12);
529 
530  std::vector<TF1 *> fitlines;
531  fitlines.push_back(fitline_residual_trackx);
532  fitlines.push_back(fitline_resslope_trackx);
533  fitlines.push_back(fitline_residual_tracky);
534  fitlines.push_back(fitline_resslope_tracky);
535  fitlines.push_back(fitline_residual_trackdxdz);
536  fitlines.push_back(fitline_resslope_trackdxdz);
537  fitlines.push_back(fitline_residual_trackdydz);
538  fitlines.push_back(fitline_resslope_trackdydz);
539 
540  double fitparameters[12] = {value(kAlignX),
541  value(kAlignY),
542  value(kAlignZ),
543  value(kAlignPhiX),
544  value(kAlignPhiY),
545  value(kAlignPhiZ),
546  mean_trackx,
547  mean_tracky,
548  mean_trackdxdz,
549  mean_trackdydz,
550  value(kAlpha),
551  mean_resslope};
553  fitparameters[10] = 0.;
554 
555  for (std::vector<TF1 *>::const_iterator itr = fitlines.begin(); itr != fitlines.end(); itr++) {
556  (*itr)->SetParameters(fitparameters);
557  (*itr)->SetLineColor(2);
558  (*itr)->SetLineWidth(2);
559  (*itr)->Write();
560  }
561 
562  for (std::vector<double *>::const_iterator resiter = residuals_begin(); resiter != residuals_end(); ++resiter) {
563  const double resid = (*resiter)[kResid];
564  const double resslope = (*resiter)[kResSlope];
565  const double positionX = (*resiter)[kPositionX];
566  const double positionY = (*resiter)[kPositionY];
567  const double angleX = (*resiter)[kAngleX];
568  const double angleY = (*resiter)[kAngleY];
569  const double redchi2 = (*resiter)[kRedChi2];
570  double weight = 1. / redchi2;
571  if (!m_weightAlignment)
572  weight = 1.;
573 
574  if (!m_weightAlignment || TMath::Prob(redchi2 * 6, 6) < 0.99) { // no spikes allowed
575  hist_alpha->Fill(1000. * resslope, 10. * resid);
576 
577  double coeff = value(kAlpha);
579  coeff = 0.;
580  double geom_resid = getResidual(value(kAlignX),
581  value(kAlignY),
582  value(kAlignZ),
583  value(kAlignPhiX),
584  value(kAlignPhiY),
585  value(kAlignPhiZ),
586  positionX,
587  positionY,
588  angleX,
589  angleY,
590  effectiveR(positionX, positionY),
591  coeff,
592  resslope);
593  hist_residual->Fill(10. * (resid - geom_resid + value(kAlignX)), weight);
594  hist_residual_trackx->Fill(positionX, 10. * resid, weight);
595  hist_residual_tracky->Fill(positionY, 10. * resid, weight);
596  hist_residual_trackdxdz->Fill(angleX, 10. * resid, weight);
597  hist_residual_trackdydz->Fill(angleY, 10. * resid, weight);
598  fit_residual_trackx->Fill(positionX, 10. * geom_resid, weight);
599  fit_residual_tracky->Fill(positionY, 10. * geom_resid, weight);
600  fit_residual_trackdxdz->Fill(angleX, 10. * geom_resid, weight);
601  fit_residual_trackdydz->Fill(angleY, 10. * geom_resid, weight);
602 
603  double geom_resslope = getResSlope(value(kAlignX),
604  value(kAlignY),
605  value(kAlignZ),
606  value(kAlignPhiX),
607  value(kAlignPhiY),
608  value(kAlignPhiZ),
609  positionX,
610  positionY,
611  angleX,
612  angleY,
613  effectiveR(positionX, positionY));
614  hist_resslope->Fill(1000. * (resslope - geom_resslope + value(kAlignPhiY)), weight);
615  hist_resslope_trackx->Fill(positionX, 1000. * resslope, weight);
616  hist_resslope_tracky->Fill(positionY, 1000. * resslope, weight);
617  hist_resslope_trackdxdz->Fill(angleX, 1000. * resslope, weight);
618  hist_resslope_trackdydz->Fill(angleY, 1000. * resslope, weight);
619  fit_resslope_trackx->Fill(positionX, 1000. * geom_resslope, weight);
620  fit_resslope_tracky->Fill(positionY, 1000. * geom_resslope, weight);
621  fit_resslope_trackdxdz->Fill(angleX, 1000. * geom_resslope, weight);
622  fit_resslope_trackdydz->Fill(angleY, 1000. * geom_resslope, weight);
623  }
624 
625  hist_residual_raw->Fill(10. * resid);
626  hist_resslope_raw->Fill(1000. * resslope);
627  if (fabs(resslope) < 0.005)
628  hist_residual_cut->Fill(10. * resid);
629  }
630 
631  double chi2 = 0.;
632  double ndof = 0.;
633  for (int i = 1; i <= hist_residual->GetNbinsX(); i++) {
634  double xi = hist_residual->GetBinCenter(i);
635  double yi = hist_residual->GetBinContent(i);
636  double yerri = hist_residual->GetBinError(i);
637  double yth = fit_residual->Eval(xi);
638  if (yerri > 0.) {
639  chi2 += pow((yth - yi) / yerri, 2);
640  ndof += 1.;
641  }
642  }
643  for (int i = 1; i <= hist_resslope->GetNbinsX(); i++) {
644  double xi = hist_resslope->GetBinCenter(i);
645  double yi = hist_resslope->GetBinContent(i);
646  double yerri = hist_resslope->GetBinError(i);
647  double yth = fit_resslope->Eval(xi);
648  if (yerri > 0.) {
649  chi2 += pow((yth - yi) / yerri, 2);
650  ndof += 1.;
651  }
652  }
653  ndof -= npar();
654 
655  return (ndof > 0. ? chi2 / ndof : -1.);
656 }
657 
659  unsigned int endcap,
660  unsigned int station,
661  unsigned int ring,
662  unsigned int chamber,
663  unsigned int preselected) {
664  TFile *f = new TFile(fname.c_str());
665  TTree *t = (TTree *)f->Get("mual_ttree");
666 
667  // Create new temporary file
668  TFile *tmpf = new TFile("small_tree_fit.root", "recreate");
669  assert(tmpf != nullptr);
670 
671  // filter the tree (temporarily lives in the new file)
672  TTree *tt = t->CopyTree(Form("!is_dt && is_plus==%d && ring_wheel==%d && station==%d && sector==%d && select==%d",
673  2 - endcap,
674  station,
675  ring,
676  chamber,
677  (bool)preselected));
678 
680  tt->SetBranchAddress("res_x", &r.res_x);
681  tt->SetBranchAddress("res_slope_x", &r.res_slope_x);
682  tt->SetBranchAddress("pos_x", &r.pos_x);
683  tt->SetBranchAddress("pos_y", &r.pos_y);
684  tt->SetBranchAddress("angle_x", &r.angle_x);
685  tt->SetBranchAddress("angle_y", &r.angle_y);
686  tt->SetBranchAddress("pz", &r.pz);
687  tt->SetBranchAddress("pt", &r.pt);
688  tt->SetBranchAddress("q", &r.q);
689 
690  Long64_t nentries = tt->GetEntries();
691  for (Long64_t i = 0; i < nentries; i++) {
692  tt->GetEntry(i);
693  double *rdata = new double[MuonResiduals6DOFrphiFitter::kNData];
694  rdata[kResid] = r.res_x;
695  rdata[kResSlope] = r.res_slope_x;
696  rdata[kPositionX] = r.pos_x;
697  rdata[kPositionY] = r.pos_y;
698  rdata[kAngleX] = r.angle_x;
699  rdata[kAngleY] = r.angle_y;
700  rdata[kRedChi2] = 0.1;
701  rdata[kPz] = r.pz;
702  rdata[kPt] = r.pt;
703  rdata[kCharge] = r.q;
704  fill(rdata);
705  }
706  delete f;
707  //delete tmpf;
708  return tt;
709 }
MuonResidualsFitter
Definition: MuonResidualsFitter.h:84
MuonResiduals6DOFrphiFitter::kAlpha
Definition: MuonResiduals6DOFrphiFitter.h:30
MuonResiduals6DOFrphiFitter::kResid
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResidualsFitter::correctBField
virtual void correctBField()=0
ApeEstimator_cff.width
width
Definition: ApeEstimator_cff.py:24
MuonResiduals6DOFrphiFitter::kPz
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResidualsFitter_pureGaussian_TF1
Double_t MuonResidualsFitter_pureGaussian_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:38
mps_fire.i
i
Definition: mps_fire.py:355
start
Definition: start.py:1
MuonResiduals6DOFrphiFitter::kRedChi2
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResiduals6DOFrphiFitter::kNData
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResidualsFitter::value
double value(int parNum)
Definition: MuonResidualsFitter.h:145
groupFilesInBlocks.tt
int tt
Definition: groupFilesInBlocks.py:144
step
step
Definition: StallMonitor.cc:94
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
MuonResidualsFitter::errorerror
double errorerror(int parNum)
Definition: MuonResidualsFitter.h:149
MuonResidualsFitter_logROOTVoigt
double MuonResidualsFitter_logROOTVoigt(double residual, double center, double sigma, double gamma)
Definition: MuonResidualsFitter.cc:120
zMuMuMuonUserData.alpha
alpha
zGenParticlesMatch = cms.InputTag(""),
Definition: zMuMuMuonUserData.py:9
relativeConstraints.station
station
Definition: relativeConstraints.py:67
mps_merge.weight
weight
Definition: mps_merge.py:88
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
MuonResidualsFitter::dofit
bool dofit(void(*fcn)(int &, double *, double &, double *, int), std::vector< int > &parNum, std::vector< std::string > &parName, std::vector< double > &start, std::vector< double > &step, std::vector< double > &low, std::vector< double > &high)
Definition: MuonResidualsFitter.cc:272
Alignable
Definition: Alignable.h:27
MuonResidualsFitter::k1100
Definition: MuonResidualsFitter.h:90
hybridSuperClusters_cfi.xi
xi
Definition: hybridSuperClusters_cfi.py:10
cms::cuda::assert
assert(be >=bs)
MuonResiduals6DOFrphiFitter::correctBField
void correctBField() override
Definition: MuonResiduals6DOFrphiFitter.cc:192
MuonResidualsFitter::kGaussPowerTails
Definition: MuonResidualsFitter.h:86
TFileDirectory
Definition: TFileDirectory.h:24
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
makeMuonMisalignmentScenario.endcap
endcap
Definition: makeMuonMisalignmentScenario.py:320
AlignableSurface::length
align::Scalar length() const
Definition: AlignableSurface.h:32
tools.TF1
TF1
Definition: tools.py:23
training_settings.idx
idx
Definition: training_settings.py:16
MuonResiduals6DOFrphiFitter::readNtuple
TTree * readNtuple(std::string fname, unsigned int endcap, unsigned int station, unsigned int ring, unsigned int chamber, unsigned int preselected=1)
Definition: MuonResiduals6DOFrphiFitter.cc:658
ndof
Definition: HIMultiTrackSelector.h:49
MuonResidualsFitter::MuonAlignmentTreeRow
Definition: MuonResidualsFitter.h:92
MuonResidualsFitter::k1111
Definition: MuonResidualsFitter.h:90
MuonResiduals6DOFrphiFitter::kAlignPhiZ
Definition: MuonResiduals6DOFrphiFitter.h:27
MuonResiduals6DOFrphiFitter::kAngleY
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResiduals6DOFrphiFitter.h
names
const std::string names[nVars_]
Definition: PhotonIDValueMapProducer.cc:122
MuonResidualsFitterFitInfo
Definition: MuonResidualsFitter.h:250
MuonResiduals6DOFrphiFitter::kResidGamma
Definition: MuonResiduals6DOFrphiFitter.h:31
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
vertices_cff.x
x
Definition: vertices_cff.py:29
MuonResidualsFitter::kPowerLawTails
Definition: MuonResidualsFitter.h:86
MuonResiduals6DOFrphiFitter_FCN
void MuonResiduals6DOFrphiFitter_FCN(int &npar, double *gin, double &fval, double *par, int iflag)
Definition: MuonResiduals6DOFrphiFitter.cc:102
MuonResidualsFitter_powerLawTails_TF1
Double_t MuonResidualsFitter_powerLawTails_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:116
MuonResiduals6DOFrphiFitter::kPositionX
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResiduals6DOFrphiFitter::sumofweights
double sumofweights() override
Definition: MuonResiduals6DOFrphiFitter.cc:196
Alignable::surface
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:132
MuonResidualsFitter::useRes
int useRes(int pattern=-1)
Definition: MuonResidualsFitter.h:119
MuonResiduals6DOFrphiFitter::kResSlope
Definition: MuonResiduals6DOFrphiFitter.h:36
OrderedSet.t
t
Definition: OrderedSet.py:90
b
double b
Definition: hdecay.h:118
MuonResiduals6DOFrphiFitter::fit
bool fit(Alignable *ali) override
Definition: MuonResiduals6DOFrphiFitter.cc:215
MuonResidualsFitter::k1010
Definition: MuonResidualsFitter.h:90
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
MuonResidualsFitter_ROOTVoigt_TF1
Double_t MuonResidualsFitter_ROOTVoigt_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:124
MuonResidualsFitter::kPureGaussian
Definition: MuonResidualsFitter.h:86
a
double a
Definition: hdecay.h:119
MuonResidualsFitter::k1110
Definition: MuonResidualsFitter.h:90
MuonResiduals6DOFrphiFitter::inform
void inform(TMinuit *tMinuit) override
Definition: MuonResiduals6DOFrphiFitter.cc:194
MuonResiduals6DOFrphiFitter::npar
int npar() override
Definition: MuonResiduals6DOFrphiFitter.h:51
MuonResidualsFitter::residuals_end
std::vector< double * >::const_iterator residuals_end() const
Definition: MuonResidualsFitter.h:191
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
MuonResiduals6DOFrphiFitter::kResSlopeSigma
Definition: MuonResiduals6DOFrphiFitter.h:29
MuonResiduals6DOFrphiFitter::kAlignPhiX
Definition: MuonResiduals6DOFrphiFitter.h:25
MuonResiduals6DOFrphiFitter::kAlignY
Definition: MuonResiduals6DOFrphiFitter.h:23
fftjetcommon_cfi.sy
sy
Definition: fftjetcommon_cfi.py:203
MuonResidualsFitter_logPowerLawTails
double MuonResidualsFitter_logPowerLawTails(double residual, double center, double sigma, double gamma)
Definition: MuonResidualsFitter.cc:85
itr
std::vector< std::pair< float, float > >::iterator itr
Definition: HGCDigitizer.cc:28
MuonResidualsFitter::fill
void fill(double *residual)
Definition: MuonResidualsFitter.cc:186
MuonResiduals6DOFrphiFitter::kPositionY
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResidualsFitter_logPureGaussian2D
double MuonResidualsFitter_logPureGaussian2D(double x, double y, double x0, double y0, double sx, double sy, double r)
Definition: MuonResidualsFitter.cc:43
MuonResidualsFitter::kROOTVoigt
Definition: MuonResidualsFitter.h:86
MuonResiduals6DOFrphiFitter::kResidSigma
Definition: MuonResiduals6DOFrphiFitter.h:28
EgammaValidation_cff.num
num
Definition: EgammaValidation_cff.py:34
alignCSCRings.r
r
Definition: alignCSCRings.py:93
MuonResidualsFitter::m_weightAlignment
bool m_weightAlignment
Definition: MuonResidualsFitter.h:227
alignmentValidation.fname
string fname
main script
Definition: alignmentValidation.py:959
LaserClient_cfi.high
high
Definition: LaserClient_cfi.py:50
MuonResiduals6DOFrphiFitter::kAlignX
Definition: MuonResiduals6DOFrphiFitter.h:22
MuonResidualsFitterFitInfo::fitter
MuonResidualsFitter * fitter()
Definition: MuonResidualsFitter.h:253
AlignableSurface::width
align::Scalar width() const
Definition: AlignableSurface.h:30
MuonResiduals6DOFrphiFitter::kAlignPhiY
Definition: MuonResiduals6DOFrphiFitter.h:26
MuonResidualsFitter::k0010
Definition: MuonResidualsFitter.h:90
relativeConstraints.ring
ring
Definition: relativeConstraints.py:68
detailsBasic3DVector::y
float float y
Definition: extBasic3DVector.h:14
MuonResidualsFitter::fix
void fix(int parNum, bool dofix=true)
Definition: MuonResidualsFitter.cc:171
relativeConstraints.chamber
chamber
Definition: relativeConstraints.py:53
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
MuonResidualsFitter_logGaussPowerTails
double MuonResidualsFitter_logGaussPowerTails(double residual, double center, double sigma)
Definition: MuonResidualsFitter.cc:128
MuonResidualsFitter::initialize_table
void initialize_table()
Definition: MuonResidualsFitter.cc:205
MuonResiduals6DOFrphiFitter::kAlignZ
Definition: MuonResiduals6DOFrphiFitter.h:24
MuonResidualsFitter_GaussPowerTails_TF1
Double_t MuonResidualsFitter_GaussPowerTails_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:142
customisers.steps
steps
Definition: customisers.py:41
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
MuonResiduals6DOFrphiFitter::kPt
Definition: MuonResiduals6DOFrphiFitter.h:36
MuonResidualsFitter::kPureGaussian2D
Definition: MuonResidualsFitter.h:86
command_line.start
start
Definition: command_line.py:167
Alignable::globalPosition
const PositionType & globalPosition() const
Return the global position of the object.
Definition: Alignable.h:135
MuonResidualsFitter::residualsModel
int residualsModel() const
Definition: MuonResidualsFitter.h:124
dttmaxenums::R
Definition: DTTMax.h:29
MuonResidualsFitter::residuals_begin
std::vector< double * >::const_iterator residuals_begin() const
Definition: MuonResidualsFitter.h:190
MuonResiduals6DOFrphiFitter::plot
double plot(std::string name, TFileDirectory *dir, Alignable *ali) override
Definition: MuonResiduals6DOFrphiFitter.cc:303
LaserClient_cfi.low
low
Definition: LaserClient_cfi.py:52
MuonResiduals6DOFrphiFitter::kAngleX
Definition: MuonResiduals6DOFrphiFitter.h:36
geometryCSVtoXML.xx
xx
Definition: geometryCSVtoXML.py:19
weight
Definition: weight.py:1
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23
MuonResidualsFitter_logPureGaussian
double MuonResidualsFitter_logPureGaussian(double residual, double center, double sigma)
Definition: MuonResidualsFitter.cc:31
MuonResiduals6DOFrphiFitter::kResSlopeGamma
Definition: MuonResiduals6DOFrphiFitter.h:32
MuonResiduals6DOFrphiFitter::kCharge
Definition: MuonResiduals6DOFrphiFitter.h:36
fftjetcommon_cfi.sx
sx
Definition: fftjetcommon_cfi.py:202