CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules 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  // NOTE: this is currently watching LS transitions, while it should watch Run transitions,
52  // even though in reality there is no Run Dependent MC (yet) in CMS
53  beamToken_ = esConsumes<SimBeamSpotObjects, SimBeamSpotObjectsRcd, edm::Transition::BeginLuminosityBlock>();
54  }
55 }
56 
58 
60  update(iEventSetup);
61 }
62 
64  if (readDB_ && parameterWatcher_.check(iEventSetup)) {
65  edm::ESHandle<SimBeamSpotObjects> beamhandle = iEventSetup.getHandle(beamToken_);
66 
67  fX0 = beamhandle->x() * cm;
68  fY0 = beamhandle->y() * cm;
69  fZ0 = beamhandle->z() * cm;
70  fSigmaZ = beamhandle->sigmaZ() * cm;
71  fTimeOffset = beamhandle->timeOffset() * ns * c_light; // HepMC time units are mm
72  fbetastar = beamhandle->betaStar() * cm;
73  femittance = beamhandle->emittance() * cm;
74  setBoost(beamhandle->alpha() * radian, beamhandle->phi() * radian);
75  }
76 }
77 
78 HepMC::FourVector BetafuncEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) const {
79  double X, Y, Z;
80 
81  double tmp_sigz = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
82  Z = tmp_sigz + fZ0;
83 
84  double tmp_sigx = BetaFunction(Z, fZ0);
85  // need sqrt(2) for beamspot width relative to single beam width
86  tmp_sigx /= sqrt(2.0);
87  X = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigx) + fX0; // + Z*fdxdz ;
88 
89  double tmp_sigy = BetaFunction(Z, fZ0);
90  // need sqrt(2) for beamspot width relative to single beam width
91  tmp_sigy /= sqrt(2.0);
92  Y = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigy) + fY0; // + Z*fdydz;
93 
94  double tmp_sigt = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
95  double T = tmp_sigt + fTimeOffset;
96 
97  return HepMC::FourVector(X, Y, Z, T);
98 }
99 
100 double BetafuncEvtVtxGenerator::BetaFunction(double z, double z0) const {
101  return sqrt(femittance * (fbetastar + (((z - z0) * (z - z0)) / fbetastar)));
102 }
103 
105  //boost_.ResizeTo(4,4);
106  //boost_ = new TMatrixD(4,4);
107  TMatrixD tmpboost(4, 4);
108 
109  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
110 
111  // Lorentz boost to frame where the collision is head-on
112  // phi is the half crossing angle in the plane ZS
113  // alpha is the angle to the S axis from the X axis in the XY plane
114 
115  tmpboost(0, 0) = 1. / cos(phi);
116  tmpboost(0, 1) = -cos(alpha) * sin(phi);
117  tmpboost(0, 2) = -tan(phi) * sin(phi);
118  tmpboost(0, 3) = -sin(alpha) * sin(phi);
119  tmpboost(1, 0) = -cos(alpha) * tan(phi);
120  tmpboost(1, 1) = 1.;
121  tmpboost(1, 2) = cos(alpha) * tan(phi);
122  tmpboost(1, 3) = 0.;
123  tmpboost(2, 0) = 0.;
124  tmpboost(2, 1) = -cos(alpha) * sin(phi);
125  tmpboost(2, 2) = cos(phi);
126  tmpboost(2, 3) = -sin(alpha) * sin(phi);
127  tmpboost(3, 0) = -sin(alpha) * tan(phi);
128  tmpboost(3, 1) = 0.;
129  tmpboost(3, 2) = sin(alpha) * tan(phi);
130  tmpboost(3, 3) = 1.;
131 
132  tmpboost.Invert();
133  boost_ = tmpboost;
134  //boost_->Print();
135 }
136 
138  if (s >= 0) {
139  fSigmaZ = s;
140  } else {
141  throw cms::Exception("LogicError") << "Error in BetafuncEvtVtxGenerator::sigmaZ: "
142  << "Illegal resolution in Z (negative)";
143  }
144 }
145 
146 TMatrixD const* BetafuncEvtVtxGenerator::GetInvLorentzBoost() const { return &boost_; }
double sigmaZ() const
get sigmaZ
void sigmaZ(double s=1.0)
set resolution in Z in cm
double BetaFunction(double z, double z0) const
beta function
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)
double alpha() const
get Alpha
float float float z
TMatrixD const * GetInvLorentzBoost() const override
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)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
double timeOffset() const
get TimeOffset
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const override
return a new event vertex
edm::ESWatcher< SimBeamSpotObjectsRcd > parameterWatcher_
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
double z() const
get Z position
double emittance() const
get Emittance
double phi() const
get Phi
double betaStar() const
get BetaStar
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
double x() const
get X position
double y() const
get Y position
long double T