CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 42 of file BSFitter.cc.

References reco::BeamSpot::Unknown.

42  {
44 }
reco::BeamSpot::BeamType fbeamtype
Definition: BSFitter.h:113
BSFitter::BSFitter ( const std::vector< BSTrkParameters > &  BSvector)

Definition at line 47 of file BSFitter.cc.

References Pi, and mathSSE::sqrt().

47  {
48 
49  ffit_type = "default";
50  ffit_variable = "default";
51 
52  fBSvector = BSvector;
53 
54  fsqrt2pi = sqrt(2.* TMath::Pi());
55 
56  fpar_name[0] = "z0 ";
57  fpar_name[1] = "SigmaZ0 ";
58  fpar_name[2] = "X0 ";
59  fpar_name[3] = "Y0 ";
60  fpar_name[4] = "dxdz ";
61  fpar_name[5] = "dydz ";
62  fpar_name[6] = "SigmaBeam ";
63 
64  //if (theGausszFcn == 0 ) {
65  thePDF = new BSpdfsFcn();
66 
67 
68 //}
69  //if (theFitter == 0 ) {
70 
71  theFitter = new VariableMetricMinimizer();
72 
73  //}
74 
75  fapplyd0cut = false;
76  fapplychi2cut = false;
77  ftmprow = 0;
78  ftmp.ResizeTo(4,1);
79  ftmp.Zero();
80  fnthite=0;
81  fMaxZ = 50.; //cm
82  fconvergence = 0.5; // stop fit when 50% of the input collection has been removed.
83  fminNtrks = 100;
84  finputBeamWidth = -1; // no input
85 
86  h1z = new TH1F("h1z","z distribution",200,-fMaxZ, fMaxZ);
87 
88 }
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:48
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 91 of file BSFitter.cc.

92 {
93  //delete fBSvector;
94  delete thePDF;
95  delete theFitter;
96 }
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 100 of file BSFitter.cc.

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

100  {
101 
102  return this->Fit(0);
103 
104 }
reco::BeamSpot Fit()
Definition: BSFitter.cc:100
reco::BeamSpot BSFitter::Fit ( double *  inipar = 0)

Definition at line 107 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().

