CMS 3D CMS Logo

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

#include <BSFitter.h>

Public Member Functions

 BSFitter ()
 
 BSFitter (const std::vector< BSTrkParameters > &BSvector)
 
void d0phi_Init ()
 
reco::BeamSpot Fit ()
 
reco::BeamSpot Fit (double *inipar)
 
reco::BeamSpot Fit_d0phi ()
 
reco::BeamSpot Fit_d_likelihood (double *inipar)
 
reco::BeamSpot Fit_d_z_likelihood (double *inipar, double *error_par)
 
reco::BeamSpot Fit_dres_z_likelihood (double *inipar)
 
reco::BeamSpot Fit_ited0phi ()
 
reco::BeamSpot Fit_z (std::string type, double *inipar)
 
reco::BeamSpot Fit_z_chi2 (double *inipar)
 
reco::BeamSpot Fit_z_likelihood (double *inipar)
 
int GetAcceptedTrks ()
 
std::vector< BSTrkParametersGetData ()
 
double GetMinimum ()
 
reco::BeamSpot::ResCovMatrix GetResMatrix ()
 
double GetResPar0 ()
 
double GetResPar0Err ()
 
double GetResPar1 ()
 
double GetResPar1Err ()
 
TH1F * GetVzHisto ()
 
double scanPDF (double *init_pars, int &tracksFailed, int option)
 
void SetChi2Cut_d0phi (double chi2cut)
 
void SetConvergence (double val)
 
void Setd0Cut_d0phi (double d0cut)
 
void SetFitType (std::string type)
 
void SetFitVariable (std::string name)
 
void SetInputBeamWidth (double val)
 
void SetMaximumZ (double z)
 
void SetMinimumNTrks (int n)
 
virtual ~BSFitter ()
 

Private Attributes

bool fapplychi2cut
 
bool fapplyd0cut
 
reco::BeamSpot::BeamType fbeamtype
 
std::vector< BSTrkParametersfBSvector
 
std::vector< BSTrkParametersfBSvectorBW
 
double fchi2cut
 
double fconvergence
 
double fd0cut
 
double ff_minimum
 
std::string ffit_type
 
std::string ffit_variable
 
double finputBeamWidth
 
double fMaxZ
 
int fminNtrks
 
int fnthite
 
std::string fpar_name [fdim]
 
double fres_c0_err
 
double fres_c1_err
 
reco::BeamSpot::ResCovMatrix fres_matrix
 
double fresolution_c0
 
double fresolution_c1
 
Double_t fsqrt2pi
 
TMatrixD ftmp
 
int ftmprow
 
bool goodfit
 
TH1F * h1z
 
ROOT::Minuit2::ModularFunctionMinimizer * theFitter
 
BSpdfsFcnthePDF
 

Static Private Attributes

static const int fdim = 7
 

Detailed Description


class: BSFitter.h package: RecoVertex/BeamSpotProducer

author: Francisco Yumiceva, Fermilab (yumic.nosp@m.eva@.nosp@m.fnal..nosp@m.gov)


Definition at line 31 of file BSFitter.h.

Constructor & Destructor Documentation

BSFitter::BSFitter ( )

Definition at line 43 of file BSFitter.cc.

References reco::BeamSpot::Unknown.

43  {
45 
46  //In order to make fitting ROOT histograms thread safe
47  // one must call this undocumented function
48  TMinuitMinimizer::UseStaticMinuit(false);
49 }
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
BSFitter::BSFitter ( const std::vector< BSTrkParameters > &  BSvector)

Definition at line 52 of file BSFitter.cc.

References Pi, and mathSSE::sqrt().

52  {
53  //In order to make fitting ROOT histograms thread safe
54  // one must call this undocumented function
55  TMinuitMinimizer::UseStaticMinuit(false);
56 
57  ffit_type = "default";
58  ffit_variable = "default";
59 
60  fBSvector = BSvector;
61 
62  fsqrt2pi = sqrt(2.* TMath::Pi());
63 
64  fpar_name[0] = "z0 ";
65  fpar_name[1] = "SigmaZ0 ";
66  fpar_name[2] = "X0 ";
67  fpar_name[3] = "Y0 ";
68  fpar_name[4] = "dxdz ";
69  fpar_name[5] = "dydz ";
70  fpar_name[6] = "SigmaBeam ";
71 
72  //if (theGausszFcn == 0 ) {
73  thePDF = new BSpdfsFcn();
74 
75 
76 //}
77  //if (theFitter == 0 ) {
78 
79  theFitter = new VariableMetricMinimizer();
80 
81  //}
82 
83  fapplyd0cut = false;
84  fapplychi2cut = false;
85  ftmprow = 0;
86  ftmp.ResizeTo(4,1);
87  ftmp.Zero();
88  fnthite=0;
89  fMaxZ = 50.; //cm
90  fconvergence = 0.5; // stop fit when 50% of the input collection has been removed.
91  fminNtrks = 100;
92  finputBeamWidth = -1; // no input
93 
94  h1z = new TH1F("h1z","z distribution",200,-fMaxZ, fMaxZ);
95 
96 }
const double Pi
TH1F * h1z
Definition: BSFitter.h:146
double fMaxZ
Definition: BSFitter.h:142
int fminNtrks
Definition: BSFitter.h:144
BSpdfsFcn * thePDF
Definition: BSFitter.h:111
std::string ffit_variable
Definition: BSFitter.h:115
T sqrt(T t)
Definition: SSEVec.h:18
bool fapplychi2cut
Definition: BSFitter.h:136
Double_t fsqrt2pi
Definition: BSFitter.h:123
std::string ffit_type
Definition: BSFitter.h:114
int ftmprow
Definition: BSFitter.h:139
TMatrixD ftmp
Definition: BSFitter.h:134
double fconvergence
Definition: BSFitter.h:143
std::string fpar_name[fdim]
Definition: BSFitter.h:121
bool fapplyd0cut
Definition: BSFitter.h:135
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
int fnthite
Definition: BSFitter.h:140
ROOT::Minuit2::ModularFunctionMinimizer * theFitter
Definition: BSFitter.h:109
double finputBeamWidth
Definition: BSFitter.h:145
BSFitter::~BSFitter ( )
virtual

Definition at line 99 of file BSFitter.cc.

100 {
101  //delete fBSvector;
102  delete thePDF;
103  delete theFitter;
104 }
BSpdfsFcn * thePDF
Definition: BSFitter.h:111
ROOT::Minuit2::ModularFunctionMinimizer * theFitter
Definition: BSFitter.h:109

Member Function Documentation

void BSFitter::d0phi_Init ( )
inline

Definition at line 68 of file BSFitter.h.

References fnthite, ftmp, ftmprow, and goodfit.

68  {
69  ftmprow = 0;
70  ftmp.ResizeTo(4,1);
71  ftmp.Zero();
72  fnthite=0;
73  goodfit=true;
74  }
bool goodfit
Definition: BSFitter.h:141
int ftmprow
Definition: BSFitter.h:139
TMatrixD ftmp
Definition: BSFitter.h:134
int fnthite
Definition: BSFitter.h:140
reco::BeamSpot BSFitter::Fit ( )

