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/GlobalSystemOfUnits.h"
26 #include "CLHEP/Units/GlobalPhysicalConstants.h"
27 #include "HepMC/SimpleVector.h"
28 
30  readDB_ = p.getParameter<bool>("readDB");
31  if (!readDB_) {
32  fX0 = p.getParameter<double>("X0") * cm;
33  fY0 = p.getParameter<double>("Y0") * cm;
34  fZ0 = p.getParameter<double>("Z0") * cm;
35  fSigmaZ = p.getParameter<double>("SigmaZ") * cm;
36  fbetastar = p.getParameter<double>("BetaStar") * cm;
37  femittance = p.getParameter<double>("Emittance") * cm; // this is not the normalized emittance
38  fTimeOffset = p.getParameter<double>("TimeOffset") * ns * c_light; // HepMC distance units are in mm
39 
40  setBoost(p.getParameter<double>("Alpha") * radian, p.getParameter<double>("Phi") * radian);
41  if (fSigmaZ <= 0) {
42  throw cms::Exception("Configuration") << "Error in BetafuncEvtVtxGenerator: "
43  << "Illegal resolution in Z (SigmaZ is negative)";
44  }
45  }
46  if (readDB_) {
47  // NOTE: this is currently watching LS transitions, while it should watch Run transitions,
48  // even though in reality there is no Run Dependent MC (yet) in CMS
49  beamToken_ = esConsumes<SimBeamSpotObjects, SimBeamSpotObjectsRcd, edm::Transition::BeginLuminosityBlock>();
50  }
51 }
52 
54  update(iEventSetup);
55 }
56 
58  if (readDB_ && parameterWatcher_.check(iEventSetup)) {
59  edm::ESHandle<SimBeamSpotObjects> beamhandle = iEventSetup.getHandle(beamToken_);
60  fX0 = beamhandle->x() * cm;
61  fY0 = beamhandle->y() * cm;
62  fZ0 = beamhandle->z() * cm;
63  fSigmaZ = beamhandle->sigmaZ() * cm;
64  fTimeOffset = beamhandle->timeOffset() * ns * c_light; // HepMC distance units are in mm
65  fbetastar = beamhandle->betaStar() * cm;
66  femittance = beamhandle->emittance() * cm;
67  setBoost(beamhandle->alpha() * radian, beamhandle->phi() * radian);
68  }
69 }
70 
71 HepMC::FourVector BetafuncEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) const {
72  double X, Y, Z;
73 
74  double tmp_sigz = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
75  Z = tmp_sigz + fZ0;
76 
77  double tmp_sigx = BetaFunction(Z, fZ0);
78  // need sqrt(2) for beamspot width relative to single beam width
79  tmp_sigx /= sqrt(2.0);
80  X = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigx) + fX0; // + Z*fdxdz ;
81 
82  double tmp_sigy = BetaFunction(Z, fZ0);
83  // need sqrt(2) for beamspot width relative to single beam width
84  tmp_sigy /= sqrt(2.0);
85  Y = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigy) + fY0; // + Z*fdydz;
86 
87  double tmp_sigt = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
88  double T = tmp_sigt + fTimeOffset;
89 
90  return HepMC::FourVector(X, Y, Z, T);
91 }
92 
93 double BetafuncEvtVtxGenerator::BetaFunction(double z, double z0) const {
94  return sqrt(femittance * (fbetastar + (((z - z0) * (z - z0)) / fbetastar)));
95 }
96 
98  //boost_.ResizeTo(4,4);
99  //boost_ = new TMatrixD(4,4);
100  TMatrixD tmpboost(4, 4);
101 
102  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
103 
104  // Lorentz boost to frame where the collision is head-on
105  // phi is the half crossing angle in the plane ZS
106  // alpha is the angle to the S axis from the X axis in the XY plane
107 
108  tmpboost(0, 0) = 1. / cos(phi);
109  tmpboost(0, 1) = -cos(alpha) * sin(phi);
110  tmpboost(0, 2) = -tan(phi) * sin(phi);
111  tmpboost(0, 3) = -sin(alpha) * sin(phi);
112  tmpboost(1, 0) = -cos(alpha) * tan(phi);
113  tmpboost(1, 1) = 1.;
114  tmpboost(1, 2) = cos(alpha) * tan(phi);
115  tmpboost(1, 3) = 0.;
116  tmpboost(2, 0) = 0.;
117  tmpboost(2, 1) = -cos(alpha) * sin(phi);
118  tmpboost(2, 2) = cos(phi);
119  tmpboost(2, 3) = -sin(alpha) * sin(phi);
120  tmpboost(3, 0) = -sin(alpha) * tan(phi);
121  tmpboost(3, 1) = 0.;
122  tmpboost(3, 2) = sin(alpha) * tan(phi);
123  tmpboost(3, 3) = 1.;
124 
125  tmpboost.Invert();
126  boost_ = tmpboost;
127  //boost_->Print();
128 }
129 
131  if (s >= 0) {
132  fSigmaZ = s;
133  } else {
134  throw cms::Exception("LogicError") << "Error in BetafuncEvtVtxGenerator::sigmaZ: "
135  << "Illegal resolution in Z (negative)";
136  }
137 }
138 
139 TMatrixD const* BetafuncEvtVtxGenerator::GetInvLorentzBoost() const { return &boost_; }
140 
143  desc.add<double>("X0", 0.0)->setComment("in cm");
144  desc.add<double>("Y0", 0.0)->setComment("in cm");
145  desc.add<double>("Z0", 0.0)->setComment("in cm");
146  desc.add<double>("SigmaZ", 0.0)->setComment("in cm");
147  desc.add<double>("BetaStar", 0.0)->setComment("in cm");
148  desc.add<double>("Emittance", 0.0)->setComment("in cm");
149  desc.add<double>("Alpha", 0.0)->setComment("in radians");
150  desc.add<double>("Phi", 0.0)->setComment("in radians");
151  desc.add<double>("TimeOffset", 0.0)->setComment("in ns");
152  desc.add<edm::InputTag>("src");
153  desc.add<bool>("readDB");
154  descriptions.add("BetafuncEvtVtxGenerator", desc);
155 }
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
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
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