CMS 3D CMS Logo

BaseProtonTransport.cc
Go to the documentation of this file.
4 #include <CLHEP/Random/RandGauss.h>
5 #include <CLHEP/Units/GlobalSystemOfUnits.h>
6 #include <cctype>
7 
9  : verbosity_(iConfig.getParameter<bool>("Verbosity")),
10  bApplyZShift(iConfig.getParameter<bool>("ApplyZShift")),
11  fPPSRegionStart_45(iConfig.getParameter<double>("PPSRegionStart_45")),
12  fPPSRegionStart_56(iConfig.getParameter<double>("PPSRegionStart_56")),
13  etaCut_(iConfig.getParameter<double>("EtaCut")),
14  momentumCut_(iConfig.getParameter<double>("MomentumCut")) {
16 }
18  TLorentzVector p_out;
19  p_out.SetPx(p->momentum().px());
20  p_out.SetPy(p->momentum().py());
21  p_out.SetPz(p->momentum().pz());
22  p_out.SetE(p->momentum().e());
23  ApplyBeamCorrection(p_out);
24  p->set_momentum(HepMC::FourVector(p_out.Px(), p_out.Py(), p_out.Pz(), p_out.E()));
25 }
26 void BaseProtonTransport::ApplyBeamCorrection(TLorentzVector& p_out) {
27  double theta = p_out.Theta();
28  double thetax = atan(p_out.Px() / fabs(p_out.Pz()));
29  double thetay = atan(p_out.Py() / fabs(p_out.Pz()));
30  double energy = p_out.E();
31  double urad = 1e-6;
32 
33  int direction = (p_out.Pz() > 0) ? 1 : -1;
34 
35  if (p_out.Pz() < 0)
36  theta = TMath::Pi() - theta;
37 
39  thetax += (p_out.Pz() > 0) ? fCrossingAngleX_45 * urad : fCrossingAngleX_56 * urad;
40 
41  double dtheta_x = (double)CLHEP::RandGauss::shoot(engine_, 0., m_sigmaSTX);
42  double dtheta_y = (double)CLHEP::RandGauss::shoot(engine_, 0., m_sigmaSTY);
43  double denergy = (double)CLHEP::RandGauss::shoot(engine_, 0., m_sig_E);
44 
45  double s_theta = sqrt(pow(thetax + dtheta_x * urad, 2) + pow(thetay + dtheta_y * urad, 2));
46  double s_phi = atan2(thetay + dtheta_y * urad, thetax + dtheta_x * urad);
47  energy += denergy;
48  double p = sqrt(pow(energy, 2) - ProtonMassSQ);
49 
50  p_out.SetPx((double)p * sin(s_theta) * cos(s_phi));
51  p_out.SetPy((double)p * sin(s_theta) * sin(s_phi));
52  p_out.SetPz((double)p * (cos(s_theta)) * direction);
53  p_out.SetE(energy);
54 }
56  NEvent++;
57  m_CorrespondenceMap.clear();
58 
59  int direction = 0;
60  HepMC::GenParticle* gpart;
61 
62  unsigned int line;
63 
64  for (auto const& it : m_beamPart) {
65  line = (it).first;
66  gpart = in_evt->barcode_to_particle(line);
67 
68  direction = (gpart->momentum().pz() > 0) ? 1 : -1;
69 
70  // Totem uses negative Z for sector 56 while Hector uses always positive distance
71  double ddd = (direction > 0) ? fPPSRegionStart_45 : fabs(fPPSRegionStart_56);
72 
73  double time = (ddd * meter - gpart->production_vertex()->position().z() * mm); // mm
74 
75  //
76  // ATTENTION: at this point, the vertex at PPS is already in mm
77  //
78  if (ddd == 0.)
79  continue;
80 
81  if (verbosity_) {
82  LogDebug("BaseProtonTransportEventProcessing")
83  << "BaseProtonTransport:: x= " << (*(m_xAtTrPoint.find(line))).second << "\n"
84  << "BaseProtonTransport:: y= " << (*(m_yAtTrPoint.find(line))).second << "\n"
85  << "BaseProtonTransport:: z= " << ddd * direction * m_to_mm << "\n"
86  << "BaseProtonTransport:: t= " << time;
87  }
88  TLorentzVector const& p_out = (it).second;
89 
90  HepMC::GenVertex* vert = new HepMC::GenVertex(HepMC::FourVector((*(m_xAtTrPoint.find(line))).second,
91  (*(m_yAtTrPoint.find(line))).second,
92  ddd * direction * m_to_mm,
93  time + time * 0.001));
94 
95  vert->add_particle_out(new HepMC::GenParticle(
96  HepMC::FourVector(p_out.Px(), p_out.Py(), p_out.Pz(), p_out.E()), gpart->pdg_id(), 1, gpart->flow()));
97  evt->add_vertex(vert);
98 
99  int ingoing = 0; //do not attach the incoming proton to this vertex to avoid duplicating data
100  int outgoing = (*vert->particles_out_const_begin())->barcode();
101 
102  LHCTransportLink theLink(ingoing, outgoing);
103  if (verbosity_)
104  LogDebug("BaseProtonTransportEventProcessing")
105  << "BaseProtonTransport:addPartToHepMC: LHCTransportLink " << theLink;
106  m_CorrespondenceMap.push_back(theLink);
107  }
108 }
110  m_beamPart.clear();
111  m_xAtTrPoint.clear();
112  m_yAtTrPoint.clear();
113  m_CorrespondenceMap.clear();
114 }
BaseProtonTransport::beamMomentum_
double beamMomentum_
Definition: BaseProtonTransport.h:64
electrons_cff.bool
bool
Definition: electrons_cff.py:366
BaseProtonTransport::m_sigmaSTY
double m_sigmaSTY
Definition: BaseProtonTransport.h:70
BaseProtonTransport::clear
void clear()
Definition: BaseProtonTransport.cc:109
BaseProtonTransport::MODE
TransportMode MODE
Definition: BaseProtonTransport.h:42
BaseProtonTransport::ApplyBeamCorrection
void ApplyBeamCorrection(HepMC::GenParticle *p)
Definition: BaseProtonTransport.cc:17
BaseProtonTransport::verbosity_
bool verbosity_
Definition: BaseProtonTransport.h:51
BaseProtonTransport::fPPSRegionStart_45
double fPPSRegionStart_45
Definition: BaseProtonTransport.h:57
protons_cff.time
time
Definition: protons_cff.py:35
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
PPSUnitConversion.h
BaseProtonTransport::m_sig_E
double m_sig_E
Definition: BaseProtonTransport.h:73
BaseProtonTransport::addPartToHepMC
void addPartToHepMC(const HepMC::GenEvent *, HepMC::GenEvent *)
Definition: BaseProtonTransport.cc:55
BaseProtonTransport::fPPSRegionStart_56
double fPPSRegionStart_56
Definition: BaseProtonTransport.h:58
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
m_to_mm
static const double m_to_mm
Definition: PPSUnitConversion.h:20
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
BaseProtonTransport::NEvent
int NEvent
Definition: BaseProtonTransport.h:44
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:125
BaseProtonTransport::fCrossingAngleX_45
double fCrossingAngleX_45
Definition: BaseProtonTransport.h:59
ProtonMassSQ
static const double ProtonMassSQ
Definition: PPSUnitConversion.h:6
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
BaseProtonTransport::m_beamPart
std::map< unsigned int, TLorentzVector > m_beamPart
Definition: BaseProtonTransport.h:47
BaseProtonTransport.h
urad
static const double urad
Definition: PPSUnitConversion.h:10
BaseProtonTransport::engine_
CLHEP::HepRandomEngine * engine_
Definition: BaseProtonTransport.h:45
BaseProtonTransport::TransportMode::TOTEM
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:69
BaseProtonTransport::beamEnergy_
double beamEnergy_
Definition: BaseProtonTransport.h:65
Exception.h
BaseProtonTransport::m_xAtTrPoint
std::map< unsigned int, double > m_xAtTrPoint
Definition: BaseProtonTransport.h:48
BaseProtonTransport::fCrossingAngleX_56
double fCrossingAngleX_56
Definition: BaseProtonTransport.h:60
Pi
const double Pi
Definition: CosmicMuonParameters.h:18
BaseProtonTransport::BaseProtonTransport
BaseProtonTransport(const edm::ParameterSet &iConfig)
Definition: BaseProtonTransport.cc:8
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
mps_splice.line
line
Definition: mps_splice.py:76
BaseProtonTransport::m_CorrespondenceMap
std::vector< LHCTransportLink > m_CorrespondenceMap
Definition: BaseProtonTransport.h:46
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37