Definition at line 108 of file BSFitter.cc.

Referenced by BeamFitter::runAllFitter(), BeamFitter::runBeamWidthFitter(), BeamFitter::runFitterNoTxt(), and SetFitVariable().

108  {
109 
110  return this->Fit(0);
111 
112 }
reco::BeamSpot Fit()
Definition: BSFitter.cc:108
reco::BeamSpot BSFitter::Fit ( double *  inipar = 0)

Definition at line 115 of file BSFitter.cc.

References reco::BeamSpot::BeamWidthX(), reco::BeamSpot::covariance(), reco::BeamSpot::dxdz(), reco::BeamSpot::dydz(), Exception, edm::isNotFinite(), makeMuonMisalignmentScenario::matrix, funct::pow(), alignCSCRings::s, reco::BeamSpot::setType(), reco::BeamSpot::sigmaZ(), reco::BeamSpot::Unknown, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().

115  {
117  if ( ffit_variable == "z" ) {
118 
119  if ( ffit_type == "chi2" ) {
120 
121  return Fit_z_chi2(inipar);
122 
123  } else if ( ffit_type == "likelihood" ) {
124 
125  return Fit_z_likelihood(inipar);
126 
127  } else if ( ffit_type == "combined" ) {
128 
129  reco::BeamSpot tmp_beamspot = Fit_z_chi2(inipar);
130  double tmp_par[2] = {tmp_beamspot.z0(), tmp_beamspot.sigmaZ()};
131  return Fit_z_likelihood(tmp_par);
132 
133  } else {
134 
135  throw cms::Exception("LogicError")
136  << "Error in BeamSpotProducer/BSFitter: "
137  << "Illegal fit type, options are chi2,likelihood or combined(ie. first chi2 then likelihood)";
138 
139  }
140  } else if ( ffit_variable == "d" ) {
141 
142  if ( ffit_type == "d0phi" ) {
143  this->d0phi_Init();
144  return Fit_d0phi();
145 
146  } else if ( ffit_type == "likelihood" ) {
147 
148  return Fit_d_likelihood(inipar);
149 
150  } else if ( ffit_type == "combined" ) {
151 
152  this->d0phi_Init();
153  reco::BeamSpot tmp_beamspot = Fit_d0phi();
154  double tmp_par[4] = {tmp_beamspot.x0(), tmp_beamspot.y0(), tmp_beamspot.dxdz(), tmp_beamspot.dydz()};
155  return Fit_d_likelihood(tmp_par);
156 
157  } else {
158  throw cms::Exception("LogicError")
159  << "Error in BeamSpotProducer/BSFitter: "
160  << "Illegal fit type, options are d0phi, likelihood or combined";
161  }
162  } else if ( ffit_variable == "d*z" || ffit_variable == "default" ) {
163 
164  if ( ffit_type == "likelihood" || ffit_type == "default" ) {
165 
167  // we are now fitting Z inside d0phi fitter
168  // first fit z distribution using a chi2 fit
169  //reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
170  //for (int j = 2 ; j < 4 ; ++j) {
171  //for(int k = j ; k < 4 ; ++k) {
172  // matrix(j,k) = tmp_z.covariance()(j,k);
173  //}
174  //}
175 
176  // use d0-phi algorithm to extract transverse position
177  this->d0phi_Init();
178  //reco::BeamSpot tmp_d0phi= Fit_d0phi(); // change to iterative procedure:
179  this->Setd0Cut_d0phi(4.0);
180  reco::BeamSpot tmp_d0phi= Fit_ited0phi();
181 
182  //for (int j = 0 ; j < 2 ; ++j) {
183  // for(int k = j ; k < 2 ; ++k) {
184  // matrix(j,k) = tmp_d0phi.covariance()(j,k);
185  //}
186  //}
187  // slopes
188  //for (int j = 4 ; j < 6 ; ++j) {
189  // for(int k = j ; k < 6 ; ++k) {
190  // matrix(j,k) = tmp_d0phi.covariance()(j,k);
191  // }
192  //}
193 
194 
195  // put everything into one object
196  reco::BeamSpot spot(reco::BeamSpot::Point(tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_d0phi.z0()),
197  tmp_d0phi.sigmaZ(),
198  tmp_d0phi.dxdz(),
199  tmp_d0phi.dydz(),
200  0.,
201  tmp_d0phi.covariance(),
202  fbeamtype );
203 
204 
205 
206  //reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
207 
208  //reco::BeamSpot tmp_d0phi = Fit_d0phi();
209 
210  // log-likelihood fit
211  if (ffit_type == "likelihood") {
212  double tmp_par[7] = {tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_d0phi.z0(),
213  tmp_d0phi.sigmaZ(), tmp_d0phi.dxdz(), tmp_d0phi.dydz(),0.0};
214 
215  double tmp_error_par[7];
216  for(int s=0;s<6;s++){ tmp_error_par[s] = pow( tmp_d0phi.covariance()(s,s),0.5);}
217  tmp_error_par[6]=0.0;
218 
219  reco::BeamSpot tmp_lh = Fit_d_z_likelihood(tmp_par,tmp_error_par);
220 
222  edm::LogWarning("BSFitter") << "BSFitter: Result is non physical. Log-Likelihood fit to extract beam width did not converge." << std::endl;
224  return tmp_lh;
225  }
226  return tmp_lh;
227 
228  } else {
229 
230  edm::LogInfo("BSFitter") << "default track-based fit does not extract beam width." << std::endl;
231  return spot;
232  }
233 
234 
235  } else if ( ffit_type == "resolution" ) {
236 
237  reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
238  this->d0phi_Init();
239  reco::BeamSpot tmp_d0phi = Fit_d0phi();
240 
241  double tmp_par[7] = {tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_z.z0(),
242  tmp_z.sigmaZ(), tmp_d0phi.dxdz(), tmp_d0phi.dydz(),0.0};
243  double tmp_error_par[7];
244  for(int s=0;s<6;s++){ tmp_error_par[s] = pow(tmp_par[s],0.5);}
245  tmp_error_par[6]=0.0;
246 
247  reco::BeamSpot tmp_beam = Fit_d_z_likelihood(tmp_par,tmp_error_par);
248 
249  double tmp_par2[7] = {tmp_beam.x0(), tmp_beam.y0(), tmp_beam.z0(),
250  tmp_beam.sigmaZ(), tmp_beam.dxdz(), tmp_beam.dydz(),
251  tmp_beam.BeamWidthX()};
252 
253  reco::BeamSpot tmp_lh = Fit_dres_z_likelihood(tmp_par2);
254 
256 
257  edm::LogWarning("BSFitter") << "Result is non physical. Log-Likelihood fit did not converge." << std::endl;
259  return tmp_lh;
260  }
261  return tmp_lh;
262 
263  } else {
264 
265  throw cms::Exception("LogicError")
266  << "Error in BeamSpotProducer/BSFitter: "
267  << "Illegal fit type, options are likelihood or resolution";
268  }
269  } else {
270 
271  throw cms::Exception("LogicError")
272  << "Error in BeamSpotProducer/BSFitter: "
273  << "Illegal variable type, options are \"z\", \"d\", or \"d*z\"";
274  }
275 
276 
277 }
reco::BeamSpot Fit_d_likelihood(double *inipar)
Definition: BSFitter.cc:657
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
double z0() const
z coordinate
Definition: BeamSpot.h:68
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
reco::BeamSpot Fit_z_likelihood(double *inipar)
Definition: BSFitter.cc:280
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
reco::BeamSpot Fit_dres_z_likelihood(double *inipar)
Definition: BSFitter.cc:888
void setType(BeamType type)
set beam type
Definition: BeamSpot.h:131
double ff_minimum
Definition: BSFitter.h:117
double dydz() const
dydz slope
Definition: BeamSpot.h:84
std::string ffit_variable
Definition: BSFitter.h:115
bool isNotFinite(T x)
Definition: isFinite.h:10
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
std::string ffit_type
Definition: BSFitter.h:114
void d0phi_Init()
Definition: BSFitter.h:68
reco::BeamSpot Fit_ited0phi()
Definition: BSFitter.cc:400
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
reco::BeamSpot Fit_z_chi2(double *inipar)
Definition: BSFitter.cc:349
void Setd0Cut_d0phi(double d0cut)
Definition: BSFitter.cc:639
double covariance(int i, int j) const
(i,j)-th element of error matrix
Definition: BeamSpot.h:112
double y0() const
y coordinate
Definition: BeamSpot.h:66
reco::BeamSpot Fit_d0phi()
Definition: BSFitter.cc:456
reco::BeamSpot Fit_d_z_likelihood(double *inipar, double *error_par)
Definition: BSFitter.cc:785
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double x0() const
x coordinate
Definition: BeamSpot.h:64
reco::BeamSpot BSFitter::Fit_d0phi ( )

