ATLAS Offline Software
Loading...
Searching...
No Matches
CrabKissingVertexPositioner.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
5// class header include
7
8// For the speed of light
9#include "GaudiKernel/PhysicalConstants.h"
10
11// CLHEP includes
12#include "CLHEP/Vector/LorentzVector.h"
13#include "CLHEP/Geometry/Point3D.h"
14#include "CLHEP/Geometry/Transform3D.h"
15#include <math.h> /* erf */
16// RandomNumber generator
18#include "CLHEP/Random/RandGaussZiggurat.h"
19#include "CLHEP/Random/RandFlat.h"
20
21namespace Simulation
22{
23
26 const std::string& n,
27 const IInterface* p )
28 : base_class(t,n,p)
29 {
31 }
32
33 void CrabKissingVertexPositioner::BunchShapeHandler(Gaudi::Details::PropertyBase&)
34 {
35 if(m_bunchShapeProp.value() == "GAUSS") m_bunchShape = BunchShape::GAUSS;
36 else if(m_bunchShapeProp.value() == "FLAT") m_bunchShape = BunchShape::FLAT;
38 }
39
42 {
43 ATH_MSG_VERBOSE("Initializing ...");
44
45 // prepare the RandonNumber generation
46 ATH_CHECK(m_rndGenSvc.retrieve());
47 ATH_CHECK(m_beamSpotKey.initialize());
48
49 m_randomEngine = m_rndGenSvc->getEngine(this, m_randomEngineName);
50 if (!m_randomEngine) {
51 ATH_MSG_ERROR("Could not get random number engine from RandomNumberService. Abort.");
52 return StatusCode::FAILURE;
53 }
54
55 ATH_MSG_DEBUG("BunchShape = " << m_bunchShapeProp.value());
56 ATH_MSG_DEBUG("BunchLength = " << m_bunchLength);
57 ATH_MSG_DEBUG("Epsilon (normalized emittance) = " << m_epsilon);
58 ATH_MSG_DEBUG("BetaStar = " << m_betaStar);
59 ATH_MSG_DEBUG("AlfaParallel (kissing angle) = " << m_alphaPar);
60 ATH_MSG_DEBUG("AlfaX (crabbing angle) = " << m_alphaX);
61 ATH_MSG_DEBUG("ThetaX (half crossing angle) = " << m_thetaX);
62
63 // everything set up properly
64 return StatusCode::SUCCESS;
65 }
66
67
70 {
71 ATH_MSG_VERBOSE("Finalizing ...");
72 return StatusCode::SUCCESS;
73 }
74
75 // beamspotFunction for rectangular bunches from Table II in S.Fartoukh Phys.Rev.ST Accel.Beams 17 (2014) no.11, 111001
76 // for z smearing: angle1=psi, angle2=phi
77 // for ct smearing: angle1=phi, angle2=psi
78 double CrabKissingVertexPositioner::beamspotFunction(double displacement, double angle1, double angle2) const
79 {
80 if ( angle1<1e-10 ) angle1 = 1e-10; // to avoid divide_by_zero errors
81 double temp(1.0-std::fabs(displacement)/m_bunchLength);
82 return 1.0/M_2_SQRTPI * std::erf(angle1*temp)/angle1 *
83 std::exp( -pow(angle2*displacement/m_bunchLength, 2) ) *
84 heaviside(temp);
85 }
86
87 double CrabKissingVertexPositioner::getDisplacement(double bunchSize, double angle1, double angle2,
88 CLHEP::HepRandomEngine* randomEngine) const
89 {
90 size_t ntries(0);
91 double yval(CLHEP::RandFlat::shoot(randomEngine, 0.0, 1.0));
92 double displ(CLHEP::RandFlat::shoot(randomEngine, -bunchSize, bunchSize));
93 while (this->beamspotFunction(displ, angle1, angle2)<yval) {
94 if(ntries>1000000) return 0.0; //just so we don't sit in this loop forever
95 yval = CLHEP::RandFlat::shoot(randomEngine, 0.0, 1.0);
96 displ = CLHEP::RandFlat::shoot(randomEngine, -bunchSize, bunchSize);
97 ++ntries;
98 }
99 return displ;
100 }
101
102 // computes the vertex displacement
103 CLHEP::HepLorentzVector *CrabKissingVertexPositioner::generate(const EventContext& ctx) const
104 {
105 // Prepare the random engine
106 m_randomEngine->setSeed( name(), ctx );
107 CLHEP::HepRandomEngine* randomEngine(m_randomEngine->getEngine(ctx));
109 // See jira issue ATLASSIM-497 for an explanation of why calling
110 // shoot outside the CLHEP::HepLorentzVector constructor is
111 // necessary/preferable.
112 double vertexX = CLHEP::RandGaussZiggurat::shoot(randomEngine)*beamSpotHandle->beamSigma(0);
113 double vertexY = CLHEP::RandGaussZiggurat::shoot(randomEngine)*beamSpotHandle->beamSigma(1);
114 double piwinski_phi = std::fabs(m_thetaX - m_alphaX) * m_bunchLength/std::sqrt(m_epsilon * m_betaStar);
115 double piwinski_psi = m_alphaPar * m_bunchLength / std::sqrt( m_epsilon * m_betaStar);
116 double vertexZ = 0;
117 double vertexT = 0;
118 // Time should be set in units of distance, the following methods generate c*t
120 double zWidth = m_bunchLength / std::sqrt(2*(1+pow(piwinski_phi,2)));
121 double timeWidth = m_bunchLength / std::sqrt(2*(1+pow(piwinski_psi,2)));
122 vertexZ = CLHEP::RandGaussZiggurat::shoot(randomEngine, 0., zWidth);
123 vertexT = CLHEP::RandGaussZiggurat::shoot(randomEngine, 0., timeWidth);
124 }
125 else if ( m_bunchShape == BunchShape::FLAT ) {
126 vertexZ = getDisplacement(m_bunchLength, piwinski_psi, piwinski_phi, randomEngine);
127 vertexT = getDisplacement(m_bunchLength, piwinski_phi, piwinski_psi, randomEngine);
128 }
129 else {
130 ATH_MSG_WARNING("Invalid BunchShape ("<<m_bunchShapeProp.value()<<"). Vertex smearing disabled");
131 return new CLHEP::HepLorentzVector(0.,0.,0.,0.);
132 }
133 ATH_MSG_VERBOSE("m_bunchLength = " << m_bunchLength <<
134 ", Zvertex = " << vertexZ << ", Tvertex = " << vertexT);
135
136 // calculate the vertexSmearing
137 CLHEP::HepLorentzVector *vertexSmearing =
138 new CLHEP::HepLorentzVector( vertexX, vertexY, vertexZ, 0. );
139
140 // (1) code from: Simulation/G4Atlas/G4AtlasUtilities/VertexPositioner.cxx
141 const double tx = tan( beamSpotHandle->beamTilt(1) );
142 const double ty = tan( beamSpotHandle->beamTilt(0) );
143
144 const double sqrt_abc = std::sqrt(1. + tx*tx + ty*ty);
145 const double sqrt_fgh = std::sqrt(1. + ty*ty);
146
147 const double a = ty/sqrt_abc;
148 const double b = tx/sqrt_abc;
149 const double c = 1./sqrt_abc;
150
151 HepGeom::Point3D<double> from1(0,0,1);
152 HepGeom::Point3D<double> to1(a,b,c);
153
154 const double f = 1./sqrt_fgh;
155 const double g = 0.;
156 const double h = -(ty)/sqrt_fgh;
157
158 HepGeom::Point3D<double> from2(1,0,0);
159 HepGeom::Point3D<double> to2(f,g,h);
160
161 // first rotation, then translation
162 HepGeom::Transform3D transform(
163 HepGeom::Rotate3D(from1, from2, to1, to2).getRotation(),
164 CLHEP::Hep3Vector( beamSpotHandle->beamPos().x(),
165 beamSpotHandle->beamPos().y(),
166 beamSpotHandle->beamPos().z() )
167 );
168
169 ATH_MSG_VERBOSE("BeamSpotSvc reported beam position as " << beamSpotHandle->beamPos());
170 ATH_MSG_VERBOSE("Width is (" << beamSpotHandle->beamSigma(0) << ", " <<
171 beamSpotHandle->beamSigma(1) << ", " << m_bunchLength << ")");
172 ATH_MSG_VERBOSE("Tilts are " << beamSpotHandle->beamTilt(0) << " and " <<
173 beamSpotHandle->beamTilt(1));
174 ATH_MSG_VERBOSE("Vertex Position before transform: " << *vertexSmearing);
175
176 // update with the tilt
177 *vertexSmearing = transform * HepGeom::Point3D<double>(*vertexSmearing);
178 vertexSmearing->setT(vertexT);
179
180 // and return it
181 return vertexSmearing;
182 }
183
184} // namespace Simulation
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static Double_t a
constexpr int pow(int base, int exp) noexcept
Header file for AthHistogramAlgorithm.
StatusCode finalize() override final
Athena algtool's Hooks.
Gaudi::Property< double > m_betaStar
parameters according to S.Fartoukh Phys.Rev.ST Accel.Beams 17 (2014) no.11, 111001 ------------------...
double beamspotFunction(double displacement, double angle1, double angle2) const
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
CrabKissingVertexPositioner(const std::string &t, const std::string &n, const IInterface *p)
Constructor with parameters.
StatusCode initialize() override final
Athena algtool's Hooks.
void BunchShapeHandler(Gaudi::Details::PropertyBase &)
double getDisplacement(double bunchSize, double angle1, double angle2, CLHEP::HepRandomEngine *rng) const
Gaudi::Property< double > m_bunchLength
Parameter in the Z distribution of the beamspot.
CLHEP::HepLorentzVector * generate(const EventContext &ctx) const override final
computes the vertex displacement