CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BSFitter.cc
Go to the documentation of this file.
1 
15 #include "Minuit2/VariableMetricMinimizer.h"
16 #include "Minuit2/FunctionMinimum.h"
17 #include "Minuit2/MnPrint.h"
18 #include "Minuit2/MnMigrad.h"
19 #include "Minuit2/MnUserParameterState.h"
20 //#include "CLHEP/config/CLHEP.h"
21 
22 // C++ standard
23 #include <vector>
24 #include <cmath>
25 
26 // CMS
31 
32 // ROOT
33 #include "TMatrixD.h"
34 #include "TMatrixDSym.h"
35 #include "TDecompBK.h"
36 #include "TH1.h"
37 #include "TF1.h"
38 
39 using namespace ROOT::Minuit2;
40 
41 
42 //_____________________________________________________________________
44  fbeamtype = reco::BeamSpot::Unknown;
45 }
46 
47 //_____________________________________________________________________
48 BSFitter::BSFitter( const std:: vector< BSTrkParameters > &BSvector ) {
49 
50  ffit_type = "default";
51  ffit_variable = "default";
52 
53  fBSvector = BSvector;
54 
55  fsqrt2pi = sqrt(2.* TMath::Pi());
56 
57  fpar_name[0] = "z0 ";
58  fpar_name[1] = "SigmaZ0 ";
59  fpar_name[2] = "X0 ";
60  fpar_name[3] = "Y0 ";
61  fpar_name[4] = "dxdz ";
62  fpar_name[5] = "dydz ";
63  fpar_name[6] = "SigmaBeam ";
64 
65  //if (theGausszFcn == 0 ) {
66  thePDF = new BSpdfsFcn();
67 
68 
69 //}
70  //if (theFitter == 0 ) {
71 
72  theFitter = new VariableMetricMinimizer();
73 
74  //}
75 
76  fapplyd0cut = false;
77  fapplychi2cut = false;
78  ftmprow = 0;
79  ftmp.ResizeTo(4,1);
80  ftmp.Zero();
81  fnthite=0;
82  fMaxZ = 50.; //cm
83  fconvergence = 0.5; // stop fit when 50% of the input collection has been removed.
84  fminNtrks = 100;
85  finputBeamWidth = -1; // no input
86 
87  h1z = new TH1F("h1z","z distribution",200,-fMaxZ, fMaxZ);
88 
89 }
90 
91 //______________________________________________________________________
93 {
94  //delete fBSvector;
95  delete thePDF;
96  delete theFitter;
97 }
98 
99 
100 //______________________________________________________________________
102 
103  return this->Fit(0);
104 
105 }
106 
107 //______________________________________________________________________
108 reco::BeamSpot BSFitter::Fit(double *inipar = 0) {
109  fbeamtype = reco::BeamSpot::Unknown;
110  if ( ffit_variable == "z" ) {
111 
112  if ( ffit_type == "chi2" ) {
113 
114  return Fit_z_chi2(inipar);
115 
116  } else if ( ffit_type == "likelihood" ) {
117 
118  return Fit_z_likelihood(inipar);
119 
120  } else if ( ffit_type == "combined" ) {
121 
122  reco::BeamSpot tmp_beamspot = Fit_z_chi2(inipar);
123  double tmp_par[2] = {tmp_beamspot.z0(), tmp_beamspot.sigmaZ()};
124  return Fit_z_likelihood(tmp_par);
125 
126  } else {
127 
128  throw cms::Exception("LogicError")
129  << "Error in BeamSpotProducer/BSFitter: "
130  << "Illegal fit type, options are chi2,likelihood or combined(ie. first chi2 then likelihood)";
131 
132  }
133  } else if ( ffit_variable == "d" ) {
134 
135  if ( ffit_type == "d0phi" ) {
136  this->d0phi_Init();
137  return Fit_d0phi();
138 
139  } else if ( ffit_type == "likelihood" ) {
140 
141  return Fit_d_likelihood(inipar);
142 
143  } else if ( ffit_type == "combined" ) {
144 
145  this->d0phi_Init();
146  reco::BeamSpot tmp_beamspot = Fit_d0phi();
147  double tmp_par[4] = {tmp_beamspot.x0(), tmp_beamspot.y0(), tmp_beamspot.dxdz(), tmp_beamspot.dydz()};
148  return Fit_d_likelihood(tmp_par);
149 
150  } else {
151  throw cms::Exception("LogicError")
152  << "Error in BeamSpotProducer/BSFitter: "
153  << "Illegal fit type, options are d0phi, likelihood or combined";
154  }
155  } else if ( ffit_variable == "d*z" || ffit_variable == "default" ) {
156 
157  if ( ffit_type == "likelihood" || ffit_type == "default" ) {
158 
160  // we are now fitting Z inside d0phi fitter
161  // first fit z distribution using a chi2 fit
162  //reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
163  //for (int j = 2 ; j < 4 ; ++j) {
164  //for(int k = j ; k < 4 ; ++k) {
165  // matrix(j,k) = tmp_z.covariance()(j,k);
166  //}
167  //}
168 
169  // use d0-phi algorithm to extract transverse position
170  this->d0phi_Init();
171  //reco::BeamSpot tmp_d0phi= Fit_d0phi(); // change to iterative procedure:
172  this->Setd0Cut_d0phi(4.0);
173  reco::BeamSpot tmp_d0phi= Fit_ited0phi();
174 
175  //for (int j = 0 ; j < 2 ; ++j) {
176  // for(int k = j ; k < 2 ; ++k) {
177  // matrix(j,k) = tmp_d0phi.covariance()(j,k);
178  //}
179  //}
180  // slopes
181  //for (int j = 4 ; j < 6 ; ++j) {
182  // for(int k = j ; k < 6 ; ++k) {
183  // matrix(j,k) = tmp_d0phi.covariance()(j,k);
184  // }
185  //}
186 
187 
188  // put everything into one object
189  reco::BeamSpot spot(reco::BeamSpot::Point(tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_d0phi.z0()),
190  tmp_d0phi.sigmaZ(),
191  tmp_d0phi.dxdz(),
192  tmp_d0phi.dydz(),
193  0.,
194  tmp_d0phi.covariance(),
195  fbeamtype );
196 
197 
198 
199  //reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
200 
201  //reco::BeamSpot tmp_d0phi = Fit_d0phi();
202 
203  // log-likelihood fit
204  if (ffit_type == "likelihood") {
205  double tmp_par[7] = {tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_d0phi.z0(),
206  tmp_d0phi.sigmaZ(), tmp_d0phi.dxdz(), tmp_d0phi.dydz(),0.0};
207 
208  double tmp_error_par[7];
209  for(int s=0;s<6;s++){ tmp_error_par[s] = pow( tmp_d0phi.covariance()(s,s),0.5);}
210  tmp_error_par[6]=0.0;
211 
212  reco::BeamSpot tmp_lh = Fit_d_z_likelihood(tmp_par,tmp_error_par);
213 
214  if (edm::isNotFinite(ff_minimum)) {
215  edm::LogWarning("BSFitter") << "BSFitter: Result is non physical. Log-Likelihood fit to extract beam width did not converge." << std::endl;
217  return tmp_lh;
218  }
219  return tmp_lh;
220 
221  } else {
222 
223  edm::LogInfo("BSFitter") << "default track-based fit does not extract beam width." << std::endl;
224  return spot;
225  }
226 
227 
228  } else if ( ffit_type == "resolution" ) {
229 
230  reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
231  this->d0phi_Init();
232  reco::BeamSpot tmp_d0phi = Fit_d0phi();
233 
234  double tmp_par[7] = {tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_z.z0(),
235  tmp_z.sigmaZ(), tmp_d0phi.dxdz(), tmp_d0phi.dydz(),0.0};
236  double tmp_error_par[7];
237  for(int s=0;s<6;s++){ tmp_error_par[s] = pow(tmp_par[s],0.5);}
238  tmp_error_par[6]=0.0;
239 
240  reco::BeamSpot tmp_beam = Fit_d_z_likelihood(tmp_par,tmp_error_par);
241 
242  double tmp_par2[7] = {tmp_beam.x0(), tmp_beam.y0(), tmp_beam.z0(),
243  tmp_beam.sigmaZ(), tmp_beam.dxdz(), tmp_beam.dydz(),
244  tmp_beam.BeamWidthX()};
245 
246  reco::BeamSpot tmp_lh = Fit_dres_z_likelihood(tmp_par2);
247 
248  if (edm::isNotFinite(ff_minimum)) {
249 
250  edm::LogWarning("BSFitter") << "Result is non physical. Log-Likelihood fit did not converge." << std::endl;
252  return tmp_lh;
253  }
254  return tmp_lh;
255 
256  } else {
257 
258  throw cms::Exception("LogicError")
259  << "Error in BeamSpotProducer/BSFitter: "
260  << "Illegal fit type, options are likelihood or resolution";
261  }
262  } else {
263 
264  throw cms::Exception("LogicError")
265  << "Error in BeamSpotProducer/BSFitter: "
266  << "Illegal variable type, options are \"z\", \"d\", or \"d*z\"";
267  }
268 
269 
270 }
271 
272 //______________________________________________________________________
274 
275  //std::cout << "Fit_z(double *) called" << std::endl;
276  //std::cout << "inipar[0]= " << inipar[0] << std::endl;
277  //std::cout << "inipar[1]= " << inipar[1] << std::endl;
278 
279  std::vector<double> par(2,0);
280  std::vector<double> err(2,0);
281 
282  par.push_back(0.0);
283  par.push_back(7.0);
284  err.push_back(0.0001);
285  err.push_back(0.0001);
286  //par[0] = 0.0; err[0] = 0.0;
287  //par[1] = 7.0; err[1] = 0.0;
288 
289  thePDF->SetPDFs("PDFGauss_z");
290  thePDF->SetData(fBSvector);
291  //std::cout << "data loaded"<< std::endl;
292 
293  //FunctionMinimum fmin = theFitter->Minimize(*theGausszFcn, par, err, 1, 500, 0.1);
294  MnUserParameters upar;
295  upar.Add("X0", 0.,0.);
296  upar.Add("Y0", 0.,0.);
297  upar.Add("Z0", inipar[0],0.001);
298  upar.Add("sigmaZ",inipar[1],0.001);
299 
300  MnMigrad migrad(*thePDF, upar);
301 
302  FunctionMinimum fmin = migrad();
303  ff_minimum = fmin.Fval();
304  //std::cout << " eval= " << ff_minimum
305  // << "/n params[0]= " << fmin.Parameters().Vec()(0) << std::endl;
306 
307  /*
308  TMinuit *gmMinuit = new TMinuit(2);
309 
310  //gmMinuit->SetFCN(z_fcn);
311  gmMinuit->SetFCN(myFitz_fcn);
312 
313 
314  int ierflg = 0;
315  double step[2] = {0.001,0.001};
316 
317  for (int i = 0; i<2; i++) {
318  gmMinuit->mnparm(i,fpar_name[i].c_str(),inipar[i],step[i],0,0,ierflg);
319  }
320  gmMinuit->Migrad();
321  */
323 
324  for (int j = 2 ; j < 4 ; ++j) {
325  for(int k = j ; k < 4 ; ++k) {
326  matrix(j,k) = fmin.Error().Matrix()(j,k);
327  }
328  }
329 
331  0.,
332  fmin.Parameters().Vec()(2)),
333  fmin.Parameters().Vec()(3),
334  0.,
335  0.,
336  0.,
337  matrix,
338  fbeamtype );
339 }
340 
341 //______________________________________________________________________
343 
344  // N.B. this fit is not performed anymore but now
345  // Z is fitted in the same track set used in the d0-phi fit after
346  // each iteration
347 
348 
349  //std::cout << "Fit_z_chi2() called" << std::endl;
350  // FIXME: include whole tracker z length for the time being
351  // ==> add protection and z0 cut
352  h1z = new TH1F("h1z","z distribution",200,-fMaxZ, fMaxZ);
353 
354  std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
355 
356  // HERE check size of track vector
357 
358  for( iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
359 
360  h1z->Fill( iparam->z0() );
361  //std::cout<<"z0="<<iparam->z0()<<"; sigZ0="<<iparam->sigz0()<<std::endl;
362  }
363 
364  h1z->Fit("gaus","QLM0");
365  //std::cout << "fitted "<< std::endl;
366 
367  TF1 *fgaus = h1z->GetFunction("gaus");
368  //std::cout << "got function" << std::endl;
369  double fpar[2] = {fgaus->GetParameter(1), fgaus->GetParameter(2) };
370  //std::cout<<"Debug fpar[2] = (" <<fpar[0]<<","<<fpar[1]<<")"<<std::endl;
372  // add matrix values.
373  matrix(2,2) = fgaus->GetParError(1) * fgaus->GetParError(1);
374  matrix(3,3) = fgaus->GetParError(2) * fgaus->GetParError(2);
375 
376  //delete h1z;
377 
379  0.,
380  fpar[0]),
381  fpar[1],
382  0.,
383  0.,
384  0.,
385  matrix,
386  fbeamtype );
387 
388 
389 }
390 
391 //______________________________________________________________________
393 
394  this->d0phi_Init();
395  edm::LogInfo("BSFitter") << "number of total input tracks: " << fBSvector.size() << std::endl;
396 
397  reco::BeamSpot theanswer;
398 
399  if ( (int)fBSvector.size() <= fminNtrks ) {
400  edm::LogWarning("BSFitter") << "need at least " << fminNtrks << " tracks to run beamline fitter." << std::endl;
401  fbeamtype = reco::BeamSpot::Fake;
402  theanswer.setType(fbeamtype);
403  return theanswer;
404  }
405 
406  theanswer = Fit_d0phi(); //get initial ftmp and ftmprow
407  if ( goodfit ) fnthite++;
408  //std::cout << "Initial tempanswer (iteration 0): " << theanswer << std::endl;
409 
410  reco::BeamSpot preanswer = theanswer;
411 
412  while ( goodfit &&
413  ftmprow > fconvergence * fBSvector.size() &&
414  ftmprow > fminNtrks ) {
415 
416  theanswer = Fit_d0phi();
417  fd0cut /= 1.5;
418  fchi2cut /= 1.5;
419  if ( goodfit &&
420  ftmprow > fconvergence * fBSvector.size() &&
421  ftmprow > fminNtrks ) {
422  preanswer = theanswer;
423  //std::cout << "Iteration " << fnthite << ": " << preanswer << std::endl;
424  fnthite++;
425  }
426  }
427  // FIXME: return fit results from previous iteration for both bad fit and for >50% tracks thrown away
428  //std::cout << "The last iteration, theanswer: " << theanswer << std::endl;
429  theanswer = preanswer;
430  //std::cout << "Use previous results from iteration #" << ( fnthite > 0 ? fnthite-1 : 0 ) << std::endl;
431  //if ( fnthite > 1 ) std::cout << theanswer << std::endl;
432 
433  edm::LogInfo("BSFitter") << "Total number of successful iterations = " << ( goodfit ? (fnthite+1) : fnthite ) << std::endl;
434  if (goodfit) {
435  fbeamtype = reco::BeamSpot::Tracker;
436  theanswer.setType(fbeamtype);
437  }
438  else {
439  edm::LogWarning("BSFitter") << "Fit doesn't converge!!!" << std::endl;
440  fbeamtype = reco::BeamSpot::Unknown;
441  theanswer.setType(fbeamtype);
442  }
443  return theanswer;
444 }
445 
446 
447 //______________________________________________________________________
449 
450  //LogDebug ("BSFitter") << " we will use " << fBSvector.size() << " tracks.";
451  if (fnthite > 0) edm::LogInfo("BSFitter") << " number of tracks used: " << ftmprow << std::endl;
452  //std::cout << " ftmp = matrix("<<ftmp.GetNrows()<<","<<ftmp.GetNcols()<<")"<<std::endl;
453  //std::cout << " ftmp(0,0)="<<ftmp(0,0)<<std::endl;
454  //std::cout << " ftmp(1,0)="<<ftmp(1,0)<<std::endl;
455  //std::cout << " ftmp(2,0)="<<ftmp(2,0)<<std::endl;
456  //std::cout << " ftmp(3,0)="<<ftmp(3,0)<<std::endl;
457 
458  h1z->Reset();
459 
460 
461  TMatrixD x_result(4,1);
462  TMatrixDSym V_result(4);
463 
464  TMatrixDSym Vint(4);
465  TMatrixD b(4,1);
466 
467  //Double_t weightsum = 0;
468 
469  Vint.Zero();
470  b.Zero();
471 
472  TMatrixD g(4,1);
473  TMatrixDSym temp(4);
474 
475  std::vector<BSTrkParameters>::iterator iparam = fBSvector.begin();
476  ftmprow=0;
477 
478 
479  //edm::LogInfo ("BSFitter") << " test";
480 
481  //std::cout << "BSFitter: fit" << std::endl;
482 
483  for( iparam = fBSvector.begin() ;
484  iparam != fBSvector.end() ; ++iparam) {
485 
486 
487  //if(i->weight2 == 0) continue;
488 
489  //if (ftmprow==0) {
490  //std::cout << "d0=" << iparam->d0() << " sigd0=" << iparam->sigd0()
491  //<< " phi0="<< iparam->phi0() << " z0=" << iparam->z0() << std::endl;
492  //std::cout << "d0phi_d0=" << iparam->d0phi_d0() << " d0phi_chi2="<<iparam->d0phi_chi2() << std::endl;
493  //}
494  g(0,0) = sin(iparam->phi0());
495  g(1,0) = -cos(iparam->phi0());
496  g(2,0) = iparam->z0() * g(0,0);
497  g(3,0) = iparam->z0() * g(1,0);
498 
499 
500  // average transverse beam width
501  double sigmabeam2 = 0.006 * 0.006;
502  if (finputBeamWidth > 0 ) sigmabeam2 = finputBeamWidth * finputBeamWidth;
503  else { edm::LogWarning("BSFitter") << "using in fit beam width = " << sqrt(sigmabeam2) << std::endl; }
504 
505  //double sigma2 = sigmabeam2 + (iparam->sigd0())* (iparam->sigd0()) / iparam->weight2;
506  // this should be 2*sigmabeam2?
507  double sigma2 = sigmabeam2 + (iparam->sigd0())* (iparam->sigd0());
508 
509  TMatrixD ftmptrans(1,4);
510  ftmptrans = ftmptrans.Transpose(ftmp);
511  TMatrixD dcor = ftmptrans * g;
512  double chi2tmp = (iparam->d0() - dcor(0,0)) * (iparam->d0() - dcor(0,0))/sigma2;
513  (*iparam) = BSTrkParameters(iparam->z0(),iparam->sigz0(),iparam->d0(),iparam->sigd0(),
514  iparam->phi0(), iparam->pt(),dcor(0,0),chi2tmp);
515 
516  bool pass = true;
517  if (fapplyd0cut && fnthite>0 ) {
518  if ( std::abs(iparam->d0() - dcor(0,0)) > fd0cut ) pass = false;
519 
520  }
521  if (fapplychi2cut && fnthite>0 ) {
522  if ( chi2tmp > fchi2cut ) pass = false;
523 
524  }
525 
526  if (pass) {
527  temp.Zero();
528  for(int j = 0 ; j < 4 ; ++j) {
529  for(int k = j ; k < 4 ; ++k) {
530  temp(j,k) += g(j,0) * g(k,0);
531  }
532  }
533 
534 
535  Vint += (temp * (1 / sigma2));
536  b += (iparam->d0() / sigma2 * g);
537  //weightsum += sqrt(i->weight2);
538  ftmprow++;
539  h1z->Fill( iparam->z0() );
540  }
541 
542 
543  }
544  Double_t determinant;
545  TDecompBK bk(Vint);
546  bk.SetTol(1e-11); //FIXME: find a better way to solve x_result
547  if (!bk.Decompose()) {
548  goodfit = false;
549  edm::LogWarning("BSFitter")
550  << "Decomposition failed, matrix singular ?" << std::endl
551  << "condition number = " << bk.Condition() << std::endl;
552  }
553  else {
554  V_result = Vint.InvertFast(&determinant);
555  x_result = V_result * b;
556  }
557  // for(int j = 0 ; j < 4 ; ++j) {
558  // for(int k = 0 ; k < 4 ; ++k) {
559  // std::cout<<"V_result("<<j<<","<<k<<")="<<V_result(j,k)<<std::endl;
560  // }
561  // }
562  // for (int j=0;j<4;++j){
563  // std::cout<<"x_result("<<j<<",0)="<<x_result(j,0)<<std::endl;
564  // }
565  //LogDebug ("BSFitter") << " d0-phi fit done.";
566  //std::cout<< " d0-phi fit done." << std::endl;
567 
568  h1z->Fit("gaus","QLM0","",h1z->GetMean() -2.*h1z->GetRMS(),h1z->GetMean() +2.*h1z->GetRMS());
569 
570  //std::cout << "fitted "<< std::endl;
571  TF1 *fgaus = h1z->GetFunction("gaus");
572  //std::cout << "got function" << std::endl;
573  if (!fgaus){
574  edm::LogError("NoBeamSpotFit")<<"gaussian fit failed. no BS d0 fit";
575  return reco::BeamSpot();
576  }
577  double fpar[2] = {fgaus->GetParameter(1), fgaus->GetParameter(2) };
578 
580  // first two parameters
581  for (int j = 0 ; j < 2 ; ++j) {
582  for(int k = j ; k < 2 ; ++k) {
583  matrix(j,k) = V_result(j,k);
584  }
585  }
586  // slope parameters
587  for (int j = 4 ; j < 6 ; ++j) {
588  for(int k = j ; k < 6 ; ++k) {
589  matrix(j,k) = V_result(j-2,k-2);
590  }
591  }
592 
593  // Z0 and sigmaZ
594  matrix(2,2) = fgaus->GetParError(1) * fgaus->GetParError(1);
595  matrix(3,3) = fgaus->GetParError(2) * fgaus->GetParError(2);
596 
597  ftmp = x_result;
598 
599  // x0 and y0 are *not* x,y at z=0, but actually at z=0
600  // to correct for this, we need to translate them to z=z0
601  // using the measured slopes
602  //
603  double x0tmp = x_result(0,0);
604  double y0tmp = x_result(1,0);
605 
606  x0tmp += x_result(2,0)*fpar[0];
607  y0tmp += x_result(3,0)*fpar[0];
608 
609 
611  y0tmp,
612  fpar[0]),
613  fpar[1],
614  x_result(2,0),
615  x_result(3,0),
616  0.,
617  matrix,
618  fbeamtype );
619 
620 }
621 
622 
623 //______________________________________________________________________
624 void BSFitter::Setd0Cut_d0phi(double d0cut) {
625 
626  fapplyd0cut = true;
627 
628  //fBSforCuts = BSfitted;
629  fd0cut = d0cut;
630 }
631 
632 //______________________________________________________________________
633 void BSFitter::SetChi2Cut_d0phi(double chi2cut) {
634 
635  fapplychi2cut = true;
636 
637  //fBSforCuts = BSfitted;
638  fchi2cut = chi2cut;
639 }
640 
641 //______________________________________________________________________
643 
644 
645  thePDF->SetPDFs("PDFGauss_d");
646  thePDF->SetData(fBSvector);
647 
648  MnUserParameters upar;
649  upar.Add("X0", inipar[0],0.001);
650  upar.Add("Y0", inipar[1],0.001);
651  upar.Add("Z0", 0.,0.001);
652  upar.Add("sigmaZ",0.,0.001);
653  upar.Add("dxdz",inipar[2],0.001);
654  upar.Add("dydz",inipar[3],0.001);
655 
656 
657  MnMigrad migrad(*thePDF, upar);
658 
659  FunctionMinimum fmin = migrad();
660  ff_minimum = fmin.Fval();
661 
663  for (int j = 0 ; j < 6 ; ++j) {
664  for(int k = j ; k < 6 ; ++k) {
665  matrix(j,k) = fmin.Error().Matrix()(j,k);
666  }
667  }
668 
669  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
670  fmin.Parameters().Vec()(1),
671  0.),
672  0.,
673  fmin.Parameters().Vec()(4),
674  fmin.Parameters().Vec()(5),
675  0.,
676  matrix,
677  fbeamtype );
678 }
679 //______________________________________________________________________
680 double BSFitter::scanPDF(double *init_pars, int & tracksfixed, int option){
681 
682  if(option==1)init_pars[6]=0.0005; //starting value for any given configuration
683 
684  //local vairables with initial values
685  double fsqrt2pi=0.0;
686  double d_sig=0.0;
687  double d_dprime=0.0;
688  double d_result=0.0;
689  double z_sig=0.0;
690  double z_result=0.0;
691  double function=0.0;
692  double tot_pdf=0.0;
693  double last_minvalue=1.0e+10;
694  double init_bw=-99.99;
695  int iters=0;
696 
697  //used to remove tracks if far away from bs by this
698  double DeltadCut=0.1000;
699  if(init_pars[6]<0.0200){DeltadCut=0.0900; } //worked for high 2.36TeV
700  if(init_pars[6]<0.0100){DeltadCut=0.0700;} //just a guesss for 7 TeV but one should scan for actual values
701 
702 
703 std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
704 
705 
706 if(option==1)iters=500;
707 if(option==2)iters=1;
708 
709 for(int p=0;p<iters;p++){
710 
711  if(iters==500)init_pars[6]+=0.0002;
712  tracksfixed=0;
713 
714 for( iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam)
715  {
716  fsqrt2pi = sqrt(2.* TMath::Pi());
717  d_sig = sqrt(init_pars[6]*init_pars[6] + (iparam->sigd0())*(iparam->sigd0()));
718  d_dprime = iparam->d0() - ( ( (init_pars[0] + iparam->z0()*(init_pars[4]))*sin(iparam->phi0()) )
719  - ( (init_pars[1] + iparam->z0()*(init_pars[5]))*cos(iparam->phi0()) ) );
720 
721  //***Remove tracks before the fit which gives low pdf values to blow up the pdf
722  if(std::abs(d_dprime)<DeltadCut && option==2){ fBSvectorBW.push_back(*iparam);}
723 
724  d_result = (exp(-(d_dprime*d_dprime)/(2.0*d_sig*d_sig)))/(d_sig*fsqrt2pi);
725  z_sig = sqrt(iparam->sigz0() * iparam->sigz0() + init_pars[3]*init_pars[3]);
726  z_result = (exp(-((iparam->z0() - init_pars[2])*(iparam->z0() - init_pars[2]))/(2.0*z_sig*z_sig)))/(z_sig*fsqrt2pi);
727  tot_pdf=z_result*d_result;
728 
729  //for those trcks which gives problems due to very tiny pdf_d values.
730  //Update: This protection will NOT be used with the dprime cut above but still kept here to get
731  // the intial value of beam width reasonably
732  //A warning will appear if there were any tracks with < 10^-5 for pdf_d so that (d-dprime) cut can be lowered
733  if(d_result < 1.0e-05){ tot_pdf=z_result*1.0e-05;
734  //if(option==2)std::cout<<"last Iter d-d' = "<<(std::abs(d_dprime))<<std::endl;
735  tracksfixed++; }
736 
737  function = function + log(tot_pdf);
738  tot_pdf=0.0;
739 
740 
741  }//loop over tracks
742 
743 
744  function= -2.0*function;
745  if(function<last_minvalue){init_bw=init_pars[6];
746  last_minvalue=function; }
747  function=0.0;
748  }//loop over beam width
749 
750  if(init_bw>0) {
751  init_bw=init_bw+(0.20*init_bw); //start with 20 % more
752 
753  }
754  else{
755 
756  if(option==1){
757  edm::LogWarning("BSFitter")
758  <<"scanPDF:====>>>> WARNING***: The initial guess value of Beam width is negative!!!!!!"<<std::endl
759  <<"scanPDF:====>>>> Assigning beam width a starting value of "<<init_bw<<" cm"<<std::endl;
760  init_bw=0.0200;
761 
762  }
763  }
764 
765 
766  return init_bw;
767 
768 }
769 
770 //________________________________________________________________________________
771 reco::BeamSpot BSFitter::Fit_d_z_likelihood(double *inipar, double *error_par) {
772 
773  int tracksFailed=0;
774 
775  //estimate first guess of beam width and tame 20% extra of it to start
776  inipar[6]=scanPDF(inipar,tracksFailed,1);
777  error_par[6]=(inipar[6])*0.20;
778 
779 
780  //Here remove the tracks which give low pdf and fill into a new vector
781  //std::cout<<"Size of Old vector = "<<(fBSvector.size())<<std::endl;
782  /* double junk= */ scanPDF(inipar,tracksFailed,2);
783  //std::cout<<"Size of New vector = "<<(fBSvectorBW.size())<<std::endl;
784 
785  //Refill the fBSVector again with new sets of tracks
786  fBSvector.clear();
787  std::vector<BSTrkParameters>::const_iterator iparamBW = fBSvectorBW.begin();
788  for( iparamBW = fBSvectorBW.begin(); iparamBW != fBSvectorBW.end(); ++iparamBW)
789  { fBSvector.push_back(*iparamBW);
790  }
791 
792 
793  thePDF->SetPDFs("PDFGauss_d*PDFGauss_z");
794  thePDF->SetData(fBSvector);
795  MnUserParameters upar;
796 
797  upar.Add("X0", inipar[0],error_par[0]);
798  upar.Add("Y0", inipar[1],error_par[1]);
799  upar.Add("Z0", inipar[2],error_par[2]);
800  upar.Add("sigmaZ",inipar[3],error_par[3]);
801  upar.Add("dxdz",inipar[4],error_par[4]);
802  upar.Add("dydz",inipar[5],error_par[5]);
803  upar.Add("BeamWidthX",inipar[6],error_par[6]);
804 
805 
806  MnMigrad migrad(*thePDF, upar);
807 
808  FunctionMinimum fmin = migrad();
809 
810  // std::cout<<"-----how the fit evoves------"<<std::endl;
811  // std::cout<<fmin<<std::endl;
812 
813  ff_minimum = fmin.Fval();
814 
815 
816  bool ff_nfcn=fmin.HasReachedCallLimit();
817  bool ff_cov=fmin.HasCovariance();
818  bool testing=fmin.IsValid();
819 
820 
821  //Print WARNINGS if minimum did not converged
822  if( ! testing )
823  {
824  edm::LogWarning("BSFitter") <<"===========>>>>>** WARNING: MINUIT DID NOT CONVERGES PROPERLY !!!!!!"<<std::endl;
825  if(ff_nfcn) edm::LogWarning("BSFitter") <<"===========>>>>>** WARNING: No. of Calls Exhausted"<<std::endl;
826  if(!ff_cov) edm::LogWarning("BSFitter") <<"===========>>>>>** WARNING: Covariance did not found"<<std::endl;
827  }
828 
829  edm::LogInfo("BSFitter") <<"The Total # Tracks used for beam width fit = "<<(fBSvectorBW.size())<<std::endl;
830 
831 
832  //Checks after fit is performed
833  double lastIter_pars[7];
834 
835  for(int ip=0;ip<7;ip++){ lastIter_pars[ip]=fmin.Parameters().Vec()(ip);
836  }
837 
838 
839 
840  tracksFailed=0;
841  /* double lastIter_scan= */ scanPDF(lastIter_pars,tracksFailed,2);
842 
843 
844  edm::LogWarning("BSFitter") <<"WARNING: # of tracks which have very low pdf value (pdf_d < 1.0e-05) are = "<<tracksFailed<<std::endl;
845 
846 
847 
848  //std::cout << " eval= " << ff_minimum
849  // << "/n params[0]= " << fmin.Parameters().Vec()(0) << std::endl;
850 
852 
853  for (int j = 0 ; j < 7 ; ++j) {
854  for(int k = j ; k < 7 ; ++k) {
855  matrix(j,k) = fmin.Error().Matrix()(j,k);
856  }
857  }
858 
859 
860  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
861  fmin.Parameters().Vec()(1),
862  fmin.Parameters().Vec()(2)),
863  fmin.Parameters().Vec()(3),
864  fmin.Parameters().Vec()(4),
865  fmin.Parameters().Vec()(5),
866  fmin.Parameters().Vec()(6),
867 
868  matrix,
869  fbeamtype );
870 }
871 
872 
873 //______________________________________________________________________
875 
876 
877  thePDF->SetPDFs("PDFGauss_d_resolution*PDFGauss_z");
878  thePDF->SetData(fBSvector);
879 
880  MnUserParameters upar;
881  upar.Add("X0", inipar[0],0.001);
882  upar.Add("Y0", inipar[1],0.001);
883  upar.Add("Z0", inipar[2],0.001);
884  upar.Add("sigmaZ",inipar[3],0.001);
885  upar.Add("dxdz",inipar[4],0.001);
886  upar.Add("dydz",inipar[5],0.001);
887  upar.Add("BeamWidthX",inipar[6],0.0001);
888  upar.Add("c0",0.0010,0.0001);
889  upar.Add("c1",0.0090,0.0001);
890 
891  // fix beam width
892  upar.Fix("BeamWidthX");
893  // number of parameters in fit are 9-1 = 8
894 
895  MnMigrad migrad(*thePDF, upar);
896 
897  FunctionMinimum fmin = migrad();
898  ff_minimum = fmin.Fval();
899 
901 
902  for (int j = 0 ; j < 6 ; ++j) {
903  for(int k = j ; k < 6 ; ++k) {
904  matrix(j,k) = fmin.Error().Matrix()(j,k);
905  }
906  }
907 
908  //std::cout << " fill resolution values" << std::endl;
909  //std::cout << " matrix size= " << fmin.Error().Matrix().size() << std::endl;
910  //std::cout << " vec(6)="<< fmin.Parameters().Vec()(6) << std::endl;
911  //std::cout << " vec(7)="<< fmin.Parameters().Vec()(7) << std::endl;
912 
913  fresolution_c0 = fmin.Parameters().Vec()(6);
914  fresolution_c1 = fmin.Parameters().Vec()(7);
915  fres_c0_err = sqrt( fmin.Error().Matrix()(6,6) );
916  fres_c1_err = sqrt( fmin.Error().Matrix()(7,7) );
917 
918  for (int j = 6 ; j < 8 ; ++j) {
919  for(int k = 6 ; k < 8 ; ++k) {
920  fres_matrix(j-6,k-6) = fmin.Error().Matrix()(j,k);
921  }
922  }
923 
924  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
925  fmin.Parameters().Vec()(1),
926  fmin.Parameters().Vec()(2)),
927  fmin.Parameters().Vec()(3),
928  fmin.Parameters().Vec()(4),
929  fmin.Parameters().Vec()(5),
930  inipar[6],
931  matrix,
932  fbeamtype );
933 }
934 
935 
936 
reco::BeamSpot Fit_d_likelihood(double *inipar)
Definition: BSFitter.cc:642
const double Pi
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:32
double z0() const
z coordinate
Definition: BeamSpot.h:69
virtual ~BSFitter()
Definition: BSFitter.cc:92
reco::BeamSpot Fit()
Definition: BSFitter.cc:101
reco::BeamSpot Fit_z_likelihood(double *inipar)
Definition: BSFitter.cc:273
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define abs(x)
Definition: mlp_lapack.h:159
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:30
double scanPDF(double *init_pars, int &tracksFailed, int option)
Definition: BSFitter.cc:680
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
reco::BeamSpot Fit_dres_z_likelihood(double *inipar)
Definition: BSFitter.cc:874
void setType(BeamType type)
set beam type
Definition: BeamSpot.h:132
double dydz() const
dydz slope
Definition: BeamSpot.h:85
bool isNotFinite(T x)
Definition: isFinite.h:10
Definition: Fit.h:34
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
BSFitter()
Definition: BSFitter.cc:43
int j
Definition: DBlmapReader.cc:9
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:87
double dxdz() const
dxdz slope
Definition: BeamSpot.h:83
int k[5][pyjets_maxn]
reco::BeamSpot Fit_ited0phi()
Definition: BSFitter.cc:392
double sigmaZ() const
sigma z
Definition: BeamSpot.h:81
double b
Definition: hdecay.h:120
reco::BeamSpot Fit_z_chi2(double *inipar)
Definition: BSFitter.cc:342
void Setd0Cut_d0phi(double d0cut)
Definition: BSFitter.cc:624
double covariance(int i, int j) const
(i,j)-th element of error matrix
Definition: BeamSpot.h:113
double y0() const
y coordinate
Definition: BeamSpot.h:67
void SetChi2Cut_d0phi(double chi2cut)
Definition: BSFitter.cc:633
reco::BeamSpot Fit_d0phi()
Definition: BSFitter.cc:448
reco::BeamSpot Fit_d_z_likelihood(double *inipar, double *error_par)
Definition: BSFitter.cc:771
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double x0() const
x coordinate
Definition: BeamSpot.h:65