Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-04 05:19:02

0001 #ifndef BeamSpotProducer_BSFitter_h
0002 #define BeamSpotProducer_BSFitter_h
0003 
0004 /**_________________________________________________________________
0005    class:   BSFitter.h
0006    package: RecoVertex/BeamSpotProducer
0007    
0008 
0009 
0010  author: Francisco Yumiceva, Fermilab (yumiceva@fnal.gov)
0011 
0012 
0013 ________________________________________________________________**/
0014 
0015 // CMS
0016 #include "RecoVertex/BeamSpotProducer/interface/BSpdfsFcn.h"
0017 #include "RecoVertex/BeamSpotProducer/interface/BSTrkParameters.h"
0018 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0019 
0020 // ROOT
0021 #include "TMatrixD.h"
0022 #include "TMath.h"
0023 #include "TH1F.h"
0024 
0025 // C++ standard
0026 #include <vector>
0027 #include <string>
0028 
0029 class BSFitter {
0030 public:
0031   //typedef std::vector <BSTrkParameters> BSTrkCollection;
0032 
0033   BSFitter();
0034   BSFitter(const std::vector<BSTrkParameters> &BSvector);
0035 
0036   virtual ~BSFitter();
0037 
0038   void SetFitType(std::string type) { ffit_type = type; }
0039 
0040   void SetFitVariable(std::string name) { ffit_variable = name; }
0041 
0042   reco::BeamSpot Fit();
0043 
0044   reco::BeamSpot Fit(double *inipar);
0045 
0046   // Fit Z distribution with a gaussian
0047   reco::BeamSpot Fit_z(std::string type, double *inipar);
0048 
0049   reco::BeamSpot Fit_z_chi2(double *inipar);
0050   reco::BeamSpot Fit_z_likelihood(double *inipar);
0051 
0052   // Fit only d0-phi distribution with a chi2
0053   reco::BeamSpot Fit_d0phi();
0054   void SetMaximumZ(double z) { fMaxZ = z; }
0055   void SetConvergence(double val) { fconvergence = val; }
0056   void SetMinimumNTrks(int n) { fminNtrks = n; }
0057   void Setd0Cut_d0phi(double d0cut);
0058   void SetChi2Cut_d0phi(double chi2cut);
0059   void SetInputBeamWidth(double val) { finputBeamWidth = val; }
0060   int GetAcceptedTrks() { return ftmprow; }
0061   void d0phi_Init() {
0062     ftmprow = 0;
0063     ftmp.ResizeTo(4, 1);
0064     ftmp.Zero();
0065     fnthite = 0;
0066     goodfit = true;
0067   }
0068   std::vector<BSTrkParameters> GetData() { return fBSvector; }
0069 
0070   reco::BeamSpot Fit_ited0phi();
0071 
0072   reco::BeamSpot Fit_d_likelihood(double *inipar);
0073   reco::BeamSpot Fit_d_z_likelihood(double *inipar, double *error_par);
0074   reco::BeamSpot Fit_dres_z_likelihood(double *inipar);
0075 
0076   double scanPDF(double *init_pars, int &tracksFailed, int option);
0077 
0078   double GetMinimum() { return ff_minimum; }
0079   double GetResPar0() { return fresolution_c0; }
0080   double GetResPar1() { return fresolution_c1; }
0081   double GetResPar0Err() { return fres_c0_err; }
0082   double GetResPar1Err() { return fres_c1_err; }
0083 
0084   reco::BeamSpot::ResCovMatrix GetResMatrix() { return fres_matrix; }
0085 
0086   TH1F *GetVzHisto() { return h1z; }
0087 
0088 private:
0089   //BSzFcn* theGausszFcn;
0090   BSpdfsFcn *thePDF;
0091 
0092   reco::BeamSpot::BeamType fbeamtype;
0093   std::string ffit_type;
0094   std::string ffit_variable;
0095 
0096   double ff_minimum;
0097 
0098   static const int fdim = 7;
0099 
0100   std::string fpar_name[fdim];
0101 
0102   Double_t fsqrt2pi;
0103 
0104   std::vector<BSTrkParameters> fBSvector;
0105   std::vector<BSTrkParameters> fBSvectorBW;
0106 
0107   double fresolution_c0;
0108   double fresolution_c1;
0109   double fres_c0_err;
0110   double fres_c1_err;
0111   reco::BeamSpot::ResCovMatrix fres_matrix;
0112   //reco::BeamSpot fBSforCuts;
0113   TMatrixD ftmp;
0114   bool fapplyd0cut;
0115   bool fapplychi2cut;
0116   double fd0cut;
0117   double fchi2cut;
0118   int ftmprow;
0119   int fnthite;
0120   bool goodfit;
0121   double fMaxZ;
0122   double fconvergence;
0123   int fminNtrks;
0124   double finputBeamWidth;
0125   TH1F *h1z;
0126 };
0127 
0128 #endif