CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ProtonReconstructionAlgorithm.h
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  * Laurent Forthomme
5  ****************************************************************************/
6 
7 #ifndef RecoPPS_ProtonReconstruction_ProtonReconstructionAlgorithm_h
8 #define RecoPPS_ProtonReconstruction_ProtonReconstructionAlgorithm_h
9 
11 
14 
18 
19 #include "TSpline.h"
20 #include "Fit/Fitter.h"
21 
22 #include <unordered_map>
23 
24 //----------------------------------------------------------------------------------------------------
25 
27 public:
28  ProtonReconstructionAlgorithm(bool fit_vtx_y,
29  bool improved_estimate,
30  const std::string &multiRPAlgorithm,
31  unsigned int verbosity);
33 
34  void init(const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions);
35  void release();
36 
39  const LHCInfo &lhcInfo,
40  std::ostream &os) const;
41 
44  const LHCInfo &lhcInfo,
45  std::ostream &os) const;
46 
47 private:
48  unsigned int verbosity_;
49  bool fitVtxY_;
53 
55  struct RPOpticsData {
57  std::shared_ptr<const TSpline3> s_x_d_vs_xi, s_L_x_vs_xi, s_xi_vs_x_d, s_y_d_vs_xi, s_v_y_vs_xi, s_L_y_vs_xi;
58  double x0;
59  double y0;
60  double ch0;
61  double ch1;
62  double la0;
63  double la1;
64  };
65 
67  std::map<unsigned int, RPOpticsData> m_rp_optics_;
68 
71  public:
72  ChiSquareCalculator() = default;
73 
74  double operator()(const double *parameters) const;
75 
77  const std::map<unsigned int, RPOpticsData> *m_rp_optics;
78  };
79 
81  std::unique_ptr<ROOT::Fit::Fitter> fitter_;
82 
84  std::unique_ptr<ChiSquareCalculator> chiSquareCalculator_;
85 
86  static void doLinearFit(const std::vector<double> &vx, const std::vector<double> &vy, double &b, double &a);
87 
88  static double newtonGoalFcn(double xi, double x_N, double x_F, const RPOpticsData &i_N, const RPOpticsData &i_F);
89 };
90 
91 #endif
std::unique_ptr< ROOT::Fit::Fitter > fitter_
fitter object
reco::ForwardProton reconstructFromSingleRP(const CTPPSLocalTrackLiteRef &track, const LHCInfo &lhcInfo, std::ostream &os) const
run proton reconstruction using single-RP strategy
auto const & tracks
cannot be loose
std::map< unsigned int, RPOpticsData > m_rp_optics_
map: RP id –&gt; optics data
const std::map< unsigned int, RPOpticsData > * m_rp_optics
void init(const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions)
std::unique_ptr< ChiSquareCalculator > chiSquareCalculator_
object to calculate chi^2
double ch0
intercept for linear approximation of
ProtonReconstructionAlgorithm(bool fit_vtx_y, bool improved_estimate, const std::string &multiRPAlgorithm, unsigned int verbosity)
enum ProtonReconstructionAlgorithm::@998 multi_rp_algorithm_
static double newtonGoalFcn(double xi, double x_N, double x_F, const RPOpticsData &i_N, const RPOpticsData &i_F)
double la1
slope for linear approximation of
reco::ForwardProton reconstructFromMultiRP(const CTPPSLocalTrackLiteRefVector &tracks, const LHCInfo &lhcInfo, std::ostream &os) const
run proton reconstruction using multiple-RP strategy
static void doLinearFit(const std::vector< double > &vx, const std::vector< double > &vy, double &b, double &a)
double la0
intercept for linear approximation of
Set of optical functions corresponding to one scoring plane along LHC, including splines for interpol...
double b
Definition: hdecay.h:118
double a
Definition: hdecay.h:119
const LHCInterpolatedOpticalFunctionsSet * optics
double ch1
slope for linear approximation of