CMS 3D CMS Logo

BetafuncEvtVtxGenerator.cc
Go to the documentation of this file.
1 /*
2 ________________________________________________________________________
3 
4  BetafuncEvtVtxGenerator
5 
6  Smear vertex according to the Beta function on the transverse plane
7  and a Gaussian on the z axis. It allows the beam to have a crossing
8  angle (slopes dxdz and dydz).
9 
10  Based on GaussEvtVtxGenerator
11  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
12 
13  FERMILAB
14  2006
15 ________________________________________________________________________
16 */
17 
23 
24 #include <CLHEP/Random/RandGaussQ.h>
25 #include <CLHEP/Units/SystemOfUnits.h>
26 #include <CLHEP/Units/GlobalPhysicalConstants.h>
27 #include "HepMC/SimpleVector.h"
28 
29 using CLHEP::cm;
30 using CLHEP::ns;
31 using CLHEP::radian;
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 distance units are in 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  update(iEventSetup);
59 }
60 
62  if (readDB_ && parameterWatcher_.check(iEventSetup)) {
63  edm::ESHandle<SimBeamSpotObjects> beamhandle = iEventSetup.getHandle(beamToken_);
64  fX0 = beamhandle->x() * cm;
65  fY0 = beamhandle->y() * cm;
66  fZ0 = beamhandle->z() * cm;
67  fSigmaZ = beamhandle->sigmaZ() * cm;
68  fTimeOffset = beamhandle->timeOffset() * ns * c_light; // HepMC distance units are in mm
69  fbetastar = beamhandle->betaStar() * cm;
70  femittance = beamhandle->emittance() * cm;
71  setBoost(beamhandle->alpha() * radian, beamhandle->phi() * radian);
72  }
73 }
74 
75 HepMC::FourVector BetafuncEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) const {
76  double X, Y, Z;
77 
78  double tmp_sigz = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
79  Z = tmp_sigz + fZ0;
80 
81  double tmp_sigx = BetaFunction(Z, fZ0);
82  // need sqrt(2) for beamspot width relative to single beam width
83  tmp_sigx /= sqrt(2.0);
84  X = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigx) + fX0; // + Z*fdxdz ;
85 
86  double tmp_sigy = BetaFunction(Z, fZ0);
87  // need sqrt(2) for beamspot width relative to single beam width
88  tmp_sigy /= sqrt(2.0);
89  Y = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigy) + fY0; // + Z*fdydz;
90 
91  double tmp_sigt = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
92  double T = tmp_sigt + fTimeOffset;
93 
94  return HepMC::FourVector(X, Y, Z, T);
95 }
96 
97 double BetafuncEvtVtxGenerator::BetaFunction(double z, double z0) const {
98  return sqrt(femittance * (fbetastar + (((z - z0) * (z - z0)) / fbetastar)));
99 }
100 
102  //boost_.ResizeTo(4,4);
103  //boost_ = new TMatrixD(4,4);
104  TMatrixD tmpboost(4, 4);
105 
106  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
107 
108  // Lorentz boost to frame where the collision is head-on
109  // phi is the half crossing angle in the plane ZS
110  // alpha is the angle to the S axis from the X axis in the XY plane
111 
112  tmpboost(0, 0) = 1. / cos(phi);
113  tmpboost(0, 1) = -cos(alpha) * sin(phi);
114  tmpboost(0, 2) = -tan(phi) * sin(phi);
115  tmpboost(0, 3) = -sin(alpha) * sin(phi);
116  tmpboost(1, 0) = -cos(alpha) * tan(phi);
117  tmpboost(1, 1) = 1.;
118  tmpboost(1, 2) = cos(alpha) * tan(phi);
119  tmpboost(1, 3) = 0.;
120  tmpboost(2, 0) = 0.;
121  tmpboost(2, 1) = -cos(alpha) * sin(phi);
122  tmpboost(2, 2) = cos(phi);
123  tmpboost(2, 3) = -sin(alpha) * sin(phi);
124  tmpboost(3, 0) = -sin(alpha) * tan(phi);
125  tmpboost(3, 1) = 0.;
126  tmpboost(3, 2) = sin(alpha) * tan(phi);
127  tmpboost(3, 3) = 1.;
128 
129  tmpboost.Invert();
130  boost_ = tmpboost;
131  //boost_->Print();
132 }
133 
135  if (s >= 0) {
136  fSigmaZ = s;
137  } else {
138  throw cms::Exception("LogicError") << "Error in BetafuncEvtVtxGenerator::sigmaZ: "
139  << "Illegal resolution in Z (negative)";
140  }
141 }
142 
143 TMatrixD const* BetafuncEvtVtxGenerator::GetInvLorentzBoost() const { return &boost_; }
144 
147  desc.add<double>("X0", 0.0)->setComment("in cm");
148  desc.add<double>("Y0", 0.0)->setComment("in cm");
149  desc.add<double>("Z0", 0.0)->setComment("in cm");
150  desc.add<double>("SigmaZ", 0.0)->setComment("in cm");
151  desc.add<double>("BetaStar", 0.0)->setComment("in cm");
152  desc.add<double>("Emittance", 0.0)->setComment("in cm");
153  desc.add<double>("Alpha", 0.0)->setComment("in radians");
154  desc.add<double>("Phi", 0.0)->setComment("in radians");
155  desc.add<double>("TimeOffset", 0.0)->setComment("in ns");
156  desc.add<edm::InputTag>("src");
157  desc.add<bool>("readDB");
158  descriptions.add("BetafuncEvtVtxGenerator", desc);
159 }
double sigmaZ() const
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
TMatrixD const * GetInvLorentzBoost() const override
T sqrt(T t)
Definition: SSEVec.h:23
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
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const override
return a new event vertex
edm::ESWatcher< SimBeamSpotObjectsRcd > parameterWatcher_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
double emittance() const
double phi() const
get Phi, Alpha and TimeOffset
double betaStar() const
get BetaStar and Emittance
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
double x() const
get X, Y, Z position
long double T