13 #include "Minuit2/FunctionMinimum.h" 14 #include "Minuit2/MnPrint.h" 15 #include "Minuit2/MnMigrad.h" 16 #include "Minuit2/MnUserParameterState.h" 31 #include "TMatrixDSym.h" 32 #include "TDecompBK.h" 35 #include "TMinuitMinimizer.h" 45 TMinuitMinimizer::UseStaticMinuit(
false);
52 TMinuitMinimizer::UseStaticMinuit(
false);
54 ffit_type =
"default";
55 ffit_variable =
"default";
62 fpar_name[1] =
"SigmaZ0 ";
65 fpar_name[4] =
"dxdz ";
66 fpar_name[5] =
"dydz ";
67 fpar_name[6] =
"SigmaBeam ";
75 fapplychi2cut =
false;
85 h1z =
new TH1F(
"h1z",
"z distribution", 200, -fMaxZ, fMaxZ);
100 if (ffit_variable ==
"z") {
101 if (ffit_type ==
"chi2") {
102 return Fit_z_chi2(inipar);
104 }
else if (ffit_type ==
"likelihood") {
105 return Fit_z_likelihood(inipar);
107 }
else if (ffit_type ==
"combined") {
109 double tmp_par[2] = {tmp_beamspot.
z0(), tmp_beamspot.
sigmaZ()};
110 return Fit_z_likelihood(tmp_par);
114 <<
"Error in BeamSpotProducer/BSFitter: " 115 <<
"Illegal fit type, options are chi2,likelihood or combined(ie. first chi2 then likelihood)";
117 }
else if (ffit_variable ==
"d") {
118 if (ffit_type ==
"d0phi") {
122 }
else if (ffit_type ==
"likelihood") {
123 return Fit_d_likelihood(inipar);
125 }
else if (ffit_type ==
"combined") {
128 double tmp_par[4] = {tmp_beamspot.
x0(), tmp_beamspot.
y0(), tmp_beamspot.
dxdz(), tmp_beamspot.
dydz()};
129 return Fit_d_likelihood(tmp_par);
132 throw cms::Exception(
"LogicError") <<
"Error in BeamSpotProducer/BSFitter: " 133 <<
"Illegal fit type, options are d0phi, likelihood or combined";
135 }
else if (ffit_variable ==
"d*z" || ffit_variable ==
"default") {
136 if (ffit_type ==
"likelihood" || ffit_type ==
"default") {
150 this->Setd0Cut_d0phi(4.0);
179 if (ffit_type ==
"likelihood") {
180 double tmp_par[7] = {
181 tmp_d0phi.
x0(), tmp_d0phi.
y0(), tmp_d0phi.
z0(), tmp_d0phi.
sigmaZ(), tmp_d0phi.
dxdz(), tmp_d0phi.
dydz(), 0.0};
183 double tmp_error_par[7];
184 for (
int s = 0;
s < 6;
s++) {
187 tmp_error_par[6] = 0.0;
189 reco::BeamSpot tmp_lh = Fit_d_z_likelihood(tmp_par, tmp_error_par);
193 <<
"BSFitter: Result is non physical. Log-Likelihood fit to extract beam width did not converge." 201 edm::LogInfo(
"BSFitter") <<
"default track-based fit does not extract beam width." << std::endl;
205 }
else if (ffit_type ==
"resolution") {
210 double tmp_par[7] = {
211 tmp_d0phi.
x0(), tmp_d0phi.
y0(), tmp_z.
z0(), tmp_z.
sigmaZ(), tmp_d0phi.
dxdz(), tmp_d0phi.
dydz(), 0.0};
212 double tmp_error_par[7];
213 for (
int s = 0;
s < 6;
s++) {
214 tmp_error_par[
s] =
pow(tmp_par[
s], 0.5);
216 tmp_error_par[6] = 0.0;
218 reco::BeamSpot tmp_beam = Fit_d_z_likelihood(tmp_par, tmp_error_par);
220 double tmp_par2[7] = {tmp_beam.
x0(),
231 edm::LogWarning(
"BSFitter") <<
"Result is non physical. Log-Likelihood fit did not converge." << std::endl;
238 throw cms::Exception(
"LogicError") <<
"Error in BeamSpotProducer/BSFitter: " 239 <<
"Illegal fit type, options are likelihood or resolution";
242 throw cms::Exception(
"LogicError") <<
"Error in BeamSpotProducer/BSFitter: " 243 <<
"Illegal variable type, options are \"z\", \"d\", or \"d*z\"";
253 std::vector<double> par(2, 0);
254 std::vector<double>
err(2, 0);
258 err.push_back(0.0001);
259 err.push_back(0.0001);
263 thePDF->SetPDFs(
"PDFGauss_z");
264 thePDF->SetData(fBSvector);
267 MnUserParameters upar;
268 upar.Add(
"X0", 0., 0.);
269 upar.Add(
"Y0", 0., 0.);
270 upar.Add(
"Z0", inipar[0], 0.001);
271 upar.Add(
"sigmaZ", inipar[1], 0.001);
273 MnMigrad migrad(*thePDF, upar);
275 FunctionMinimum fmin = migrad();
276 ff_minimum = fmin.Fval();
297 for (
int j = 2;
j < 4; ++
j) {
298 for (
int k =
j;
k < 4; ++
k) {
304 fmin.Parameters().Vec()(3),
321 h1z =
new TH1F(
"h1z",
"z distribution", 200, -fMaxZ, fMaxZ);
323 std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
327 for (iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
328 h1z->Fill(iparam->z0());
334 TF1 fgaus(
"fgaus",
"gaus", 0., 1., TF1::EAddToList::kNo);
335 h1z->Fit(&fgaus,
"QLMN0 SERIAL");
339 double fpar[2] = {fgaus.GetParameter(1), fgaus.GetParameter(2)};
343 matrix(2, 2) = fgaus.GetParError(1) * fgaus.GetParError(1);
344 matrix(3, 3) = fgaus.GetParError(2) * fgaus.GetParError(2);
354 edm::LogInfo(
"BSFitter") <<
"number of total input tracks: " << fBSvector.size() << std::endl;
358 if ((
int)fBSvector.size() <= fminNtrks) {
359 edm::LogWarning(
"BSFitter") <<
"need at least " << fminNtrks <<
" tracks to run beamline fitter." << std::endl;
365 theanswer = Fit_d0phi();
372 while (goodfit && ftmprow > fconvergence * fBSvector.size() && ftmprow > fminNtrks) {
373 theanswer = Fit_d0phi();
376 if (goodfit && ftmprow > fconvergence * fBSvector.size() && ftmprow > fminNtrks) {
377 preanswer = theanswer;
384 theanswer = preanswer;
388 edm::LogInfo(
"BSFitter") <<
"Total number of successful iterations = " << (goodfit ? (fnthite + 1) : fnthite)
405 edm::LogInfo(
"BSFitter") <<
" number of tracks used: " << ftmprow << std::endl;
414 TMatrixD x_result(4, 1);
415 TMatrixDSym V_result(4);
428 std::vector<BSTrkParameters>::iterator iparam = fBSvector.begin();
435 for (iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
443 g(0, 0) =
sin(iparam->phi0());
444 g(1, 0) = -
cos(iparam->phi0());
445 g(2, 0) = iparam->z0() *
g(0, 0);
446 g(3, 0) = iparam->z0() *
g(1, 0);
449 double sigmabeam2 = 0.006 * 0.006;
450 if (finputBeamWidth > 0)
451 sigmabeam2 = finputBeamWidth * finputBeamWidth;
458 double sigma2 = sigmabeam2 + (iparam->sigd0()) * (iparam->sigd0());
460 TMatrixD ftmptrans(1, 4);
461 ftmptrans = ftmptrans.Transpose(ftmp);
462 TMatrixD dcor = ftmptrans *
g;
463 double chi2tmp = (iparam->d0() - dcor(0, 0)) * (iparam->d0() - dcor(0, 0)) / sigma2;
465 iparam->z0(), iparam->sigz0(), iparam->d0(), iparam->sigd0(), iparam->phi0(), iparam->pt(), dcor(0, 0), chi2tmp);
468 if (fapplyd0cut && fnthite > 0) {
469 if (
std::abs(iparam->d0() - dcor(0, 0)) > fd0cut)
472 if (fapplychi2cut && fnthite > 0) {
473 if (chi2tmp > fchi2cut)
479 for (
int j = 0;
j < 4; ++
j) {
480 for (
int k =
j;
k < 4; ++
k) {
485 Vint += (
temp * (1 / sigma2));
486 b += (iparam->d0() / sigma2 *
g);
489 h1z->Fill(iparam->z0());
492 Double_t determinant;
495 if (!
bk.Decompose()) {
497 edm::LogWarning(
"BSFitter") <<
"Decomposition failed, matrix singular ?" << std::endl
498 <<
"condition number = " <<
bk.Condition() << std::endl;
500 V_result = Vint.InvertFast(&determinant);
501 x_result = V_result *
b;
516 TF1 fgaus(
"fgaus",
"gaus", 0., 1., TF1::EAddToList::kNo);
520 h1z->Fit(&fgaus,
"QLN0 SERIAL",
"", h1z->GetMean() - 2. * h1z->GetRMS(), h1z->GetMean() + 2. * h1z->GetRMS());
531 double fpar[2] = {fgaus.GetParameter(1), fgaus.GetParameter(2)};
535 for (
int j = 0;
j < 2; ++
j) {
536 for (
int k =
j;
k < 2; ++
k) {
541 for (
int j = 4;
j < 6; ++
j) {
542 for (
int k =
j;
k < 6; ++
k) {
548 matrix(2, 2) = fgaus.GetParError(1) * fgaus.GetParError(1);
549 matrix(3, 3) = fgaus.GetParError(2) * fgaus.GetParError(2);
557 double x0tmp = x_result(0, 0);
558 double y0tmp = x_result(1, 0);
560 x0tmp += x_result(2, 0) * fpar[0];
561 y0tmp += x_result(3, 0) * fpar[0];
564 reco::BeamSpot::Point(x0tmp, y0tmp, fpar[0]), fpar[1], x_result(2, 0), x_result(3, 0), 0.,
matrix, fbeamtype);
577 fapplychi2cut =
true;
585 thePDF->SetPDFs(
"PDFGauss_d");
586 thePDF->SetData(fBSvector);
588 MnUserParameters upar;
589 upar.Add(
"X0", inipar[0], 0.001);
590 upar.Add(
"Y0", inipar[1], 0.001);
591 upar.Add(
"Z0", 0., 0.001);
592 upar.Add(
"sigmaZ", 0., 0.001);
593 upar.Add(
"dxdz", inipar[2], 0.001);
594 upar.Add(
"dydz", inipar[3], 0.001);
596 MnMigrad migrad(*thePDF, upar);
598 FunctionMinimum fmin = migrad();
599 ff_minimum = fmin.Fval();
602 for (
int j = 0;
j < 6; ++
j) {
603 for (
int k =
j;
k < 6; ++
k) {
610 fmin.Parameters().Vec()(4),
611 fmin.Parameters().Vec()(5),
619 init_pars[6] = 0.0005;
622 double fsqrt2pi = 0.0;
624 double d_dprime = 0.0;
625 double d_result = 0.0;
627 double z_result = 0.0;
628 double function = 0.0;
629 double tot_pdf = 0.0;
630 double last_minvalue = 1.0e+10;
631 double init_bw = -99.99;
635 double DeltadCut = 0.1000;
636 if (init_pars[6] < 0.0200) {
639 if (init_pars[6] < 0.0100) {
643 std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
650 for (
int p = 0;
p < iters;
p++) {
652 init_pars[6] += 0.0002;
655 for (iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
657 d_sig =
sqrt(init_pars[6] * init_pars[6] + (iparam->sigd0()) * (iparam->sigd0()));
658 d_dprime = iparam->d0() - (((init_pars[0] + iparam->z0() * (init_pars[4])) *
sin(iparam->phi0())) -
659 ((init_pars[1] + iparam->z0() * (init_pars[5])) *
cos(iparam->phi0())));
663 fBSvectorBW.push_back(*iparam);
666 d_result = (
exp(-(d_dprime * d_dprime) / (2.0 * d_sig * d_sig))) / (d_sig * fsqrt2pi);
667 z_sig =
sqrt(iparam->sigz0() * iparam->sigz0() + init_pars[3] * init_pars[3]);
668 z_result = (
exp(-((iparam->z0() - init_pars[2]) * (iparam->z0() - init_pars[2])) / (2.0 * z_sig * z_sig))) /
670 tot_pdf = z_result * d_result;
676 if (d_result < 1.0
e-05) {
677 tot_pdf = z_result * 1.0e-05;
682 function =
function +
log(tot_pdf);
686 function = -2.0 *
function;
687 if (
function < last_minvalue) {
688 init_bw = init_pars[6];
689 last_minvalue =
function;
695 init_bw = init_bw + (0.20 * init_bw);
700 <<
"scanPDF:====>>>> WARNING***: The initial guess value of Beam width is negative!!!!!!" << std::endl
701 <<
"scanPDF:====>>>> Assigning beam width a starting value of " << init_bw <<
" cm" << std::endl;
711 int tracksFailed = 0;
714 inipar[6] = scanPDF(inipar, tracksFailed, 1);
715 error_par[6] = (inipar[6]) * 0.20;
719 scanPDF(inipar, tracksFailed, 2);
724 std::vector<BSTrkParameters>::const_iterator iparamBW = fBSvectorBW.begin();
725 for (iparamBW = fBSvectorBW.begin(); iparamBW != fBSvectorBW.end(); ++iparamBW) {
726 fBSvector.push_back(*iparamBW);
729 thePDF->SetPDFs(
"PDFGauss_d*PDFGauss_z");
730 thePDF->SetData(fBSvector);
731 MnUserParameters upar;
733 upar.Add(
"X0", inipar[0], error_par[0]);
734 upar.Add(
"Y0", inipar[1], error_par[1]);
735 upar.Add(
"Z0", inipar[2], error_par[2]);
736 upar.Add(
"sigmaZ", inipar[3], error_par[3]);
737 upar.Add(
"dxdz", inipar[4], error_par[4]);
738 upar.Add(
"dydz", inipar[5], error_par[5]);
739 upar.Add(
"BeamWidthX", inipar[6], error_par[6]);
741 MnMigrad migrad(*thePDF, upar);
743 FunctionMinimum fmin = migrad();
748 ff_minimum = fmin.Fval();
750 bool ff_nfcn = fmin.HasReachedCallLimit();
751 bool ff_cov = fmin.HasCovariance();
756 edm::LogWarning(
"BSFitter") <<
"===========>>>>>** WARNING: MINUIT DID NOT CONVERGES PROPERLY !!!!!!" << std::endl;
758 edm::LogWarning(
"BSFitter") <<
"===========>>>>>** WARNING: No. of Calls Exhausted" << std::endl;
760 edm::LogWarning(
"BSFitter") <<
"===========>>>>>** WARNING: Covariance did not found" << std::endl;
763 edm::LogInfo(
"BSFitter") <<
"The Total # Tracks used for beam width fit = " << (fBSvectorBW.size()) << std::endl;
766 double lastIter_pars[7];
768 for (
int ip = 0; ip < 7; ip++) {
769 lastIter_pars[ip] = fmin.Parameters().Vec()(ip);
773 scanPDF(lastIter_pars, tracksFailed, 2);
775 edm::LogWarning(
"BSFitter") <<
"WARNING: # of tracks which have very low pdf value (pdf_d < 1.0e-05) are = " 776 << tracksFailed << std::endl;
783 for (
int j = 0;
j < 7; ++
j) {
784 for (
int k =
j;
k < 7; ++
k) {
790 reco::BeamSpot::Point(fmin.Parameters().Vec()(0), fmin.Parameters().Vec()(1), fmin.Parameters().Vec()(2)),
791 fmin.Parameters().Vec()(3),
792 fmin.Parameters().Vec()(4),
793 fmin.Parameters().Vec()(5),
794 fmin.Parameters().Vec()(6),
802 thePDF->SetPDFs(
"PDFGauss_d_resolution*PDFGauss_z");
803 thePDF->SetData(fBSvector);
805 MnUserParameters upar;
806 upar.Add(
"X0", inipar[0], 0.001);
807 upar.Add(
"Y0", inipar[1], 0.001);
808 upar.Add(
"Z0", inipar[2], 0.001);
809 upar.Add(
"sigmaZ", inipar[3], 0.001);
810 upar.Add(
"dxdz", inipar[4], 0.001);
811 upar.Add(
"dydz", inipar[5], 0.001);
812 upar.Add(
"BeamWidthX", inipar[6], 0.0001);
813 upar.Add(
"c0", 0.0010, 0.0001);
814 upar.Add(
"c1", 0.0090, 0.0001);
817 upar.Fix(
"BeamWidthX");
820 MnMigrad migrad(*thePDF, upar);
822 FunctionMinimum fmin = migrad();
823 ff_minimum = fmin.Fval();
827 for (
int j = 0;
j < 6; ++
j) {
828 for (
int k =
j;
k < 6; ++
k) {
838 fresolution_c0 = fmin.Parameters().Vec()(6);
839 fresolution_c1 = fmin.Parameters().Vec()(7);
840 fres_c0_err =
sqrt(fmin.Error().Matrix()(6, 6));
841 fres_c1_err =
sqrt(fmin.Error().Matrix()(7, 7));
843 for (
int j = 6;
j < 8; ++
j) {
844 for (
int k = 6;
k < 8; ++
k) {
845 fres_matrix(
j - 6,
k - 6) = fmin.Error().Matrix()(
j,
k);
850 reco::BeamSpot::Point(fmin.Parameters().Vec()(0), fmin.Parameters().Vec()(1), fmin.Parameters().Vec()(2)),
851 fmin.Parameters().Vec()(3),
852 fmin.Parameters().Vec()(4),
853 fmin.Parameters().Vec()(5),
reco::BeamSpot Fit_d_likelihood(double *inipar)
math::Error< dimension >::type CovarianceMatrix
double BeamWidthX() const
beam width X
constexpr bool isNotFinite(T x)
reco::BeamSpot Fit_z_likelihood(double *inipar)
Sin< T >::type sin(const T &t)
math::XYZPoint Point
point in the space
double dydz() const
dydz slope
double scanPDF(double *init_pars, int &tracksFailed, int option)
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
reco::BeamSpot Fit_dres_z_likelihood(double *inipar)
void setType(BeamType type)
set beam type
double x0() const
x coordinate
Cos< T >::type cos(const T &t)
Abs< T >::type abs(const T &t)
double y0() const
y coordinate
Log< level::Info, false > LogInfo
double sigmaZ() const
sigma z
reco::BeamSpot Fit_ited0phi()
double dxdz() const
dxdz slope
reco::BeamSpot Fit_z_chi2(double *inipar)
void Setd0Cut_d0phi(double d0cut)
double covariance(int i, int j) const
(i,j)-th element of error matrix
double z0() const
z coordinate
void SetChi2Cut_d0phi(double chi2cut)
Log< level::Warning, false > LogWarning
reco::BeamSpot Fit_d0phi()
reco::BeamSpot Fit_d_z_likelihood(double *inipar, double *error_par)
Power< A, B >::type pow(const A &a, const B &b)