107  {
109  if ( ffit_variable == "z" ) {
110 
111  if ( ffit_type == "chi2" ) {
112 
113  return Fit_z_chi2(inipar);
114 
115  } else if ( ffit_type == "likelihood" ) {
116 
117  return Fit_z_likelihood(inipar);
118 
119  } else if ( ffit_type == "combined" ) {
120 
121  reco::BeamSpot tmp_beamspot = Fit_z_chi2(inipar);
122  double tmp_par[2] = {tmp_beamspot.z0(), tmp_beamspot.sigmaZ()};
123  return Fit_z_likelihood(tmp_par);
124 
125  } else {
126 
127  throw cms::Exception("LogicError")
128  << "Error in BeamSpotProducer/BSFitter: "
129  << "Illegal fit type, options are chi2,likelihood or combined(ie. first chi2 then likelihood)";
130 
131  }
132  } else if ( ffit_variable == "d" ) {
133 
134  if ( ffit_type == "d0phi" ) {
135  this->d0phi_Init();
136  return Fit_d0phi();
137 
138  } else if ( ffit_type == "likelihood" ) {
139 
140  return Fit_d_likelihood(inipar);
141 
142  } else if ( ffit_type == "combined" ) {
143 
144  this->d0phi_Init();
145  reco::BeamSpot tmp_beamspot = Fit_d0phi();
146  double tmp_par[4] = {tmp_beamspot.x0(), tmp_beamspot.y0(), tmp_beamspot.dxdz(), tmp_beamspot.dydz()};
147  return Fit_d_likelihood(tmp_par);
148 
149  } else {
150  throw cms::Exception("LogicError")
151  << "Error in BeamSpotProducer/BSFitter: "
152  << "Illegal fit type, options are d0phi, likelihood or combined";
153  }
154  } else if ( ffit_variable == "d*z" || ffit_variable == "default" ) {
155 
156  if ( ffit_type == "likelihood" || ffit_type == "default" ) {
157 
159  // we are now fitting Z inside d0phi fitter
160  // first fit z distribution using a chi2 fit
161  //reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
162  //for (int j = 2 ; j < 4 ; ++j) {
163  //for(int k = j ; k < 4 ; ++k) {
164  // matrix(j,k) = tmp_z.covariance()(j,k);
165  //}
166  //}
167 
168  // use d0-phi algorithm to extract transverse position
169  this->d0phi_Init();
170  //reco::BeamSpot tmp_d0phi= Fit_d0phi(); // change to iterative procedure:
171  this->Setd0Cut_d0phi(4.0);
172  reco::BeamSpot tmp_d0phi= Fit_ited0phi();
173 
174  //for (int j = 0 ; j < 2 ; ++j) {
175  // for(int k = j ; k < 2 ; ++k) {
176  // matrix(j,k) = tmp_d0phi.covariance()(j,k);
177  //}
178  //}
179  // slopes
180  //for (int j = 4 ; j < 6 ; ++j) {
181  // for(int k = j ; k < 6 ; ++k) {
182  // matrix(j,k) = tmp_d0phi.covariance()(j,k);
183  // }
184  //}
185 
186 
187  // put everything into one object
188  reco::BeamSpot spot(reco::BeamSpot::Point(tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_d0phi.z0()),
189  tmp_d0phi.sigmaZ(),
190  tmp_d0phi.dxdz(),
191  tmp_d0phi.dydz(),
192  0.,
193  tmp_d0phi.covariance(),
194  fbeamtype );
195 
196 
197 
198  //reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
199 
200  //reco::BeamSpot tmp_d0phi = Fit_d0phi();
201 
202  // log-likelihood fit
203  if (ffit_type == "likelihood") {
204  double tmp_par[7] = {tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_d0phi.z0(),
205  tmp_d0phi.sigmaZ(), tmp_d0phi.dxdz(), tmp_d0phi.dydz(),0.0};
206 
207  double tmp_error_par[7];
208  for(int s=0;s<6;s++){ tmp_error_par[s] = pow( tmp_d0phi.covariance()(s,s),0.5);}
209  tmp_error_par[6]=0.0;
210 
211  reco::BeamSpot tmp_lh = Fit_d_z_likelihood(tmp_par,tmp_error_par);
212 
214  edm::LogWarning("BSFitter") << "BSFitter: Result is non physical. Log-Likelihood fit to extract beam width did not converge." << std::endl;
216  return tmp_lh;
217  }
218  return tmp_lh;
219 
220  } else {
221 
222  edm::LogInfo("BSFitter") << "default track-based fit does not extract beam width." << std::endl;
223  return spot;
224  }
225 
226 
227  } else if ( ffit_type == "resolution" ) {
228 
229  reco::BeamSpot tmp_z = Fit_z_chi2(inipar);
230  this->d0phi_Init();
231  reco::BeamSpot tmp_d0phi = Fit_d0phi();
232 
233  double tmp_par[7] = {tmp_d0phi.x0(), tmp_d0phi.y0(), tmp_z.z0(),
234  tmp_z.sigmaZ(), tmp_d0phi.dxdz(), tmp_d0phi.dydz(),0.0};
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;
238 
239  reco::BeamSpot tmp_beam = Fit_d_z_likelihood(tmp_par,tmp_error_par);
240 
241  double tmp_par2[7] = {tmp_beam.x0(), tmp_beam.y0(), tmp_beam.z0(),
242  tmp_beam.sigmaZ(), tmp_beam.dxdz(), tmp_beam.dydz(),
243  tmp_beam.BeamWidthX()};
244 
245  reco::BeamSpot tmp_lh = Fit_dres_z_likelihood(tmp_par2);
246 
248 
249  edm::LogWarning("BSFitter") << "Result is non physical. Log-Likelihood fit did not converge." << std::endl;
251  return tmp_lh;
252  }
253  return tmp_lh;
254 
255  } else {
256 
257  throw cms::Exception("LogicError")
258  << "Error in BeamSpotProducer/BSFitter: "
259  << "Illegal fit type, options are likelihood or resolution";
260  }
261  } else {
262 
263  throw cms::Exception("LogicError")
264  << "Error in BeamSpotProducer/BSFitter: "
265  << "Illegal variable type, options are \"z\", \"d\", or \"d*z\"";
266  }
267 
268 
269 }
reco::BeamSpot Fit_d_likelihood(double *inipar)
Definition: BSFitter.cc:649
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:272
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
reco::BeamSpot Fit_dres_z_likelihood(double *inipar)
Definition: BSFitter.cc:881
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:392
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
reco::BeamSpot Fit_z_chi2(double *inipar)
Definition: BSFitter.cc:341
void Setd0Cut_d0phi(double d0cut)
Definition: BSFitter.cc:631
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:448
reco::BeamSpot Fit_d_z_likelihood(double *inipar, double *error_par)
Definition: BSFitter.cc:778
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 448 of file BSFitter.cc.