Definition at line 456 of file BSFitter.cc.

References funct::abs(), b, align::BeamSpot, funct::cos(), MillePedeFileConverter_cfg::e, g, gen::k, makeMuonMisalignmentScenario::matrix, funct::sin(), mps_update::status, and groupFilesInBlocks::temp.

Referenced by BeamFitter::runAllFitter(), and SetFitVariable().

456  {
457 
458  //LogDebug ("BSFitter") << " we will use " << fBSvector.size() << " tracks.";
459  if (fnthite > 0) edm::LogInfo("BSFitter") << " number of tracks used: " << ftmprow << std::endl;
460  //std::cout << " ftmp = matrix("<<ftmp.GetNrows()<<","<<ftmp.GetNcols()<<")"<<std::endl;
461  //std::cout << " ftmp(0,0)="<<ftmp(0,0)<<std::endl;
462  //std::cout << " ftmp(1,0)="<<ftmp(1,0)<<std::endl;
463  //std::cout << " ftmp(2,0)="<<ftmp(2,0)<<std::endl;
464  //std::cout << " ftmp(3,0)="<<ftmp(3,0)<<std::endl;
465 
466  h1z->Reset();
467 
468 
469  TMatrixD x_result(4,1);
470  TMatrixDSym V_result(4);
471 
472  TMatrixDSym Vint(4);
473  TMatrixD b(4,1);
474 
475  //Double_t weightsum = 0;
476 
477  Vint.Zero();
478  b.Zero();
479 
480  TMatrixD g(4,1);
481  TMatrixDSym temp(4);
482 
483  std::vector<BSTrkParameters>::iterator iparam = fBSvector.begin();
484  ftmprow=0;
485 
486 
487  //edm::LogInfo ("BSFitter") << " test";
488 
489  //std::cout << "BSFitter: fit" << std::endl;
490 
491  for( iparam = fBSvector.begin() ;
492  iparam != fBSvector.end() ; ++iparam) {
493 
494 
495  //if(i->weight2 == 0) continue;
496 
497  //if (ftmprow==0) {
498  //std::cout << "d0=" << iparam->d0() << " sigd0=" << iparam->sigd0()
499  //<< " phi0="<< iparam->phi0() << " z0=" << iparam->z0() << std::endl;
500  //std::cout << "d0phi_d0=" << iparam->d0phi_d0() << " d0phi_chi2="<<iparam->d0phi_chi2() << std::endl;
501  //}
502  g(0,0) = sin(iparam->phi0());
503  g(1,0) = -cos(iparam->phi0());
504  g(2,0) = iparam->z0() * g(0,0);
505  g(3,0) = iparam->z0() * g(1,0);
506 
507 
508  // average transverse beam width
509  double sigmabeam2 = 0.006 * 0.006;
510  if (finputBeamWidth > 0 ) sigmabeam2 = finputBeamWidth * finputBeamWidth;
511  else {
512  //edm::LogWarning("BSFitter") << "using in fit beam width = " << sqrt(sigmabeam2) << std::endl;
513  }
514 
515  //double sigma2 = sigmabeam2 + (iparam->sigd0())* (iparam->sigd0()) / iparam->weight2;
516  // this should be 2*sigmabeam2?
517  double sigma2 = sigmabeam2 + (iparam->sigd0())* (iparam->sigd0());
518 
519  TMatrixD ftmptrans(1,4);
520  ftmptrans = ftmptrans.Transpose(ftmp);
521  TMatrixD dcor = ftmptrans * g;
522  double chi2tmp = (iparam->d0() - dcor(0,0)) * (iparam->d0() - dcor(0,0))/sigma2;
523  (*iparam) = BSTrkParameters(iparam->z0(),iparam->sigz0(),iparam->d0(),iparam->sigd0(),
524  iparam->phi0(), iparam->pt(),dcor(0,0),chi2tmp);
525 
526  bool pass = true;
527  if (fapplyd0cut && fnthite>0 ) {
528  if ( std::abs(iparam->d0() - dcor(0,0)) > fd0cut ) pass = false;
529 
530  }
531  if (fapplychi2cut && fnthite>0 ) {
532  if ( chi2tmp > fchi2cut ) pass = false;
533 
534  }
535 
536  if (pass) {
537  temp.Zero();
538  for(int j = 0 ; j < 4 ; ++j) {
539  for(int k = j ; k < 4 ; ++k) {
540  temp(j,k) += g(j,0) * g(k,0);
541  }
542  }
543 
544 
545  Vint += (temp * (1 / sigma2));
546  b += (iparam->d0() / sigma2 * g);
547  //weightsum += sqrt(i->weight2);
548  ftmprow++;
549  h1z->Fill( iparam->z0() );
550  }
551 
552 
553  }
554  Double_t determinant;
555  TDecompBK bk(Vint);
556  bk.SetTol(1e-11); //FIXME: find a better way to solve x_result
557  if (!bk.Decompose()) {
558  goodfit = false;
559  edm::LogWarning("BSFitter")
560  << "Decomposition failed, matrix singular ?" << std::endl
561  << "condition number = " << bk.Condition() << std::endl;
562  }
563  else {
564  V_result = Vint.InvertFast(&determinant);
565  x_result = V_result * b;
566  }
567  // for(int j = 0 ; j < 4 ; ++j) {
568  // for(int k = 0 ; k < 4 ; ++k) {
569  // std::cout<<"V_result("<<j<<","<<k<<")="<<V_result(j,k)<<std::endl;
570  // }
571  // }
572  // for (int j=0;j<4;++j){
573  // std::cout<<"x_result("<<j<<",0)="<<x_result(j,0)<<std::endl;
574  // }
575  //LogDebug ("BSFitter") << " d0-phi fit done.";
576  //std::cout<< " d0-phi fit done." << std::endl;
577 
578  //Use our own copy for thread safety
579  TF1 fgaus("fgaus","gaus");
580  //returns 0 if OK
581  //auto status = h1z->Fit(&fgaus,"QLM0","",h1z->GetMean() -2.*h1z->GetRMS(),h1z->GetMean() +2.*h1z->GetRMS());
582  auto status = h1z->Fit(&fgaus,"QLN0","",h1z->GetMean() -2.*h1z->GetRMS(),h1z->GetMean() +2.*h1z->GetRMS());
583 
584  //std::cout << "fitted "<< std::endl;
585 
586  //std::cout << "got function" << std::endl;
587  if (status){
588  //edm::LogError("NoBeamSpotFit")<<"gaussian fit failed. no BS d0 fit";
589 
590  return reco::BeamSpot();
591  }
592  double fpar[2] = {fgaus.GetParameter(1), fgaus.GetParameter(2) };
593 
595  // first two parameters
596  for (int j = 0 ; j < 2 ; ++j) {
597  for(int k = j ; k < 2 ; ++k) {
598  matrix(j,k) = V_result(j,k);
599  }
600  }
601  // slope parameters
602  for (int j = 4 ; j < 6 ; ++j) {
603  for(int k = j ; k < 6 ; ++k) {
604  matrix(j,k) = V_result(j-2,k-2);
605  }
606  }
607 
608  // Z0 and sigmaZ
609  matrix(2,2) = fgaus.GetParError(1) * fgaus.GetParError(1);
610  matrix(3,3) = fgaus.GetParError(2) * fgaus.GetParError(2);
611 
612  ftmp = x_result;
613 
614  // x0 and y0 are *not* x,y at z=0, but actually at z=0
615  // to correct for this, we need to translate them to z=z0
616  // using the measured slopes
617  //
618  double x0tmp = x_result(0,0);
619  double y0tmp = x_result(1,0);
620 
621  x0tmp += x_result(2,0)*fpar[0];
622  y0tmp += x_result(3,0)*fpar[0];
623 
624 
626  y0tmp,
627  fpar[0]),
628  fpar[1],
629  x_result(2,0),
630  x_result(3,0),
631  0.,
632  matrix,
633  fbeamtype );
634 
635 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
double fchi2cut
Definition: BSFitter.h:138
TH1F * h1z
Definition: BSFitter.h:146
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
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
bool goodfit
Definition: BSFitter.h:141
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool fapplychi2cut
Definition: BSFitter.h:136
int k[5][pyjets_maxn]
int ftmprow
Definition: BSFitter.h:139
TMatrixD ftmp
Definition: BSFitter.h:134
double b
Definition: hdecay.h:120
bool fapplyd0cut
Definition: BSFitter.h:135
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
int fnthite
Definition: BSFitter.h:140
double fd0cut
Definition: BSFitter.h:137
double finputBeamWidth
Definition: BSFitter.h:145
reco::BeamSpot BSFitter::Fit_d_likelihood ( double *  inipar)

