CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
HectorTransport Class Reference

#include <HectorTransport.h>

Inheritance diagram for HectorTransport:
ProtonTransport

Public Member Functions

 HectorTransport ()
 
 HectorTransport (const edm::ParameterSet &ps, bool verbosity)
 
void process (const HepMC::GenEvent *, const edm::EventSetup &, CLHEP::HepRandomEngine *) override
 
 ~HectorTransport () override
 
- Public Member Functions inherited from ProtonTransport
void addPartToHepMC (HepMC::GenEvent *)
 
void ApplyBeamCorrection (HepMC::GenParticle *p)
 
void ApplyBeamCorrection (TLorentzVector &p)
 
void clear ()
 
std::vector< LHCTransportLink > & getCorrespondenceMap ()
 
 ProtonTransport ()
 
virtual ~ProtonTransport ()
 

Private Member Functions

void genProtonsLoop (const HepMC::GenEvent *, const edm::EventSetup &)
 Clears BeamParticle, prepares PPSHector for a next Aperture check or/and a next event. More...
 
double getBeamEnergy ()
 
double getBeamMomentum ()
 
void setBeamEnergy (double e)
 
bool setBeamLine ()
 
bool transportProton (const HepMC::GenParticle *)
 propagate the particles through a beamline to PPS More...
 

Private Attributes

double m_b_ctpps_b
 
std::string m_beam1filename
 
std::string m_beam2filename
 
std::unique_ptr< H_BeamLine > m_beamline45
 
std::unique_ptr< H_BeamLine > m_beamline56
 
double m_f_ctpps_f
 
double m_fEtacut
 
double m_fMomentumMin
 
double m_lengthctpps
 
edm::ESHandle< ParticleDataTablem_pdt
 

Additional Inherited Members

- Protected Types inherited from ProtonTransport
enum  TransportMode { TransportMode::HECTOR, TransportMode::TOTEM }
 
- Protected Attributes inherited from ProtonTransport
bool bApplyZShift
 
CLHEP::HepRandomEngine * engine
 
double fBeamEnergy
 
double fBeamMomentum
 
double fBeamXatIP
 
double fBeamYatIP
 
double fCrossingAngle_45
 
double fCrossingAngle_56
 
double fPPSRegionStart_45
 
double fPPSRegionStart_56
 
double fVtxMeanX
 
double fVtxMeanY
 
double fVtxMeanZ
 
std::map< unsigned int, TLorentzVector > m_beamPart
 
std::vector< LHCTransportLinkm_CorrespondenceMap
 
double m_sig_E
 
double m_sigmaSTX
 
double m_sigmaSTY
 
double m_sigmaSX
 
double m_sigmaSY
 
bool m_verbosity
 
std::map< unsigned int, double > m_xAtTrPoint
 
std::map< unsigned int, double > m_yAtTrPoint
 
TransportMode MODE
 
int NEvent
 

Detailed Description

Definition at line 37 of file HectorTransport.h.

Constructor & Destructor Documentation

HectorTransport::HectorTransport ( )

Definition at line 10 of file HectorTransport.cc.

10 {};
HectorTransport::HectorTransport ( const edm::ParameterSet ps,
bool  verbosity 
)

Definition at line 12 of file HectorTransport.cc.

References ProtonTransport::bApplyZShift, ProtonTransport::fBeamEnergy, ProtonTransport::fBeamMomentum, ProtonTransport::fBeamXatIP, ProtonTransport::fBeamYatIP, ProtonTransport::fCrossingAngle_45, ProtonTransport::fCrossingAngle_56, objects.autophobj::float, ProtonTransport::fPPSRegionStart_45, ProtonTransport::fPPSRegionStart_56, ProtonTransport::fVtxMeanX, ProtonTransport::fVtxMeanY, ProtonTransport::fVtxMeanZ, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), ProtonTransport::HECTOR, m_b_ctpps_b, m_beam1filename, m_beam2filename, m_f_ctpps_f, m_fEtacut, m_fMomentumMin, m_lengthctpps, ProtonTransport::m_sig_E, ProtonTransport::m_sigmaSTX, ProtonTransport::m_sigmaSTY, ProtonTransport::m_sigmaSX, ProtonTransport::m_sigmaSY, ProtonTransport::m_verbosity, ProtonTransport::MODE, PPSTools::ProtonMassSQ, setBeamLine(), mathSSE::sqrt(), and HIPAlignmentAlgorithm_cfi::verbosity.