References funct::abs(), b, align::BeamSpot, funct::cos(), alignCSCRings::e, g, j, relval_steps::k, makeMuonMisalignmentScenario::matrix, funct::sin(), ntuplemaker::status, and groupFilesInBlocks::temp.

Referenced by BeamFitter::runAllFitter().

448  {
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 {
504  //edm::LogWarning("BSFitter") << "using in fit beam width = " << sqrt(sigmabeam2) << std::endl;
505  }
506 
507  //double sigma2 = sigmabeam2 + (iparam->sigd0())* (iparam->sigd0()) / iparam->weight2;
508  // this should be 2*sigmabeam2?
509  double sigma2 = sigmabeam2 + (iparam->sigd0())* (iparam->sigd0());
510 
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);
517 
518  bool pass = true;
519  if (fapplyd0cut && fnthite>0 ) {
520  if ( std::abs(iparam->d0() - dcor(0,0)) > fd0cut ) pass = false;
521 
522  }
523  if (fapplychi2cut && fnthite>0 ) {
524  if ( chi2tmp > fchi2cut ) pass = false;
525 
526  }
527 
528  if (pass) {
529  temp.Zero();
530  for(int j = 0 ; j < 4 ; ++j) {
531  for(int k = j ; k < 4 ; ++k) {
532  temp(j,k) += g(j,0) * g(k,0);
533  }
534  }
535 
536 
537  Vint += (temp * (1 / sigma2));
538  b += (iparam->d0() / sigma2 * g);
539  //weightsum += sqrt(i->weight2);
540  ftmprow++;
541  h1z->Fill( iparam->z0() );
542  }
543 
544 
545  }
546  Double_t determinant;
547  TDecompBK bk(Vint);
548  bk.SetTol(1e-11); //FIXME: find a better way to solve x_result
549  if (!bk.Decompose()) {
550  goodfit = false;
551  edm::LogWarning("BSFitter")
552  << "Decomposition failed, matrix singular ?" << std::endl
553  << "condition number = " << bk.Condition() << std::endl;
554  }
555  else {
556  V_result = Vint.InvertFast(&determinant);
557  x_result = V_result * b;
558  }
559  // for(int j = 0 ; j < 4 ; ++j) {
560  // for(int k = 0 ; k < 4 ; ++k) {
561  // std::cout<<"V_result("<<j<<","<<k<<")="<<V_result(j,k)<<std::endl;
562  // }
563  // }
564  // for (int j=0;j<4;++j){
565  // std::cout<<"x_result("<<j<<",0)="<<x_result(j,0)<<std::endl;
566  // }
567  //LogDebug ("BSFitter") << " d0-phi fit done.";
568  //std::cout<< " d0-phi fit done." << std::endl;
569 
570  //Use our own copy for thread safety
571  TF1 fgaus("fgaus","gaus");
572  //returns 0 if OK
573  //auto status = h1z->Fit(&fgaus,"QLM0","",h1z->GetMean() -2.*h1z->GetRMS(),h1z->GetMean() +2.*h1z->GetRMS());
574  auto status = h1z->Fit(&fgaus,"QL0","",h1z->GetMean() -2.*h1z->GetRMS(),h1z->GetMean() +2.*h1z->GetRMS());
575 
576  //std::cout << "fitted "<< std::endl;
577 
578  //std::cout << "got function" << std::endl;
579  if (status){
580  //edm::LogError("NoBeamSpotFit")<<"gaussian fit failed. no BS d0 fit";
581 
582  return reco::BeamSpot();
583  }
584  double fpar[2] = {fgaus.GetParameter(1), fgaus.GetParameter(2) };
585 
587  // first two parameters
588  for (int j = 0 ; j < 2 ; ++j) {
589  for(int k = j ; k < 2 ; ++k) {
590  matrix(j,k) = V_result(j,k);
591  }
592  }
593  // slope parameters
594  for (int j = 4 ; j < 6 ; ++j) {
595  for(int k = j ; k < 6 ; ++k) {
596  matrix(j,k) = V_result(j-2,k-2);
597  }
598  }
599 
600  // Z0 and sigmaZ
601  matrix(2,2) = fgaus.GetParError(1) * fgaus.GetParError(1);
602  matrix(3,3) = fgaus.GetParError(2) * fgaus.GetParError(2);
603 
604  ftmp = x_result;
605 
606  // x0 and y0 are *not* x,y at z=0, but actually at z=0
607  // to correct for this, we need to translate them to z=z0
608  // using the measured slopes
609  //
610  double x0tmp = x_result(0,0);
611  double y0tmp = x_result(1,0);
612 
613  x0tmp += x_result(2,0)*fpar[0];
614  y0tmp += x_result(3,0)*fpar[0];
615 
616 
618  y0tmp,
619  fpar[0]),
620  fpar[1],
621  x_result(2,0),
622  x_result(3,0),
623  0.,
624  matrix,
625  fbeamtype );
626 
627 }
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
int j
Definition: DBlmapReader.cc:9
bool fapplychi2cut
Definition: BSFitter.h:136
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
tuple status
Definition: ntuplemaker.py:245
double finputBeamWidth
Definition: BSFitter.h:145
reco::BeamSpot BSFitter::Fit_d_likelihood ( double *  inipar)

