CMS 3D CMS Logo

HectorTransport.cc
Go to the documentation of this file.
3 #include <CLHEP/Random/RandGauss.h>
4 #include "TLorentzVector.h"
5 //Hector headers
6 #include "H_BeamLine.h"
7 #include "H_BeamParticle.h"
8 #include <string>
9 
11 
13  : BaseProtonTransport(iConfig), m_beamline45(nullptr), m_beamline56(nullptr) {
14  // Create LHC beam line
15  MODE = TransportMode::HECTOR; // defines the MODE for the transport
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  //PPS
33  edm::LogInfo("ProtonTransport") << "=============================================================================\n"
34  << " Bulding LHC Proton transporter based on HECTOR model\n"
35  << "=============================================================================\n";
36  setBeamLine();
37 }
38 //
39 // this method is the same for all propagator, but since transportProton is different for each derived class
40 // it needes to be overriden
41 //
43  const edm::EventSetup& iSetup,
44  CLHEP::HepRandomEngine* _engine) {
45  this->clear();
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 
61  if (gpart->pdg_id() != 2212)
62  continue; // only transport stable protons
63  if (gpart->status() != 1)
64  continue;
65  if (m_beamPart.find(line) != m_beamPart.end())
66  continue; // assures this protons has not been already propagated
67 
68  transportProton(gpart);
69  }
70 }
72  edm::LogInfo("ProtonTransport") << "Starting proton transport using HECTOR method\n";
73 
74  double px, py, pz, e;
75  unsigned int line = (gpart)->barcode();
76 
77  double mass = gpart->generatedMass();
78  double charge = 1.;
79 
80  px = gpart->momentum().px();
81  py = gpart->momentum().py();
82  pz = gpart->momentum().pz();
83  e = gpart->momentum().e();
84 
85  e = sqrt(pow(mass, 2) + pow(px, 2) + pow(py, 2) + pow(pz, 2));
86 
87  int direction = (pz > 0) ? 1 : -1;
88 
89  // Apply Beam and Crossing Angle Corrections
90  TLorentzVector p_out(px, py, pz, e);
92  "LAB",
93  {fCrossingAngle_56, //Beam1
94  fCrossingAngle_45, //Beam2
96  beamEnergy_});
97 
98  ApplyBeamCorrection(p_out);
99 
100  // from mm to cm
101  double XforPosition = gpart->production_vertex()->position().x() / cm; //cm
102  double YforPosition = gpart->production_vertex()->position().y() / cm; //cm
103  double ZforPosition = gpart->production_vertex()->position().z() / cm; //cm
104 
105  H_BeamParticle h_p(mass, charge);
106  h_p.set4Momentum(-direction * p_out.Px(), p_out.Py(), fabs(p_out.Pz()), p_out.E());
107  // shift the beam position to the given beam position at IP (in cm)
108  XforPosition = (XforPosition - fVtxMeanX) + fBeamXatIP * mm_to_cm;
109  YforPosition = (YforPosition - fVtxMeanY) + fBeamYatIP * mm_to_cm;
110  //ZforPosition stays the same, move the closed orbit only in the X,Y plane
111  //
112  // shift the starting position of the track to Z=0 if configured so (all the variables below are still in cm)
113  if (bApplyZShift) {
114  double fCrossingAngle = (p_out.Pz() > 0) ? fCrossingAngle_45 : -fCrossingAngle_56;
115  XforPosition = XforPosition +
116  (tan((long double)fCrossingAngle * urad) - ((long double)p_out.Px()) / ((long double)p_out.Pz())) *
117  ZforPosition;
118  YforPosition = YforPosition - ((long double)p_out.Py()) / ((long double)p_out.Pz()) * ZforPosition;
119  ZforPosition = 0.;
120  }
121 
122  //
123  // set position, but do not invert the coordinate for the angles (TX,TY) because it is done by set4Momentum
124  h_p.setPosition(-direction * XforPosition * cm_to_um,
125  YforPosition * cm_to_um,
126  h_p.getTX(),
127  h_p.getTY(),
128  -direction * ZforPosition * cm_to_m);
129  bool is_stop;
130  float x1_ctpps;
131  float y1_ctpps;
132 
133  H_BeamLine* _beamline = nullptr;
134  double _targetZ = 0;
135  switch (direction) {
136  case -1:
137  _beamline = &*m_beamline56; // negative side propagation
138  _targetZ = fPPSRegionStart_56;
139  break;
140  case 1:
141  _beamline = &*m_beamline45;
142  _targetZ = fPPSRegionStart_45;
143  break;
144  }
145  // insert protection for NULL beamlines here
146  h_p.computePath(&*_beamline);
147  is_stop = h_p.stopped(&*_beamline);
148  if (verbosity_)
149  LogDebug("HectorTransportEventProcessing")
150  << "HectorTransport:filterPPS: barcode = " << line << " is_stop= " << is_stop;
151  if (is_stop) {
152  return false;
153  }
154  //
155  //propagating
156  //
157  h_p.propagate(_targetZ);
158 
159  p_out = PPSTools::HectorParticle2LorentzVector(h_p, direction);
160 
161  p_out.SetPx(direction * p_out.Px());
162  x1_ctpps = direction * h_p.getX() * um_to_mm;
163  y1_ctpps = h_p.getY() * um_to_mm;
164 
165  if (verbosity_)
166  LogDebug("HectorTransportEventProcessing")
167  << "HectorTransport:filterPPS: barcode = " << line << " x= " << x1_ctpps << " y= " << y1_ctpps;
168 
169  m_beamPart[line] = p_out;
170  m_xAtTrPoint[line] = x1_ctpps;
171  m_yAtTrPoint[line] = y1_ctpps;
172  return true;
173 }
177 
178  // construct beam line for PPS (forward 1 backward 2):
179  if (fPPSBeamLineLength_ > 0.) {
180  m_beamline45 = std::unique_ptr<H_BeamLine>(new H_BeamLine(-1, fPPSBeamLineLength_ + 0.1)); // (direction, length)
181  m_beamline45->fill(b2.fullPath(), 1, "IP5");
182  m_beamline56 = std::unique_ptr<H_BeamLine>(new H_BeamLine(1, fPPSBeamLineLength_ + 0.1)); //
183  m_beamline56->fill(b1.fullPath(), 1, "IP5");
184  } else {
185  if (verbosity_)
186  LogDebug("HectorTransportSetup") << "HectorTransport: WARNING: lengthctpps= " << fPPSBeamLineLength_;
187  return false;
188  }
189  if (verbosity_) {
190  edm::LogInfo("HectorTransportSetup") << "===================================================================\n"
191  << " * * * * * * * * * * * * * * * * * * * * * * * * * * * * \n"
192  << " * * \n"
193  << " * --<--<-- A fast simulator --<--<-- * \n"
194  << " * | --<--<-- of particle --<--<-- * \n"
195  << " * ----HECTOR----< * \n"
196  << " * | -->-->-- transport through-->-->-- * \n"
197  << " * -->-->-- generic beamlines -->-->-- * \n"
198  << " * * \n"
199  << " * JINST 2:P09005 (2007) * \n"
200  << " * X Rouby, J de Favereau, K Piotrzkowski (CP3) * \n"
201  << " * http://www.fynu.ucl.ac.be/hector.html * \n"
202  << " * * \n"
203  << " * Center for Cosmology, Particle Physics and Phenomenology * \n"
204  << " * Universite catholique de Louvain * \n"
205  << " * Louvain-la-Neuve, Belgium * \n"
206  << " * * \n"
207  << " * * * * * * * * * * * * * * * * * * * * * * * * * * * * \n"
208  << " HectorTransport configuration: \n"
209  << " Beam line length = " << fPPSBeamLineLength_ << "\n"
210  << " PPS Region Start 44 = " << fPPSRegionStart_45 << "\n"
211  << " PPS Region Start 56 = " << fPPSRegionStart_56 << "\n"
212  << "===================================================================\n";
213  }
214  if (verbosity_) {
215  edm::LogInfo("HectorTransportSetup") << "====================================================================\n"
216  << " Forward beam line elements \n";
217  m_beamline45->showElements();
218  edm::LogInfo("HectorTransportSetup") << "====================================================================\n"
219  << " Backward beam line elements \n";
220  m_beamline56->showElements();
221  }
222  return true;
223 }
BaseProtonTransport::beam2Filename_
std::string beam2Filename_
Definition: BaseProtonTransport.h:55
BaseProtonTransport::fVtxMeanY
double fVtxMeanY
Definition: BaseProtonTransport.h:67
BaseProtonTransport::beamMomentum_
double beamMomentum_
Definition: BaseProtonTransport.h:62
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
HectorTransport::transportProton
bool transportProton(const HepMC::GenParticle *)
propagate the particles through a beamline to PPS
Definition: HectorTransport.cc:71
BaseProtonTransport::verbosity_
bool verbosity_
Definition: BaseProtonTransport.h:51
edm::LogInfo
Definition: MessageLogger.h:254
BaseProtonTransport::fPPSRegionStart_45
double fPPSRegionStart_45
Definition: BaseProtonTransport.h:57
PPSTools::HectorParticle2LorentzVector
TLorentzVector HectorParticle2LorentzVector(H_BeamParticle hp, int)
Definition: PPSUtilities.cc:6
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
BaseProtonTransport::TransportMode::HECTOR
b2
static constexpr float b2
Definition: L1EGammaCrystalsEmulatorProducer.cc:83
um_to_mm
static const double um_to_mm
Definition: PPSUnitConversion.h:11
BaseProtonTransport::momentumCut_
double momentumCut_
Definition: BaseProtonTransport.h:65
BaseProtonTransport::m_sig_E
double m_sig_E
Definition: BaseProtonTransport.h:74
PPSTools::LorentzBoost
void LorentzBoost(H_BeamParticle &h_p, int dir, const std::string &frame, FullBeamInfo const &bi)
BaseProtonTransport::m_sigmaSY
double m_sigmaSY
Definition: BaseProtonTransport.h:73
PPSUtilities.h
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
b1
static constexpr float b1
Definition: L1EGammaCrystalsEmulatorProducer.cc:83
edm::FileInPath
Definition: FileInPath.h:64
BaseProtonTransport::fVtxMeanX
double fVtxMeanX
Definition: BaseProtonTransport.h:66
BaseProtonTransport::m_sigmaSX
double m_sigmaSX
Definition: BaseProtonTransport.h:72
HectorTransport::fPPSBeamLineLength_
static constexpr double fPPSBeamLineLength_
Definition: HectorTransport.h:44
cm_to_m
static const double cm_to_m
Definition: PPSUnitConversion.h:18
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
HectorTransport::process
void process(const HepMC::GenEvent *ev, const edm::EventSetup &es, CLHEP::HepRandomEngine *engine) override
Definition: HectorTransport.cc:42
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
HectorTransport::m_beamline56
std::unique_ptr< H_BeamLine > m_beamline56
Definition: HectorTransport.h:59
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
edm::ParameterSet
Definition: ParameterSet.h:36
BaseProtonTransport::etaCut_
double etaCut_
Definition: BaseProtonTransport.h:64
BaseProtonTransport::fBeamXatIP
double fBeamXatIP
Definition: BaseProtonTransport.h:75
HectorTransport::HectorTransport
HectorTransport(const edm::ParameterSet &ps)
Definition: HectorTransport.cc:12
BaseProtonTransport::bApplyZShift
bool bApplyZShift
Definition: BaseProtonTransport.h:52
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
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
cm_to_um
static const double cm_to_um
Definition: PPSUnitConversion.h:16
HectorTransport.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
HectorTransport::m_beamline45
std::unique_ptr< H_BeamLine > m_beamline45
Definition: HectorTransport.h:58
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
GenParticle.GenParticle
GenParticle
Definition: GenParticle.py:18
BaseProtonTransport::m_yAtTrPoint
std::map< unsigned int, double > m_yAtTrPoint
Definition: BaseProtonTransport.h:49
BaseProtonTransport::m_sigmaSTX
double m_sigmaSTX
Definition: BaseProtonTransport.h:70
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
BaseProtonTransport::beamEnergy_
double beamEnergy_
Definition: BaseProtonTransport.h:63
BaseProtonTransport::m_xAtTrPoint
std::map< unsigned int, double > m_xAtTrPoint
Definition: BaseProtonTransport.h:48
BaseProtonTransport::fCrossingAngle_45
double fCrossingAngle_45
Definition: BaseProtonTransport.h:59
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
BaseProtonTransport::fBeamYatIP
double fBeamYatIP
Definition: BaseProtonTransport.h:76
mps_splice.line
line
Definition: mps_splice.py:76
HectorTransport::~HectorTransport
~HectorTransport() override
Definition: HectorTransport.cc:10
BaseProtonTransport::fCrossingAngle_56
double fCrossingAngle_56
Definition: BaseProtonTransport.h:60
HectorTransport::setBeamLine
bool setBeamLine()
Definition: HectorTransport.cc:174
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
mm_to_cm
static const double mm_to_cm
Definition: PPSUnitConversion.h:14