CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Friends
LHCOpticsApproximator Class Reference

Class finds the parametrisation of MADX proton transport and transports the protons according to it 5 phase space variables are taken in to configuration: x, y, theta_x, theta_y, xi xi < 0 for momentum losses (that is for diffractive protons) More...

#include <LHCOpticsApproximator.h>

Inheritance diagram for LHCOpticsApproximator:
LHCApertureApproximator

Public Types

enum  beam_type { lhcb1, lhcb2 }
 
enum  polynomials_selection { AUTOMATIC, PREDEFINED }
 

Public Member Functions

void AddRectEllipseAperture (const LHCOpticsApproximator &in, double rect_x, double rect_y, double r_el_x, double r_el_y)
 
bool CheckInputRange (const double *in, bool invert_beam_coord_sytems=true) const
 
beam_type GetBeamType () const
 
double GetDx (double mad_init_x, double mad_init_thx, double mad_init_y, double mad_init_thy, double mad_init_xi, double d_mad_xi=0.001)
 
double GetDxds (double mad_init_x, double mad_init_thx, double mad_init_y, double mad_init_thy, double mad_init_xi, double d_mad_xi=0.001)
 
void GetLinearApproximation (double atPoint[], double &Cx, double &Lx, double &vx, double &Cy, double &Ly, double &vy, double &D, double ep=1E-5)
 
void GetLineariasedTransportMatrixX (double mad_init_x, double mad_init_thx, double mad_init_y, double mad_init_thy, double mad_init_xi, TMatrixD &tr_matrix, double d_mad_x=10e-6, double d_mad_thx=10e-6)
 returns linearised transport matrix for x projection | dx_out/dx_in dx_out/dthx_in | | dthx_out/dx_in dthx_out/dthx_in | More...
 
void GetLineariasedTransportMatrixY (double mad_init_x, double mad_init_thx, double mad_init_y, double mad_init_thy, double mad_init_xi, TMatrixD &tr_matrix, double d_mad_y=10e-6, double d_mad_thy=10e-6)
 returns linearised transport matrix for y projection | dy_out/dy_in dy_out/dthy_in | | dthy_out/dy_in dthy_out/dthy_in | More...
 
 LHCOpticsApproximator ()
 
 LHCOpticsApproximator (std::string name, std::string title, TMultiDimFet::EMDFPolyType polynom_type, std::string beam_direction, double nominal_beam_momentum)
 begin and end position along the beam of the particle to transport, training_tree, prefix of data branch in the tree More...
 
 LHCOpticsApproximator (const LHCOpticsApproximator &org)
 
const LHCOpticsApproximatoroperator= (const LHCOpticsApproximator &org)
 
double ParameterOutOfRangePenalty (double par_m[], bool invert_beam_coord_sytems=true) const
 
void PrintCoordinateOpticalFunctions (TMultiDimFet &parametrization, const std::string &coord_name, const std::vector< std::string > &input_vars)
 
void PrintInputRange ()
 
void PrintOpticalFunctions ()
 
void Test (TTree *inp_tree, TFile *f_out, std::string data_prefix=std::string("def"), std::string base_out_dir=std::string(""))
 
void TestAperture (TTree *in_tree, TTree *out_tree)
 x, theta_x, y, theta_y, ksi, mad_accepted, parametriz_accepted More...
 
void Train (TTree *inp_tree, std::string data_prefix=std::string("def"), polynomials_selection mode=PREDEFINED, int max_degree_x=10, int max_degree_tx=10, int max_degree_y=10, int max_degree_ty=10, bool common_terms=false, double *prec=nullptr)
 
