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 // $Id: BetafuncEvtVtxGenerator.cc,v 1.14 2012/01/18 16:22:25 vlimant Exp $
3 /*
4 ________________________________________________________________________
5 
6  BetafuncEvtVtxGenerator
7 
8  Smear vertex according to the Beta function on the transverse plane
9  and a Gaussian on the z axis. It allows the beam to have a crossing
10  angle (slopes dxdz and dydz).
11 
12  Based on GaussEvtVtxGenerator
13  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
14 
15  FERMILAB
16  2006
17 ________________________________________________________________________
18 */
19 
20 
21 
27 
28 #include "CLHEP/Random/RandGaussQ.h"
29 #include "CLHEP/Units/GlobalSystemOfUnits.h"
30 #include "CLHEP/Units/GlobalPhysicalConstants.h"
31 //#include "CLHEP/Vector/ThreeVector.h"
32 #include "HepMC/SimpleVector.h"
33 
36 
37 #include <iostream>
38 
39 
40 
43 {
44 
45  fRandom = new CLHEP::RandGaussQ(getEngine());
46 
47  readDB_=p.getParameter<bool>("readDB");
48  if (!readDB_){
49  fX0 = p.getParameter<double>("X0")*cm;
50  fY0 = p.getParameter<double>("Y0")*cm;
51  fZ0 = p.getParameter<double>("Z0")*cm;
52  fSigmaZ = p.getParameter<double>("SigmaZ")*cm;
53  alpha_ = p.getParameter<double>("Alpha")*radian;
54  phi_ = p.getParameter<double>("Phi")*radian;
55  fbetastar = p.getParameter<double>("BetaStar")*cm;
56  femittance = p.getParameter<double>("Emittance")*cm; // this is not the normalized emittance
57  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light; // HepMC time units are mm
58 
59  if (fSigmaZ <= 0) {
60  throw cms::Exception("Configuration")
61  << "Error in BetafuncEvtVtxGenerator: "
62  << "Illegal resolution in Z (SigmaZ is negative)";
63  }
64  }
65 
66 
67 }
68 
70 {
71  delete fRandom;
72 }
73 
75  update(iEventSetup);
76 }
78  update(iEventSetup);
79 }
80 
82  if (readDB_ && parameterWatcher_.check(iEventSetup)){
84  iEventSetup.get<SimBeamSpotObjectsRcd>().get(beamhandle);
85 
86  fX0=beamhandle->fX0;
87  fY0=beamhandle->fY0;
88  fZ0=beamhandle->fZ0;
89  // falpha=beamhandle->fAlpha;
90  alpha_=beamhandle->fAlpha;
91  phi_=beamhandle->fPhi;
92  fSigmaZ=beamhandle->fSigmaZ;
93  fTimeOffset=beamhandle->fTimeOffset;
94  fbetastar=beamhandle->fbetastar;
95  femittance=beamhandle->femittance;
96 
97  //re-initialize the boost matrix
98  delete boost_;
99  boost_=0;
100  }
101 }
102 
103 //Hep3Vector* BetafuncEvtVtxGenerator::newVertex() {
105 
106 
107  double X,Y,Z;
108 
109  double tmp_sigz = fRandom->fire(0., fSigmaZ);
110  Z = tmp_sigz + fZ0;
111 
112  double tmp_sigx = BetaFunction(Z,fZ0);
113  // need sqrt(2) for beamspot width relative to single beam width
114  tmp_sigx /= sqrt(2.0);
115  X = fRandom->fire(0.,tmp_sigx) + fX0; // + Z*fdxdz ;
116 
117  double tmp_sigy = BetaFunction(Z,fZ0);
118  // need sqrt(2) for beamspot width relative to single beam width
119  tmp_sigy /= sqrt(2.0);
120  Y = fRandom->fire(0.,tmp_sigy) + fY0; // + Z*fdydz;
121 
122  double tmp_sigt = fRandom->fire(0., fSigmaZ);
123  double T = tmp_sigt + fTimeOffset;
124 
125  if ( fVertex == 0 ) fVertex = new HepMC::FourVector();
126  fVertex->set(X,Y,Z,T);
127 
128  return fVertex;
129 }
130 
131 double BetafuncEvtVtxGenerator::BetaFunction(double z, double z0)
132 {
133  return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
134 
135 }
136 
137 
139 {
140  if (s>=0 ) {
141  fSigmaZ=s;
142  }
143  else {
144  throw cms::Exception("LogicError")
145  << "Error in BetafuncEvtVtxGenerator::sigmaZ: "
146  << "Illegal resolution in Z (negative)";
147  }
148 }
149 
151 
152  //alpha_ = 0;
153  //phi_ = 142.e-6;
154 
155  if (boost_ != 0 ) return boost_;
156 
157  //boost_.ResizeTo(4,4);
158  //boost_ = new TMatrixD(4,4);
159  TMatrixD tmpboost(4,4);
160 
161  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
162 
163  // Lorentz boost to frame where the collision is head-on
164  // phi is the half crossing angle in the plane ZS
165  // alpha is the angle to the S axis from the X axis in the XY plane
166 
167  tmpboost(0,0) = 1./cos(phi_);
168  tmpboost(0,1) = - cos(alpha_)*sin(phi_);
169  tmpboost(0,2) = - tan(phi_)*sin(phi_);
170  tmpboost(0,3) = - sin(alpha_)*sin(phi_);
171  tmpboost(1,0) = - cos(alpha_)*tan(phi_);
172  tmpboost(1,1) = 1.;
173  tmpboost(1,2) = cos(alpha_)*tan(phi_);
174  tmpboost(1,3) = 0.;
175  tmpboost(2,0) = 0.;
176  tmpboost(2,1) = - cos(alpha_)*sin(phi_);
177  tmpboost(2,2) = cos(phi_);
178  tmpboost(2,3) = - sin(alpha_)*sin(phi_);
179  tmpboost(3,0) = - sin(alpha_)*tan(phi_);
180  tmpboost(3,1) = 0.;
181  tmpboost(3,2) = sin(alpha_)*tan(phi_);
182  tmpboost(3,3) = 1.;
183 
184  tmpboost.Invert();
185  boost_ = new TMatrixD(tmpboost);
186  //boost_->Print();
187 
188  return boost_;
189 }
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:49
virtual void beginRun(edm::Run &, const edm::EventSetup &)
virtual void beginLuminosityBlock(edm::LuminosityBlock &, edm::EventSetup const &)
double double double z
double BetaFunction(double z, double z0)
beta function
virtual TMatrixD * GetInvLorentzBoost()
T sqrt(T t)
Definition: SSEVec.h:46
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:59
long double T
Definition: Run.h:33