CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
19 
20 
26 
27 #include "CLHEP/Random/RandGaussQ.h"
28 #include "CLHEP/Units/GlobalSystemOfUnits.h"
29 #include "CLHEP/Units/GlobalPhysicalConstants.h"
30 //#include "CLHEP/Vector/ThreeVector.h"
31 #include "HepMC/SimpleVector.h"
32 
35 
36 #include <iostream>
37 
38 
39 
42 {
43  readDB_=p.getParameter<bool>("readDB");
44  if (!readDB_){
45  fX0 = p.getParameter<double>("X0")*cm;
46  fY0 = p.getParameter<double>("Y0")*cm;
47  fZ0 = p.getParameter<double>("Z0")*cm;
48  fSigmaZ = p.getParameter<double>("SigmaZ")*cm;
49  alpha_ = p.getParameter<double>("Alpha")*radian;
50  phi_ = p.getParameter<double>("Phi")*radian;
51  fbetastar = p.getParameter<double>("BetaStar")*cm;
52  femittance = p.getParameter<double>("Emittance")*cm; // this is not the normalized emittance
53  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light; // HepMC time units are mm
54 
55  if (fSigmaZ <= 0) {
56  throw cms::Exception("Configuration")
57  << "Error in BetafuncEvtVtxGenerator: "
58  << "Illegal resolution in Z (SigmaZ is negative)";
59  }
60  }
61 }
62 
64 {
65 }
66 
68  update(iEventSetup);
69 }
70 void BetafuncEvtVtxGenerator::beginRun(const edm::Run & , const edm::EventSetup& iEventSetup){
71  update(iEventSetup);
72 }
73 
75  if (readDB_ && parameterWatcher_.check(iEventSetup)){
77  iEventSetup.get<SimBeamSpotObjectsRcd>().get(beamhandle);
78 
79  fX0=beamhandle->fX0;
80  fY0=beamhandle->fY0;
81  fZ0=beamhandle->fZ0;
82  // falpha=beamhandle->fAlpha;
83  alpha_=beamhandle->fAlpha;
84  phi_=beamhandle->fPhi;
85  fSigmaZ=beamhandle->fSigmaZ;
86  fTimeOffset=beamhandle->fTimeOffset;
87  fbetastar=beamhandle->fbetastar;
88  femittance=beamhandle->femittance;
89 
90  //re-initialize the boost matrix
91  delete boost_;
92  boost_=0;
93  }
94 }
95 
96 //Hep3Vector* BetafuncEvtVtxGenerator::newVertex() {
97 HepMC::FourVector* BetafuncEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) {
98 
99 
100  double X,Y,Z;
101 
102  double tmp_sigz = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
103  Z = tmp_sigz + fZ0;
104 
105  double tmp_sigx = BetaFunction(Z,fZ0);
106  // need sqrt(2) for beamspot width relative to single beam width
107  tmp_sigx /= sqrt(2.0);
108  X = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigx) + fX0; // + Z*fdxdz ;
109 
110  double tmp_sigy = BetaFunction(Z,fZ0);
111  // need sqrt(2) for beamspot width relative to single beam width
112  tmp_sigy /= sqrt(2.0);
113  Y = CLHEP::RandGaussQ::shoot(engine, 0., tmp_sigy) + fY0; // + Z*fdydz;
114 
115  double tmp_sigt = CLHEP::RandGaussQ::shoot(engine, 0., fSigmaZ);
116  double T = tmp_sigt + fTimeOffset;
117 
118  if ( fVertex == 0 ) fVertex = new HepMC::FourVector();
119  fVertex->set(X,Y,Z,T);
120 
121  return fVertex;
122 }
123 
124 double BetafuncEvtVtxGenerator::BetaFunction(double z, double z0)
125 {
126  return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
127 
128 }
129 
130 
132 {
133  if (s>=0 ) {
134  fSigmaZ=s;
135  }
136  else {
137  throw cms::Exception("LogicError")
138  << "Error in BetafuncEvtVtxGenerator::sigmaZ: "
139  << "Illegal resolution in Z (negative)";
140  }
141 }
142 
144 
145  //alpha_ = 0;
146  //phi_ = 142.e-6;
147 
148  if (boost_ != 0 ) return boost_;
149 
150  //boost_.ResizeTo(4,4);
151  //boost_ = new TMatrixD(4,4);
152  TMatrixD tmpboost(4,4);
153 
154  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
155 
156  // Lorentz boost to frame where the collision is head-on
157  // phi is the half crossing angle in the plane ZS
158  // alpha is the angle to the S axis from the X axis in the XY plane
159 
160  tmpboost(0,0) = 1./cos(phi_);
161  tmpboost(0,1) = - cos(alpha_)*sin(phi_);
162  tmpboost(0,2) = - tan(phi_)*sin(phi_);
163  tmpboost(0,3) = - sin(alpha_)*sin(phi_);
164  tmpboost(1,0) = - cos(alpha_)*tan(phi_);
165  tmpboost(1,1) = 1.;
166  tmpboost(1,2) = cos(alpha_)*tan(phi_);
167  tmpboost(1,3) = 0.;
168  tmpboost(2,0) = 0.;
169  tmpboost(2,1) = - cos(alpha_)*sin(phi_);
170  tmpboost(2,2) = cos(phi_);
171  tmpboost(2,3) = - sin(alpha_)*sin(phi_);
172  tmpboost(3,0) = - sin(alpha_)*tan(phi_);
173  tmpboost(3,1) = 0.;
174  tmpboost(3,2) = sin(alpha_)*tan(phi_);
175  tmpboost(3,3) = 1.;
176 
177  tmpboost.Invert();
178  boost_ = new TMatrixD(tmpboost);
179  //boost_->Print();
180 
181  return boost_;
182 }
const double Z[kNumberCalorimeter]
T getParameter(std::string const &) const
void sigmaZ(double s=1.0)
set resolution in Z in cm
HepMC::FourVector * fVertex
BetafuncEvtVtxGenerator(const edm::ParameterSet &p)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define X(str)
Definition: MuonsGrabber.cc:48
float float float z
double BetaFunction(double z, double z0)
beta function
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
void update(const edm::EventSetup &iEventSetup)
edm::ESWatcher< SimBeamSpotObjectsRcd > parameterWatcher_
virtual TMatrixD * GetInvLorentzBoost() override
const T & get() const
Definition: EventSetup.h:56
virtual HepMC::FourVector * newVertex(CLHEP::HepRandomEngine *) override
return a new event vertex
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
virtual void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
long double T
Definition: Run.h:42