CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BetaFuncPrimaryVertexGenerator.cc
Go to the documentation of this file.
1 //Framework Headers
3 
4 //Famos Headers
7 
9  const edm::ParameterSet& vtx, const RandomEngine* engine) :
10  PrimaryVertexGenerator(engine),
11  fX0(vtx.getParameter<double>("X0")),
12  fY0(vtx.getParameter<double>("Y0")),
13  fZ0(vtx.getParameter<double>("Z0")),
14  fSigmaZ(vtx.getParameter<double>("SigmaZ")),
15  alpha_(vtx.getParameter<double>("Alpha")),
16  phi_(vtx.getParameter<double>("Phi")),
17  fbetastar(vtx.getParameter<double>("BetaStar")),
18  femittance(vtx.getParameter<double>("Emittance"))
19 {
20 
23 
24 }
25 
26 
28 
29  double tmp_sigz = random->gaussShoot(0., fSigmaZ);
30  this->SetZ(tmp_sigz + fZ0);
31 
32  double tmp_sigx = BetaFunction(tmp_sigz,fZ0);
33  // need to divide by sqrt(2) for beamspot width relative to single beam width
34  tmp_sigx *= 0.707107;
35  this->SetX(random->gaussShoot(fX0,tmp_sigx));
36 
37  double tmp_sigy = BetaFunction(tmp_sigz,fZ0);
38  // need to divide by sqrt(2) for beamspot width relative to single beam width
39  tmp_sigy *= 0.707107;
40  this->SetY(random->gaussShoot(fY0,tmp_sigy));
41 
42 }
43 
45 {
46  return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
47 
48 }
49 
50 
51 TMatrixD*
53 
54  TMatrixD* aBoost = 0;
55  if ( fabs(alpha_) < 1E-12 && fabs(phi_) < 1E-12 ) return aBoost;
56 
57  TMatrixD tmpboost(4,4);
58 
59  // Lorentz boost to frame where the collision is head-on
60  // phi is the half crossing angle in the plane ZS
61  // alpha is the angle to the S axis from the X axis in the XY plane
62 
63  double calpha = std::cos(alpha_);
64  double salpha = std::sin(alpha_);
65  double cphi = std::cos(phi_);
66  double sphi = std::sin(phi_);
67  double tphi = sphi/cphi;
68  tmpboost(0,0) = 1./cphi;
69  tmpboost(0,1) = - calpha*sphi;
70  tmpboost(0,2) = - tphi*sphi;
71  tmpboost(0,3) = - salpha*sphi;
72  tmpboost(1,0) = - calpha*tphi;
73  tmpboost(1,1) = 1.;
74  tmpboost(1,2) = calpha*tphi;
75  tmpboost(1,3) = 0.;
76  tmpboost(2,0) = 0.;
77  tmpboost(2,1) = -calpha*sphi;
78  tmpboost(2,2) = cphi;
79  tmpboost(2,3) = - salpha*sphi;
80  tmpboost(3,0) = - salpha*tphi;
81  tmpboost(3,1) = 0.;
82  tmpboost(3,2) = salpha*tphi;
83  tmpboost(3,3) = 1.;
84 
85  tmpboost.Invert();
86  aBoost = new TMatrixD(tmpboost);
87 
88  return aBoost;
89 
90 }
BetaFuncPrimaryVertexGenerator(const edm::ParameterSet &vtx, const RandomEngine *engine)
Default constructor.
const RandomEngine * random
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual void generate()
Generation process (to be implemented)
double gaussShoot(double mean=0.0, double sigma=1.0) const
Definition: RandomEngine.h:37
Definition: DDAxes.h:10
T sqrt(T t)
Definition: SSEVec.h:28
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13