13 {
14  // Create LHC beam line
15  MODE=TransportMode::HECTOR; // defines the MODE for the transport
17  edm::ParameterSet hector_par = param.getParameter<edm::ParameterSet>("PPSHector");
18 
19  m_f_ctpps_f = (float) hector_par.getParameter<double>("PPSf");
20  m_b_ctpps_b = (float) hector_par.getParameter<double>("PPSb");
21  fCrossingAngle_56 = hector_par.getParameter<double>("CrossingAngleBeam1");
22  fCrossingAngle_45 = hector_par.getParameter<double>("CrossingAngleBeam2");
23  fBeamEnergy = hector_par.getParameter<double>("BeamEnergy"); // beam energy in GeV
24  m_fEtacut = hector_par.getParameter<double>("EtaCutForHector");
25  m_fMomentumMin = hector_par.getParameter<double>("MomentumMin");
26  fBeamMomentum = sqrt(fBeamEnergy*fBeamEnergy - PPSTools::ProtonMassSQ);
27 
28  // User definitons
29  m_lengthctpps = hector_par.getParameter<double>("BeamLineLengthPPS" );
30 
31  m_beam1filename = hector_par.getParameter<string>("Beam1");
32  m_beam2filename = hector_par.getParameter<string>("Beam2");
33  fVtxMeanX = param.getParameter<double>("VtxMeanX");
34  fVtxMeanY = param.getParameter<double>("VtxMeanY");
35  fVtxMeanZ = param.getParameter<double>("VtxMeanZ");
36  m_sigmaSTX = hector_par.getParameter<double>("sigmaSTX" );
37  m_sigmaSTY = hector_par.getParameter<double>("sigmaSTY" );
38  m_sigmaSX = hector_par.getParameter<double>("sigmaSX");
39  m_sigmaSY = hector_par.getParameter<double>("sigmaSY");
40  m_sig_E = hector_par.getParameter<double>("sigmaEnergy");
41  fBeamXatIP = hector_par.getUntrackedParameter<double>("BeamXatIP",fVtxMeanX);
42  fBeamYatIP = hector_par.getUntrackedParameter<double>("BeamYatIP",fVtxMeanY);
43  bApplyZShift = hector_par.getParameter<bool>("ApplyZShift");
44  //PPS
45  edm::LogInfo("ProtonTransport")
46  <<"=============================================================================\n"
47  <<" Bulding LHC Proton transporter based on HECTOR model\n"
48  <<"=============================================================================\n";
49  setBeamLine();
52 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
double fPPSRegionStart_45
TransportMode MODE
std::string m_beam2filename
double fCrossingAngle_45
double fCrossingAngle_56
T sqrt(T t)
Definition: SSEVec.h:18
double fPPSRegionStart_56
const double ProtonMassSQ
Definition: PPSUtilities.h:33
std::string m_beam1filename
HectorTransport::~HectorTransport ( )
override

Definition at line 53 of file HectorTransport.cc.

References ProtonTransport::clear().

54 {
55  this->clear();
56 
57 }

Member Function Documentation

void HectorTransport::genProtonsLoop ( const HepMC::GenEvent evt,
const edm::EventSetup iSetup 
)
private

Clears BeamParticle, prepares PPSHector for a next Aperture check or/and a next event.

Adds the stable protons from the event ev to a beamline

Definition at line 156 of file HectorTransport.cc.

References GenParticle::GenParticle, mps_splice::line, ProtonTransport::m_beamPart, m_fEtacut, m_fMomentumMin, and transportProton().

Referenced by process().

157 {
158 /*
159  Loop over genVertex looking for transportable protons
160 */
161  unsigned int line;
162 
163  for (HepMC::GenEvent::particle_const_iterator eventParticle =evt->particles_begin(); eventParticle != evt->particles_end(); ++eventParticle ) {
164  if (!((*eventParticle)->status() == 1 && (*eventParticle)->pdg_id()==2212 )) continue;
165  if (!(fabs((*eventParticle)->momentum().eta())>m_fEtacut && fabs((*eventParticle)->momentum().pz())>m_fMomentumMin)) continue;
166  line = (*eventParticle)->barcode();
167 
168  if ( m_beamPart.find(line) != m_beamPart.end() ) continue;
169 
170  HepMC::GenParticle * gpart = (*eventParticle);
171 
172  transportProton(gpart);
173  }
174 }
std::map< unsigned int, TLorentzVector > m_beamPart
bool transportProton(const HepMC::GenParticle *)
propagate the particles through a beamline to PPS
double HectorTransport::getBeamEnergy ( )
inlineprivate

Definition at line 57 of file HectorTransport.h.

57 {return fBeamEnergy;}
double HectorTransport::getBeamMomentum ( )
inlineprivate

Definition at line 59 of file HectorTransport.h.

59 {return fBeamMomentum;}
void HectorTransport::process ( const HepMC::GenEvent ev,
const edm::EventSetup iSetup,
CLHEP::HepRandomEngine *  _engine 
)
overridevirtual

Implements ProtonTransport.

Definition at line 58 of file HectorTransport.cc.

References ProtonTransport::addPartToHepMC(), ProtonTransport::engine, genProtonsLoop(), edm::EventSetup::getData(), and m_pdt.

59 {
60  engine = _engine;
61  iSetup.getData( m_pdt );
62  genProtonsLoop(ev,iSetup);
63  addPartToHepMC(const_cast<HepMC::GenEvent*>(ev));
64 }
void genProtonsLoop(const HepMC::GenEvent *, const edm::EventSetup &)
Clears BeamParticle, prepares PPSHector for a next Aperture check or/and a next event.
bool getData(T &iHolder) const
Definition: EventSetup.h:111
edm::ESHandle< ParticleDataTable > m_pdt
void addPartToHepMC(HepMC::GenEvent *)
CLHEP::HepRandomEngine * engine
void HectorTransport::setBeamEnergy ( double  e)
inlineprivate
bool HectorTransport::setBeamLine ( )
private

Definition at line 175 of file HectorTransport.cc.

References LogDebug, m_b_ctpps_b, m_beam1filename, m_beam2filename, m_beamline45, m_beamline56, m_f_ctpps_f, m_lengthctpps, and ProtonTransport::m_verbosity.

Referenced by HectorTransport().

176 {
177 
178  m_beamline45=nullptr;
179  m_beamline56=nullptr;
180  edm::FileInPath b1(m_beam1filename.c_str());
181  edm::FileInPath b2(m_beam2filename.c_str());
182  if(m_verbosity) {
183  edm::LogInfo("HectorTransportSetup") << "===================================================================\n"
184  << " * * * * * * * * * * * * * * * * * * * * * * * * * * * * \n"
185  << " * * \n"
186  << " * --<--<-- A fast simulator --<--<-- * \n"
187  << " * | --<--<-- of particle --<--<-- * \n"
188  << " * ----HECTOR----< * \n"
189  << " * | -->-->-- transport through-->-->-- * \n"
190  << " * -->-->-- generic beamlines -->-->-- * \n"
191  << " * * \n"
192  << " * JINST 2:P09005 (2007) * \n"
193  << " * X Rouby, J de Favereau, K Piotrzkowski (CP3) * \n"
194  << " * http://www.fynu.ucl.ac.be/hector.html * \n"
195  << " * * \n"
196  << " * Center for Cosmology, Particle Physics and Phenomenology * \n"
197  << " * Universite catholique de Louvain * \n"
198  << " * Louvain-la-Neuve, Belgium * \n"
199  << " * * \n"
200  << " * * * * * * * * * * * * * * * * * * * * * * * * * * * * \n"
201  << " HectorTransport configuration: \n"
202  << " lengthctpps = " << m_lengthctpps << "\n"
203  << " m_f_ctpps_f = " << m_f_ctpps_f << "\n"
204  << " m_b_ctpps_b = " << m_b_ctpps_b << "\n"
205  << "===================================================================\n";
206  }
207 
208  // construct beam line for PPS (forward 1 backward 2):
209  if(m_lengthctpps>0. ) {
210  m_beamline45 = std::unique_ptr<H_BeamLine>(new H_BeamLine(-1, m_lengthctpps + 0.1 )); // (direction, length)
211  m_beamline45->fill( b2.fullPath(), 1, "IP5");
212  m_beamline56 = std::unique_ptr<H_BeamLine>(new H_BeamLine( 1, m_lengthctpps + 0.1 )); //
213  m_beamline56->fill( b1.fullPath(), 1, "IP5");
214  }
215  else {
216  if ( m_verbosity ) LogDebug("HectorTransportSetup") << "HectorTransport: WARNING: lengthctpps= " << m_lengthctpps;
217  return false;
218  }
219  if (m_verbosity) {
220  edm::LogInfo("HectorTransportSetup")
221  << "====================================================================\n"
222  << " Forward beam line elements \n";
223  m_beamline45->showElements();
224  edm::LogInfo("HectorTransportSetup")
225  << "====================================================================\n"
226  << " Backward beam line elements \n";
227  m_beamline56->showElements();
228  }
229 
230  return true;
231 }
#define LogDebug(id)
std::string m_beam2filename
std::string m_beam1filename
std::unique_ptr< H_BeamLine > m_beamline56
std::unique_ptr< H_BeamLine > m_beamline45
bool HectorTransport::transportProton ( const HepMC::GenParticle *  gpart)
private

propagate the particles through a beamline to PPS

Definition at line 65 of file HectorTransport.cc.

References ProtonTransport::ApplyBeamCorrection(), ProtonTransport::bApplyZShift, ALCARECOTkAlJpsiMuMu_cff::charge, cm_to_m, cm_to_um, MillePedeFileConverter_cfg::e, ProtonTransport::fBeamEnergy, ProtonTransport::fBeamMomentum, ProtonTransport::fBeamXatIP, ProtonTransport::fBeamYatIP, ProtonTransport::fCrossingAngle_45, ProtonTransport::fCrossingAngle_56, ProtonTransport::fVtxMeanX, ProtonTransport::fVtxMeanY, PPSTools::HectorParticle2LorentzVector(), mps_splice::line, LogDebug, PPSTools::LorentzBoost(), m_b_ctpps_b, m_beamline45, m_beamline56, ProtonTransport::m_beamPart, m_f_ctpps_f, ProtonTransport::m_verbosity, ProtonTransport::m_xAtTrPoint, ProtonTransport::m_yAtTrPoint, ResonanceBuilder::mass, mm_to_cm, funct::pow(), mathSSE::sqrt(), funct::tan(), um_to_mm, and urad.

Referenced by genProtonsLoop().

66 {
67  edm::LogInfo("ProtonTransport")<<"Starting proton transport using HECTOR method\n";
68 
69  double px,py,pz,e;
70  unsigned int line = (gpart)->barcode();
71 
72  double mass = gpart->generatedMass();
73  double charge = 1.;
74 
75  px = gpart->momentum().px();
76  py = gpart->momentum().py();
77  pz = gpart->momentum().pz();
78  e = gpart->momentum().e();
79 
80  e = sqrt(pow(mass,2)+pow(px,2)+pow(py,2)+pow(pz,2));
81 
82  int direction = (pz>0)?1:-1;
83 
84  // Apply Beam and Crossing Angle Corrections
85  TLorentzVector p_out(px,py,pz,e);
86  PPSTools::LorentzBoost(p_out,"LAB", {fCrossingAngle_56 ,//Beam1
87  fCrossingAngle_45, //Beam2
89 
90  ApplyBeamCorrection(p_out);
91 
92  // from mm to cm
93  double XforPosition = gpart->production_vertex()->position().x()/cm;//cm
94  double YforPosition = gpart->production_vertex()->position().y()/cm;//cm
95  double ZforPosition = gpart->production_vertex()->position().z()/cm;//cm
96 
97  H_BeamParticle h_p(mass,charge);
98  h_p.set4Momentum(-direction*p_out.Px(), p_out.Py(), fabs(p_out.Pz()), p_out.E());
99 // shift the beam position to the given beam position at IP (in cm)
100  XforPosition=(XforPosition-fVtxMeanX)+fBeamXatIP*mm_to_cm;
101  YforPosition=(YforPosition-fVtxMeanY)+fBeamYatIP*mm_to_cm;
102  //ZforPosition stays the same, move the closed orbit only in the X,Y plane
103 //
104 // shift the starting position of the track to Z=0 if configured so (all the variables below are still in cm)
105  if (bApplyZShift) {
106  double fCrossingAngle = (p_out.Pz()>0)?fCrossingAngle_45:-fCrossingAngle_56;
107  XforPosition = XforPosition+(tan((long double)fCrossingAngle*urad)-((long double)p_out.Px())/((long double)p_out.Pz()))*ZforPosition;
108  YforPosition = YforPosition-((long double)p_out.Py())/((long double)p_out.Pz())*ZforPosition;
109  ZforPosition = 0.;
110  }
111 
112 //
113 // set position, but do not invert the coordinate for the angles (TX,TY) because it is done by set4Momentum
114  h_p.setPosition(-direction*XforPosition*cm_to_um,YforPosition*cm_to_um,h_p.getTX(),h_p.getTY(),-direction*ZforPosition*cm_to_m);
115  bool is_stop;
116  float x1_ctpps;
117  float y1_ctpps;
118 
119  H_BeamLine* _beamline = nullptr;
120  double _targetZ=0;
121  switch (direction) {
122  case -1: _beamline = &*m_beamline56;// negative side propagation
123  _targetZ = m_b_ctpps_b;
124  break;
125  case 1: _beamline = &*m_beamline45;
126  _targetZ = m_f_ctpps_f;
127  break;
128  }
129  // insert protection for NULL beamlines here
130  h_p.computePath(&*_beamline);
131  is_stop = h_p.stopped(&*_beamline);
132  if(m_verbosity) LogDebug("HectorTransportEventProcessing") << "HectorTransport:filterPPS: barcode = "
133  << line << " is_stop= "<< is_stop;
134  if (is_stop) {
135  return false;
136  }
137  //
138  //propagating
139  //
140  h_p.propagate( _targetZ );
141 
142  p_out = PPSTools::HectorParticle2LorentzVector(h_p,direction);
143 
144  p_out.SetPx(direction*p_out.Px());
145  x1_ctpps = direction*h_p.getX()*um_to_mm;
146  y1_ctpps = h_p.getY()*um_to_mm;
147 
148  if(m_verbosity) LogDebug("HectorTransportEventProcessing") <<
149  "HectorTransport:filterPPS: barcode = " << line << " x= "<< x1_ctpps <<" y= " << y1_ctpps;
150 
151  m_beamPart[line] = p_out;
152  m_xAtTrPoint[line] = x1_ctpps;
153  m_yAtTrPoint[line] = y1_ctpps;
154  return true;
155 }
#define LogDebug(id)
void LorentzBoost(H_BeamParticle &h_p, int dir, const std::string &frame, FullBeamInfo const &bi)
TLorentzVector HectorParticle2LorentzVector(H_BeamParticle hp, int)
Definition: PPSUtilities.cc:6
void ApplyBeamCorrection(HepMC::GenParticle *p)
std::map< unsigned int, TLorentzVector > m_beamPart
std::map< unsigned int, double > m_yAtTrPoint
static const double urad
static const double cm_to_um
double fCrossingAngle_45
std::map< unsigned int, double > m_xAtTrPoint
double fCrossingAngle_56
T sqrt(T t)
Definition: SSEVec.h:18
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
static const double um_to_mm
static const double mm_to_cm
static const double cm_to_m
std::unique_ptr< H_BeamLine > m_beamline56
std::unique_ptr< H_BeamLine > m_beamline45
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40

Member Data Documentation

double HectorTransport::m_b_ctpps_b
private

Definition at line 75 of file HectorTransport.h.

Referenced by HectorTransport(), setBeamLine(), and transportProton().

std::string HectorTransport::m_beam1filename
private

Definition at line 82 of file HectorTransport.h.

Referenced by HectorTransport(), and setBeamLine().

std::string HectorTransport::m_beam2filename
private

Definition at line 83 of file HectorTransport.h.

Referenced by HectorTransport(), and setBeamLine().

std::unique_ptr<H_BeamLine> HectorTransport::m_beamline45
private

Definition at line 79 of file HectorTransport.h.

Referenced by setBeamLine(), and transportProton().

std::unique_ptr<H_BeamLine> HectorTransport::m_beamline56
private

Definition at line 80 of file HectorTransport.h.

Referenced by setBeamLine(), and transportProton().

double HectorTransport::m_f_ctpps_f
private

Definition at line 74 of file HectorTransport.h.

Referenced by HectorTransport(), setBeamLine(), and transportProton().

double HectorTransport::m_fEtacut
private

Definition at line 70 of file HectorTransport.h.

Referenced by genProtonsLoop(), and HectorTransport().

double HectorTransport::m_fMomentumMin
private

Definition at line 71 of file HectorTransport.h.

Referenced by genProtonsLoop(), and HectorTransport().

double HectorTransport::m_lengthctpps
private

Definition at line 73 of file HectorTransport.h.

Referenced by HectorTransport(), and setBeamLine().

edm::ESHandle< ParticleDataTable > HectorTransport::m_pdt
private

Definition at line 68 of file HectorTransport.h.

Referenced by process().