CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BetafuncEvtVtxGenerator.cc
Go to the documentation of this file.
1 
2 /*
3 ________________________________________________________________________
4 
5  BetafuncEvtVtxGenerator
6 
7  Smear vertex according to the Beta function on the transverse plane
8  and a Gaussian on the z axis. It allows the beam to have a crossing
9  angle (slopes dxdz and dydz).
10 
11  Based on GaussEvtVtxGenerator
12  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
13 
14  FERMILAB
15  2006
16 ________________________________________________________________________
17 */
18 
24 
25 #include "CLHEP/Random/RandGaussQ.h"
26 #include "CLHEP/Units/GlobalSystemOfUnits.h"
27 #include "CLHEP/Units/GlobalPhysicalConstants.h"
28 //#include "CLHEP/Vector/ThreeVector.h"
29 #include "HepMC/SimpleVector.h"
30 
31 #include <iostream>
32 
34  readDB_ = p.getParameter<bool>("readDB");
35  if (!readDB_) {
36  fX0 = p.getParameter<double>("X0") * cm;
37  fY0 = p.getParameter<double>("Y0") * cm;
38  fZ0 = p.getParameter<double>("Z0") * cm;
39  fSigmaZ = p.getParameter<double>("SigmaZ") * cm;
40  fbetastar = p.getParameter<double>("BetaStar") * cm;
41  femittance = p.getParameter<double>("Emittance") * cm; // this is not the normalized emittance
42  fTimeOffset = p.getParameter<double>("TimeOffset") * ns * c_light; // HepMC time units are mm
43 
44  setBoost(p.getParameter<double>("Alpha") * radian, p.getParameter<double>("Phi") * radian);
45  if (fSigmaZ <= 0) {
46  throw cms::Exception("Configuration") << "Error in BetafuncEvtVtxGenerator: "
47  << "Illegal resolution in Z (SigmaZ is negative)";
48  }
49  }
50  if (readDB_) {
51  beamToken_ = esConsumes<SimBeamSpotObjects, SimBeamSpotObjectsRcd, edm::Transition::BeginLuminosityBlock>();
52  }
53 }
54 
56 
58  update(iEventSetup);
59 }
60 
62  if (readDB_ && parameterWatcher_.check(iEventSetup)) {
63  edm::ESHandle<SimBeamSpotObjects> beamhandle = iEventSetup.getHandle(beamToken_);
64 
65  fX0 = beamhandle->fX0;
66  fY0 = beamhandle->fY0;
67  fZ0 = beamhandle->fZ0;
68  // falpha=beamhandle->fAlpha;
69  fSigmaZ = beamhandle->fSigmaZ;
70  fTimeOffset = beamhandle->fTimeOffset;
71  fbetastar = beamhandle->fbetastar;
72  femittance = beamhandle->femittance;
73  setBoost(beamhandle->fAlpha, beamhandle->fPhi);
74  }
75 }
76 
77 HepMC::FourVector BetafuncEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) const {
78  double X, Y, Z;
79 
80  double tmp_sigz = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
81  Z = tmp_sigz + fZ0;
82 
83  double tmp_sigx = BetaFunction(Z, fZ0);
84  // need sqrt(2) for beamspot width relative to single beam width
85  tmp_sigx /= sqrt(2.0);
86  X = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigx) + fX0; // + Z*fdxdz ;
87 
88  double tmp_sigy = BetaFunction(Z, fZ0);
89  // need sqrt(2) for beamspot width relative to single beam width
90  tmp_sigy /= sqrt(2.0);
91  Y = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigy) + fY0; // + Z*fdydz;
92 
93  double tmp_sigt = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
94  double T = tmp_sigt + fTimeOffset;
95 
96  return HepMC::FourVector(X, Y, Z, T);
97 }
98 
99 double BetafuncEvtVtxGenerator::BetaFunction(double z, double z0) const {
100  return sqrt(femittance * (fbetastar + (((z - z0) * (z - z0)) / fbetastar)));
101 }
102 
104  //boost_.ResizeTo(4,4);
105  //boost_ = new TMatrixD(4,4);
106  TMatrixD tmpboost(4, 4);
107 
108  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
109 
110  // Lorentz boost to frame where the collision is head-on
111  // phi is the half crossing angle in the plane ZS
112  // alpha is the angle to the S axis from the X axis in the XY plane
113 
114  tmpboost(0, 0) = 1. / cos(phi);
115  tmpboost(0, 1) = -cos(alpha) * sin(phi);
116  tmpboost(0, 2) = -tan(phi) * sin(phi);
117  tmpboost(0, 3) = -sin(alpha) * sin(phi);
118  tmpboost(1, 0) = -cos(alpha) * tan(phi);
119  tmpboost(1, 1) = 1.;
120  tmpboost(1, 2) = cos(alpha) * tan(phi);
121  tmpboost(1, 3) = 0.;
122  tmpboost(2, 0) = 0.;
123  tmpboost(2, 1) = -cos(alpha) * sin(phi);
124  tmpboost(2, 2) = cos(phi);
125  tmpboost(2, 3) = -sin(alpha) * sin(phi);
126  tmpboost(3, 0) = -sin(alpha) * tan(phi);
127  tmpboost(3, 1) = 0.;
128  tmpboost(3, 2) = sin(alpha) * tan(phi);
129  tmpboost(3, 3) = 1.;
130 
131  tmpboost.Invert();
132  boost_ = tmpboost;
133  //boost_->Print();
134 }
135 
137  if (s >= 0) {
138  fSigmaZ = s;
139  } else {
140  throw cms::Exception("LogicError") << "Error in BetafuncEvtVtxGenerator::sigmaZ: "
141  << "Illegal resolution in Z (negative)";
142  }
143 }
144 
145 TMatrixD const* BetafuncEvtVtxGenerator::GetInvLorentzBoost() const { return &boost_; }
float alpha
Definition: AMPTWrapper.h:105
void sigmaZ(double s=1.0)
set resolution in Z in cm
TMatrixD const * GetInvLorentzBoost() const override
BetafuncEvtVtxGenerator(const edm::ParameterSet &p)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define X(str)
Definition: MuonsGrabber.cc:38
void setBoost(double alpha, double phi)
float float float z
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
edm::ESGetToken< SimBeamSpotObjects, SimBeamSpotObjectsRcd > beamToken_
void update(const edm::EventSetup &iEventSetup)
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const override
return a new event vertex
edm::ESWatcher< SimBeamSpotObjectsRcd > parameterWatcher_
double BetaFunction(double z, double z0) const
beta function
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
long double T