bool Transport (const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const
 
bool Transport (const MadKinematicDescriptor *in, MadKinematicDescriptor *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const
 
bool Transport2D (const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const
 
bool Transport_m_GeV (double in_pos[3], double in_momentum[3], double out_pos[3], double out_momentum[3], bool check_apertures, double z2_z1_dist) const
 pos, momentum: x,y,z; pos in m, momentum in GeV/c More...
 

Private Types

enum  VariableType { VariableType::X, VariableType::THETA_X, VariableType::Y, VariableType::THETA_Y }
 

Private Member Functions

void AllocateErrorHists (TH1D *err_hists[4])
 
void AllocateErrorInputCorHists (TH2D *err_inp_cor_hists[4][5])
 
void AllocateErrorOutputCorHists (TH2D *err_out_cor_hists[4][5])
 
void DeleteErrorCorHistograms (TH2D *err_cor_hists[4][5])
 
void DeleteErrorHists (TH1D *err_hists[4])
 
void FillErrorDataCorHistograms (double errors[4], double var[5], TH2D *err_cor_hists[4][5])
 
void FillErrorHistograms (double errors[4], TH1D *err_hists[4])
 
void Init ()
 
void InitializeApproximators (polynomials_selection mode, int max_degree_x, int max_degree_tx, int max_degree_y, int max_degree_ty, bool common_terms)
 
void SetDefaultAproximatorSettings (TMultiDimFet &approximator, VariableType var_type, int max_degree)
 
void SetTermsManually (TMultiDimFet &approximator, VariableType variable, int max_degree, bool common_terms)
 
void WriteHistograms (TH1D *err_hists[4], TH2D *err_inp_cor_hists[4][5], TH2D *err_out_cor_hists[4][5], TFile *f_out, std::string base_out_dir)
 

Private Attributes

std::vector< LHCApertureApproximatorapertures_
 apertures on the way More...
 
beam_type beam
 
std::vector< std::string > coord_names
 pointers to polynomials More...
 
double nominal_beam_energy_
 GeV. More...
 
double nominal_beam_momentum_
 GeV/c. More...
 
std::vector< TMultiDimFet * > out_polynomials
 
double s_begin_
 begin of transport along the reference orbit More...
 
double s_end_
 end of transport along the reference orbit More...
 
TMultiDimFet theta_x_parametrisation
 polynomial approximation for theta_x More...
 
TMultiDimFet theta_y_parametrisation
 polynomial approximation for theta_y More...
 
bool trained_
 trained polynomials More...
 
TMultiDimFet x_parametrisation
 polynomial approximation for x More...
 
TMultiDimFet y_parametrisation
 polynomial approximation for y More...
 

Friends

class ProtonTransportFunctionsESSource
 

Detailed Description

Class finds the parametrisation of MADX proton transport and transports the protons according to it 5 phase space variables are taken in to configuration: x, y, theta_x, theta_y, xi xi < 0 for momentum losses (that is for diffractive protons)

Definition at line 29 of file LHCOpticsApproximator.h.

Member Enumeration Documentation

◆ beam_type

◆ polynomials_selection

◆ VariableType

Enumerator
THETA_X 
THETA_Y 

Definition at line 182 of file LHCOpticsApproximator.h.

182 { X, THETA_X, Y, THETA_Y };
#define X(str)
Definition: MuonsGrabber.cc:38

Constructor & Destructor Documentation

◆ LHCOpticsApproximator() [1/3]

LHCOpticsApproximator::LHCOpticsApproximator ( )

◆ LHCOpticsApproximator() [2/3]

LHCOpticsApproximator::LHCOpticsApproximator ( std::string  name,
std::string  title,
TMultiDimFet::EMDFPolyType  polynom_type,
std::string  beam_direction,
double  nominal_beam_momentum 
)

begin and end position along the beam of the particle to transport, training_tree, prefix of data branch in the tree

Definition at line 35 of file LHCOpticsApproximator.cc.

References beam, Init(), lhcb1, lhcb2, mergeVDriftHistosByStation::name, nominal_beam_energy_, nominal_beam_momentum_, and runGCPTkAlMap::title.

40  : x_parametrisation(5, polynom_type, "k"),
41  theta_x_parametrisation(5, polynom_type, "k"),
42  y_parametrisation(5, polynom_type, "k"),
43  theta_y_parametrisation(5, polynom_type, "k") {
44  this->SetName(name.c_str());
45  this->SetTitle(title.c_str());
46  Init();
47 
48  if (beam_direction == "lhcb1")
49  beam = lhcb1;
50  else if (beam_direction == "lhcb2")
51  beam = lhcb2;
52  else
53  beam = lhcb1;
54 
55  nominal_beam_momentum_ = nominal_beam_momentum;
56  nominal_beam_energy_ = TMath::Sqrt(nominal_beam_momentum_ * nominal_beam_momentum_ + 0.938272029 * 0.938272029);
57 }
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
double nominal_beam_momentum_
GeV/c.
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ LHCOpticsApproximator() [3/3]

LHCOpticsApproximator::LHCOpticsApproximator ( const LHCOpticsApproximator org)

Definition at line 249 of file LHCOpticsApproximator.cc.

References apertures_, beam, Init(), nominal_beam_energy_, nominal_beam_momentum_, s_begin_, s_end_, and trained_.

250  : TNamed(org),
255  Init();
256  s_begin_ = org.s_begin_;
257  s_end_ = org.s_end_;
258  trained_ = org.trained_;
259  apertures_ = org.apertures_;
260  beam = org.beam;
263 }
std::vector< LHCApertureApproximator > apertures_
apertures on the way
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
double nominal_beam_momentum_
GeV/c.
double s_begin_
begin of transport along the reference orbit
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
bool trained_
trained polynomials
double s_end_
end of transport along the reference orbit
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

Member Function Documentation

◆ AddRectEllipseAperture()

void LHCOpticsApproximator::AddRectEllipseAperture ( const LHCOpticsApproximator in,
double  rect_x,
double  rect_y,
double  r_el_x,
double  r_el_y 
)

◆ AllocateErrorHists()

void LHCOpticsApproximator::AllocateErrorHists ( TH1D *  err_hists[4])
private

Definition at line 724 of file LHCOpticsApproximator.cc.

References mps_fire::i.

Referenced by Test().

724  {
725  std::vector<std::string> error_labels;
726  error_labels.push_back("x error");
727  error_labels.push_back("theta_x error");
728  error_labels.push_back("y error");
729  error_labels.push_back("theta_y error");
730 
731  for (int i = 0; i < 4; ++i) {
732  err_hists[i] = new TH1D(error_labels[i].c_str(), error_labels[i].c_str(), 100, -0.0000000001, 0.0000000001);
733  err_hists[i]->SetXTitle(error_labels[i].c_str());
734  err_hists[i]->SetYTitle("counts");
735  err_hists[i]->SetDirectory(nullptr);
736  err_hists[i]->SetCanExtend(TH1::kAllAxes);
737  }
738 }

◆ AllocateErrorInputCorHists()

void LHCOpticsApproximator::AllocateErrorInputCorHists ( TH2D *  err_inp_cor_hists[4][5])
private

Definition at line 783 of file LHCOpticsApproximator.cc.

References mergeVDriftHistosByStation::name, AlCaHLTBitMon_QueryRunRegistry::string, and runGCPTkAlMap::title.

Referenced by Test().

783  {
784  std::vector<std::string> error_labels;
785  std::vector<std::string> data_labels;
786 
787  error_labels.push_back("x error");
788  error_labels.push_back("theta_x error");
789  error_labels.push_back("y error");
790  error_labels.push_back("theta_y error");
791 
792  data_labels.push_back("x input");
793  data_labels.push_back("theta_x input");
794  data_labels.push_back("y input");
795  data_labels.push_back("theta_y input");
796  data_labels.push_back("ksi input");
797 
798  for (int eri = 0; eri < 4; ++eri) {
799  for (int dati = 0; dati < 5; ++dati) {
800  std::string name = error_labels[eri] + " vs. " + data_labels[dati];
801  const std::string &title = name;
802  err_inp_cor_hists[eri][dati] =
803  new TH2D(name.c_str(), title.c_str(), 100, -0.0000000001, 0.0000000001, 100, -0.0000000001, 0.0000000001);
804  err_inp_cor_hists[eri][dati]->SetXTitle(error_labels[eri].c_str());
805  err_inp_cor_hists[eri][dati]->SetYTitle(data_labels[dati].c_str());
806  err_inp_cor_hists[eri][dati]->SetDirectory(nullptr);
807  err_inp_cor_hists[eri][dati]->SetCanExtend(TH1::kAllAxes);
808  }
809  }
810 }

◆ AllocateErrorOutputCorHists()

void LHCOpticsApproximator::AllocateErrorOutputCorHists ( TH2D *  err_out_cor_hists[4][5])
private

Definition at line 812 of file LHCOpticsApproximator.cc.

References mergeVDriftHistosByStation::name, AlCaHLTBitMon_QueryRunRegistry::string, and runGCPTkAlMap::title.

Referenced by Test().

812  {
813  std::vector<std::string> error_labels;
814  std::vector<std::string> data_labels;
815 
816  error_labels.push_back("x error");
817  error_labels.push_back("theta_x error");
818  error_labels.push_back("y error");
819  error_labels.push_back("theta_y error");
820 
821  data_labels.push_back("x output");
822  data_labels.push_back("theta_x output");
823  data_labels.push_back("y output");
824  data_labels.push_back("theta_y output");
825  data_labels.push_back("ksi output");
826 
827  for (int eri = 0; eri < 4; ++eri) {
828  for (int dati = 0; dati < 5; ++dati) {
829  std::string name = error_labels[eri] + " vs. " + data_labels[dati];
830  const std::string &title = name;
831  err_out_cor_hists[eri][dati] =
832  new TH2D(name.c_str(), title.c_str(), 100, -0.0000000001, 0.0000000001, 100, -0.0000000001, 0.0000000001);
833  err_out_cor_hists[eri][dati]->SetXTitle(error_labels[eri].c_str());
834  err_out_cor_hists[eri][dati]->SetYTitle(data_labels[dati].c_str());
835  err_out_cor_hists[eri][dati]->SetDirectory(nullptr);
836  err_out_cor_hists[eri][dati]->SetCanExtend(TH1::kAllAxes);
837  }
838  }
839 }

◆ CheckInputRange()

bool LHCOpticsApproximator::CheckInputRange ( const double *  in,
bool  invert_beam_coord_sytems = true 
) const

Definition at line 935 of file LHCOpticsApproximator.cc.

References EcalCondDBWriter_cfi::beam, mps_fire::i, and recoMuon::in.

Referenced by Transport().

937 {
938  double in_corrected[5];
939  if (beam == lhcb1 || !invert_beam_coord_sytems) {
940  in_corrected[0] = in[0];
941  in_corrected[1] = in[1];
942  in_corrected[2] = in[2];
943  in_corrected[3] = in[3];
944  in_corrected[4] = in[4];
945  } else {
946  in_corrected[0] = -in[0];
947  in_corrected[1] = -in[1];
948  in_corrected[2] = in[2];
949  in_corrected[3] = in[3];
950  in_corrected[4] = in[4];
951  }
952 
953  const TVectorD *min_var = x_parametrisation.GetMinVariables();
954  const TVectorD *max_var = x_parametrisation.GetMaxVariables();
955  bool res = true;
956 
957  for (int i = 0; i < 5; i++) {
958  res = res && in_corrected[i] >= (*min_var)(i) && in_corrected[i] <= (*max_var)(i);
959  }
960  return res;
961 }
Definition: Electron.h:6
const TVectorD * GetMinVariables() const
Definition: TMultiDimFet.h:166
const TVectorD * GetMaxVariables() const
Definition: TMultiDimFet.h:160
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ DeleteErrorCorHistograms()

void LHCOpticsApproximator::DeleteErrorCorHistograms ( TH2D *  err_cor_hists[4][5])
private

Definition at line 861 of file LHCOpticsApproximator.cc.

Referenced by Test().

861  {
862  for (int eri = 0; eri < 4; ++eri) {
863  for (int dati = 0; dati < 5; ++dati) {
864  delete err_cor_hists[eri][dati];
865  }
866  }
867 }

◆ DeleteErrorHists()

void LHCOpticsApproximator::DeleteErrorHists ( TH1D *  err_hists[4])
private

Definition at line 855 of file LHCOpticsApproximator.cc.

References mps_fire::i.

Referenced by Test().

855  {
856  for (int i = 0; i < 4; ++i) {
857  delete err_hists[i];
858  }
859 }

◆ FillErrorDataCorHistograms()

void LHCOpticsApproximator::FillErrorDataCorHistograms ( double  errors[4],
double  var[5],
TH2D *  err_cor_hists[4][5] 
)
private

Definition at line 847 of file LHCOpticsApproximator.cc.

References ALCARECOEcalPhiSym_cff::var.

Referenced by Test().

847  {
848  for (int eri = 0; eri < 4; ++eri) {
849  for (int dati = 0; dati < 5; ++dati) {
850  err_cor_hists[eri][dati]->Fill(errors[eri], var[dati]);
851  }
852  }
853 }
Definition: errors.py:1

◆ FillErrorHistograms()

void LHCOpticsApproximator::FillErrorHistograms ( double  errors[4],
TH1D *  err_hists[4] 
)
private

Definition at line 841 of file LHCOpticsApproximator.cc.

References mps_fire::i.

Referenced by Test().

841  {
842  for (int i = 0; i < 4; ++i) {
843  err_hists[i]->Fill(errors[i]);
844  }
845 }
Definition: errors.py:1

◆ GetBeamType()

beam_type LHCOpticsApproximator::GetBeamType ( ) const
inline

Definition at line 145 of file LHCOpticsApproximator.h.

References beam.

145 { return beam; }

◆ GetDx()

double LHCOpticsApproximator::GetDx ( double  mad_init_x,
double  mad_init_thx,
double  mad_init_y,
double  mad_init_thy,
double  mad_init_xi,
double  d_mad_xi = 0.001 
)

Definition at line 1173 of file LHCOpticsApproximator.cc.

References recoMuon::in, MillePedeFileConverter_cfg::out, Transport(), and testProducerWithPsetDescEmpty_cfi::x1.

1178  {
1179  double in[5];
1180  in[0] = mad_init_x;
1181  in[1] = mad_init_thx;
1182  in[2] = mad_init_y;
1183  in[3] = mad_init_thy;
1184  in[4] = mad_init_xi;
1185 
1186  double out[5];
1187 
1188  if (!Transport(in, out)) {
1189  return 0.0;
1190  }
1191  double x1 = out[0];
1192 
1193  in[4] = mad_init_xi + d_mad_xi;
1194  if (!Transport(in, out)) {
1195  return 0.0;
1196  }
1197  double x2_dxi = out[0];
1198  double dispersion = (x2_dxi - x1) / d_mad_xi;
1199 
1200  return dispersion;
1201 }
bool Transport(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const

◆ GetDxds()

double LHCOpticsApproximator::GetDxds ( double  mad_init_x,
double  mad_init_thx,
double  mad_init_y,
double  mad_init_thy,
double  mad_init_xi,
double  d_mad_xi = 0.001 
)

Definition at line 1205 of file LHCOpticsApproximator.cc.

References recoMuon::in, MillePedeFileConverter_cfg::out, and Transport().

1210  {
1211  double MADX_momentum_correction_factor = 1.0 + mad_init_xi;
1212  double in[5];
1213  in[0] = mad_init_x;
1214  in[1] = mad_init_thx;
1215  in[2] = mad_init_y;
1216  in[3] = mad_init_thy;
1217  in[4] = mad_init_xi;
1218 
1219  double out[5];
1220 
1221  if (!Transport(in, out)) {
1222  return 0.0;
1223  }
1224  double thx1 = out[1] / MADX_momentum_correction_factor;
1225 
1226  in[4] = mad_init_xi + d_mad_xi;
1227  if (!Transport(in, out)) {
1228  return 0.0;
1229  }
1230  double thx2_dxi = out[1] / MADX_momentum_correction_factor;
1231  double dispersion = (thx2_dxi - thx1) / d_mad_xi;
1232 
1233  return dispersion;
1234 }
bool Transport(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const

◆ GetLinearApproximation()

void LHCOpticsApproximator::GetLinearApproximation ( double  atPoint[],
double &  Cx,
double &  Lx,
double &  vx,
double &  Cy,
double &  Ly,
double &  vy,
double &  D,
double  ep = 1E-5 
)

returns linear approximation of the transport parameterization takes numerical derivatives (see parameter ep) around point ‘atPoint’ (this array has the same structure as ‘in’ parameter in Transport method) the linearized transport: x = Cx + Lx*theta_x + vx*x_star

Definition at line 1037 of file LHCOpticsApproximator.cc.

References ecalHexDisplay_cfg::ep, mps_fire::i, MillePedeFileConverter_cfg::out, Transport2D(), btvMC_cff::vx, and btvMC_cff::vy.

1038  {
1039  double out[2];
1040  if (!Transport2D(atPoint, out)) {
1041  return;
1042  }
1043  Cx = out[0];
1044  Cy = out[1];
1045 
1046  for (int i = 0; i < 5; ++i) {
1047  atPoint[i] += ep;
1048  if (!Transport2D(atPoint, out)) {
1049  continue;
1050  }
1051  switch (i) {
1052  case 0:
1053  vx = (out[0] - Cx) / ep;
1054  break;
1055  case 1:
1056  Lx = (out[0] - Cx) / ep;
1057  break;
1058  case 2:
1059  vy = (out[1] - Cy) / ep;
1060  break;
1061  case 3:
1062  Ly = (out[1] - Cy) / ep;
1063  break;
1064  case 4:
1065  D = (out[0] - Cx) / ep;
1066  break;
1067  }
1068  atPoint[i] -= ep;
1069  }
1070 }
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
bool Transport2D(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const

◆ GetLineariasedTransportMatrixX()

void LHCOpticsApproximator::GetLineariasedTransportMatrixX ( double  mad_init_x,
double  mad_init_thx,
double  mad_init_y,
double  mad_init_thy,
double  mad_init_xi,
TMatrixD &  tr_matrix,
double  d_mad_x = 10e-6,
double  d_mad_thx = 10e-6 
)

returns linearised transport matrix for x projection | dx_out/dx_in dx_out/dthx_in | | dthx_out/dx_in dthx_out/dthx_in |

input: [m], [rad], xi:-1...0

Definition at line 1073 of file LHCOpticsApproximator.cc.

References recoMuon::in, MillePedeFileConverter_cfg::out, Transport(), and testProducerWithPsetDescEmpty_cfi::x1.

1080  {
1081  double MADX_momentum_correction_factor = 1.0 + mad_init_xi;
1082  transp_matrix.ResizeTo(2, 2);
1083  double in[5];
1084  in[0] = mad_init_x;
1085  in[1] = mad_init_thx;
1086  in[2] = mad_init_y;
1087  in[3] = mad_init_thy;
1088  in[4] = mad_init_xi;
1089 
1090  double out[5];
1091 
1092  if (!Transport(in, out)) {
1093  return;
1094  };
1095  double x1 = out[0];
1096  double thx1 = out[1];
1097 
1098  in[0] = mad_init_x + d_mad_x;
1099  if (!Transport(in, out)) {
1100  return;
1101  }
1102  double x2_dx = out[0];
1103  double thx2_dx = out[1];
1104 
1105  in[0] = mad_init_x;
1106  in[1] = mad_init_thx + d_mad_thx; //?
1107  if (!Transport(in, out)) {
1108  return;
1109  }
1110  double x2_dthx = out[0];
1111  double thx2_dthx = out[1];
1112 
1113  // | dx/dx, dx/dthx |
1114  // | dthx/dx, dtchx/dthx |
1115 
1116  transp_matrix(0, 0) = (x2_dx - x1) / d_mad_x;
1117  transp_matrix(1, 0) = (thx2_dx - thx1) / (d_mad_x * MADX_momentum_correction_factor);
1118  transp_matrix(0, 1) = MADX_momentum_correction_factor * (x2_dthx - x1) / d_mad_thx;
1119  transp_matrix(1, 1) = (thx2_dthx - thx1) / d_mad_thx;
1120 }
bool Transport(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const

◆ GetLineariasedTransportMatrixY()

void LHCOpticsApproximator::GetLineariasedTransportMatrixY ( double  mad_init_x,
double  mad_init_thx,
double  mad_init_y,
double  mad_init_thy,
double  mad_init_xi,
TMatrixD &  tr_matrix,
double  d_mad_y = 10e-6,
double  d_mad_thy = 10e-6 
)

returns linearised transport matrix for y projection | dy_out/dy_in dy_out/dthy_in | | dthy_out/dy_in dthy_out/dthy_in |

input: [m], [rad], xi:-1...0

Definition at line 1123 of file LHCOpticsApproximator.cc.

References recoMuon::in, MillePedeFileConverter_cfg::out, Transport(), and testProducerWithPsetDescEmpty_cfi::y1.

1130  {
1131  double MADX_momentum_correction_factor = 1.0 + mad_init_xi;
1132  transp_matrix.ResizeTo(2, 2);
1133  double in[5];
1134  in[0] = mad_init_x;
1135  in[1] = mad_init_thx;
1136  in[2] = mad_init_y;
1137  in[3] = mad_init_thy;
1138  in[4] = mad_init_xi;
1139 
1140  double out[5];
1141 
1142  if (!Transport(in, out)) {
1143  return;
1144  };
1145  double y1 = out[2];
1146  double thy1 = out[3];
1147 
1148  in[2] = mad_init_y + d_mad_y;
1149  if (!Transport(in, out)) {
1150  return;
1151  }
1152  double y2_dy = out[2];
1153  double thy2_dy = out[3];
1154 
1155  in[2] = mad_init_y;
1156  in[3] = mad_init_thy + d_mad_thy; //?
1157  if (!Transport(in, out)) {
1158  return;
1159  }
1160  double y2_dthy = out[2];
1161  double thy2_dthy = out[3];
1162 
1163  // | dy/dy, dy/dthy |
1164  // | dthy/dy, dtchy/dthy |
1165 
1166  transp_matrix(0, 0) = (y2_dy - y1) / d_mad_y;
1167  transp_matrix(1, 0) = (thy2_dy - thy1) / (d_mad_y * MADX_momentum_correction_factor);
1168  transp_matrix(0, 1) = MADX_momentum_correction_factor * (y2_dthy - y1) / d_mad_thy;
1169  transp_matrix(1, 1) = (thy2_dthy - thy1) / d_mad_thy;
1170 }
bool Transport(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const

◆ Init()

void LHCOpticsApproximator::Init ( void  )
private

Definition at line 15 of file LHCOpticsApproximator.cc.

References apertures_, coord_names, out_polynomials, s_begin_, s_end_, theta_x_parametrisation, theta_y_parametrisation, trained_, x_parametrisation, and y_parametrisation.

Referenced by LHCOpticsApproximator(), and operator=().

15  {
16  out_polynomials.clear();
17  apertures_.clear();
22 
23  coord_names.clear();
24  coord_names.push_back("x");
25  coord_names.push_back("theta_x");
26  coord_names.push_back("y");
27  coord_names.push_back("theta_y");
28  coord_names.push_back("ksi");
29 
30  s_begin_ = 0.0;
31  s_end_ = 0.0;
32  trained_ = false;
33 }
std::vector< LHCApertureApproximator > apertures_
apertures on the way
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
double s_begin_
begin of transport along the reference orbit
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
bool trained_
trained polynomials
std::vector< TMultiDimFet * > out_polynomials
double s_end_
end of transport along the reference orbit
std::vector< std::string > coord_names
pointers to polynomials
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ InitializeApproximators()

void LHCOpticsApproximator::InitializeApproximators ( polynomials_selection  mode,
int  max_degree_x,
int  max_degree_tx,
int  max_degree_y,
int  max_degree_ty,
bool  common_terms 
)
private

Definition at line 395 of file LHCOpticsApproximator.cc.

References ALCARECOPromptCalibProdSiPixelAli0T_cff::mode, PREDEFINED, SetDefaultAproximatorSettings(), SetTermsManually(), THETA_X, theta_x_parametrisation, THETA_Y, theta_y_parametrisation, X, x_parametrisation, Y, and y_parametrisation.

Referenced by Train().

400  {
405 
406  if (mode == PREDEFINED) {
407  SetTermsManually(x_parametrisation, VariableType::X, max_degree_x, common_terms);
408  SetTermsManually(theta_x_parametrisation, VariableType::THETA_X, max_degree_tx, common_terms);
409  SetTermsManually(y_parametrisation, VariableType::Y, max_degree_y, common_terms);
410  SetTermsManually(theta_y_parametrisation, VariableType::THETA_Y, max_degree_ty, common_terms);
411  }
412 }
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
void SetTermsManually(TMultiDimFet &approximator, VariableType variable, int max_degree, bool common_terms)
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
void SetDefaultAproximatorSettings(TMultiDimFet &approximator, VariableType var_type, int max_degree)
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ operator=()

const LHCOpticsApproximator & LHCOpticsApproximator::operator= ( const LHCOpticsApproximator org)

Definition at line 265 of file LHCOpticsApproximator.cc.

References apertures_, beam, Init(), nominal_beam_energy_, nominal_beam_momentum_, operator=(), s_begin_, s_end_, theta_x_parametrisation, theta_y_parametrisation, trained_, x_parametrisation, and y_parametrisation.

265  {
266  if (this != &org) {
271  Init();
272 
273  TNamed::operator=(org);
274  s_begin_ = org.s_begin_;
275  s_end_ = org.s_end_;
276  trained_ = org.trained_;
277 
278  apertures_ = org.apertures_;
279  beam = org.beam;
282  }
283  return org;
284 }
Basic3DVector & operator=(const Basic3DVector &)=default
Assignment operator.
std::vector< LHCApertureApproximator > apertures_
apertures on the way
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
double nominal_beam_momentum_
GeV/c.
double s_begin_
begin of transport along the reference orbit
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
bool trained_
trained polynomials
double s_end_
end of transport along the reference orbit
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ ParameterOutOfRangePenalty()

double LHCOpticsApproximator::ParameterOutOfRangePenalty ( double  par_m[],
bool  invert_beam_coord_sytems = true 
) const

Definition at line 68 of file LHCOpticsApproximator.cc.

References beam, TMultiDimFet::GetMaxVariables(), TMultiDimFet::GetMinVariables(), mps_fire::i, recoMuon::in, lhcb1, and x_parametrisation.

68  {
69  double in_corrected[5];
70  if (beam == lhcb1 || !invert_beam_coord_sytems) {
71  in_corrected[0] = in[0];
72  in_corrected[1] = in[1];
73  in_corrected[2] = in[2];
74  in_corrected[3] = in[3];
75  in_corrected[4] = in[4];
76  } else {
77  in_corrected[0] = -in[0];
78  in_corrected[1] = -in[1];
79  in_corrected[2] = in[2];
80  in_corrected[3] = in[3];
81  in_corrected[4] = in[4];
82  }
83 
84  const TVectorD *min_var = x_parametrisation.GetMinVariables();
85  const TVectorD *max_var = x_parametrisation.GetMaxVariables();
86  double res = 0.;
87 
88  for (int i = 0; i < 5; i++) {
89  if (in_corrected[i] < (*min_var)(i)) {
90  double dist = TMath::Abs(((*min_var)(i)-in_corrected[i]) / ((*max_var)(i) - (*min_var)(i)));
91  res += 8 * (TMath::Exp(dist) - 1.0);
92  in_corrected[i] = (*min_var)(i);
93  } else if (in_corrected[i] > (*max_var)(i)) {
94  double dist = TMath::Abs((in_corrected[i] - (*max_var)(i)) / ((*max_var)(i) - (*min_var)(i)));
95  res += 8 * (TMath::Exp(dist) - 1.0);
96  in_corrected[i] = (*max_var)(i);
97  }
98  }
99  return res;
100 }
Definition: Electron.h:6
const TVectorD * GetMinVariables() const
Definition: TMultiDimFet.h:166
const TVectorD * GetMaxVariables() const
Definition: TMultiDimFet.h:160
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ PrintCoordinateOpticalFunctions()

void LHCOpticsApproximator::PrintCoordinateOpticalFunctions ( TMultiDimFet parametrization,
const std::string &  coord_name,
const std::vector< std::string > &  input_vars 
)

Definition at line 1006 of file LHCOpticsApproximator.cc.

References MillePedeFileConverter_cfg::e, mps_fire::i, recoMuon::in, dqmiolumiharvest::j, and stringResolutionProvider_cfi::parametrization.

Referenced by PrintOpticalFunctions().

1008  {
1009  double in[5];
1010  // double out;
1011  double d_out_d_in[5];
1012  double d_par = 1e-5;
1013  double bias = 0;
1014 
1015  for (int j = 0; j < 5; j++)
1016  in[j] = 0.0;
1017 
1018  bias = parametrization.Eval(in);
1019 
1020  for (int i = 0; i < 5; i++) {
1021  for (int j = 0; j < 5; j++)
1022  in[j] = 0.0;
1023 
1024  in[i] = d_par;
1025  d_out_d_in[i] = parametrization.Eval(in);
1026  in[i] = 0.0;
1027  d_out_d_in[i] = d_out_d_in[i] - parametrization.Eval(in);
1028  d_out_d_in[i] = d_out_d_in[i] / d_par;
1029  }
1030  edm::LogInfo("LHCOpticsApproximator") << coord_name << " = " << bias;
1031  for (int i = 0; i < 5; i++) {
1032  edm::LogInfo("LHCOpticsApproximator") << " + " << d_out_d_in[i] << "*" << input_vars[i];
1033  }
1034  edm::LogInfo("LHCOpticsApproximator") << "\n";
1035 }
parametrization
specify parametrization (see SWGuidePATKinematicResolutions for more details)
Log< level::Info, false > LogInfo

◆ PrintInputRange()

void LHCOpticsApproximator::PrintInputRange ( )

Definition at line 923 of file LHCOpticsApproximator.cc.

References coord_names, TMultiDimFet::GetMaxVariables(), TMultiDimFet::GetMinVariables(), mps_fire::i, and x_parametrisation.

Referenced by Train().

923  {
924  const TVectorD *min_var = x_parametrisation.GetMinVariables();
925  const TVectorD *max_var = x_parametrisation.GetMaxVariables();
926 
927  edm::LogInfo("LHCOpticsApproximator") << "Covered input parameters range:"
928  << "\n";
929  for (int i = 0; i < 5; i++) {
930  edm::LogInfo("LHCOpticsApproximator") << (*min_var)(i) << " < " << coord_names[i] << " < " << (*max_var)(i) << "\n";
931  }
932  edm::LogInfo("LHCOpticsApproximator") << "\n";
933 }
Log< level::Info, false > LogInfo
const TVectorD * GetMinVariables() const
Definition: TMultiDimFet.h:166
const TVectorD * GetMaxVariables() const
Definition: TMultiDimFet.h:160
std::vector< std::string > coord_names
pointers to polynomials
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ PrintOpticalFunctions()

void LHCOpticsApproximator::PrintOpticalFunctions ( )

Definition at line 997 of file LHCOpticsApproximator.cc.

References coord_names, mps_fire::i, out_polynomials, and PrintCoordinateOpticalFunctions().

997  {
998  edm::LogInfo("LHCOpticsApproximator") << std::endl
999  << "Linear terms of optical functions:"
1000  << "\n";
1001  for (int i = 0; i < 4; i++) {
1003  }
1004 }
void PrintCoordinateOpticalFunctions(TMultiDimFet &parametrization, const std::string &coord_name, const std::vector< std::string > &input_vars)
std::vector< TMultiDimFet * > out_polynomials
Log< level::Info, false > LogInfo
std::vector< std::string > coord_names
pointers to polynomials

◆ SetDefaultAproximatorSettings()

void LHCOpticsApproximator::SetDefaultAproximatorSettings ( TMultiDimFet approximator,
VariableType  var_type,
int  max_degree 
)
private

Definition at line 414 of file LHCOpticsApproximator.cc.

References MillePedeFileConverter_cfg::e, TMultiDimFet::SetMaxAngle(), TMultiDimFet::SetMaxFunctions(), TMultiDimFet::SetMaxPowers(), TMultiDimFet::SetMaxStudy(), TMultiDimFet::SetMaxTerms(), TMultiDimFet::SetMinAngle(), TMultiDimFet::SetMinRelativeError(), TMultiDimFet::SetPowerLimit(), THETA_X, THETA_Y, X, and Y.

Referenced by InitializeApproximators().

416  {
417  if (max_degree < 1 || max_degree > 20)
418  max_degree = 10;
419 
420  if (var_type == VariableType::X || var_type == VariableType::THETA_X) {
421  Int_t mPowers[] = {1, 1, 0, 0, max_degree};
422  approximator.SetMaxPowers(mPowers);
423  approximator.SetMaxFunctions(900);
424  approximator.SetMaxStudy(1000);
425  approximator.SetMaxTerms(900);
426  approximator.SetPowerLimit(1.50);
427  approximator.SetMinAngle(10);
428  approximator.SetMaxAngle(10);
429  approximator.SetMinRelativeError(1e-13);
430  }
431 
432  if (var_type == VariableType::Y || var_type == VariableType::THETA_Y) {
433  Int_t mPowers[] = {0, 0, 1, 1, max_degree};
434  approximator.SetMaxPowers(mPowers);
435  approximator.SetMaxFunctions(900);
436  approximator.SetMaxStudy(1000);
437  approximator.SetMaxTerms(900);
438  approximator.SetPowerLimit(1.50);
439  approximator.SetMinAngle(10);
440  approximator.SetMaxAngle(10);
441  approximator.SetMinRelativeError(1e-13);
442  }
443 }
void SetPowerLimit(Double_t limit=1e-3)
void SetMaxPowers(const Int_t *powers)
void SetMaxTerms(Int_t terms)
Definition: TMultiDimFet.h:206
void SetMinRelativeError(Double_t error)
void SetMaxStudy(Int_t n)
Definition: TMultiDimFet.h:205
void SetMaxFunctions(Int_t n)
Definition: TMultiDimFet.h:203
void SetMinAngle(Double_t angle=1)
void SetMaxAngle(Double_t angle=0)

◆ SetTermsManually()

void LHCOpticsApproximator::SetTermsManually ( TMultiDimFet approximator,
VariableType  variable,
int  max_degree,
bool  common_terms 
)
private

Definition at line 445 of file LHCOpticsApproximator.cc.

References mps_fire::i, TMultiDimFet::SetPowers(), THETA_X, THETA_Y, HPSPFTaus_cff::variable, X, and Y.

Referenced by InitializeApproximators().

448  {
449  if (max_degree < 1 || max_degree > 20)
450  max_degree = 10;
451 
452  //put terms of shape:
453  //1,0,0,0,t 0,1,0,0,t 0,2,0,0,t 0,3,0,0,t 0,0,0,0,t
454  //t: 0,1,...,max_degree
455 
456  std::vector<Int_t> term_literals;
457  term_literals.reserve(5000);
458 
460  //1,0,0,0,t
461  for (int i = 0; i <= max_degree; ++i) {
462  term_literals.push_back(1);
463  term_literals.push_back(0);
464  term_literals.push_back(0);
465  term_literals.push_back(0);
466  term_literals.push_back(i);
467  }
468  //0,1,0,0,t
469  for (int i = 0; i <= max_degree; ++i) {
470  term_literals.push_back(0);
471  term_literals.push_back(1);
472  term_literals.push_back(0);
473  term_literals.push_back(0);
474  term_literals.push_back(i);
475  }
476  //0,2,0,0,t
477  for (int i = 0; i <= max_degree; ++i) {
478  term_literals.push_back(0);
479  term_literals.push_back(2);
480  term_literals.push_back(0);
481  term_literals.push_back(0);
482  term_literals.push_back(i);
483  }
484  //0,3,0,0,t
485  for (int i = 0; i <= max_degree; ++i) {
486  term_literals.push_back(0);
487  term_literals.push_back(3);
488  term_literals.push_back(0);
489  term_literals.push_back(0);
490  term_literals.push_back(i);
491  }
492  //0,0,0,0,t
493  for (int i = 0; i <= max_degree; ++i) {
494  term_literals.push_back(0);
495  term_literals.push_back(0);
496  term_literals.push_back(0);
497  term_literals.push_back(0);
498  term_literals.push_back(i);
499  }
500  }
501 
503  //0,0,1,0,t
504  for (int i = 0; i <= max_degree; ++i) {
505  term_literals.push_back(0);
506  term_literals.push_back(0);
507  term_literals.push_back(1);
508  term_literals.push_back(0);
509  term_literals.push_back(i);
510  }
511  //0,0,0,1,t
512  for (int i = 0; i <= max_degree; ++i) {
513  term_literals.push_back(0);
514  term_literals.push_back(0);
515  term_literals.push_back(0);
516  term_literals.push_back(1);
517  term_literals.push_back(i);
518  }
519  //0,0,0,2,t
520  for (int i = 0; i <= max_degree; ++i) {
521  term_literals.push_back(0);
522  term_literals.push_back(0);
523  term_literals.push_back(0);
524  term_literals.push_back(2);
525  term_literals.push_back(i);
526  }
527  //0,0,0,3,t
528  for (int i = 0; i <= max_degree; ++i) {
529  term_literals.push_back(0);
530  term_literals.push_back(0);
531  term_literals.push_back(0);
532  term_literals.push_back(3);
533  term_literals.push_back(i);
534  }
535  //0,0,0,0,t
536  for (int i = 0; i <= max_degree; ++i) {
537  term_literals.push_back(0);
538  term_literals.push_back(0);
539  term_literals.push_back(0);
540  term_literals.push_back(0);
541  term_literals.push_back(i);
542  }
543  }
544 
545  //push common terms
546  if (common_terms) {
547  term_literals.push_back(1), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0),
548  term_literals.push_back(0);
549  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0),
550  term_literals.push_back(0);
551  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1),
552  term_literals.push_back(0);
553  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(1), term_literals.push_back(0),
554  term_literals.push_back(0);
555  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(1),
556  term_literals.push_back(0);
557  term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(1),
558  term_literals.push_back(0);
559 
560  term_literals.push_back(1), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0),
561  term_literals.push_back(1);
562  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0),
563  term_literals.push_back(1);
564  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1),
565  term_literals.push_back(1);
566  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(1), term_literals.push_back(0),
567  term_literals.push_back(1);
568  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(1),
569  term_literals.push_back(1);
570  term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(1),
571  term_literals.push_back(1);
572 
573  term_literals.push_back(1), term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(0),
574  term_literals.push_back(0);
575  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(0),
576  term_literals.push_back(0);
577  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(2),
578  term_literals.push_back(0);
579  term_literals.push_back(2), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0),
580  term_literals.push_back(0);
581  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(2), term_literals.push_back(0),
582  term_literals.push_back(0);
583  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(2),
584  term_literals.push_back(0);
585  term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0),
586  term_literals.push_back(0);
587  term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(1), term_literals.push_back(0),
588  term_literals.push_back(0);
589  term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(2),
590  term_literals.push_back(0);
591  term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1),
592  term_literals.push_back(0);
593  term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(1),
594  term_literals.push_back(0);
595  term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(1),
596  term_literals.push_back(0);
597 
598  term_literals.push_back(1), term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(0),
599  term_literals.push_back(1);
600  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(0),
601  term_literals.push_back(1);
602  term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(2),
603  term_literals.push_back(1);
604  term_literals.push_back(2), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(0),
605  term_literals.push_back(1);
606  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(2), term_literals.push_back(0),
607  term_literals.push_back(1);
608  term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0), term_literals.push_back(2),
609  term_literals.push_back(1);
610  term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(0),
611  term_literals.push_back(1);
612  term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(1), term_literals.push_back(0),
613  term_literals.push_back(1);
614  term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1), term_literals.push_back(2),
615  term_literals.push_back(1);
616  term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(1),
617  term_literals.push_back(1);
618  term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(0), term_literals.push_back(1),
619  term_literals.push_back(1);
620  term_literals.push_back(0), term_literals.push_back(0), term_literals.push_back(2), term_literals.push_back(1),
621  term_literals.push_back(1);
622  }
623 
624  std::vector<Int_t> powers;
625  powers.resize(term_literals.size());
626 
627  for (unsigned int i = 0; i < term_literals.size(); ++i) {
628  powers[i] = term_literals[i];
629  }
630  approximator.SetPowers(&powers[0], term_literals.size() / 5);
631 }
virtual void SetPowers(const Int_t *powers, Int_t terms)