Definition at line 657 of file BSFitter.cc.

References align::BeamSpot, gen::k, and makeMuonMisalignmentScenario::matrix.

Referenced by GetData().

657  {
658 
659 
660  thePDF->SetPDFs("PDFGauss_d");
662 
663  MnUserParameters upar;
664  upar.Add("X0", inipar[0],0.001);
665  upar.Add("Y0", inipar[1],0.001);
666  upar.Add("Z0", 0.,0.001);
667  upar.Add("sigmaZ",0.,0.001);
668  upar.Add("dxdz",inipar[2],0.001);
669  upar.Add("dydz",inipar[3],0.001);
670 
671 
672  MnMigrad migrad(*thePDF, upar);
673 
674  FunctionMinimum fmin = migrad();
675  ff_minimum = fmin.Fval();
676 
678  for (int j = 0 ; j < 6 ; ++j) {
679  for(int k = j ; k < 6 ; ++k) {
680  matrix(j,k) = fmin.Error().Matrix()(j,k);
681  }
682  }
683 
684  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
685  fmin.Parameters().Vec()(1),
686  0.),
687  0.,
688  fmin.Parameters().Vec()(4),
689  fmin.Parameters().Vec()(5),
690  0.,
691  matrix,
692  fbeamtype );
693 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
void SetPDFs(std::string usepdfs)
Definition: BSpdfsFcn.h:32
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
BSpdfsFcn * thePDF
Definition: BSFitter.h:111
double ff_minimum
Definition: BSFitter.h:117
int k[5][pyjets_maxn]
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
void SetData(const std::vector< BSTrkParameters > &a_BSvector)
Definition: BSpdfsFcn.h:26
reco::BeamSpot BSFitter::Fit_d_z_likelihood ( double *  inipar,
double *  error_par 
)

Definition at line 785 of file BSFitter.cc.

References align::BeamSpot, gen::k, and makeMuonMisalignmentScenario::matrix.

Referenced by GetData().

