CMS 3D CMS Logo

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