Definition at line 649 of file BSFitter.cc.

References align::BeamSpot, j, relval_steps::k, and makeMuonMisalignmentScenario::matrix.

649  {
650 
651 
652  thePDF->SetPDFs("PDFGauss_d");
654 
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);
662 
663 
664  MnMigrad migrad(*thePDF, upar);
665 
666  FunctionMinimum fmin = migrad();
667  ff_minimum = fmin.Fval();
668 
670  for (int j = 0 ; j < 6 ; ++j) {
671  for(int k = j ; k < 6 ; ++k) {
672  matrix(j,k) = fmin.Error().Matrix()(j,k);
673  }
674  }
675 
676  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
677  fmin.Parameters().Vec()(1),
678  0.),
679  0.,
680  fmin.Parameters().Vec()(4),
681  fmin.Parameters().Vec()(5),
682  0.,
683  matrix,
684  fbeamtype );
685 }
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 j
Definition: DBlmapReader.cc:9
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 778 of file BSFitter.cc.

References align::BeamSpot, j, relval_steps::k, and makeMuonMisalignmentScenario::matrix.

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

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

881  {
882 
883 
884  thePDF->SetPDFs("PDFGauss_d_resolution*PDFGauss_z");
886 
887  MnUserParameters upar;
888  upar.Add("X0", inipar[0],0.001);
889  upar.Add("Y0", inipar[1],0.001);
890  upar.Add("Z0", inipar[2],0.001);
891  upar.Add("sigmaZ",inipar[3],0.001);
892  upar.Add("dxdz",inipar[4],0.001);
893  upar.Add("dydz",inipar[5],0.001);
894  upar.Add("BeamWidthX",inipar[6],0.0001);
895  upar.Add("c0",0.0010,0.0001);
896  upar.Add("c1",0.0090,0.0001);
897 
898  // fix beam width
899  upar.Fix("BeamWidthX");
900  // number of parameters in fit are 9-1 = 8
901 
902  MnMigrad migrad(*thePDF, upar);
903 
904  FunctionMinimum fmin = migrad();
905  ff_minimum = fmin.Fval();
906 
908 
909  for (int j = 0 ; j < 6 ; ++j) {
910  for(int k = j ; k < 6 ; ++k) {
911  matrix(j,k) = fmin.Error().Matrix()(j,k);
912  }
913  }
914 
915  //std::cout << " fill resolution values" << std::endl;
916  //std::cout << " matrix size= " << fmin.Error().Matrix().size() << std::endl;
917  //std::cout << " vec(6)="<< fmin.Parameters().Vec()(6) << std::endl;
918  //std::cout << " vec(7)="<< fmin.Parameters().Vec()(7) << std::endl;
919 
920  fresolution_c0 = fmin.Parameters().Vec()(6);
921  fresolution_c1 = fmin.Parameters().Vec()(7);
922  fres_c0_err = sqrt( fmin.Error().Matrix()(6,6) );
923  fres_c1_err = sqrt( fmin.Error().Matrix()(7,7) );
924 
925  for (int j = 6 ; j < 8 ; ++j) {
926  for(int k = 6 ; k < 8 ; ++k) {
927  fres_matrix(j-6,k-6) = fmin.Error().Matrix()(j,k);
928  }
929  }
930 
931  return reco::BeamSpot( reco::BeamSpot::Point(fmin.Parameters().Vec()(0),
932  fmin.Parameters().Vec()(1),
933  fmin.Parameters().Vec()(2)),
934  fmin.Parameters().Vec()(3),
935  fmin.Parameters().Vec()(4),
936  fmin.Parameters().Vec()(5),
937  inipar[6],
938  matrix,
939  fbeamtype );
940 }
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:48
double fres_c0_err
Definition: BSFitter.h:130
int j
Definition: DBlmapReader.cc:9
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 392 of file BSFitter.cc.

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

