CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TotemTransport.cc
Go to the documentation of this file.
3 #include <CLHEP/Units/GlobalSystemOfUnits.h>
4 #include <CLHEP/Random/RandGauss.h>
5 #include "TLorentzVector.h"
6 #include "TFile.h"
7 
8 #include <cmath>
9 
11  : BaseProtonTransport(iConfig),
12  m_model_ip_150_r_name(iConfig.getParameter<std::string>("Model_IP_150_R_Name")),
13  m_model_ip_150_l_name(iConfig.getParameter<std::string>("Model_IP_150_L_Name")),
14  m_beampipe_aperture_radius(iConfig.getParameter<double>("BeampipeApertureRadius")) {
16  beam1Filename_ = iConfig.getParameter<std::string>("Beam1Filename");
17  beam2Filename_ = iConfig.getParameter<std::string>("Beam2Filename");
18  fCrossingAngleX_45 = iConfig.getParameter<double>("halfCrossingAngleSector45");
19  fCrossingAngleX_56 = iConfig.getParameter<double>("halfCrossingAngleSector56");
20  beamEnergy_ = iConfig.getParameter<double>("BeamEnergy");
21  m_sigmaSTX = iConfig.getParameter<double>("BeamDivergenceX");
22  m_sigmaSTY = iConfig.getParameter<double>("BeamDivergenceY");
23  m_sigmaSX = iConfig.getParameter<double>("BeamSigmaX");
24  m_sigmaSY = iConfig.getParameter<double>("BeamSigmaY");
25  m_sig_E = iConfig.getParameter<double>("BeamEnergyDispersion");
26 
27  if (fPPSRegionStart_56 > 0)
28  fPPSRegionStart_56 *= -1; // make sure sector 56 has negative position, as TOTEM convention
29 
30  edm::LogInfo("TotemTransport") << "=============================================================================\n"
31  << " Bulding LHC Proton transporter based on TOTEM model\n"
32  << "=============================================================================\n";
33 
36 
37  if (m_aprox_ip_150_r == nullptr || m_aprox_ip_150_l == nullptr) {
38  edm::LogError("TotemTransport") << "Parameterisation " << m_model_ip_150_r_name << " or " << m_model_ip_150_l_name
39  << " missing in file. Cannot proceed. ";
40  exit(1);
41  }
42  edm::LogInfo("TotemTransport") << "Parameterizations read from file, pointers:" << m_aprox_ip_150_r << " "
43  << m_aprox_ip_150_l << " ";
44 }
45 //
46 // this method is the same for all propagator, but since transportProton is different for each derived class
47 // it needes to be overriden
48 //
50  const edm::EventSetup& iSetup,
51  CLHEP::HepRandomEngine* _engine) {
52  this->clear();
53 
54  engine_ = _engine; // the engine needs to be updated for each event
55 
56  for (HepMC::GenEvent::particle_const_iterator eventParticle = evt->particles_begin();
57  eventParticle != evt->particles_end();
58  ++eventParticle) {
59  if (!((*eventParticle)->status() == 1 && (*eventParticle)->pdg_id() == 2212))
60  continue;
61 
62  if (!(fabs((*eventParticle)->momentum().eta()) > etaCut_ && fabs((*eventParticle)->momentum().pz()) > momentumCut_))
63  continue; // discard protons outside kinematic acceptance
64 
65  unsigned int line = (*eventParticle)->barcode();
66  HepMC::GenParticle* gpart = (*eventParticle);
67  if (gpart->pdg_id() != 2212)
68  continue; // only transport stable protons
69  if (gpart->status() != 1)
70  continue;
71  if (m_beamPart.find(line) != m_beamPart.end()) // assures this protons has not been already propagated
72  continue;
73 
74  transportProton(gpart);
75  }
76 }
77 //
78 //
79 // here comes the real thing
80 //
81 //
83  //
84  //
85  edm::LogInfo("TotemTransport") << "Starting proton transport using TOTEM method\n";
86  //
87  ApplyBeamCorrection(const_cast<HepMC::GenParticle*>(in_trk));
88 
89  const HepMC::GenVertex* in_pos = in_trk->production_vertex();
90  const HepMC::FourVector in_mom = in_trk->momentum();
91  //
92  // ATTENTION: HepMC uses mm, vertex config of CMS uses cm and SimTransport uses mm
93  //
94  double in_position[3] = {in_pos->position().x(), in_pos->position().y(), in_pos->position().z()}; //in LHC ref. frame
95 
96  double fCrossingAngleX = (in_mom.z() > 0) ? fCrossingAngleX_45 : fCrossingAngleX_56;
97 
98  // Move the position to z=0. Do it in the CMS ref frame. Totem parameterization does the rotation internatlly
99  in_position[0] =
100  in_position[0] - in_position[2] * (in_mom.x() / in_mom.z() - fCrossingAngleX * urad); // in CMS ref. frame
101  in_position[1] = in_position[1] - in_position[2] * (in_mom.y() / (in_mom.z()));
102  in_position[2] = 0.;
103  double in_momentum[3] = {in_mom.x(), in_mom.y(), in_mom.z()};
104  double out_position[3];
105  double out_momentum[3];
106  edm::LogInfo("TotemTransport") << "before transport ->"
107  << " position: " << in_position[0] << ", " << in_position[1] << ", " << in_position[2]
108  << " momentum: " << in_momentum[0] << ", " << in_momentum[1] << ", " << in_momentum[2];
109 
110  LHCOpticsApproximator* approximator_ = nullptr;
111  double m_Zin_;
112  double m_Zout_;
113  if (in_mom.z() > 0) {
114  approximator_ = m_aprox_ip_150_l;
115  m_Zin_ = 0.0; // Totem propagations assumes the starting point at 0 (zero)
116  m_Zout_ = fPPSRegionStart_45;
117  } else {
118  approximator_ = m_aprox_ip_150_r;
119  m_Zin_ = 0.0; // Totem propagations assumes the starting point at 0 (zero)
120  m_Zout_ = fPPSRegionStart_56;
121  }
122 
123  bool invert_beam_coord_system =
124  true; // it doesn't matter the option here, it is hard coded as TRUE inside LHCOpticsApproximator!
125 
126  bool tracked = approximator_->Transport_m_GeV(
127  in_position, in_momentum, out_position, out_momentum, invert_beam_coord_system, m_Zout_ - m_Zin_);
128 
129  if (!tracked)
130  return false;
131 
132  edm::LogInfo("TotemTransport") << "after transport -> "
133  << "position: " << out_position[0] << ", " << out_position[1] << ", "
134  << out_position[2] << "momentum: " << out_momentum[0] << ", " << out_momentum[1]
135  << ", " << out_momentum[2];
136 
137  if (out_position[0] * out_position[0] + out_position[1] * out_position[1] >
139  edm::LogInfo("TotemTransport") << "Proton ouside beampipe";
140  edm::LogInfo("TotemTransport") << "===== END Transport "
141  << "====================";
142  return false;
143  }
144 
145  TVector3 out_pos(out_position[0] * meter, out_position[1] * meter, out_position[2] * meter);
146  TVector3 out_mom(out_momentum[0], out_momentum[1], out_momentum[2]);
147 
148  if (verbosity_) {
149  LogDebug("TotemTransport") << "output -> position: ";
150  out_pos.Print();
151  LogDebug("TotemTransport") << " momentum: ";
152  out_mom.Print();
153  }
154 
155  double px = -out_momentum[0]; // tote calculates px by means of TH_X, which is in the LHC ref. frame.
156  double py = out_momentum[1]; // this need to be checked again, since it seems an invertion is occuring in the prop.
157  double pz =
158  out_momentum[2]; // totem calculates output pz already in the CMS ref. frame, it doesn't need to be converted
159  double e = sqrt(px * px + py * py + pz * pz + ProtonMassSQ);
160  TLorentzVector p_out(px, py, pz, e);
161  double x1_ctpps = -out_position[0] * meter; // Totem parameterization uses meter, one need it in millimeter
162  double y1_ctpps = out_position[1] * meter;
163 
164  unsigned int line = in_trk->barcode();
165 
166  if (verbosity_)
167  LogDebug("TotemTransport") << "TotemTransport:filterPPS: barcode = " << line << " x= " << x1_ctpps
168  << " y= " << y1_ctpps;
169 
170  m_beamPart[line] = p_out;
171  m_xAtTrPoint[line] = x1_ctpps;
172  m_yAtTrPoint[line] = y1_ctpps;
173  return true;
174 }
175 //
177  const std::string& rootfile) {
178  edm::FileInPath fileName(rootfile.c_str());
179  TFile* f = TFile::Open(fileName.fullPath().c_str(), "read");
180  if (!f) {
181  edm::LogError("TotemTransport") << "File " << fileName << " not found. Exiting.";
182  return nullptr;
183  }
184  edm::LogInfo("TotemTransport") << "Root file opened, pointer:" << f;
185 
186  // read parametrization
187  LHCOpticsApproximator* aprox = (LHCOpticsApproximator*)f->Get(m_model_name.c_str());
188  f->Close();
189  return aprox;
190 }
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
static const double urad
std::string m_model_ip_150_r_name
Log< level::Error, false > LogError
LHCOpticsApproximator * m_aprox_ip_150_l
std::string m_model_ip_150_l_name
TotemTransport(const edm::ParameterSet &ps)
T sqrt(T t)
Definition: SSEVec.h:19
Class finds the parametrisation of MADX proton transport and transports the protons according to it 5...
CLHEP::HepRandomEngine * engine_
Log< level::Info, false > LogInfo
bool transportProton(const HepMC::GenParticle *)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::map< unsigned int, double > m_xAtTrPoint
void ApplyBeamCorrection(HepMC::GenParticle *p)
double m_beampipe_aperture_radius
std::map< unsigned int, TLorentzVector > m_beamPart
LHCOpticsApproximator * ReadParameterization(const std::string &, const std::string &)
LHCOpticsApproximator * m_aprox_ip_150_r
void process(const HepMC::GenEvent *ev, const edm::EventSetup &es, CLHEP::HepRandomEngine *engine) override
static const double ProtonMassSQ
std::map< unsigned int, double > m_yAtTrPoint
#define LogDebug(id)