785  {
786 
787  int tracksFailed=0;
788 
789  //estimate first guess of beam width and tame 20% extra of it to start
790  inipar[6]=scanPDF(inipar,tracksFailed,1);
791  error_par[6]=(inipar[6])*0.20;
792 
793 
794  //Here remove the tracks which give low pdf and fill into a new vector
795  //std::cout<<"Size of Old vector = "<<(fBSvector.size())<<std::endl;
796  /* double junk= */ scanPDF(inipar,tracksFailed,2);
797  //std::cout<<"Size of New vector = "<<(fBSvectorBW.size())<<std::endl;
798 
799  //Refill the fBSVector again with new sets of tracks
800  fBSvector.clear();
801  std::vector<BSTrkParameters>::const_iterator iparamBW = fBSvectorBW.begin();
802  for( iparamBW = fBSvectorBW.begin(); iparamBW != fBSvectorBW.end(); ++iparamBW)
803  { fBSvector.push_back(*iparamBW);
804  }
805 
806 
807  thePDF->SetPDFs("PDFGauss_d*PDFGauss_z");
809  MnUserParameters upar;
810 
811  upar.Add("X0", inipar[0],error_par[0]);
812  upar.Add("Y0", inipar[1],error_par[1]);
813  upar.Add("Z0", inipar[2],error_par[2]);
814  upar.Add("sigmaZ",inipar[3],error_par[3]);
815  upar.Add("dxdz",inipar[4],error_par[4]);
816  upar.Add("dydz",inipar[5],error_par[5]);
817  upar.Add("BeamWidthX",inipar[6],error_par[6]);
818 
819 
820  MnMigrad migrad(*thePDF, upar);
821 
822  FunctionMinimum fmin = migrad();
823 
824  // std::cout<<"-----how the fit evoves------"<<std::endl;
825  // std::cout<<fmin<<std::endl;
826 
827  ff_minimum = fmin.Fval();
828 
829 
830  bool ff_nfcn=fmin.HasReachedCallLimit();
831  bool ff_cov=fmin.HasCovariance();
832  bool testing=fmin.IsValid();
833 
834 
835  //Print WARNINGS if minimum did not converged
836  if( ! testing )
837  {
838  edm::LogWarning("BSFitter") <<"===========>>>>>** WARNING: MINUIT DID NOT CONVERGES PROPERLY !!!!!!"<<std::endl;
839  if(ff_nfcn) edm::LogWarning("BSFitter") <<"===========>>>>>** WARNING: No. of Calls Exhausted"<<std::endl;
840  if(!ff_cov) edm::LogWarning("BSFitter") <<"===========>>>>>** WARNING: Covariance did not found"<<std::endl;
841  }
842 
843  edm::LogInfo("BSFitter") <<"The Total # Tracks used for beam width fit = "<<(fBSvectorBW.size())<<std::endl;
844 
845 
846  //Checks after fit is performed
847  double lastIter_pars[7];
848 
849  for(int ip=0;ip<7;ip++){ lastIter_pars[ip]=fmin.Parameters().Vec()(ip);
850  }
851 
852 
853 
854  tracksFailed=0;
855  /* double lastIter_scan= */ scanPDF(lastIter_pars,tracksFailed,2);
856 
857 
858  edm::LogWarning("BSFitter") <<"WARNING: # of tracks which have very low pdf value (pdf_d < 1.0e-05) are = "<<tracksFailed<<std::endl;
859 
860 
861 
862  //std::cout << " eval= " << ff_minimum
863  // << "/n params[0]= " << fmin.Parameters().Vec()(0) << std::endl;
864 
866 
867  for (int j = 0 ; j < 7 ; ++j) {
868  for(int k = j ; k < 7 ; ++k) {
869  matrix(j,k) = fmin.Error().Matrix()(j,k);
870  }
871  }
872 
873 
874  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
875  fmin.Parameters().Vec()(1),
876  fmin.Parameters().Vec()(2)),
877  fmin.Parameters().Vec()(3),
878  fmin.Parameters().Vec()(4),
879  fmin.Parameters().Vec()(5),
880  fmin.Parameters().Vec()(6),
881 
882  matrix,
883  fbeamtype );
884 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
std::vector< BSTrkParameters > fBSvectorBW
Definition: BSFitter.h:126
void SetPDFs(std::string usepdfs)
Definition: BSpdfsFcn.h:32
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
double scanPDF(double *init_pars, int &tracksFailed, int option)
Definition: BSFitter.cc:695
BSpdfsFcn * thePDF
Definition: BSFitter.h:111
double ff_minimum
Definition: BSFitter.h:117
int k[5][pyjets_maxn]
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
void SetData(const std::vector< BSTrkParameters > &a_BSvector)
Definition: BSpdfsFcn.h:26
reco::BeamSpot BSFitter::Fit_dres_z_likelihood ( double *  inipar)

Definition at line 888 of file BSFitter.cc.

References align::BeamSpot, gen::k, makeMuonMisalignmentScenario::matrix, and mathSSE::sqrt().

Referenced by GetData().

888  {
889 
890 
891  thePDF->SetPDFs("PDFGauss_d_resolution*PDFGauss_z");
893 
894  MnUserParameters upar;
895  upar.Add("X0", inipar[0],0.001);
896  upar.Add("Y0", inipar[1],0.001);
897  upar.Add("Z0", inipar[2],0.001);
898  upar.Add("sigmaZ",inipar[3],0.001);
899  upar.Add("dxdz",inipar[4],0.001);
900  upar.Add("dydz",inipar[5],0.001);
901  upar.Add("BeamWidthX",inipar[6],0.0001);
902  upar.Add("c0",0.0010,0.0001);
903  upar.Add("c1",0.0090,0.0001);
904 
905  // fix beam width
906  upar.Fix("BeamWidthX");
907  // number of parameters in fit are 9-1 = 8
908 
909  MnMigrad migrad(*thePDF, upar);
910 
911  FunctionMinimum fmin = migrad();
912  ff_minimum = fmin.Fval();
913 
915 
916  for (int j = 0 ; j < 6 ; ++j) {
917  for(int k = j ; k < 6 ; ++k) {
918  matrix(j,k) = fmin.Error().Matrix()(j,k);
919  }
920  }
921 
922  //std::cout << " fill resolution values" << std::endl;
923  //std::cout << " matrix size= " << fmin.Error().Matrix().size() << std::endl;
924  //std::cout << " vec(6)="<< fmin.Parameters().Vec()(6) << std::endl;
925  //std::cout << " vec(7)="<< fmin.Parameters().Vec()(7) << std::endl;
926 
927  fresolution_c0 = fmin.Parameters().Vec()(6);
928  fresolution_c1 = fmin.Parameters().Vec()(7);
929  fres_c0_err = sqrt( fmin.Error().Matrix()(6,6) );
930  fres_c1_err = sqrt( fmin.Error().Matrix()(7,7) );
931 
932  for (int j = 6 ; j < 8 ; ++j) {
933  for(int k = 6 ; k < 8 ; ++k) {
934  fres_matrix(j-6,k-6) = fmin.Error().Matrix()(j,k);
935  }
936  }
937 
938  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
939  fmin.Parameters().Vec()(1),
940  fmin.Parameters().Vec()(2)),
941  fmin.Parameters().Vec()(3),
942  fmin.Parameters().Vec()(4),
943  fmin.Parameters().Vec()(5),
944  inipar[6],
945  matrix,
946  fbeamtype );
947 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
void SetPDFs(std::string usepdfs)
Definition: BSpdfsFcn.h:32
double fresolution_c1
Definition: BSFitter.h:129
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
BSpdfsFcn * thePDF
Definition: BSFitter.h:111
double fresolution_c0
Definition: BSFitter.h:128
double ff_minimum
Definition: BSFitter.h:117
T sqrt(T t)
Definition: SSEVec.h:18
double fres_c0_err
Definition: BSFitter.h:130
int k[5][pyjets_maxn]
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
reco::BeamSpot::ResCovMatrix fres_matrix
Definition: BSFitter.h:132
double fres_c1_err
Definition: BSFitter.h:131
void SetData(const std::vector< BSTrkParameters > &a_BSvector)
Definition: BSpdfsFcn.h:26
reco::BeamSpot BSFitter::Fit_ited0phi ( )

