File indexing completed on 2024-04-06 12:31:09
0001 #ifndef TotemRPProtonTransportParametrization_LHC_OPTICS_APPROXIMATOR_H
0002 #define TotemRPProtonTransportParametrization_LHC_OPTICS_APPROXIMATOR_H
0003
0004 #include <string>
0005 #include "TNamed.h"
0006 #include "TTree.h"
0007 #include "TH1D.h"
0008 #include "TH2D.h"
0009 #include <memory>
0010 #include "TMatrixD.h"
0011
0012 #include "SimTransport/TotemRPProtonTransportParametrization/interface/TMultiDimFet.h"
0013
0014 struct MadKinematicDescriptor {
0015 double x;
0016 double theta_x;
0017 double y;
0018 double theta_y;
0019 double ksi;
0020 };
0021
0022 class LHCApertureApproximator;
0023
0024
0025
0026
0027
0028
0029 class LHCOpticsApproximator : public TNamed {
0030 public:
0031 LHCOpticsApproximator();
0032
0033 LHCOpticsApproximator(std::string name,
0034 std::string title,
0035 TMultiDimFet::EMDFPolyType polynom_type,
0036 std::string beam_direction,
0037 double nominal_beam_momentum);
0038 LHCOpticsApproximator(const LHCOpticsApproximator &org);
0039 const LHCOpticsApproximator &operator=(const LHCOpticsApproximator &org);
0040
0041 enum polynomials_selection { AUTOMATIC, PREDEFINED };
0042 enum beam_type { lhcb1, lhcb2 };
0043 void Train(TTree *inp_tree,
0044 std::string data_prefix = std::string("def"),
0045 polynomials_selection mode = PREDEFINED,
0046 int max_degree_x = 10,
0047 int max_degree_tx = 10,
0048 int max_degree_y = 10,
0049 int max_degree_ty = 10,
0050 bool common_terms = false,
0051 double *prec = nullptr);
0052 void Test(TTree *inp_tree,
0053 TFile *f_out,
0054 std::string data_prefix = std::string("def"),
0055 std::string base_out_dir = std::string(""));
0056 void TestAperture(TTree *in_tree,
0057 TTree *out_tree);
0058
0059 double ParameterOutOfRangePenalty(double par_m[], bool invert_beam_coord_sytems = true) const;
0060
0061
0062
0063
0064
0065
0066 bool Transport(const double *in,
0067 double *out,
0068 bool check_apertures = false,
0069 bool invert_beam_coord_sytems = true) const;
0070 bool Transport(const MadKinematicDescriptor *in,
0071 MadKinematicDescriptor *out,
0072 bool check_apertures = false,
0073 bool invert_beam_coord_sytems = true) const;
0074
0075
0076
0077
0078
0079
0080 bool Transport2D(const double *in,
0081 double *out,
0082 bool check_apertures = false,
0083 bool invert_beam_coord_sytems = true) const;
0084
0085 bool Transport_m_GeV(double in_pos[3],
0086 double in_momentum[3],
0087 double out_pos[3],
0088 double out_momentum[3],
0089 bool check_apertures,
0090 double z2_z1_dist) const;
0091
0092 void PrintInputRange();
0093 bool CheckInputRange(const double *in, bool invert_beam_coord_sytems = true) const;
0094 void AddRectEllipseAperture(
0095 const LHCOpticsApproximator &in, double rect_x, double rect_y, double r_el_x, double r_el_y);
0096 void PrintOpticalFunctions();
0097 void PrintCoordinateOpticalFunctions(TMultiDimFet ¶metrization,
0098 const std::string &coord_name,
0099 const std::vector<std::string> &input_vars);
0100
0101
0102
0103
0104
0105
0106
0107
0108 void GetLineariasedTransportMatrixX(double mad_init_x,
0109 double mad_init_thx,
0110 double mad_init_y,
0111 double mad_init_thy,
0112 double mad_init_xi,
0113 TMatrixD &tr_matrix,
0114 double d_mad_x = 10e-6,
0115 double d_mad_thx = 10e-6);
0116
0117
0118
0119
0120
0121
0122
0123
0124 void GetLineariasedTransportMatrixY(double mad_init_x,
0125 double mad_init_thx,
0126 double mad_init_y,
0127 double mad_init_thy,
0128 double mad_init_xi,
0129 TMatrixD &tr_matrix,
0130 double d_mad_y = 10e-6,
0131 double d_mad_thy = 10e-6);
0132
0133 double GetDx(double mad_init_x,
0134 double mad_init_thx,
0135 double mad_init_y,
0136 double mad_init_thy,
0137 double mad_init_xi,
0138 double d_mad_xi = 0.001);
0139 double GetDxds(double mad_init_x,
0140 double mad_init_thx,
0141 double mad_init_y,
0142 double mad_init_thy,
0143 double mad_init_xi,
0144 double d_mad_xi = 0.001);
0145 inline beam_type GetBeamType() const { return beam; }
0146
0147
0148
0149
0150 void GetLinearApproximation(double atPoint[],
0151 double &Cx,
0152 double &Lx,
0153 double &vx,
0154 double &Cy,
0155 double &Ly,
0156 double &vy,
0157 double &D,
0158 double ep = 1E-5);
0159
0160 private:
0161 void Init();
0162 double s_begin_;
0163 double s_end_;
0164 beam_type beam;
0165 double nominal_beam_energy_;
0166 double nominal_beam_momentum_;
0167 bool trained_;
0168 std::vector<TMultiDimFet *> out_polynomials;
0169 std::vector<std::string> coord_names;
0170 std::vector<LHCApertureApproximator> apertures_;
0171
0172 #ifndef __CINT_
0173 friend class ProtonTransportFunctionsESSource;
0174 #endif
0175
0176 TMultiDimFet x_parametrisation;
0177 TMultiDimFet theta_x_parametrisation;
0178 TMultiDimFet y_parametrisation;
0179 TMultiDimFet theta_y_parametrisation;
0180
0181
0182 enum class VariableType { X, THETA_X, Y, THETA_Y };
0183
0184 void InitializeApproximators(polynomials_selection mode,
0185 int max_degree_x,
0186 int max_degree_tx,
0187 int max_degree_y,
0188 int max_degree_ty,
0189 bool common_terms);
0190 void SetDefaultAproximatorSettings(TMultiDimFet &approximator, VariableType var_type, int max_degree);
0191 void SetTermsManually(TMultiDimFet &approximator, VariableType variable, int max_degree, bool common_terms);
0192
0193 void AllocateErrorHists(TH1D *err_hists[4]);
0194 void AllocateErrorInputCorHists(TH2D *err_inp_cor_hists[4][5]);
0195 void AllocateErrorOutputCorHists(TH2D *err_out_cor_hists[4][5]);
0196
0197 void DeleteErrorHists(TH1D *err_hists[4]);
0198 void DeleteErrorCorHistograms(TH2D *err_cor_hists[4][5]);
0199
0200 void FillErrorHistograms(double errors[4], TH1D *err_hists[4]);
0201 void FillErrorDataCorHistograms(double errors[4], double var[5], TH2D *err_cor_hists[4][5]);
0202
0203 void WriteHistograms(TH1D *err_hists[4],
0204 TH2D *err_inp_cor_hists[4][5],
0205 TH2D *err_out_cor_hists[4][5],
0206 TFile *f_out,
0207 std::string base_out_dir);
0208
0209 ClassDef(LHCOpticsApproximator, 1)
0210 };
0211
0212 class LHCApertureApproximator : public LHCOpticsApproximator {
0213 public:
0214 enum class ApertureType { NO_APERTURE, RECTELLIPSE };
0215
0216 LHCApertureApproximator();
0217 LHCApertureApproximator(const LHCOpticsApproximator &in,
0218 double rect_x,
0219 double rect_y,
0220 double r_el_x,
0221 double r_el_y,
0222 ApertureType type = ApertureType::RECTELLIPSE);
0223
0224 bool CheckAperture(const double *in, bool invert_beam_coord_sytems = true) const;
0225
0226 private:
0227 double rect_x_, rect_y_, r_el_x_, r_el_y_;
0228 ApertureType ap_type_;
0229
0230 ClassDef(LHCApertureApproximator, 1)
0231 };
0232 #endif