CMS 3D CMS Logo

OpticalFunctionsTransport.cc
Go to the documentation of this file.
3 
5  : BaseProtonTransport(iConfig),
6  lhcInfoToken_(iC.esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("lhcInfoLabel")))),
7  beamParametersToken_(iC.esConsumes()),
8  opticsToken_(iC.esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("opticsLabel")))),
9  beamspotToken_(iC.esConsumes()),
10  useEmpiricalApertures_(iConfig.getParameter<bool>("useEmpiricalApertures")),
11  empiricalAperture45_xi0_int_(iConfig.getParameter<double>("empiricalAperture45_xi0_int")),
12  empiricalAperture45_xi0_slp_(iConfig.getParameter<double>("empiricalAperture45_xi0_slp")),
13  empiricalAperture45_a_int_(iConfig.getParameter<double>("empiricalAperture45_a_int")),
14  empiricalAperture45_a_slp_(iConfig.getParameter<double>("empiricalAperture45_a_slp")),
15  empiricalAperture56_xi0_int_(iConfig.getParameter<double>("empiricalAperture56_xi0_int")),
16  empiricalAperture56_xi0_slp_(iConfig.getParameter<double>("empiricalAperture56_xi0_slp")),
17  empiricalAperture56_a_int_(iConfig.getParameter<double>("empiricalAperture56_a_int")),
18  empiricalAperture56_a_slp_(iConfig.getParameter<double>("empiricalAperture56_a_slp")) {
20 }
22 
24  const edm::EventSetup& iSetup,
25  CLHEP::HepRandomEngine* engine) {
26  clear();
27 
28  lhcInfo_ = &iSetup.getData(lhcInfoToken_);
31  beamspot_ = &iSetup.getData(beamspotToken_);
32 
33  // Choose the optical function corresponding to the first station ono each side (it is in lhc ref. frame)
34  optFunctionId45_ = 0;
35  optFunctionId56_ = 0;
36  for (const auto& ofp : (*opticalFunctions_)) {
37  if (ofp.second.getScoringPlaneZ() < 0) {
38  if (optFunctionId45_ == 0)
39  optFunctionId45_ = ofp.first;
40  if (opticalFunctions_->at(optFunctionId45_).getScoringPlaneZ() < ofp.second.getScoringPlaneZ())
41  optFunctionId45_ = ofp.first;
42  }
43  if (ofp.second.getScoringPlaneZ() > 0) {
44  if (optFunctionId56_ == 0)
45  optFunctionId56_ = ofp.first;
46  if (opticalFunctions_->at(optFunctionId56_).getScoringPlaneZ() > ofp.second.getScoringPlaneZ())
47  optFunctionId56_ = ofp.first;
48  }
49  }
50  //
51  engine_ = engine; // the engine needs to be updated for each event
52 
53  for (HepMC::GenEvent::particle_const_iterator eventParticle = evt->particles_begin();
54  eventParticle != evt->particles_end();
55  ++eventParticle) {
56  if (!((*eventParticle)->status() == 1 && (*eventParticle)->pdg_id() == 2212))
57  continue;
58 
59  if (!(fabs((*eventParticle)->momentum().eta()) > etaCut_ && fabs((*eventParticle)->momentum().pz()) > momentumCut_))
60  continue; // discard protons outside kinematic acceptance
61 
62  unsigned int line = (*eventParticle)->barcode();
63  HepMC::GenParticle* gpart = (*eventParticle);
64  if (gpart->pdg_id() != 2212)
65  continue; // only transport stable protons
66  if (gpart->status() != 1)
67  continue;
68  if (m_beamPart.find(line) != m_beamPart.end()) // assures this protons has not been already propagated
69  continue;
70 
71  transportProton(gpart);
72  }
73 }
74 
76  const HepMC::FourVector& vtx_cms = in_trk->production_vertex()->position(); // in mm
77  const HepMC::FourVector& mom_cms = in_trk->momentum();
78 
79  // transformation to LHC/TOTEM convention
80  HepMC::FourVector vtx_lhc(-vtx_cms.x(), vtx_cms.y(), -vtx_cms.z(), vtx_cms.t());
81  HepMC::FourVector mom_lhc(-mom_cms.x(), mom_cms.y(), -mom_cms.z(), mom_cms.t());
82 
83  // determine the LHC arm and related parameters
84  double urad = 1.e-6;
85  double beamMomentum = 0.;
86  double xangle = 0.;
87  double empiricalAperture_xi0_int, empiricalAperture_xi0_slp;
88  double empiricalAperture_a_int, empiricalAperture_a_slp;
89  unsigned int optFunctionId;
90  // get the beam position at the IP in mm and in the LHC ref. frame
91  double vtxXoffset;
92  double vtxYoffset;
94  vtxXoffset = -beamParameters_->getVtxOffsetX45() * cm_to_mm;
95  vtxYoffset = beamParameters_->getVtxOffsetY45() * cm_to_mm;
96  } else {
97  vtxXoffset = -beamspot_->x() * cm_to_mm;
98  vtxYoffset = beamspot_->y() * cm_to_mm;
99  }
100 
101  if (mom_lhc.z() < 0) // sector 45
102  {
103  optFunctionId = optFunctionId45_;
106  empiricalAperture_xi0_int = empiricalAperture45_xi0_int_;
107  empiricalAperture_xi0_slp = empiricalAperture45_xi0_slp_;
108  empiricalAperture_a_int = empiricalAperture45_a_int_;
109  empiricalAperture_a_slp = empiricalAperture45_a_slp_;
110  } else { // sector 56
111  optFunctionId = optFunctionId56_;
114  empiricalAperture_xi0_int = empiricalAperture56_xi0_int_;
115  empiricalAperture_xi0_slp = empiricalAperture56_xi0_slp_;
116  empiricalAperture_a_int = empiricalAperture56_a_int_;
117  empiricalAperture_a_slp = empiricalAperture56_a_slp_;
118  }
119  if (xangle > 1.0)
120  xangle *= urad;
121  // calculate kinematics for optics parametrisation, avoid the aproximation for small angles xangle -> tan(xangle)
122  const double p = mom_lhc.rho();
123  const double xi = 1. - p / beamMomentum;
124  const double th_x_phys = mom_lhc.x() / abs(mom_lhc.z()) - tan(xangle); //"-" in the LHC ref. frame
125  const double th_y_phys = mom_lhc.y() / abs(mom_lhc.z());
126  const double vtx_lhc_eff_x = vtx_lhc.x() - vtx_lhc.z() * (mom_lhc.x() / mom_lhc.z() + tan(xangle)) - (vtxXoffset);
127  const double vtx_lhc_eff_y = vtx_lhc.y() - vtx_lhc.z() * (mom_lhc.y() / mom_lhc.z()) - (vtxYoffset);
128 
129  if (verbosity_) {
130  LogDebug("OpticalFunctionsTransport")
131  << "simu: xi = " << xi << ", th_x_phys = " << th_x_phys << ", th_y_phys = " << th_y_phys
132  << ", vtx_lhc_eff_x = " << vtx_lhc_eff_x << ", vtx_lhc_eff_y = " << vtx_lhc_eff_y;
133  }
134 
135  // check empirical aperture
137  const auto& xangle =
139  const double xi_th = (empiricalAperture_xi0_int + xangle * empiricalAperture_xi0_slp) +
140  (empiricalAperture_a_int + xangle * empiricalAperture_a_slp) * th_x_phys;
141 
142  if (xi > xi_th) {
143  if (verbosity_) {
144  LogDebug("OpticalFunctionsTransport") << "stop because of empirical appertures";
145  }
146  return false;
147  }
148  }
149 
150  // transport the proton into pot/scoring plane
151  auto ofp = opticalFunctions_->at(optFunctionId);
152  CTPPSDetId rpId(optFunctionId);
153  const unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
154 
155  if (verbosity_)
156  LogDebug("OpticalFunctionsTransport") << " RP " << rpDecId << std::endl;
157 
158  // transport proton
160  vtx_lhc_eff_x * cm_to_mm, th_x_phys, vtx_lhc_eff_y * cm_to_mm, th_y_phys, xi}; // conversions: mm -> cm
161 
163  ofp.transport(k_in, k_out, true);
164 
165  // Original code uses mm, but CMS uses cm, so keep it in cm
166  double b_x = k_out.x * cm_to_mm, b_y = k_out.y * cm_to_mm; // conversions: cm -> mm
167  double a_x = k_out.th_x, a_y = k_out.th_y;
168 
169  // if needed, subtract beam position and angle
171  // determine beam position
172  LHCInterpolatedOpticalFunctionsSet::Kinematics k_be_in = {0., -tan(xangle), 0., 0., 0.};
174  ofp.transport(k_be_in, k_be_out, true);
175 
176  a_x -= k_be_out.th_x;
177  a_y -= k_be_out.th_y;
178  b_x -= k_be_out.x * cm_to_mm;
179  b_y -= k_be_out.y * cm_to_mm; // conversions: cm -> mm
180  }
181 
182  const double z_scoringPlane = ofp.getScoringPlaneZ() * cm_to_mm; // conversion: cm --> mm
183 
184  if (verbosity_) {
185  LogDebug("OpticalFunctionsTransport")
186  << " proton transported: a_x = " << a_x << " rad, a_y = " << a_y << " rad, b_x = " << b_x
187  << " mm, b_y = " << b_y << " mm, z = " << z_scoringPlane << " mm";
188  }
189  //
190  // Project the track back to the starting of PPS region in mm
191  b_x -= (abs(z_scoringPlane) - ((z_scoringPlane < 0) ? fPPSRegionStart_45_ : fPPSRegionStart_56_) * 1e3) * a_x;
192  b_y -= (abs(z_scoringPlane) - ((z_scoringPlane < 0) ? fPPSRegionStart_45_ : fPPSRegionStart_56_) * 1e3) * a_y;
193 
194  unsigned int line = in_trk->barcode();
195  double px = -p * a_x;
196  double py = p * a_y;
197  double pz = std::copysign(sqrt(p * p - px * px - py * py), mom_cms.z());
198  double e = sqrt(px * px + py * py + pz * pz + ProtonMassSQ);
199  TLorentzVector p_out(px, py, pz, e);
200  m_beamPart[line] = p_out;
201  m_xAtTrPoint[line] = -b_x;
202  m_yAtTrPoint[line] = b_y;
203  return true;
204 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
uint32_t arm() const
Definition: CTPPSDetId.h:51
double getHalfXangleX56() const
edm::ESGetToken< LHCInfo, LHCInfoRcd > lhcInfoToken_
OpticalFunctionsTransport(const edm::ParameterSet &ps, edm::ConsumesCollector iC)
edm::ESGetToken< LHCInterpolatedOpticalFunctionsSetCollection, CTPPSInterpolatedOpticsRcd > opticsToken_
static const double urad
double getVtxOffsetX45() const
void process(const HepMC::GenEvent *ev, const edm::EventSetup &es, CLHEP::HepRandomEngine *engine) override
T sqrt(T t)
Definition: SSEVec.h:19
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::ESGetToken< CTPPSBeamParameters, CTPPSBeamParametersRcd > beamParametersToken_
double x() const
get X beam position
float const crossingAngle() const
Definition: LHCInfo.cc:182
double getVtxOffsetY45() const
const LHCInterpolatedOpticalFunctionsSetCollection * opticalFunctions_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const BeamSpotObjects * beamspot_
CLHEP::HepRandomEngine * engine_
double getBeamMom45() const
double y() const
get Y beam position
uint32_t rp() const
Definition: CTPPSDetId.h:65
uint32_t station() const
Definition: CTPPSDetId.h:58
std::map< unsigned int, double > m_xAtTrPoint
edm::ESGetToken< BeamSpotObjects, BeamSpotObjectsRcd > beamspotToken_
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
std::map< unsigned int, TLorentzVector > m_beamPart
HLT enums.
double getBeamMom56() const
const CTPPSBeamParameters * beamParameters_
bool transportProton(const HepMC::GenParticle *)
static const double ProtonMassSQ
std::map< unsigned int, double > m_yAtTrPoint
double getHalfXangleX45() const
static const double cm_to_mm
#define LogDebug(id)