Referenced by BeamFitter::runAllFitter().

392  {
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;
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) {
436  theanswer.setType(fbeamtype);
437  }
438  else {
439  edm::LogWarning("BSFitter") << "Fit doesn't converge!!!" << std::endl;
441  theanswer.setType(fbeamtype);
442  }
443  return theanswer;
444 }
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:448
reco::BeamSpot BSFitter::Fit_z ( std::string  type,
double *  inipar 
)
reco::BeamSpot BSFitter::Fit_z_chi2 ( double *  inipar)

Definition at line 341 of file BSFitter.cc.

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

341  {
342 
343  // N.B. this fit is not performed anymore but now
344  // Z is fitted in the same track set used in the d0-phi fit after
345  // each iteration
346 
347 
348  //std::cout << "Fit_z_chi2() called" << std::endl;
349  // FIXME: include whole tracker z length for the time being
350  // ==> add protection and z0 cut
351  h1z = new TH1F("h1z","z distribution",200,-fMaxZ, fMaxZ);
352 
353  std::vector<BSTrkParameters>::const_iterator iparam = fBSvector.begin();
354 
355  // HERE check size of track vector
356 
357  for( iparam = fBSvector.begin(); iparam != fBSvector.end(); ++iparam) {
358 
359  h1z->Fill( iparam->z0() );
360  //std::cout<<"z0="<<iparam->z0()<<"; sigZ0="<<iparam->sigz0()<<std::endl;
361  }
362 
363  //Use our own copy for thread safety
364  TF1 fgaus("fgaus","gaus");
365  h1z->Fit(&fgaus,"QLM0");
366  //std::cout << "fitted "<< std::endl;
367 
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 }
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 272 of file BSFitter.cc.