Definition at line 400 of file BSFitter.cc.

References reco::BeamSpot::Fake, reco::BeamSpot::setType(), reco::BeamSpot::Tracker, and reco::BeamSpot::Unknown.

Referenced by GetData(), and BeamFitter::runAllFitter().

400  {
401 
402  this->d0phi_Init();
403  edm::LogInfo("BSFitter") << "number of total input tracks: " << fBSvector.size() << std::endl;
404 
405  reco::BeamSpot theanswer;
406 
407  if ( (int)fBSvector.size() <= fminNtrks ) {
408  edm::LogWarning("BSFitter") << "need at least " << fminNtrks << " tracks to run beamline fitter." << std::endl;
410  theanswer.setType(fbeamtype);
411  return theanswer;
412  }
413 
414  theanswer = Fit_d0phi(); //get initial ftmp and ftmprow
415  if ( goodfit ) fnthite++;
416  //std::cout << "Initial tempanswer (iteration 0): " << theanswer << std::endl;
417 
418  reco::BeamSpot preanswer = theanswer;
419 
420  while ( goodfit &&
421  ftmprow > fconvergence * fBSvector.size() &&
422  ftmprow > fminNtrks ) {
423 
424  theanswer = Fit_d0phi();
425  fd0cut /= 1.5;
426  fchi2cut /= 1.5;
427  if ( goodfit &&
428  ftmprow > fconvergence * fBSvector.size() &&
429  ftmprow > fminNtrks ) {
430  preanswer = theanswer;
431  //std::cout << "Iteration " << fnthite << ": " << preanswer << std::endl;
432  fnthite++;
433  }
434  }
435  // FIXME: return fit results from previous iteration for both bad fit and for >50% tracks thrown away
436  //std::cout << "The last iteration, theanswer: " << theanswer << std::endl;
437  theanswer = preanswer;
438  //std::cout << "Use previous results from iteration #" << ( fnthite > 0 ? fnthite-1 : 0 ) << std::endl;
439  //if ( fnthite > 1 ) std::cout << theanswer << std::endl;
440 
441  edm::LogInfo("BSFitter") << "Total number of successful iterations = " << ( goodfit ? (fnthite+1) : fnthite ) << std::endl;
442  if (goodfit) {
444  theanswer.setType(fbeamtype);
445  }
446  else {
447  edm::LogWarning("BSFitter") << "Fit doesn't converge!!!" << std::endl;
449  theanswer.setType(fbeamtype);
450  }
451  return theanswer;
452 }
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
double fchi2cut
Definition: BSFitter.h:138
int fminNtrks
Definition: BSFitter.h:144
bool goodfit
Definition: BSFitter.h:141
void setType(BeamType type)
set beam type
Definition: BeamSpot.h:131
void d0phi_Init()
Definition: BSFitter.h:68
int ftmprow
Definition: BSFitter.h:139
double fconvergence
Definition: BSFitter.h:143
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
int fnthite
Definition: BSFitter.h:140
double fd0cut
Definition: BSFitter.h:137
reco::BeamSpot Fit_d0phi()
Definition: BSFitter.cc:456
reco::BeamSpot BSFitter::Fit_z ( std::string  type,
double *  inipar 
)

Referenced by SetFitVariable().

reco::BeamSpot BSFitter::Fit_z_chi2 ( double *  inipar)

Definition at line 349 of file BSFitter.cc.

References align::BeamSpot, and makeMuonMisalignmentScenario::matrix.

Referenced by SetFitVariable().

349  {
350 
351  // N.B. this fit is not performed anymore but now
352  // Z is fitted in the same track set used in the d0-phi fit after
353  // each iteration
354 
355 
356  //std::cout << "Fit_z_chi2() called" << std::endl;
357  // FIXME: include whole tracker z length for the time being
358  // ==> add protection and z0 cut
359  h1z = new TH1F("h1z","z distribution",200,-fMaxZ, fMaxZ);
360 
361  std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
362 
363  // HERE check size of track vector
364 
365  for( iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
366 
367  h1z->Fill( iparam->z0() );
368  //std::cout<<"z0="<<iparam->z0()<<"; sigZ0="<<iparam->sigz0()<<std::endl;
369  }
370 
371  //Use our own copy for thread safety
372  TF1 fgaus("fgaus","gaus");
373  h1z->Fit(&fgaus,"QLMN0");
374  //std::cout << "fitted "<< std::endl;
375 
376  //std::cout << "got function" << std::endl;
377  double fpar[2] = {fgaus.GetParameter(1), fgaus.GetParameter(2) };
378  //std::cout<<"Debug fpar[2] = (" <<fpar[0]<<","<<fpar[1]<<")"<<std::endl;
380  // add matrix values.
381  matrix(2,2) = fgaus.GetParError(1) * fgaus.GetParError(1);
382  matrix(3,3) = fgaus.GetParError(2) * fgaus.GetParError(2);
383 
384  //delete h1z;
385 
387  0.,
388  fpar[0]),
389  fpar[1],
390  0.,
391  0.,
392  0.,
393  matrix,
394  fbeamtype );
395 
396 
397 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
TH1F * h1z
Definition: BSFitter.h:146
double fMaxZ
Definition: BSFitter.h:142
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
reco::BeamSpot BSFitter::Fit_z_likelihood ( double *  inipar)

Definition at line 280 of file BSFitter.cc.

References align::BeamSpot, gen::k, and makeMuonMisalignmentScenario::matrix.

Referenced by SetFitVariable().

280  {
281 
282  //std::cout << "Fit_z(double *) called" << std::endl;
283  //std::cout << "inipar[0]= " << inipar[0] << std::endl;
284  //std::cout << "inipar[1]= " << inipar[1] << std::endl;
285 
286  std::vector<double> par(2,0);
287  std::vector<double> err(2,0);
288 
289  par.push_back(0.0);
290  par.push_back(7.0);
291  err.push_back(0.0001);
292  err.push_back(0.0001);
293  //par[0] = 0.0; err[0] = 0.0;
294  //par[1] = 7.0; err[1] = 0.0;
295 
296  thePDF->SetPDFs("PDFGauss_z");
298  //std::cout << "data loaded"<< std::endl;
299 
300  //FunctionMinimum fmin = theFitter->Minimize(*theGausszFcn, par, err, 1, 500, 0.1);
301  MnUserParameters upar;
302  upar.Add("X0", 0.,0.);
303  upar.Add("Y0", 0.,0.);
304  upar.Add("Z0", inipar[0],0.001);
305  upar.Add("sigmaZ",inipar[1],0.001);
306 
307  MnMigrad migrad(*thePDF, upar);
308 
309  FunctionMinimum fmin = migrad();
310  ff_minimum = fmin.Fval();
311  //std::cout << " eval= " << ff_minimum
312  // << "/n params[0]= " << fmin.Parameters().Vec()(0) << std::endl;
313 
314  /*
315  TMinuit *gmMinuit = new TMinuit(2);
316 
317  //gmMinuit->SetFCN(z_fcn);
318  gmMinuit->SetFCN(myFitz_fcn);
319 
320 
321  int ierflg = 0;
322  double step[2] = {0.001,0.001};
323 
324  for (int i = 0; i<2; i++) {
325  gmMinuit->mnparm(i,fpar_name[i].c_str(),inipar[i],step[i],0,0,ierflg);
326  }
327  gmMinuit->Migrad();
328  */
330 
331  for (int j = 2 ; j < 4 ; ++j) {
332  for(int k = j ; k < 4 ; ++k) {
333  matrix(j,k) = fmin.Error().Matrix()(j,k);
334  }
335  }
336 
338  0.,
339  fmin.Parameters().Vec()(2)),
340  fmin.Parameters().Vec()(3),
341  0.,
342  0.,
343  0.,
344  matrix,
345  fbeamtype );
346 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
void SetPDFs(std::string usepdfs)
Definition: BSpdfsFcn.h:32
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
BSpdfsFcn * thePDF
Definition: BSFitter.h:111
double ff_minimum
Definition: BSFitter.h:117
int k[5][pyjets_maxn]
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
void SetData(const std::vector< BSTrkParameters > &a_BSvector)
Definition: BSpdfsFcn.h:26
int BSFitter::GetAcceptedTrks ( )
inline

