CMS 3D CMS Logo

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 
17 
18 #include "TSpline.h"
19 #include "Fit/Fitter.h"
20 
21 #include <unordered_map>
22 
23 //----------------------------------------------------------------------------------------------------
24 
26 public:
27  ProtonReconstructionAlgorithm(bool fit_vtx_y,
28  bool improved_estimate,
29  const std::string &multiRPAlgorithm,
30  unsigned int verbosity);
32 
34  void release();
35 
38  const float energy,
39  std::ostream &os) const;
40 
43  const float energy,
44  std::ostream &os) const;
45 
46 private:
47  unsigned int verbosity_;
48  bool fitVtxY_;
52 
54  struct RPOpticsData {
56  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;
57  double x0;
58  double y0;
59  double ch0;
60  double ch1;
61  double la0;
62  double la1;
63  };
64 
66  std::map<unsigned int, RPOpticsData> m_rp_optics_;
67 
70  public:
71  ChiSquareCalculator() = default;
72 
73  double operator()(const double *parameters) const;
74 
76  const std::map<unsigned int, RPOpticsData> *m_rp_optics;
77  };
78 
80  std::unique_ptr<ROOT::Fit::Fitter> fitter_;
81 
83  std::unique_ptr<ChiSquareCalculator> chiSquareCalculator_;
84 
85  static void doLinearFit(const std::vector<double> &vx, const std::vector<double> &vy, double &b, double &a);
86 
87  static double newtonGoalFcn(double xi, double x_N, double x_F, const RPOpticsData &i_N, const RPOpticsData &i_F);
88 };
89 
90 #endif
std::unique_ptr< ROOT::Fit::Fitter > fitter_
fitter object
reco::ForwardProton reconstructFromSingleRP(const CTPPSLocalTrackLiteRef &track, const float energy, std::ostream &os) const
run proton reconstruction using single-RP strategy
std::map< unsigned int, RPOpticsData > m_rp_optics_
map: RP id –> 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)
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
static void doLinearFit(const std::vector< double > &vx, const std::vector< double > &vy, double &b, double &a)
reco::ForwardProton reconstructFromMultiRP(const CTPPSLocalTrackLiteRefVector &tracks, const float energy, std::ostream &os) const
run proton reconstruction using multiple-RP strategy
double la0
intercept for linear approximation of
const int verbosity
Set of optical functions corresponding to one scoring plane along LHC, including splines for interpol...
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
const LHCInterpolatedOpticalFunctionsSet * optics
double ch1
slope for linear approximation of
enum ProtonReconstructionAlgorithm::@991 multi_rp_algorithm_