3 #include <CLHEP/Random/RandGauss.h>
4 #include "TLorentzVector.h"
6 #include "H_BeamLine.h"
7 #include "H_BeamParticle.h"
16 m_beamline45(nullptr),
17 m_beamline56(nullptr),
36 edm::LogInfo(
"ProtonTransport") <<
"=============================================================================\n"
37 <<
" Bulding LHC Proton transporter based on HECTOR model\n"
38 <<
"=============================================================================\n";
47 CLHEP::HepRandomEngine* _engine) {
55 for (HepMC::GenEvent::particle_const_iterator eventParticle = evt->particles_begin();
56 eventParticle != evt->particles_end();
58 if (!((*eventParticle)->status() == 1 && (*eventParticle)->pdg_id() == 2212))
61 if (!(fabs((*eventParticle)->momentum().eta()) >
etaCut_ && fabs((*eventParticle)->momentum().pz()) >
momentumCut_))
64 unsigned int line = (*eventParticle)->barcode();
67 if (gpart->pdg_id() != 2212)
69 if (gpart->status() != 1)
78 edm::LogInfo(
"ProtonTransport") <<
"Starting proton transport using HECTOR method\n";
81 unsigned int line = (gpart)->barcode();
83 double mass = gpart->generatedMass();
101 px = -gpart->momentum().px();
102 py = gpart->momentum().py();
103 pz = -gpart->momentum().pz();
104 e = gpart->momentum().e();
106 int direction = (pz > 0) ? 1 : -1;
108 double XforPosition = -gpart->production_vertex()->position().x();
109 double YforPosition = gpart->production_vertex()->position().y();
110 double ZforPosition = -gpart->production_vertex()->position().z();
114 H_BeamParticle h_p(mass, charge);
115 h_p.set4Momentum(px, py, pz, e);
118 XforPosition = XforPosition - ZforPosition * (px / pz + fCrossingAngleX *
urad);
119 YforPosition = YforPosition - ZforPosition * (py / pz + fCrossingAngleY *
urad);
120 XforPosition -= (-vtxXoffset_);
121 YforPosition -= vtxYoffset_;
125 h_p.setPosition(XforPosition *
mm_to_um, YforPosition * mm_to_um, h_p.getTX(), h_p.getTY(), 0.);
131 H_BeamLine* _beamline =
nullptr;
144 h_p.computePath(&*_beamline);
145 is_stop = h_p.stopped(&*_beamline);
147 LogDebug(
"HectorTransportEventProcessing")
148 <<
"HectorTransport:filterPPS: barcode = " << line <<
" is_stop= " << is_stop;
155 h_p.propagate(_targetZ);
156 x1_ctpps = h_p.getX();
157 y1_ctpps = h_p.getY();
159 double thx = h_p.getTX();
160 double thy = h_p.getTY();
163 H_BeamParticle p_beam(mass, charge);
165 p_beam.setPosition(0., 0., fCrossingAngleX *
urad, fCrossingAngleY * urad, 0.);
166 p_beam.computePath(&*_beamline);
167 thx -= p_beam.getTX();
168 thy -= p_beam.getTY();
169 x1_ctpps -= p_beam.getX();
170 y1_ctpps -= p_beam.getY();
172 <<
"Shifting the hit relative to beam : X = " << p_beam.getX() <<
"(mm) Y = " << p_beam.getY() <<
"(mm)";
179 TLorentzVector p_out(-
tan(thx *
urad) * partP *
cos(theta),
180 tan(thy * urad) * partP *
cos(theta),
181 -direction * partP *
cos(theta),
188 LogDebug(
"HectorTransportEventProcessing")
189 <<
"HectorTransport:filterPPS: barcode = " << line <<
" x= " << x1_ctpps <<
" y= " << y1_ctpps;
209 edm::LogInfo(
"HectorTransportSetup") <<
"===================================================================\n"
210 <<
" * * * * * * * * * * * * * * * * * * * * * * * * * * * * \n"
212 <<
" * --<--<-- A fast simulator --<--<-- * \n"
213 <<
" * | --<--<-- of particle --<--<-- * \n"
214 <<
" * ----HECTOR----< * \n"
215 <<
" * | -->-->-- transport through-->-->-- * \n"
216 <<
" * -->-->-- generic beamlines -->-->-- * \n"
218 <<
" * JINST 2:P09005 (2007) * \n"
219 <<
" * X Rouby, J de Favereau, K Piotrzkowski (CP3) * \n"
220 <<
" * http://www.fynu.ucl.ac.be/hector.html * \n"
222 <<
" * Center for Cosmology, Particle Physics and Phenomenology * \n"
223 <<
" * Universite catholique de Louvain * \n"
224 <<
" * Louvain-la-Neuve, Belgium * \n"
226 <<
" * * * * * * * * * * * * * * * * * * * * * * * * * * * * \n"
227 <<
" HectorTransport configuration: \n"
231 <<
"===================================================================\n";
234 edm::LogInfo(
"HectorTransportSetup") <<
"====================================================================\n"
235 <<
" Forward beam line elements \n";
237 edm::LogInfo(
"HectorTransportSetup") <<
"====================================================================\n"
238 <<
" Backward beam line elements \n";
double fPPSRegionStart_45
std::string beam1Filename_
double GetY() const
get Y beam position
~HectorTransport() override
std::string beam2Filename_
bool produceHitsRelativeToBeam_
bool useBeamPositionFromLHCInfo_
Geom::Theta< T > theta() const
double getVtxOffsetX45() const
bool transportProton(const HepMC::GenParticle *)
propagate the particles through a beamline to PPS
double fCrossingAngleY_56
bool getData(T &iHolder) const
double fCrossingAngleX_56
edm::ESGetToken< CTPPSBeamParameters, CTPPSBeamParametersRcd > beamParametersToken_
const CTPPSBeamParameters * beamParameters_
Cos< T >::type cos(const T &t)
Tan< T >::type tan(const T &t)
double GetZ() const
get Z beam position
double fPPSRegionStart_56
const BeamSpotObjects * beamspot_
CLHEP::HepRandomEngine * engine_
static const double um_to_mm
HectorTransport(const edm::ParameterSet &ps, edm::ConsumesCollector iC)
Log< level::Info, false > LogInfo
double GetX() const
get X beam position
void process(const HepMC::GenEvent *ev, const edm::EventSetup &es, CLHEP::HepRandomEngine *engine) override
T getParameter(std::string const &) const
std::map< unsigned int, double > m_xAtTrPoint
double getVtxOffsetZ45() const
static const double mm_to_um
std::map< unsigned int, TLorentzVector > m_beamPart
double getVtxOffsetY45() const
static constexpr double fPPSBeamLineLength_
double fCrossingAngleY_45
std::unique_ptr< H_BeamLine > m_beamline56
std::unique_ptr< H_BeamLine > m_beamline45
static constexpr float b2
edm::ESGetToken< BeamSpotObjects, BeamSpotObjectsRcd > beamspotToken_
static const double ProtonMassSQ
std::map< unsigned int, double > m_yAtTrPoint
Power< A, B >::type pow(const A &a, const B &b)
static constexpr float b1
static const double cm_to_mm
double fCrossingAngleX_45