Definition at line 67 of file BSFitter.h.

References ftmprow.

67 { return ftmprow; }
int ftmprow
Definition: BSFitter.h:139
std::vector< BSTrkParameters > BSFitter::GetData ( )
inline

Definition at line 75 of file BSFitter.h.

References fBSvector, Fit_d_likelihood(), Fit_d_z_likelihood(), Fit_dres_z_likelihood(), Fit_ited0phi(), TSGForRoadSearch_cfi::option, and scanPDF().

75 { return fBSvector; }
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
double BSFitter::GetMinimum ( )
inline

Definition at line 85 of file BSFitter.h.

References ff_minimum.

85  {
86  return ff_minimum;
87  }
double ff_minimum
Definition: BSFitter.h:117
reco::BeamSpot::ResCovMatrix BSFitter::GetResMatrix ( )
inline

Definition at line 101 of file BSFitter.h.

References fres_matrix.

101  {
102  return fres_matrix;
103  }
reco::BeamSpot::ResCovMatrix fres_matrix
Definition: BSFitter.h:132
double BSFitter::GetResPar0 ( )
inline

Definition at line 88 of file BSFitter.h.

References fresolution_c0.

88  {
89  return fresolution_c0;
90  }
double fresolution_c0
Definition: BSFitter.h:128
double BSFitter::GetResPar0Err ( )
inline

Definition at line 94 of file BSFitter.h.

References fres_c0_err.

94  {
95  return fres_c0_err;
96  }
double fres_c0_err
Definition: BSFitter.h:130
double BSFitter::GetResPar1 ( )
inline

Definition at line 91 of file BSFitter.h.

References fresolution_c1.

91  {
92  return fresolution_c1;
93  }
double fresolution_c1
Definition: BSFitter.h:129
double BSFitter::GetResPar1Err ( )
inline

Definition at line 97 of file BSFitter.h.

References fres_c1_err.

97  {
98  return fres_c1_err;
99  }
double fres_c1_err
Definition: BSFitter.h:131
TH1F* BSFitter::GetVzHisto ( )
inline

Definition at line 105 of file BSFitter.h.

References h1z.

Referenced by BeamFitter::runFitterNoTxt().

105 { return h1z; }
TH1F * h1z
Definition: BSFitter.h:146
double BSFitter::scanPDF ( double *  init_pars,
int &  tracksFailed,
int  option 
)

Definition at line 695 of file BSFitter.cc.

References funct::abs(), funct::cos(), MillePedeFileConverter_cfg::e, JetChargeProducer_cfi::exp, cmsBatch::log, AlCaHLTBitMon_ParallelJobs::p, Pi, funct::sin(), and mathSSE::sqrt().

Referenced by GetData().

695  {
696 
697  if(option==1)init_pars[6]=0.0005; //starting value for any given configuration
698 
699  //local vairables with initial values
700  double fsqrt2pi=0.0;
701  double d_sig=0.0;
702  double d_dprime=0.0;
703  double d_result=0.0;
704  double z_sig=0.0;
705  double z_result=0.0;
706  double function=0.0;
707  double tot_pdf=0.0;
708  double last_minvalue=1.0e+10;
709  double init_bw=-99.99;
710  int iters=0;
711 
712  //used to remove tracks if far away from bs by this
713  double DeltadCut=0.1000;
714  if(init_pars[6]<0.0200){DeltadCut=0.0900; } //worked for high 2.36TeV
715  if(init_pars[6]<0.0100){DeltadCut=0.0700;} //just a guesss for 7 TeV but one should scan for actual values
716 
717 
718 std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
719 
720 
721 if(option==1)iters=500;
722 if(option==2)iters=1;
723 
724 for(int p=0;p<iters;p++){
725 
726  if(iters==500)init_pars[6]+=0.0002;
727  tracksfixed=0;
728 
729 for( iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam)
730  {
731  fsqrt2pi = sqrt(2.* TMath::Pi());
732  d_sig = sqrt(init_pars[6]*init_pars[6] + (iparam->sigd0())*(iparam->sigd0()));
733  d_dprime = iparam->d0() - ( ( (init_pars[0] + iparam->z0()*(init_pars[4]))*sin(iparam->phi0()) )
734  - ( (init_pars[1] + iparam->z0()*(init_pars[5]))*cos(iparam->phi0()) ) );
735 
736  //***Remove tracks before the fit which gives low pdf values to blow up the pdf
737  if(std::abs(d_dprime)<DeltadCut && option==2){ fBSvectorBW.push_back(*iparam);}
738 
739  d_result = (exp(-(d_dprime*d_dprime)/(2.0*d_sig*d_sig)))/(d_sig*fsqrt2pi);
740  z_sig = sqrt(iparam->sigz0() * iparam->sigz0() + init_pars[3]*init_pars[3]);
741  z_result = (exp(-((iparam->z0() - init_pars[2])*(iparam->z0() - init_pars[2]))/(2.0*z_sig*z_sig)))/(z_sig*fsqrt2pi);
742  tot_pdf=z_result*d_result;
743 
744  //for those trcks which gives problems due to very tiny pdf_d values.
745  //Update: This protection will NOT be used with the dprime cut above but still kept here to get
746  // the intial value of beam width reasonably
747  //A warning will appear if there were any tracks with < 10^-5 for pdf_d so that (d-dprime) cut can be lowered
748  if(d_result < 1.0e-05){ tot_pdf=z_result*1.0e-05;
749  //if(option==2)std::cout<<"last Iter d-d' = "<<(std::abs(d_dprime))<<std::endl;
750  tracksfixed++; }
751 
752  function = function + log(tot_pdf);
753 
754 
755  }//loop over tracks
756 
757 
758  function= -2.0*function;
759  if(function<last_minvalue){init_bw=init_pars[6];
760  last_minvalue=function; }
761  function=0.0;
762  }//loop over beam width
763 
764  if(init_bw>0) {
765  init_bw=init_bw+(0.20*init_bw); //start with 20 % more
766 
767  }
768  else{
769 
770  if(option==1){
771  edm::LogWarning("BSFitter")
772  <<"scanPDF:====>>>> WARNING***: The initial guess value of Beam width is negative!!!!!!"<<std::endl
773  <<"scanPDF:====>>>> Assigning beam width a starting value of "<<init_bw<<" cm"<<std::endl;
774  init_bw=0.0200;
775 
776  }
777  }
778 
779 
780  return init_bw;
781 
782 }
const double Pi
std::vector< BSTrkParameters > fBSvectorBW
Definition: BSFitter.h:126
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Double_t fsqrt2pi
Definition: BSFitter.h:123
std::vector< BSTrkParameters > fBSvector
Definition: BSFitter.h:125
void BSFitter::SetChi2Cut_d0phi ( double  chi2cut)