◆ Test()

void LHCOpticsApproximator::Test ( TTree *  inp_tree,
TFile *  f_out,
std::string  data_prefix = std::string("def"),
std::string  base_out_dir = std::string("") 
)

Definition at line 633 of file LHCOpticsApproximator.cc.

References AllocateErrorHists(), AllocateErrorInputCorHists(), AllocateErrorOutputCorHists(), DeleteErrorCorHistograms(), DeleteErrorHists(), TMultiDimFet::Eval(), FillErrorDataCorHistograms(), FillErrorHistograms(), mps_fire::i, AlCaHLTBitMon_QueryRunRegistry::string, theta_x_parametrisation, theta_y_parametrisation, WriteHistograms(), x_parametrisation, and y_parametrisation.

633  {
634  if (inp_tree == nullptr || f_out == nullptr)
635  return;
636 
637  //in-variables
638  //x_in, theta_x_in, y_in, theta_y_in, ksi_in, s_in
639  double in_var[6];
640 
641  //out-variables
642  //x_out, theta_x_out, y_out, theta_y_out, ksi_out, s_out, valid_out;
643  double out_var[7];
644 
645  //in- out-lables
646  std::string x_in_lab = "x_in";
647  std::string theta_x_in_lab = "theta_x_in";
648  std::string y_in_lab = "y_in";
649  std::string theta_y_in_lab = "theta_y_in";
650  std::string ksi_in_lab = "ksi_in";
651  std::string s_in_lab = "s_in";
652 
653  std::string x_out_lab = data_prefix + "_x_out";
654  std::string theta_x_out_lab = data_prefix + "_theta_x_out";
655  std::string y_out_lab = data_prefix + "_y_out";
656  std::string theta_y_out_lab = data_prefix + "_theta_y_out";
657  std::string ksi_out_lab = data_prefix + "_ksi_out";
658  std::string s_out_lab = data_prefix + "_s_out";
659  std::string valid_out_lab = data_prefix + "_valid_out";
660 
661  //disable not needed branches to speed up the readin
662  inp_tree->SetBranchStatus("*", false); //disable all branches
663  inp_tree->SetBranchStatus(x_in_lab.c_str(), true);
664  inp_tree->SetBranchStatus(theta_x_in_lab.c_str(), true);
665  inp_tree->SetBranchStatus(y_in_lab.c_str(), true);
666  inp_tree->SetBranchStatus(theta_y_in_lab.c_str(), true);
667  inp_tree->SetBranchStatus(ksi_in_lab.c_str(), true);
668  inp_tree->SetBranchStatus(x_out_lab.c_str(), true);
669  inp_tree->SetBranchStatus(theta_x_out_lab.c_str(), true);
670  inp_tree->SetBranchStatus(y_out_lab.c_str(), true);
671  inp_tree->SetBranchStatus(theta_y_out_lab.c_str(), true);
672  inp_tree->SetBranchStatus(ksi_out_lab.c_str(), true);
673  inp_tree->SetBranchStatus(valid_out_lab.c_str(), true);
674 
675  //set input data adresses
676  inp_tree->SetBranchAddress(x_in_lab.c_str(), &(in_var[0]));
677  inp_tree->SetBranchAddress(theta_x_in_lab.c_str(), &(in_var[1]));
678  inp_tree->SetBranchAddress(y_in_lab.c_str(), &(in_var[2]));
679  inp_tree->SetBranchAddress(theta_y_in_lab.c_str(), &(in_var[3]));
680  inp_tree->SetBranchAddress(ksi_in_lab.c_str(), &(in_var[4]));
681  inp_tree->SetBranchAddress(s_in_lab.c_str(), &(in_var[5]));
682 
683  //set output data adresses
684  inp_tree->SetBranchAddress(x_out_lab.c_str(), &(out_var[0]));
685  inp_tree->SetBranchAddress(theta_x_out_lab.c_str(), &(out_var[1]));
686  inp_tree->SetBranchAddress(y_out_lab.c_str(), &(out_var[2]));
687  inp_tree->SetBranchAddress(theta_y_out_lab.c_str(), &(out_var[3]));
688  inp_tree->SetBranchAddress(ksi_out_lab.c_str(), &(out_var[4]));
689  inp_tree->SetBranchAddress(s_out_lab.c_str(), &(out_var[5]));
690  inp_tree->SetBranchAddress(valid_out_lab.c_str(), &(out_var[6]));
691 
692  //test histogramms
693  TH1D *err_hists[4];
694  TH2D *err_inp_cor_hists[4][5];
695  TH2D *err_out_cor_hists[4][5];
696 
697  AllocateErrorHists(err_hists);
698  AllocateErrorInputCorHists(err_inp_cor_hists);
699  AllocateErrorOutputCorHists(err_out_cor_hists);
700 
701  Long64_t entries = inp_tree->GetEntries();
702  //set input and output variables for fitting
703  for (Long64_t i = 0; i < entries; ++i) {
704  double errors[4];
705  inp_tree->GetEntry(i);
706 
707  errors[0] = out_var[0] - x_parametrisation.Eval(in_var);
708  errors[1] = out_var[1] - theta_x_parametrisation.Eval(in_var);
709  errors[2] = out_var[2] - y_parametrisation.Eval(in_var);
710  errors[3] = out_var[3] - theta_y_parametrisation.Eval(in_var);
711 
712  FillErrorHistograms(errors, err_hists);
713  FillErrorDataCorHistograms(errors, in_var, err_inp_cor_hists);
714  FillErrorDataCorHistograms(errors, out_var, err_out_cor_hists);
715  }
716 
717  WriteHistograms(err_hists, err_inp_cor_hists, err_out_cor_hists, f_out, base_out_dir);
718 
719  DeleteErrorHists(err_hists);
720  DeleteErrorCorHistograms(err_inp_cor_hists);
721  DeleteErrorCorHistograms(err_out_cor_hists);
722 }
void FillErrorDataCorHistograms(double errors[4], double var[5], TH2D *err_cor_hists[4][5])
void DeleteErrorHists(TH1D *err_hists[4])
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
void AllocateErrorOutputCorHists(TH2D *err_out_cor_hists[4][5])
virtual Double_t Eval(const Double_t *x, const Double_t *coeff=nullptr) const
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
void AllocateErrorHists(TH1D *err_hists[4])
void DeleteErrorCorHistograms(TH2D *err_cor_hists[4][5])
void FillErrorHistograms(double errors[4], TH1D *err_hists[4])
void AllocateErrorInputCorHists(TH2D *err_inp_cor_hists[4][5])
TMultiDimFet y_parametrisation
polynomial approximation for y
Definition: errors.py:1
TMultiDimFet x_parametrisation
polynomial approximation for x
void WriteHistograms(TH1D *err_hists[4], TH2D *err_inp_cor_hists[4][5], TH2D *err_out_cor_hists[4][5], TFile *f_out, std::string base_out_dir)

