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