14 #include "Minuit2/VariableMetricMinimizer.h"
15 #include "Minuit2/FunctionMinimum.h"
16 #include "Minuit2/MnPrint.h"
17 #include "Minuit2/MnMigrad.h"
18 #include "Minuit2/MnUserParameterState.h"
33 #include "TMatrixDSym.h"
34 #include "TDecompBK.h"
38 using namespace ROOT::Minuit2;
49 ffit_type =
"default";
50 ffit_variable =
"default";
57 fpar_name[1] =
"SigmaZ0 ";
60 fpar_name[4] =
"dxdz ";
61 fpar_name[5] =
"dydz ";
62 fpar_name[6] =
"SigmaBeam ";
71 theFitter =
new VariableMetricMinimizer();
76 fapplychi2cut =
false;
86 h1z =
new TH1F(
"h1z",
"z distribution",200,-fMaxZ, fMaxZ);
109 if ( ffit_variable ==
"z" ) {
111 if ( ffit_type ==
"chi2" ) {
113 return Fit_z_chi2(inipar);
115 }
else if ( ffit_type ==
"likelihood" ) {
117 return Fit_z_likelihood(inipar);
119 }
else if ( ffit_type ==
"combined" ) {
122 double tmp_par[2] = {tmp_beamspot.
z0(), tmp_beamspot.
sigmaZ()};
123 return Fit_z_likelihood(tmp_par);
128 <<
"Error in BeamSpotProducer/BSFitter: "
129 <<
"Illegal fit type, options are chi2,likelihood or combined(ie. first chi2 then likelihood)";
132 }
else if ( ffit_variable ==
"d" ) {
134 if ( ffit_type ==
"d0phi" ) {
138 }
else if ( ffit_type ==
"likelihood" ) {
140 return Fit_d_likelihood(inipar);
142 }
else if ( ffit_type ==
"combined" ) {
146 double tmp_par[4] = {tmp_beamspot.
x0(), tmp_beamspot.
y0(), tmp_beamspot.
dxdz(), tmp_beamspot.
dydz()};
147 return Fit_d_likelihood(tmp_par);
151 <<
"Error in BeamSpotProducer/BSFitter: "
152 <<
"Illegal fit type, options are d0phi, likelihood or combined";
154 }
else if ( ffit_variable ==
"d*z" || ffit_variable ==
"default" ) {
156 if ( ffit_type ==
"likelihood" || ffit_type ==
"default" ) {
171 this->Setd0Cut_d0phi(4.0);
203 if (ffit_type ==
"likelihood") {
204 double tmp_par[7] = {tmp_d0phi.
x0(), tmp_d0phi.
y0(), tmp_d0phi.
z0(),
207 double tmp_error_par[7];
209 tmp_error_par[6]=0.0;
211 reco::BeamSpot tmp_lh = Fit_d_z_likelihood(tmp_par,tmp_error_par);
214 edm::LogWarning(
"BSFitter") <<
"BSFitter: Result is non physical. Log-Likelihood fit to extract beam width did not converge." << std::endl;
222 edm::LogInfo(
"BSFitter") <<
"default track-based fit does not extract beam width." << std::endl;
227 }
else if ( ffit_type ==
"resolution" ) {
233 double tmp_par[7] = {tmp_d0phi.
x0(), tmp_d0phi.
y0(), tmp_z.
z0(),
235 double tmp_error_par[7];
236 for(
int s=0;
s<6;
s++){ tmp_error_par[
s] =
pow(tmp_par[
s],0.5);}
237 tmp_error_par[6]=0.0;
239 reco::BeamSpot tmp_beam = Fit_d_z_likelihood(tmp_par,tmp_error_par);
241 double tmp_par2[7] = {tmp_beam.
x0(), tmp_beam.
y0(), tmp_beam.
z0(),
249 edm::LogWarning(
"BSFitter") <<
"Result is non physical. Log-Likelihood fit did not converge." << std::endl;
258 <<
"Error in BeamSpotProducer/BSFitter: "
259 <<
"Illegal fit type, options are likelihood or resolution";
264 <<
"Error in BeamSpotProducer/BSFitter: "
265 <<
"Illegal variable type, options are \"z\", \"d\", or \"d*z\"";
278 std::vector<double> par(2,0);
279 std::vector<double> err(2,0);
283 err.push_back(0.0001);
284 err.push_back(0.0001);
288 thePDF->SetPDFs(
"PDFGauss_z");
289 thePDF->SetData(fBSvector);
293 MnUserParameters upar;
294 upar.Add(
"X0", 0.,0.);
295 upar.Add(
"Y0", 0.,0.);
296 upar.Add(
"Z0", inipar[0],0.001);
297 upar.Add(
"sigmaZ",inipar[1],0.001);
299 MnMigrad migrad(*thePDF, upar);
301 FunctionMinimum fmin = migrad();
302 ff_minimum = fmin.Fval();
323 for (
int j = 2 ;
j < 4 ; ++
j) {
324 for(
int k =
j ;
k < 4 ; ++
k) {
331 fmin.Parameters().Vec()(2)),
332 fmin.Parameters().Vec()(3),
351 h1z =
new TH1F(
"h1z",
"z distribution",200,-fMaxZ, fMaxZ);
353 std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
357 for( iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
359 h1z->Fill( iparam->z0() );
364 TF1 fgaus(
"fgaus",
"gaus");
365 h1z->Fit(&fgaus,
"QLM0");
369 double fpar[2] = {fgaus.GetParameter(1), fgaus.GetParameter(2) };
373 matrix(2,2) = fgaus.GetParError(1) * fgaus.GetParError(1);
374 matrix(3,3) = fgaus.GetParError(2) * fgaus.GetParError(2);
395 edm::LogInfo(
"BSFitter") <<
"number of total input tracks: " << fBSvector.size() << std::endl;
399 if ( (
int)fBSvector.size() <= fminNtrks ) {
400 edm::LogWarning(
"BSFitter") <<
"need at least " << fminNtrks <<
" tracks to run beamline fitter." << std::endl;
406 theanswer = Fit_d0phi();
407 if ( goodfit ) fnthite++;
413 ftmprow > fconvergence * fBSvector.size() &&
414 ftmprow > fminNtrks ) {
416 theanswer = Fit_d0phi();
420 ftmprow > fconvergence * fBSvector.size() &&
421 ftmprow > fminNtrks ) {
422 preanswer = theanswer;
429 theanswer = preanswer;
433 edm::LogInfo(
"BSFitter") <<
"Total number of successful iterations = " << ( goodfit ? (fnthite+1) : fnthite ) << std::endl;
451 if (fnthite > 0)
edm::LogInfo(
"BSFitter") <<
" number of tracks used: " << ftmprow << std::endl;
461 TMatrixD x_result(4,1);
462 TMatrixDSym V_result(4);
475 std::vector<BSTrkParameters>::iterator iparam = fBSvector.begin();
483 for( iparam = fBSvector.begin() ;
484 iparam != fBSvector.end() ; ++iparam) {
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);
501 double sigmabeam2 = 0.006 * 0.006;
502 if (finputBeamWidth > 0 ) sigmabeam2 = finputBeamWidth * finputBeamWidth;
509 double sigma2 = sigmabeam2 + (iparam->sigd0())* (iparam->sigd0());
511 TMatrixD ftmptrans(1,4);
512 ftmptrans = ftmptrans.Transpose(ftmp);
513 TMatrixD dcor = ftmptrans *
g;
514 double chi2tmp = (iparam->d0() - dcor(0,0)) * (iparam->d0() - dcor(0,0))/sigma2;
515 (*iparam) =
BSTrkParameters(iparam->z0(),iparam->sigz0(),iparam->d0(),iparam->sigd0(),
516 iparam->phi0(), iparam->pt(),dcor(0,0),chi2tmp);
519 if (fapplyd0cut && fnthite>0 ) {
520 if (
std::abs(iparam->d0() - dcor(0,0)) > fd0cut ) pass =
false;
523 if (fapplychi2cut && fnthite>0 ) {
524 if ( chi2tmp > fchi2cut ) pass =
false;
530 for(
int j = 0 ;
j < 4 ; ++
j) {
531 for(
int k =
j ;
k < 4 ; ++
k) {
537 Vint += (temp * (1 / sigma2));
538 b += (iparam->d0() / sigma2 *
g);
541 h1z->Fill( iparam->z0() );
546 Double_t determinant;
549 if (!bk.Decompose()) {
552 <<
"Decomposition failed, matrix singular ?" << std::endl
553 <<
"condition number = " << bk.Condition() << std::endl;
556 V_result = Vint.InvertFast(&determinant);
557 x_result = V_result *
b;
571 TF1 fgaus(
"fgaus",
"gaus");
574 auto status = h1z->Fit(&fgaus,
"QL0",
"",h1z->GetMean() -2.*h1z->GetRMS(),h1z->GetMean() +2.*h1z->GetRMS());
584 double fpar[2] = {fgaus.GetParameter(1), fgaus.GetParameter(2) };
588 for (
int j = 0 ;
j < 2 ; ++
j) {
589 for(
int k =
j ;
k < 2 ; ++
k) {
594 for (
int j = 4 ;
j < 6 ; ++
j) {
595 for(
int k =
j ;
k < 6 ; ++
k) {
601 matrix(2,2) = fgaus.GetParError(1) * fgaus.GetParError(1);
602 matrix(3,3) = fgaus.GetParError(2) * fgaus.GetParError(2);
610 double x0tmp = x_result(0,0);
611 double y0tmp = x_result(1,0);
613 x0tmp += x_result(2,0)*fpar[0];
614 y0tmp += x_result(3,0)*fpar[0];
642 fapplychi2cut =
true;
652 thePDF->SetPDFs(
"PDFGauss_d");
653 thePDF->SetData(fBSvector);
655 MnUserParameters upar;
656 upar.Add(
"X0", inipar[0],0.001);
657 upar.Add(
"Y0", inipar[1],0.001);
658 upar.Add(
"Z0", 0.,0.001);
659 upar.Add(
"sigmaZ",0.,0.001);
660 upar.Add(
"dxdz",inipar[2],0.001);
661 upar.Add(
"dydz",inipar[3],0.001);
664 MnMigrad migrad(*thePDF, upar);
666 FunctionMinimum fmin = migrad();
667 ff_minimum = fmin.Fval();
670 for (
int j = 0 ;
j < 6 ; ++
j) {
671 for(
int k =
j ;
k < 6 ; ++
k) {
677 fmin.Parameters().Vec()(1),
680 fmin.Parameters().Vec()(4),
681 fmin.Parameters().Vec()(5),
689 if(option==1)init_pars[6]=0.0005;
700 double last_minvalue=1.0e+10;
701 double init_bw=-99.99;
705 double DeltadCut=0.1000;
706 if(init_pars[6]<0.0200){DeltadCut=0.0900; }
707 if(init_pars[6]<0.0100){DeltadCut=0.0700;}
710 std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
713 if(option==1)iters=500;
714 if(option==2)iters=1;
716 for(
int p=0;
p<iters;
p++){
718 if(iters==500)init_pars[6]+=0.0002;
721 for( iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam)
724 d_sig =
sqrt(init_pars[6]*init_pars[6] + (iparam->sigd0())*(iparam->sigd0()));
725 d_dprime = iparam->d0() - ( ( (init_pars[0] + iparam->z0()*(init_pars[4]))*
sin(iparam->phi0()) )
726 - ( (init_pars[1] + iparam->z0()*(init_pars[5]))*
cos(iparam->phi0()) ) );
729 if(
std::abs(d_dprime)<DeltadCut && option==2){ fBSvectorBW.push_back(*iparam);}
731 d_result = (
exp(-(d_dprime*d_dprime)/(2.0*d_sig*d_sig)))/(d_sig*fsqrt2pi);
732 z_sig =
sqrt(iparam->sigz0() * iparam->sigz0() + init_pars[3]*init_pars[3]);
733 z_result = (
exp(-((iparam->z0() - init_pars[2])*(iparam->z0() - init_pars[2]))/(2.0*z_sig*z_sig)))/(z_sig*fsqrt2pi);
734 tot_pdf=z_result*d_result;
740 if(d_result < 1.0
e-05){ tot_pdf=z_result*1.0e-05;
744 function =
function +
log(tot_pdf);
750 function= -2.0*
function;
751 if(
function<last_minvalue){init_bw=init_pars[6];
752 last_minvalue=
function; }
757 init_bw=init_bw+(0.20*init_bw);
764 <<
"scanPDF:====>>>> WARNING***: The initial guess value of Beam width is negative!!!!!!"<<std::endl
765 <<
"scanPDF:====>>>> Assigning beam width a starting value of "<<init_bw<<
" cm"<<std::endl;
782 inipar[6]=scanPDF(inipar,tracksFailed,1);
783 error_par[6]=(inipar[6])*0.20;
788 scanPDF(inipar,tracksFailed,2);
793 std::vector<BSTrkParameters>::const_iterator iparamBW = fBSvectorBW.begin();
794 for( iparamBW = fBSvectorBW.begin(); iparamBW != fBSvectorBW.end(); ++iparamBW)
795 { fBSvector.push_back(*iparamBW);
799 thePDF->SetPDFs(
"PDFGauss_d*PDFGauss_z");
800 thePDF->SetData(fBSvector);
801 MnUserParameters upar;
803 upar.Add(
"X0", inipar[0],error_par[0]);
804 upar.Add(
"Y0", inipar[1],error_par[1]);
805 upar.Add(
"Z0", inipar[2],error_par[2]);
806 upar.Add(
"sigmaZ",inipar[3],error_par[3]);
807 upar.Add(
"dxdz",inipar[4],error_par[4]);
808 upar.Add(
"dydz",inipar[5],error_par[5]);
809 upar.Add(
"BeamWidthX",inipar[6],error_par[6]);
812 MnMigrad migrad(*thePDF, upar);
814 FunctionMinimum fmin = migrad();
819 ff_minimum = fmin.Fval();
822 bool ff_nfcn=fmin.HasReachedCallLimit();
823 bool ff_cov=fmin.HasCovariance();
824 bool testing=fmin.IsValid();
830 edm::LogWarning(
"BSFitter") <<
"===========>>>>>** WARNING: MINUIT DID NOT CONVERGES PROPERLY !!!!!!"<<std::endl;
831 if(ff_nfcn)
edm::LogWarning(
"BSFitter") <<
"===========>>>>>** WARNING: No. of Calls Exhausted"<<std::endl;
832 if(!ff_cov)
edm::LogWarning(
"BSFitter") <<
"===========>>>>>** WARNING: Covariance did not found"<<std::endl;
835 edm::LogInfo(
"BSFitter") <<
"The Total # Tracks used for beam width fit = "<<(fBSvectorBW.size())<<std::endl;
839 double lastIter_pars[7];
841 for(
int ip=0;ip<7;ip++){ lastIter_pars[ip]=fmin.Parameters().Vec()(ip);
847 scanPDF(lastIter_pars,tracksFailed,2);
850 edm::LogWarning(
"BSFitter") <<
"WARNING: # of tracks which have very low pdf value (pdf_d < 1.0e-05) are = "<<tracksFailed<<std::endl;
859 for (
int j = 0 ;
j < 7 ; ++
j) {
860 for(
int k =
j ;
k < 7 ; ++
k) {
867 fmin.Parameters().Vec()(1),
868 fmin.Parameters().Vec()(2)),
869 fmin.Parameters().Vec()(3),
870 fmin.Parameters().Vec()(4),
871 fmin.Parameters().Vec()(5),
872 fmin.Parameters().Vec()(6),
883 thePDF->SetPDFs(
"PDFGauss_d_resolution*PDFGauss_z");
884 thePDF->SetData(fBSvector);
886 MnUserParameters upar;
887 upar.Add(
"X0", inipar[0],0.001);
888 upar.Add(
"Y0", inipar[1],0.001);
889 upar.Add(
"Z0", inipar[2],0.001);
890 upar.Add(
"sigmaZ",inipar[3],0.001);
891 upar.Add(
"dxdz",inipar[4],0.001);
892 upar.Add(
"dydz",inipar[5],0.001);
893 upar.Add(
"BeamWidthX",inipar[6],0.0001);
894 upar.Add(
"c0",0.0010,0.0001);
895 upar.Add(
"c1",0.0090,0.0001);
898 upar.Fix(
"BeamWidthX");
901 MnMigrad migrad(*thePDF, upar);
903 FunctionMinimum fmin = migrad();
904 ff_minimum = fmin.Fval();
908 for (
int j = 0 ;
j < 6 ; ++
j) {
909 for(
int k =
j ;
k < 6 ; ++
k) {
919 fresolution_c0 = fmin.Parameters().Vec()(6);
920 fresolution_c1 = fmin.Parameters().Vec()(7);
921 fres_c0_err =
sqrt( fmin.Error().Matrix()(6,6) );
922 fres_c1_err =
sqrt( fmin.Error().Matrix()(7,7) );
924 for (
int j = 6 ;
j < 8 ; ++
j) {
925 for(
int k = 6 ;
k < 8 ; ++
k) {
926 fres_matrix(
j-6,
k-6) = fmin.Error().Matrix()(
j,
k);
931 fmin.Parameters().Vec()(1),
932 fmin.Parameters().Vec()(2)),
933 fmin.Parameters().Vec()(3),
934 fmin.Parameters().Vec()(4),
935 fmin.Parameters().Vec()(5),
reco::BeamSpot Fit_d_likelihood(double *inipar)
math::Error< dimension >::type CovarianceMatrix
double z0() const
z coordinate
reco::BeamSpot Fit_z_likelihood(double *inipar)
Sin< T >::type sin(const T &t)
math::XYZPoint Point
point in the space
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 dydz() const
dydz slope
Cos< T >::type cos(const T &t)
Abs< T >::type abs(const T &t)
double BeamWidthX() const
beam width X
double dxdz() const
dxdz slope
reco::BeamSpot Fit_ited0phi()
double sigmaZ() const
sigma z
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 y0() const
y coordinate
void SetChi2Cut_d0phi(double chi2cut)
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)
double x0() const
x coordinate