◆ TestAperture()

void LHCOpticsApproximator::TestAperture ( TTree *  in_tree,
TTree *  out_tree 
)

x, theta_x, y, theta_y, ksi, mad_accepted, parametriz_accepted

Definition at line 740 of file LHCOpticsApproximator.cc.

References mps_splice::entry, mps_fire::i, and Transport().

742 {
743  if (inp_tree == nullptr || out_tree == nullptr)
744  return;
745 
746  Long64_t entries = inp_tree->GetEntries();
747  double entry[7];
748  double parametrization_out[5];
749 
750  inp_tree->SetBranchAddress("x", &(entry[0]));
751  inp_tree->SetBranchAddress("theta_x", &(entry[1]));
752  inp_tree->SetBranchAddress("y", &(entry[2]));
753  inp_tree->SetBranchAddress("theta_y", &(entry[3]));
754  inp_tree->SetBranchAddress("ksi", &(entry[4]));
755  inp_tree->SetBranchAddress("mad_accept", &(entry[5]));
756  inp_tree->SetBranchAddress("par_accept", &(entry[6]));
757 
758  out_tree->SetBranchAddress("x", &(entry[0]));
759  out_tree->SetBranchAddress("theta_x", &(entry[1]));
760  out_tree->SetBranchAddress("y", &(entry[2]));
761  out_tree->SetBranchAddress("theta_y", &(entry[3]));
762  out_tree->SetBranchAddress("ksi", &(entry[4]));
763  out_tree->SetBranchAddress("mad_accept", &(entry[5]));
764  out_tree->SetBranchAddress("par_accept", &(entry[6]));
765 
766  // int ind=0;
767  for (Long64_t i = 0; i < entries; ++i) {
768  inp_tree->GetEntry(i);
769 
770  //Don't invert the coordinate systems, appertures are defined in the
771  //coordinate system of the beam - perhaps to be changed
772  bool res = Transport(entry, parametrization_out, true, false);
773 
774  if (res)
775  entry[6] = 1.0;
776  else
777  entry[6] = 0.0;
778 
779  out_tree->Fill();
780  }
781 }
Definition: Electron.h:6
bool Transport(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const

◆ Train()

void LHCOpticsApproximator::Train ( TTree *  inp_tree,
std::string  data_prefix = std::string("def"),
polynomials_selection  mode = PREDEFINED,
int  max_degree_x = 10,
int  max_degree_tx = 10,
int  max_degree_y = 10,
int  max_degree_ty = 10,
bool  common_terms = false,
double *  prec = nullptr 
)

Definition at line 286 of file LHCOpticsApproximator.cc.

References TMultiDimFet::AddRow(), coord_names, mps_fire::i, InitializeApproximators(), ALCARECOPromptCalibProdSiPixelAli0T_cff::mode, out_polynomials, PrintInputRange(), s_begin_, s_end_, AlCaHLTBitMon_QueryRunRegistry::string, theta_x_parametrisation, theta_y_parametrisation, trained_, x_parametrisation, and y_parametrisation.

294  {
295  if (inp_tree == nullptr)
296  return;
297 
298  InitializeApproximators(mode, max_degree_x, max_degree_tx, max_degree_y, max_degree_ty, common_terms);
299 
300  //in-variables
301  //x_in, theta_x_in, y_in, theta_y_in, ksi_in, s_in
302  double in_var[6];
303 
304  //out-variables
305  //x_out, theta_x_out, y_out, theta_y_out, ksi_out, s_out, valid_out;
306  double out_var[7];
307 
308  //in- out-lables
309  std::string x_in_lab = "x_in";
310  std::string theta_x_in_lab = "theta_x_in";
311  std::string y_in_lab = "y_in";
312  std::string theta_y_in_lab = "theta_y_in";
313  std::string ksi_in_lab = "ksi_in";
314  std::string s_in_lab = "s_in";
315 
316  std::string x_out_lab = data_prefix + "_x_out";
317  std::string theta_x_out_lab = data_prefix + "_theta_x_out";
318  std::string y_out_lab = data_prefix + "_y_out";
319  std::string theta_y_out_lab = data_prefix + "_theta_y_out";
320  std::string ksi_out_lab = data_prefix + "_ksi_out";
321  std::string s_out_lab = data_prefix + "_s_out";
322  std::string valid_out_lab = data_prefix + "_valid_out";
323 
324  //disable not needed branches to speed up the readin
325  inp_tree->SetBranchStatus("*", false); //disable all branches
326  inp_tree->SetBranchStatus(x_in_lab.c_str(), true);
327  inp_tree->SetBranchStatus(theta_x_in_lab.c_str(), true);
328  inp_tree->SetBranchStatus(y_in_lab.c_str(), true);
329  inp_tree->SetBranchStatus(theta_y_in_lab.c_str(), true);
330  inp_tree->SetBranchStatus(ksi_in_lab.c_str(), true);
331  inp_tree->SetBranchStatus(x_out_lab.c_str(), true);
332  inp_tree->SetBranchStatus(theta_x_out_lab.c_str(), true);
333  inp_tree->SetBranchStatus(y_out_lab.c_str(), true);
334  inp_tree->SetBranchStatus(theta_y_out_lab.c_str(), true);
335  inp_tree->SetBranchStatus(ksi_out_lab.c_str(), true);
336  inp_tree->SetBranchStatus(valid_out_lab.c_str(), true);
337 
338  //set input data adresses
339  inp_tree->SetBranchAddress(x_in_lab.c_str(), &(in_var[0]));
340  inp_tree->SetBranchAddress(theta_x_in_lab.c_str(), &(in_var[1]));
341  inp_tree->SetBranchAddress(y_in_lab.c_str(), &(in_var[2]));
342  inp_tree->SetBranchAddress(theta_y_in_lab.c_str(), &(in_var[3]));
343  inp_tree->SetBranchAddress(ksi_in_lab.c_str(), &(in_var[4]));
344  inp_tree->SetBranchAddress(s_in_lab.c_str(), &(in_var[5]));
345 
346  //set output data adresses
347  inp_tree->SetBranchAddress(x_out_lab.c_str(), &(out_var[0]));
348  inp_tree->SetBranchAddress(theta_x_out_lab.c_str(), &(out_var[1]));
349  inp_tree->SetBranchAddress(y_out_lab.c_str(), &(out_var[2]));
350  inp_tree->SetBranchAddress(theta_y_out_lab.c_str(), &(out_var[3]));
351  inp_tree->SetBranchAddress(ksi_out_lab.c_str(), &(out_var[4]));
352  inp_tree->SetBranchAddress(s_out_lab.c_str(), &(out_var[5]));
353  inp_tree->SetBranchAddress(valid_out_lab.c_str(), &(out_var[6]));
354 
355  Long64_t entries = inp_tree->GetEntries();
356  if (entries > 0) {
357  inp_tree->SetBranchStatus(s_in_lab.c_str(), true);
358  inp_tree->SetBranchStatus(s_out_lab.c_str(), true);
359  inp_tree->GetEntry(0);
360  s_begin_ = in_var[5];
361  s_end_ = out_var[5];
362  inp_tree->SetBranchStatus(s_in_lab.c_str(), false);
363  inp_tree->SetBranchStatus(s_out_lab.c_str(), false);
364  }
365 
366  //set input and output variables for fitting
367  for (Long64_t i = 0; i < entries; ++i) {
368  inp_tree->GetEntry(i);
369  if (out_var[6] != 0) //if out data valid
370  {
371  x_parametrisation.AddRow(in_var, out_var[0], 0);
372  theta_x_parametrisation.AddRow(in_var, out_var[1], 0);
373  y_parametrisation.AddRow(in_var, out_var[2], 0);
374  theta_y_parametrisation.AddRow(in_var, out_var[3], 0);
375  }
376  }
377 
378  edm::LogInfo("LHCOpticsApproximator") << "Optical functions parametrizations from " << s_begin_ << " to " << s_end_
379  << "\n";
380  PrintInputRange();
381  for (int i = 0; i < 4; i++) {
382  double best_precision = 0.0;
383  if (prec)
384  best_precision = prec[i];
385  out_polynomials[i]->FindParameterization(best_precision);
386  edm::LogInfo("LHCOpticsApproximator") << "Out variable " << coord_names[i] << " polynomial"
387  << "\n";
388  out_polynomials[i]->PrintPolynomialsSpecial("M");
389  edm::LogInfo("LHCOpticsApproximator") << "\n";
390  }
391 
392  trained_ = true;
393 }
void InitializeApproximators(polynomials_selection mode, int max_degree_x, int max_degree_tx, int max_degree_y, int max_degree_ty, bool common_terms)
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
double s_begin_
begin of transport along the reference orbit
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
bool trained_
trained polynomials
std::vector< TMultiDimFet * > out_polynomials
virtual void AddRow(const Double_t *x, Double_t D, Double_t E=0)
Log< level::Info, false > LogInfo
double s_end_
end of transport along the reference orbit
std::vector< std::string > coord_names
pointers to polynomials
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ Transport() [1/2]

bool LHCOpticsApproximator::Transport ( const double *  in,
double *  out,
bool  check_apertures = false,
bool  invert_beam_coord_sytems = true 
) const

Basic 3D transport method MADX canonical variables IN/OUT: (x, theta_x, y, theta_y, xi) [m, rad, m, rad, 1] returns true if transport possible if theta is calculated from momentum p, use theta_x = p.x() / p.mag() and theta_y = p.y() / p.mag()

Definition at line 102 of file LHCOpticsApproximator.cc.

References apertures_, beam, CheckInputRange(), TMultiDimFet::Eval(), mps_fire::i, recoMuon::in, lhcb1, MillePedeFileConverter_cfg::out, theta_x_parametrisation, theta_y_parametrisation, trained_, x_parametrisation, and y_parametrisation.

Referenced by GetDx(), GetDxds(), GetLineariasedTransportMatrixX(), GetLineariasedTransportMatrixY(), TestAperture(), Transport(), and Transport_m_GeV().

105  {
106  if (in == nullptr || out == nullptr || !trained_)
107  return false;
108 
109  bool res = CheckInputRange(in);
110  double in_corrected[5];
111 
112  if (beam == lhcb1 || !invert_beam_coord_sytems) {
113  in_corrected[0] = in[0];
114  in_corrected[1] = in[1];
115  in_corrected[2] = in[2];
116  in_corrected[3] = in[3];
117  in_corrected[4] = in[4];
118  out[0] = x_parametrisation.Eval(in_corrected);
119  out[1] = theta_x_parametrisation.Eval(in_corrected);
120  out[2] = y_parametrisation.Eval(in_corrected);
121  out[3] = theta_y_parametrisation.Eval(in_corrected);
122  out[4] = in[4];
123  } else {
124  in_corrected[0] = -in[0];
125  in_corrected[1] = -in[1];
126  in_corrected[2] = in[2];
127  in_corrected[3] = in[3];
128  in_corrected[4] = in[4];
129  out[0] = -x_parametrisation.Eval(in_corrected);
130  out[1] = -theta_x_parametrisation.Eval(in_corrected);
131  out[2] = y_parametrisation.Eval(in_corrected);
132  out[3] = theta_y_parametrisation.Eval(in_corrected);
133  out[4] = in[4];
134  }
135 
136  if (check_apertures) {
137  for (unsigned int i = 0; i < apertures_.size(); i++) {
138  res = res && apertures_[i].CheckAperture(in);
139  }
140  }
141  return res;
142 }
std::vector< LHCApertureApproximator > apertures_
apertures on the way
TMultiDimFet theta_x_parametrisation
polynomial approximation for theta_x
virtual Double_t Eval(const Double_t *x, const Double_t *coeff=nullptr) const
Definition: Electron.h:6
TMultiDimFet theta_y_parametrisation
polynomial approximation for theta_y
bool trained_
trained polynomials
bool CheckInputRange(const double *in, bool invert_beam_coord_sytems=true) const
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ Transport() [2/2]

bool LHCOpticsApproximator::Transport ( const MadKinematicDescriptor in,
MadKinematicDescriptor out,
bool  check_apertures = false,
bool  invert_beam_coord_sytems = true 
) const

Definition at line 220 of file LHCOpticsApproximator.cc.

References recoMuon::in, input, MillePedeFileConverter_cfg::out, trained_, and Transport().

223  {
224  if (in == nullptr || out == nullptr || !trained_)
225  return false;
226 
227  Double_t input[5];
228  Double_t output[5];
229  input[0] = in->x;
230  input[1] = in->theta_x;
231  input[2] = in->y;
232  input[3] = in->theta_y;
233  input[4] = in->ksi;
234 
235  //transport inverts the coordinate systems
236  if (!Transport(input, output, check_apertures, invert_beam_coord_sytems)) {
237  return false;
238  }
239 
240  out->x = output[0];
241  out->theta_x = output[1];
242  out->y = output[2];
243  out->theta_y = output[3];
244  out->ksi = output[4];
245 
246  return true;
247 }
static std::string const input
Definition: EdmProvDump.cc:50
bool trained_
trained polynomials
bool Transport(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const
Definition: output.py:1

◆ Transport2D()

bool LHCOpticsApproximator::Transport2D ( const double *  in,
double *  out,
bool  check_apertures = false,
bool  invert_beam_coord_sytems = true 
) const

Basic 2D transport method MADX canonical variables IN : (x, theta_x, y, theta_y, xi) [m, rad, m, rad, 1] OUT : (x, y) [m, m] returns true if transport possible

Definition at line 144 of file LHCOpticsApproximator.cc.

References EcalCondDBWriter_cfi::beam, mps_fire::i, recoMuon::in, and MillePedeFileConverter_cfg::out.

Referenced by GetLinearApproximation().

149 {
150  if (in == nullptr || out == nullptr || !trained_)
151  return false;
152 
153  bool res = CheckInputRange(in);
154  double in_corrected[5];
155 
156  if (beam == lhcb1 || !invert_beam_coord_sytems) {
157  in_corrected[0] = in[0];
158  in_corrected[1] = in[1];
159  in_corrected[2] = in[2];
160  in_corrected[3] = in[3];
161  in_corrected[4] = in[4];
162  out[0] = x_parametrisation.Eval(in_corrected);
163  out[1] = y_parametrisation.Eval(in_corrected);
164  } else {
165  in_corrected[0] = -in[0];
166  in_corrected[1] = -in[1];
167  in_corrected[2] = in[2];
168  in_corrected[3] = in[3];
169  in_corrected[4] = in[4];
170  out[0] = -x_parametrisation.Eval(in_corrected);
171  out[1] = y_parametrisation.Eval(in_corrected);
172  }
173 
174  if (check_apertures) {
175  for (unsigned int i = 0; i < apertures_.size(); i++) {
176  res = res && apertures_[i].CheckAperture(in);
177  }
178  }
179  return res;
180 }
std::vector< LHCApertureApproximator > apertures_
apertures on the way
virtual Double_t Eval(const Double_t *x, const Double_t *coeff=nullptr) const
Definition: Electron.h:6
bool trained_
trained polynomials
bool CheckInputRange(const double *in, bool invert_beam_coord_sytems=true) const
TMultiDimFet y_parametrisation
polynomial approximation for y
TMultiDimFet x_parametrisation
polynomial approximation for x

◆ Transport_m_GeV()

bool LHCOpticsApproximator::Transport_m_GeV ( double  in_pos[3],
double  in_momentum[3],
double  out_pos[3],
double  out_momentum[3],
bool  check_apertures,
double  z2_z1_dist 
) const

pos, momentum: x,y,z; pos in m, momentum in GeV/c

Definition at line 182 of file LHCOpticsApproximator.cc.

References mps_fire::i, recoMuon::in, nominal_beam_momentum_, MillePedeFileConverter_cfg::out, mathSSE::sqrt(), and Transport().

Referenced by TotemTransport::transportProton().

187  {
188  double in[5];
189  double out[5];
190  double part_mom = 0.0;
191  for (int i = 0; i < 3; ++i)
192  part_mom += in_momentum[i] * in_momentum[i];
193 
194  part_mom = std::sqrt(part_mom);
195 
196  in[0] = in_pos[0];
197  in[1] = in_momentum[0] / nominal_beam_momentum_;
198  in[2] = in_pos[1];
199  in[3] = in_momentum[1] / nominal_beam_momentum_;
201 
202  if (!Transport(in, out, check_apertures, true)) {
203  return false;
204  }
205 
206  out_pos[0] = out[0];
207  out_pos[1] = out[2];
208  out_pos[2] = in_pos[2] + z2_z1_dist;
209 
210  out_momentum[0] = out[1] * nominal_beam_momentum_;
211  out_momentum[1] = out[3] * nominal_beam_momentum_;
212  double part_out_total_mom = (out[4] + 1) * nominal_beam_momentum_;
213  out_momentum[2] = std::sqrt(part_out_total_mom * part_out_total_mom - out_momentum[0] * out_momentum[0] -
214  out_momentum[1] * out_momentum[1]);
215  out_momentum[2] = TMath::Sign(out_momentum[2], in_momentum[2]);
216 
217  return true;
218 }
double nominal_beam_momentum_
GeV/c.
T sqrt(T t)
Definition: SSEVec.h:23
bool Transport(const double *in, double *out, bool check_apertures=false, bool invert_beam_coord_sytems=true) const

◆ WriteHistograms()

void LHCOpticsApproximator::WriteHistograms ( TH1D *  err_hists[4],
TH2D *  err_inp_cor_hists[4][5],
TH2D *  err_out_cor_hists[4][5],
TFile *  f_out,
std::string  base_out_dir 
)
private

Definition at line 869 of file LHCOpticsApproximator.cc.

References mps_fire::i.

Referenced by Test().

873  {
874  if (f_out == nullptr)
875  return;
876 
877  f_out->cd();
878  if (!gDirectory->cd(base_out_dir.c_str()))
879  gDirectory->mkdir(base_out_dir.c_str());
880 
881  gDirectory->cd(base_out_dir.c_str());
882  gDirectory->mkdir(this->GetName());
883  gDirectory->cd(this->GetName());
884  gDirectory->mkdir("x");
885  gDirectory->mkdir("theta_x");
886  gDirectory->mkdir("y");
887  gDirectory->mkdir("theta_y");
888 
889  gDirectory->cd("x");
890  err_hists[0]->Write("", TObject::kWriteDelete);
891  for (int i = 0; i < 5; i++) {
892  err_inp_cor_hists[0][i]->Write("", TObject::kWriteDelete);
893  err_out_cor_hists[0][i]->Write("", TObject::kWriteDelete);
894  }
895 
896  gDirectory->cd("..");
897  gDirectory->cd("theta_x");
898  err_hists[1]->Write("", TObject::kWriteDelete);
899  for (int i = 0; i < 5; i++) {
900  err_inp_cor_hists[1][i]->Write("", TObject::kWriteDelete);
901  err_out_cor_hists[1][i]->Write("", TObject::kWriteDelete);
902  }
903 
904  gDirectory->cd("..");
905  gDirectory->cd("y");
906  err_hists[2]->Write("", TObject::kWriteDelete);
907  for (int i = 0; i < 5; i++) {
908  err_inp_cor_hists[2][i]->Write("", TObject::kWriteDelete);
909  err_out_cor_hists[2][i]->Write("", TObject::kWriteDelete);
910  }
911 
912  gDirectory->cd("..");
913  gDirectory->cd("theta_y");
914  err_hists[3]->Write("", TObject::kWriteDelete);
915  for (int i = 0; i < 5; i++) {
916  err_inp_cor_hists[3][i]->Write("", TObject::kWriteDelete);
917  err_out_cor_hists[3][i]->Write("", TObject::kWriteDelete);
918  }
919  gDirectory->cd("..");
920  gDirectory->cd("..");
921 }

Friends And Related Function Documentation

◆ ProtonTransportFunctionsESSource

friend class ProtonTransportFunctionsESSource
friend

Definition at line 173 of file LHCOpticsApproximator.h.

Member Data Documentation

◆ apertures_

std::vector<LHCApertureApproximator> LHCOpticsApproximator::apertures_
private

apertures on the way

Definition at line 170 of file LHCOpticsApproximator.h.

Referenced by AddRectEllipseAperture(), Init(), LHCOpticsApproximator(), operator=(), and Transport().

◆ beam

beam_type LHCOpticsApproximator::beam
private

◆ coord_names

std::vector<std::string> LHCOpticsApproximator::coord_names
private

pointers to polynomials

Definition at line 169 of file LHCOpticsApproximator.h.

Referenced by Init(), PrintInputRange(), PrintOpticalFunctions(), and Train().

◆ nominal_beam_energy_

double LHCOpticsApproximator::nominal_beam_energy_
private

GeV.

Definition at line 165 of file LHCOpticsApproximator.h.

Referenced by LHCOpticsApproximator(), and operator=().

◆ nominal_beam_momentum_

double LHCOpticsApproximator::nominal_beam_momentum_
private

GeV/c.

Definition at line 166 of file LHCOpticsApproximator.h.

Referenced by LHCOpticsApproximator(), operator=(), and Transport_m_GeV().

◆ out_polynomials

std::vector<TMultiDimFet *> LHCOpticsApproximator::out_polynomials
private

Definition at line 168 of file LHCOpticsApproximator.h.

Referenced by Init(), PrintOpticalFunctions(), and Train().

◆ s_begin_

double LHCOpticsApproximator::s_begin_
private

begin of transport along the reference orbit

Definition at line 162 of file LHCOpticsApproximator.h.

Referenced by Init(), LHCOpticsApproximator(), operator=(), and Train().

◆ s_end_

double LHCOpticsApproximator::s_end_
private

end of transport along the reference orbit

Definition at line 163 of file LHCOpticsApproximator.h.

Referenced by Init(), LHCOpticsApproximator(), operator=(), and Train().

◆ theta_x_parametrisation

TMultiDimFet LHCOpticsApproximator::theta_x_parametrisation
private

polynomial approximation for theta_x

Definition at line 177 of file LHCOpticsApproximator.h.

Referenced by Init(), InitializeApproximators(), operator=(), Test(), Train(), and Transport().

◆ theta_y_parametrisation

TMultiDimFet LHCOpticsApproximator::theta_y_parametrisation
private

polynomial approximation for theta_y

Definition at line 179 of file LHCOpticsApproximator.h.

Referenced by Init(), InitializeApproximators(), operator=(), Test(), Train(), and Transport().

◆ trained_

bool LHCOpticsApproximator::trained_
private

trained polynomials

Definition at line 167 of file LHCOpticsApproximator.h.

Referenced by Init(), LHCOpticsApproximator(), operator=(), Train(), and Transport().

◆ x_parametrisation

TMultiDimFet LHCOpticsApproximator::x_parametrisation
private

polynomial approximation for x

Definition at line 176 of file LHCOpticsApproximator.h.

Referenced by Init(), InitializeApproximators(), operator=(), ParameterOutOfRangePenalty(), PrintInputRange(), Test(), Train(), and Transport().

◆ y_parametrisation

TMultiDimFet LHCOpticsApproximator::y_parametrisation
private

polynomial approximation for y

Definition at line 178 of file LHCOpticsApproximator.h.

Referenced by Init(), InitializeApproximators(), operator=(), Test(), Train(), and Transport().