CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
RecoProducerFP420.cc
Go to the documentation of this file.
3 
6 
7 #include <cmath>
8 
9 //#include <iostream>
10 
12  // Create LHC beam line
13  // double length = param.getParameter<double>("BeamLineLength");
14  // std::cout << " BeamLineLength = " << length << std::endl;
15 
16  length = conf_.getParameter<double>("BeamLineLength"); //m
17  verbosity = conf_.getUntrackedParameter<int>("VerbosityLevel");
20 
21  edm::LogInfo("RecoProducerFP420") << "RecoProducerFP420 parameters: \n"
22  << " Verbosity: " << verbosity << "\n"
23  << " length: " << length << "\n";
24  if (verbosity > 1) {
25  std::cout << " RecoProducerFP420: constructor " << std::endl;
26  std::cout << " BeamLineLength: " << length << std::endl;
27  }
28 
29  // edm::FileInPath b1("SimTransport/HectorData/twiss_ip5_b1_v6.5.txt");
30  // edm::FileInPath b2("SimTransport/HectorData/twiss_ip5_b2_v6.5.txt");
31 
32  edm::FileInPath b1(beam1filename.c_str());
33  edm::FileInPath b2(beam2filename.c_str());
34 
35  m_beamline1 = new H_BeamLine(1, length + 0.1); // (direction, length)
36  m_beamline2 = new H_BeamLine(-1, length + 0.1); //
37 
38  try {
39  m_beamline1->fill(b1.fullPath(), 1, "IP5");
40  m_beamline2->fill(b2.fullPath(), -1, "IP5");
41  } catch (const edm::Exception& e) {
42  std::string msg = e.what();
43  msg += " caught in RecoProducerFP420... \nERROR: Could not locate SimTransport/HectorData data files.";
44  edm::LogError("DataNotFound") << msg;
45  }
46 
47  m_beamline1->offsetElements(120, -0.097);
48  m_beamline2->offsetElements(120, +0.097);
49 
50  m_beamline1->calcMatrix();
51  m_beamline2->calcMatrix();
52 
53  edm::LogInfo("RecoProducerFP420") << "==============================\n";
54 }
55 
57 
58 std::vector<RecoFP420> RecoProducerFP420::reconstruct(
59  int direction, double x1_420, double y1_420, double x2_420, double y2_420, double z1_420, double z2_420) {
60  // ==================================================================================
61  // ==================================================================================
62  std::vector<RecoFP420> rhits;
63  int restracks = 10; // max # tracks
64  rhits.reserve(restracks);
65  rhits.clear();
66  // ==================================================================================
67  // trivial (TM), angle compensation (AM) and position compensation (PM) methods
68  // #define TM 1 #define AM 2 #define PM 3
69  m_tx0 = -100., m_ty0 = -100., m_x0 = -100., m_y0 = -100.;
70  if (direction == 1) {
71  m_rp420_f = new H_RecRPObject(z1_420 * 0.001, z2_420 * 0.001, *m_beamline1); // m
72  if (m_rp420_f) {
73  if (verbosity > 1) {
74  std::cout << " RecoProducerFP420: input coord. in um " << std::endl;
75  std::cout << " x1_420: " << x1_420 << " y1_420: " << y1_420 << std::endl;
76  std::cout << " x2_420: " << x2_420 << " y2_420: " << y2_420 << std::endl;
77  }
78  m_rp420_f->setPositions(x1_420, y1_420, x2_420, y2_420); //input coord. in um
79  m_e = m_rp420_f->getE(AM); // GeV
80  // std::cout << " m_e1: " << m_rp420_f->getE( TM ) << std::endl;
81  // std::cout << " m_e2: " << m_rp420_f->getE( AM ) << std::endl;
82  m_tx0 = m_rp420_f->getTXIP(); // urad
83  m_ty0 = m_rp420_f->getTYIP(); // urad
84  m_x0 = m_rp420_f->getX0(); // um
85  m_y0 = m_rp420_f->getY0(); // um
86  m_q2 = m_rp420_f->getQ2(); // GeV^2
87  } // if ( m_rp420_f
88  } // if ( dire
89  else if (direction == 2) {
90  m_rp420_b = new H_RecRPObject(z1_420 * 0.001, z2_420 * 0.001, *m_beamline2); // m
91  if (m_rp420_b) {
92  m_rp420_b->setPositions(x1_420, y1_420, x2_420, y2_420); // input coord. in um
93  m_e = m_rp420_b->getE(AM); // GeV
94  m_tx0 = m_rp420_b->getTXIP(); // urad
95  m_ty0 = m_rp420_b->getTYIP(); // urad
96  m_x0 = m_rp420_b->getX0(); // um
97  m_y0 = m_rp420_b->getY0(); // um
98  m_q2 = m_rp420_b->getQ2(); // GeV^2
99  } // if ( m_rp420_b
100  } else {
101  return rhits;
102  }
103 
104  // ==============================
105  if (verbosity > 1) {
106  std::cout << " RecoProducerFP420: rhits.push_back " << std::endl;
107  std::cout << " m_e: " << m_e << std::endl;
108  std::cout << " m_x0: " << m_x0 << std::endl;
109  std::cout << " m_y0: " << m_y0 << std::endl;
110  std::cout << " m_tx0: " << m_tx0 << std::endl;
111  std::cout << " m_ty0: " << m_ty0 << std::endl;
112  std::cout << " m_q2: " << m_q2 << std::endl;
113  std::cout << " direction: " << direction << std::endl;
114  }
115  rhits.push_back(RecoFP420(m_e, m_x0, m_y0, m_tx0, m_ty0, m_q2, direction));
116  // ==============================
118  return rhits;
119  //============
120 }
T getUntrackedParameter(std::string const &, T const &) const
std::vector< RecoFP420 > rhits
std::string beam1filename
H_BeamLine * m_beamline2
Log< level::Error, false > LogError
edm::ParameterSet conf_
RecoProducerFP420(const edm::ParameterSet &conf)
char const * what() const noexceptoverride
Definition: Exception.cc:103
H_RecRPObject * m_rp420_f
Log< level::Info, false > LogInfo
std::string beam2filename
std::vector< RecoFP420 > reconstruct(int, double, double, double, double, double, double)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
tuple msg
Definition: mps_check.py:285
H_RecRPObject * m_rp420_b
tuple cout
Definition: gather_cfg.py:144
H_BeamLine * m_beamline1
static constexpr float b2
virtual ~RecoProducerFP420()
static constexpr float b1