References align::BeamSpot, j, relval_steps::k, and makeMuonMisalignmentScenario::matrix.

272  {
273 
274  //std::cout << "Fit_z(double *) called" << std::endl;
275  //std::cout << "inipar[0]= " << inipar[0] << std::endl;
276  //std::cout << "inipar[1]= " << inipar[1] << std::endl;
277 
278  std::vector<double> par(2,0);
279  std::vector<double> err(2,0);
280 
281  par.push_back(0.0);
282  par.push_back(7.0);
283  err.push_back(0.0001);
284  err.push_back(0.0001);
285  //par[0] = 0.0; err[0] = 0.0;
286  //par[1] = 7.0; err[1] = 0.0;
287 
288  thePDF->SetPDFs("PDFGauss_z");
290  //std::cout << "data loaded"<< std::endl;
291 
292  //FunctionMinimum fmin = theFitter->Minimize(*theGausszFcn, par, err, 1, 500, 0.1);
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);
298 
299  MnMigrad migrad(*thePDF, upar);
300 
301  FunctionMinimum fmin = migrad();
302  ff_minimum = fmin.Fval();
303  //std::cout << " eval= " << ff_minimum
304  // << "/n params[0]= " << fmin.Parameters().Vec()(0) << std::endl;
305 
306  /*
307  TMinuit *gmMinuit = new TMinuit(2);
308 
309  //gmMinuit->SetFCN(z_fcn);
310  gmMinuit->SetFCN(myFitz_fcn);
311 
312 
313  int ierflg = 0;
314  double step[2] = {0.001,0.001};
315 
316  for (int i = 0; i<2; i++) {
317  gmMinuit->mnparm(i,fpar_name[i].c_str(),inipar[i],step[i],0,0,ierflg);
318  }
319  gmMinuit->Migrad();
320  */
322 
323  for (int j = 2 ; j < 4 ; ++j) {
324  for(int k = j ; k < 4 ; ++k) {
325  matrix(j,k) = fmin.Error().Matrix()(j,k);
326  }
327  }
328 
330  0.,
331  fmin.Parameters().Vec()(2)),
332  fmin.Parameters().Vec()(3),
333  0.,
334  0.,
335  0.,
336  matrix,
337  fbeamtype );
338 }
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 j
Definition: DBlmapReader.cc:9
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.

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 687 of file BSFitter.cc.

References funct::abs(), funct::cos(), alignCSCRings::e, create_public_lumi_plots::exp, create_public_lumi_plots::log, AlCaHLTBitMon_ParallelJobs::p, Pi, funct::sin(), and mathSSE::sqrt().

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

640  {
641 
642  fapplychi2cut = true;
643 
644  //fBSforCuts = BSfitted;
645  fchi2cut = chi2cut;
646 }
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.

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 631 of file BSFitter.cc.

Referenced by BeamFitter::runAllFitter().

631  {
632 
633  fapplyd0cut = true;
634 
635  //fBSforCuts = BSfitted;
636  fd0cut = d0cut;
637 }
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

Definition at line 45 of file BSFitter.h.

References ffit_variable, and mergeVDriftHistosByStation::name.

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

45  {
47  }
std::string ffit_variable
Definition: BSFitter.h:115
void BSFitter::SetInputBeamWidth ( double  val)
inline

Definition at line 66 of file BSFitter.h.

References finputBeamWidth.

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

Definition at line 63 of file BSFitter.h.

References fminNtrks, and gen::n.

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

63 { fminNtrks = n; }
int fminNtrks
Definition: BSFitter.h:144

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.