CMS 3D CMS Logo

OpticalFunctionsTransport.cc
Go to the documentation of this file.
3 
5  : BaseProtonTransport(iConfig),
6  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
7  opticsLabel_(iConfig.getParameter<std::string>("opticsLabel")),
8  useEmpiricalApertures_(iConfig.getParameter<bool>("useEmpiricalApertures")),
9  empiricalAperture45_xi0_int_(iConfig.getParameter<double>("empiricalAperture45_xi0_int")),
10  empiricalAperture45_xi0_slp_(iConfig.getParameter<double>("empiricalAperture45_xi0_slp")),
11  empiricalAperture45_a_int_(iConfig.getParameter<double>("empiricalAperture45_a_int")),
12  empiricalAperture45_a_slp_(iConfig.getParameter<double>("empiricalAperture45_a_slp")),
13  empiricalAperture56_xi0_int_(iConfig.getParameter<double>("empiricalAperture56_xi0_int")),
14  empiricalAperture56_xi0_slp_(iConfig.getParameter<double>("empiricalAperture56_xi0_slp")),
15  empiricalAperture56_a_int_(iConfig.getParameter<double>("empiricalAperture56_a_int")),
16  empiricalAperture56_a_slp_(iConfig.getParameter<double>("empiricalAperture56_a_slp")),
17  produceHitsRelativeToBeam_(iConfig.getParameter<bool>("produceHitsRelativeToBeam")) {
19 }
21  const edm::EventSetup& iSetup,
22  CLHEP::HepRandomEngine* _engine) {
23  this->clear();
24 
28 
29  // Choose the optical function corresponding to the first station ono each side (it is in lhc ref. frame)
30  optFunctionId45_ = 0;
31  optFunctionId56_ = 0;
32  for (const auto& ofp : (*opticalFunctions_)) {
33  if (ofp.second.getScoringPlaneZ() < 0) {
34  if (optFunctionId45_ == 0)
35  optFunctionId45_ = ofp.first;
36  if (opticalFunctions_->at(optFunctionId45_).getScoringPlaneZ() < ofp.second.getScoringPlaneZ())
37  optFunctionId45_ = ofp.first;
38  }
39  if (ofp.second.getScoringPlaneZ() > 0) {
40  if (optFunctionId56_ == 0)
41  optFunctionId56_ = ofp.first;
42  if (opticalFunctions_->at(optFunctionId56_).getScoringPlaneZ() > ofp.second.getScoringPlaneZ())
43  optFunctionId56_ = ofp.first;
44  }
45  }
46  //
47  engine_ = _engine; // the engine needs to be updated for each event
48 
49  for (HepMC::GenEvent::particle_const_iterator eventParticle = evt->particles_begin();
50  eventParticle != evt->particles_end();
51  ++eventParticle) {
52  if (!((*eventParticle)->status() == 1 && (*eventParticle)->pdg_id() == 2212))
53  continue;
54 
55  if (!(fabs((*eventParticle)->momentum().eta()) > etaCut_ && fabs((*eventParticle)->momentum().pz()) > momentumCut_))
56  continue; // discard protons outside kinematic acceptance
57 
58  unsigned int line = (*eventParticle)->barcode();
59  HepMC::GenParticle* gpart = (*eventParticle);
60  if (gpart->pdg_id() != 2212)
61  continue; // only transport stable protons
62  if (gpart->status() != 1)
63  continue;
64  if (m_beamPart.find(line) != m_beamPart.end()) // assures this protons has not been already propagated
65  continue;
66 
67  transportProton(gpart);
68  }
69 }
70 
72  const HepMC::FourVector& vtx_cms = in_trk->production_vertex()->position(); // in mm
73  const HepMC::FourVector& mom_cms = in_trk->momentum();
74 
75  // transformation to LHC/TOTEM convention
76  HepMC::FourVector vtx_lhc(-vtx_cms.x(), vtx_cms.y(), -vtx_cms.z(), vtx_cms.t());
77  HepMC::FourVector mom_lhc(-mom_cms.x(), mom_cms.y(), -mom_cms.z(), mom_cms.t());
78 
79  // determine the LHC arm and related parameters
80  double urad = 1.e-6;
81  double beamMomentum = 0.;
82  double xangle = 0.;
83  double empiricalAperture_xi0_int, empiricalAperture_xi0_slp;
84  double empiricalAperture_a_int, empiricalAperture_a_slp;
85  unsigned int optFunctionId_;
86 
87  if (mom_lhc.z() < 0) // sector 45
88  {
89  optFunctionId_ = optFunctionId45_;
92  empiricalAperture_xi0_int = empiricalAperture45_xi0_int_;
93  empiricalAperture_xi0_slp = empiricalAperture45_xi0_slp_;
94  empiricalAperture_a_int = empiricalAperture45_a_int_;
95  empiricalAperture_a_slp = empiricalAperture45_a_slp_;
96  } else { // sector 56
97  optFunctionId_ = optFunctionId56_;
100  empiricalAperture_xi0_int = empiricalAperture56_xi0_int_;
101  empiricalAperture_xi0_slp = empiricalAperture56_xi0_slp_;
102  empiricalAperture_a_int = empiricalAperture56_a_int_;
103  empiricalAperture_a_slp = empiricalAperture56_a_slp_;
104  }
105  if (xangle > 1.0)
106  xangle *= urad;
107  // calculate kinematics for optics parametrisation
108  const double p = mom_lhc.rho();
109  const double xi = 1. - p / beamMomentum;
110  const double th_x_phys = mom_lhc.x() / p;
111  const double th_y_phys = mom_lhc.y() / p;
112  const double vtx_lhc_eff_x = vtx_lhc.x() - vtx_lhc.z() * (mom_lhc.x() / mom_lhc.z() + xangle);
113  const double vtx_lhc_eff_y = vtx_lhc.y() - vtx_lhc.z() * (mom_lhc.y() / mom_lhc.z());
114 
115  if (verbosity_) {
116  LogDebug("OpticalFunctionsTransport")
117  << "simu: xi = " << xi << ", th_x_phys = " << th_x_phys << ", th_y_phys = " << th_y_phys
118  << ", vtx_lhc_eff_x = " << vtx_lhc_eff_x << ", vtx_lhc_eff_y = " << vtx_lhc_eff_y << std::endl;
119  }
120 
121  // check empirical aperture
123  const auto& xangle =
125  const double xi_th = (empiricalAperture_xi0_int + xangle * empiricalAperture_xi0_slp) +
126  (empiricalAperture_a_int + xangle * empiricalAperture_a_slp) * th_x_phys;
127 
128  if (xi > xi_th) {
129  if (verbosity_) {
130  LogDebug("OpticalFunctionsTransport") << "stop because of empirical appertures";
131  }
132  return false;
133  }
134  }
135 
136  // transport the proton into pot/scoring plane
137  auto ofp = opticalFunctions_->at(optFunctionId_);
138  CTPPSDetId rpId(optFunctionId_);
139  const unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
140 
141  if (verbosity_)
142  LogDebug("OpticalFunctionsTransport") << " RP " << rpDecId << std::endl;
143 
144  // transport proton
146  vtx_lhc_eff_x * 1E-1, th_x_phys, vtx_lhc_eff_y * 1E-1, th_y_phys, xi}; // conversions: mm -> cm
147 
149  ofp.transport(k_in, k_out, true);
150 
151  // Original code uses mm, but CMS uses cm, so keep it in cm
152  double b_x = k_out.x * 1E1, b_y = k_out.y * 1E1; // conversions: cm -> mm
153  double a_x = k_out.th_x, a_y = k_out.th_y;
154 
155  // if needed, subtract beam position and angle
157  // determine beam position
158  LHCInterpolatedOpticalFunctionsSet::Kinematics k_be_in = {0., 0., 0., 0., 0.};
160  ofp.transport(k_be_in, k_be_out, true);
161 
162  a_x -= k_be_out.th_x;
163  a_y -= k_be_out.th_y;
164  b_x -= k_be_out.x * 1E1;
165  b_y -= k_be_out.y * 1E1; // conversions: cm -> mm
166  }
167 
168  const double z_scoringPlane = ofp.getScoringPlaneZ() * 1E1; // conversion: cm --> mm
169 
170  if (verbosity_) {
171  LogDebug("OpticalFunctionsTransport")
172  << " proton transported: a_x = " << a_x << " rad, a_y = " << a_y << " rad, b_x = " << b_x
173  << " mm, b_y = " << b_y << " mm, z = " << z_scoringPlane << " mm" << std::endl;
174  }
175 
176  unsigned int line = in_trk->barcode();
177  double px = -p * a_x;
178  double py = p * a_y;
179  double pz = std::copysign(sqrt(p * p - px * px - py * py), mom_cms.z());
180  double e = sqrt(px * px + py * py + pz * pz + ProtonMassSQ);
181  TLorentzVector p_out(px, py, pz, e);
182  m_beamPart[line] = p_out;
183  m_xAtTrPoint[line] = -b_x;
184  m_yAtTrPoint[line] = b_y;
185  return true;
186 }
electrons_cff.bool
bool
Definition: electrons_cff.py:393
BaseProtonTransport::TransportMode::OPTICALFUNCTIONS
OpticalFunctionsTransport.h
OpticalFunctionsTransport::opticsLabel_
std::string opticsLabel_
Definition: OpticalFunctionsTransport.h:33
LHCInterpolatedOpticalFunctionsSet::Kinematics
proton kinematics description
Definition: LHCInterpolatedOpticalFunctionsSet.h:28
OpticalFunctionsTransport::empiricalAperture45_a_int_
double empiricalAperture45_a_int_
Definition: OpticalFunctionsTransport.h:42
OpticalFunctionsTransport::empiricalAperture45_xi0_int_
double empiricalAperture45_xi0_int_
Definition: OpticalFunctionsTransport.h:42
OpticalFunctionsTransport::lhcInfoLabel_
std::string lhcInfoLabel_
Definition: OpticalFunctionsTransport.h:32
BaseProtonTransport::clear
void clear()
Definition: BaseProtonTransport.cc:109
BaseProtonTransport::MODE
TransportMode MODE
Definition: BaseProtonTransport.h:42
OpticalFunctionsTransport::empiricalAperture45_a_slp_
double empiricalAperture45_a_slp_
Definition: OpticalFunctionsTransport.h:42
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
OpticalFunctionsTransport::transportProton
bool transportProton(const HepMC::GenParticle *)
Definition: OpticalFunctionsTransport.cc:71
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
CTPPSBeamParameters::getHalfXangleX56
double getHalfXangleX56() const
Definition: CTPPSBeamParameters.cc:70
BaseProtonTransport::verbosity_
bool verbosity_
Definition: BaseProtonTransport.h:51
hybridSuperClusters_cfi.xi
xi
Definition: hybridSuperClusters_cfi.py:10
OpticalFunctionsTransport::process
void process(const HepMC::GenEvent *ev, const edm::EventSetup &es, CLHEP::HepRandomEngine *engine) override
Definition: OpticalFunctionsTransport.cc:20
CTPPSBeamParameters::getBeamMom56
double getBeamMom56() const
Definition: CTPPSBeamParameters.cc:56
CTPPSBeamParameters::getHalfXangleX45
double getHalfXangleX45() const
Definition: CTPPSBeamParameters.cc:68
OpticalFunctionsTransport::optFunctionId45_
unsigned int optFunctionId45_
Definition: OpticalFunctionsTransport.h:38
OpticalFunctionsTransport::optFunctionId56_
unsigned int optFunctionId56_
Definition: OpticalFunctionsTransport.h:39
year_2016_postTS2_cff.rpId
rpId
Definition: year_2016_postTS2_cff.py:23
BaseProtonTransport::momentumCut_
double momentumCut_
Definition: BaseProtonTransport.h:67
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
OpticalFunctionsTransport::lhcInfo_
edm::ESHandle< LHCInfo > lhcInfo_
Definition: OpticalFunctionsTransport.h:35
OpticalFunctionsTransport::empiricalAperture56_a_slp_
double empiricalAperture56_a_slp_
Definition: OpticalFunctionsTransport.h:44
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
LHCInfo::crossingAngle
const float crossingAngle() const
Definition: LHCInfo.cc:182
LHCInterpolatedOpticalFunctionsSet::Kinematics::y
double y
Definition: LHCInterpolatedOpticalFunctionsSet.h:31
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
OpticalFunctionsTransport::empiricalAperture56_a_int_
double empiricalAperture56_a_int_
Definition: OpticalFunctionsTransport.h:44
OpticalFunctionsTransport::empiricalAperture56_xi0_int_
double empiricalAperture56_xi0_int_
Definition: OpticalFunctionsTransport.h:44
LHCInterpolatedOpticalFunctionsSet::Kinematics::x
double x
Definition: LHCInterpolatedOpticalFunctionsSet.h:29
OpticalFunctionsTransport::produceHitsRelativeToBeam_
bool produceHitsRelativeToBeam_
Definition: OpticalFunctionsTransport.h:47
LHCInfoRcd
Definition: LHCInfoRcd.h:24
OpticalFunctionsTransport::beamParameters_
edm::ESHandle< CTPPSBeamParameters > beamParameters_
Definition: OpticalFunctionsTransport.h:36
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
OpticalFunctionsConfig_cfi.xangle
xangle
Definition: OpticalFunctionsConfig_cfi.py:17
LHCInterpolatedOpticalFunctionsSet::Kinematics::th_y
double th_y
Definition: LHCInterpolatedOpticalFunctionsSet.h:32
ProtonMassSQ
static const double ProtonMassSQ
Definition: PPSUnitConversion.h:6
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
edm::ParameterSet
Definition: ParameterSet.h:47
BaseProtonTransport::etaCut_
double etaCut_
Definition: BaseProtonTransport.h:66
BaseProtonTransport::m_beamPart
std::map< unsigned int, TLorentzVector > m_beamPart
Definition: BaseProtonTransport.h:47
CTPPSDetId
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:31
urad
static const double urad
Definition: PPSUnitConversion.h:10
OpticalFunctionsTransport::opticalFunctions_
edm::ESHandle< LHCInterpolatedOpticalFunctionsSetCollection > opticalFunctions_
Definition: OpticalFunctionsTransport.h:37
BaseProtonTransport
Definition: BaseProtonTransport.h:17
edm::EventSetup
Definition: EventSetup.h:57
CTPPSBeamParameters::getBeamMom45
double getBeamMom45() const
Definition: CTPPSBeamParameters.cc:55
BaseProtonTransport::engine_
CLHEP::HepRandomEngine * engine_
Definition: BaseProtonTransport.h:45
get
#define get
CTPPSBeamParametersRcd
Definition: CTPPSBeamParametersRcd.h:14
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
GenParticle.GenParticle
GenParticle
Definition: GenParticle.py:18
std
Definition: JetResolutionObject.h:76
OpticalFunctionsTransport::empiricalAperture56_xi0_slp_
double empiricalAperture56_xi0_slp_
Definition: OpticalFunctionsTransport.h:44
BaseProtonTransport::m_yAtTrPoint
std::map< unsigned int, double > m_yAtTrPoint
Definition: BaseProtonTransport.h:49
OpticalFunctionsTransport::empiricalAperture45_xi0_slp_
double empiricalAperture45_xi0_slp_
Definition: OpticalFunctionsTransport.h:42
OpticalFunctionsTransport::OpticalFunctionsTransport
OpticalFunctionsTransport(const edm::ParameterSet &ps)
Definition: OpticalFunctionsTransport.cc:4
CTPPSDetId.h
BaseProtonTransport::m_xAtTrPoint
std::map< unsigned int, double > m_xAtTrPoint
Definition: BaseProtonTransport.h:48
BaseProtonTransport::beamMomentum
double beamMomentum()
Definition: BaseProtonTransport.h:38
LHCInterpolatedOpticalFunctionsSet::Kinematics::th_x
double th_x
Definition: LHCInterpolatedOpticalFunctionsSet.h:30
OpticalFunctionsTransport::useEmpiricalApertures_
bool useEmpiricalApertures_
Definition: OpticalFunctionsTransport.h:41
mps_splice.line
line
Definition: mps_splice.py:76
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
CTPPSInterpolatedOpticsRcd
Definition: CTPPSInterpolatedOpticsRcd.h:13