Definition at line 648 of file BSFitter.cc.

References trackAssociatorByChi2_cfi::chi2cut.

Referenced by SetMinimumNTrks().

648  {
649 
650  fapplychi2cut = true;
651 
652  //fBSforCuts = BSfitted;
653  fchi2cut = chi2cut;
654 }
double fchi2cut
Definition: BSFitter.h:138
bool fapplychi2cut
Definition: BSFitter.h:136
void BSFitter::SetConvergence ( double  val)
inline

Definition at line 62 of file BSFitter.h.

References fconvergence, and heppy_batch::val.

Referenced by BeamFitter::runBeamWidthFitter(), and BeamFitter::runFitterNoTxt().

62 { fconvergence = val; }
double fconvergence
Definition: BSFitter.h:143
void BSFitter::Setd0Cut_d0phi ( double  d0cut)

Definition at line 639 of file BSFitter.cc.

Referenced by BeamFitter::runAllFitter(), and SetMinimumNTrks().

639  {
640 
641  fapplyd0cut = true;
642 
643  //fBSforCuts = BSfitted;
644  fd0cut = d0cut;
645 }
bool fapplyd0cut
Definition: BSFitter.h:135
double fd0cut
Definition: BSFitter.h:137
void BSFitter::SetFitType ( std::string  type)
inline

Definition at line 41 of file BSFitter.h.

References ffit_type.

Referenced by BeamFitter::runAllFitter(), and BeamFitter::runBeamWidthFitter().

41  {
42  ffit_type = type;
43  }
type
Definition: HCALResponse.h:21
std::string ffit_type
Definition: BSFitter.h:114
void BSFitter::SetFitVariable ( std::string  name)
inline
void BSFitter::SetInputBeamWidth ( double  val)
inline

Definition at line 66 of file BSFitter.h.

References finputBeamWidth, and heppy_batch::val.

Referenced by BeamFitter::runBeamWidthFitter(), and BeamFitter::runFitterNoTxt().

66 { finputBeamWidth = val; }
double finputBeamWidth
Definition: BSFitter.h:145
void BSFitter::SetMaximumZ ( double  z)
inline

Definition at line 61 of file BSFitter.h.

References fMaxZ, and z.

Referenced by BeamFitter::runBeamWidthFitter(), and BeamFitter::runFitterNoTxt().

61 { fMaxZ = z; }
double fMaxZ
Definition: BSFitter.h:142
void BSFitter::SetMinimumNTrks ( int  n)
inline

Member Data Documentation

bool BSFitter::fapplychi2cut
private

Definition at line 136 of file BSFitter.h.

bool BSFitter::fapplyd0cut
private

Definition at line 135 of file BSFitter.h.

reco::BeamSpot::BeamType BSFitter::fbeamtype
private

Definition at line 113 of file BSFitter.h.

std::vector< BSTrkParameters > BSFitter::fBSvector
private

Definition at line 125 of file BSFitter.h.

Referenced by GetData().

std::vector< BSTrkParameters > BSFitter::fBSvectorBW
private

Definition at line 126 of file BSFitter.h.

double BSFitter::fchi2cut
private

Definition at line 138 of file BSFitter.h.

double BSFitter::fconvergence
private

Definition at line 143 of file BSFitter.h.

Referenced by SetConvergence().

double BSFitter::fd0cut
private

Definition at line 137 of file BSFitter.h.

const int BSFitter::fdim = 7
staticprivate

Definition at line 119 of file BSFitter.h.

double BSFitter::ff_minimum
private

Definition at line 117 of file BSFitter.h.

Referenced by GetMinimum().

std::string BSFitter::ffit_type
private

Definition at line 114 of file BSFitter.h.

Referenced by SetFitType().

std::string BSFitter::ffit_variable
private

Definition at line 115 of file BSFitter.h.

Referenced by SetFitVariable().

double BSFitter::finputBeamWidth
private

Definition at line 145 of file BSFitter.h.

Referenced by SetInputBeamWidth().

double BSFitter::fMaxZ
private

Definition at line 142 of file BSFitter.h.

Referenced by SetMaximumZ().

int BSFitter::fminNtrks
private

Definition at line 144 of file BSFitter.h.

Referenced by SetMinimumNTrks().

int BSFitter::fnthite
private

Definition at line 140 of file BSFitter.h.

Referenced by d0phi_Init().

std::string BSFitter::fpar_name[fdim]
private

Definition at line 121 of file BSFitter.h.

double BSFitter::fres_c0_err
private

Definition at line 130 of file BSFitter.h.

Referenced by GetResPar0Err().

double BSFitter::fres_c1_err
private

Definition at line 131 of file BSFitter.h.

Referenced by GetResPar1Err().

reco::BeamSpot::ResCovMatrix BSFitter::fres_matrix
private

Definition at line 132 of file BSFitter.h.

Referenced by GetResMatrix().

double BSFitter::fresolution_c0
private

Definition at line 128 of file BSFitter.h.

Referenced by GetResPar0().

double BSFitter::fresolution_c1
private

Definition at line 129 of file BSFitter.h.

Referenced by GetResPar1().

Double_t BSFitter::fsqrt2pi
private

Definition at line 123 of file BSFitter.h.

TMatrixD BSFitter::ftmp
private

Definition at line 134 of file BSFitter.h.

Referenced by d0phi_Init().

int BSFitter::ftmprow
private

Definition at line 139 of file BSFitter.h.

Referenced by d0phi_Init(), and GetAcceptedTrks().

bool BSFitter::goodfit
private

Definition at line 141 of file BSFitter.h.

Referenced by d0phi_Init().

TH1F* BSFitter::h1z
private

Definition at line 146 of file BSFitter.h.

Referenced by GetVzHisto().

ROOT::Minuit2::ModularFunctionMinimizer* BSFitter::theFitter
private

Definition at line 109 of file BSFitter.h.

BSpdfsFcn* BSFitter::thePDF
private

Definition at line 111 of file BSFitter.h.