![]() |
![]() |
#include <FastSimulation/Event/interface/BetaFuncPrimaryVertexGenerator.h>
Public Member Functions | |
BetaFuncPrimaryVertexGenerator (const edm::ParameterSet &vtx, const RandomEngine *engine) | |
Default constructor. | |
double | BetaFunction (double z, double z0) |
set resolution in Z in cm set mean in X in cm beta function | |
virtual void | generate () |
Generation process (to be implemented). | |
~BetaFuncPrimaryVertexGenerator () | |
Destructor. | |
Private Member Functions | |
TMatrixD * | inverseLorentzBoost () |
Private Attributes | |
double | alpha_ |
double | fbetastar |
double | femittance |
double | fSigmaZ |
double | fX0 |
double | fY0 |
double | fZ0 |
double | phi_ |
Definition at line 15 of file BetaFuncPrimaryVertexGenerator.h.
BetaFuncPrimaryVertexGenerator::BetaFuncPrimaryVertexGenerator | ( | const edm::ParameterSet & | vtx, | |
const RandomEngine * | engine | |||
) |
Default constructor.
Definition at line 8 of file BetaFuncPrimaryVertexGenerator.cc.
References PrimaryVertexGenerator::beamSpot_, fX0, fY0, fZ0, inverseLorentzBoost(), and PrimaryVertexGenerator::setBoost().
00009 : 00010 PrimaryVertexGenerator(engine), 00011 fX0(vtx.getParameter<double>("X0")), 00012 fY0(vtx.getParameter<double>("Y0")), 00013 fZ0(vtx.getParameter<double>("Z0")), 00014 fSigmaZ(vtx.getParameter<double>("SigmaZ")), 00015 alpha_(vtx.getParameter<double>("Alpha")), 00016 phi_(vtx.getParameter<double>("Phi")), 00017 fbetastar(vtx.getParameter<double>("BetaStar")), 00018 femittance(vtx.getParameter<double>("Emittance")) 00019 { 00020 00021 this->setBoost(inverseLorentzBoost()); 00022 beamSpot_ = math::XYZPoint(fX0,fY0,fZ0); 00023 00024 }
BetaFuncPrimaryVertexGenerator::~BetaFuncPrimaryVertexGenerator | ( | ) | [inline] |
double BetaFuncPrimaryVertexGenerator::BetaFunction | ( | double | z, | |
double | z0 | |||
) |
set resolution in Z in cm set mean in X in cm beta function
Definition at line 40 of file BetaFuncPrimaryVertexGenerator.cc.
References fbetastar, femittance, and funct::sqrt().
Referenced by generate().
void BetaFuncPrimaryVertexGenerator::generate | ( | ) | [virtual] |
Generation process (to be implemented).
Implements PrimaryVertexGenerator.
Definition at line 27 of file BetaFuncPrimaryVertexGenerator.cc.
References BetaFunction(), fSigmaZ, fX0, fY0, fZ0, RandomEngine::gaussShoot(), and PrimaryVertexGenerator::random.
00027 { 00028 00029 double tmp_sigz = random->gaussShoot(0., fSigmaZ); 00030 this->SetZ(tmp_sigz + fZ0); 00031 00032 double tmp_sigx = BetaFunction(tmp_sigz,fZ0); 00033 this->SetX(random->gaussShoot(fX0,tmp_sigx)); 00034 00035 double tmp_sigy = BetaFunction(tmp_sigz,fZ0); 00036 this->SetY(random->gaussShoot(fY0,tmp_sigy)); 00037 00038 }
TMatrixD * BetaFuncPrimaryVertexGenerator::inverseLorentzBoost | ( | ) | [private] |
Definition at line 48 of file BetaFuncPrimaryVertexGenerator.cc.
References alpha_, funct::cos(), phi_, and funct::sin().
Referenced by BetaFuncPrimaryVertexGenerator().
00048 { 00049 00050 TMatrixD* aBoost = 0; 00051 if ( fabs(alpha_) < 1E-12 && fabs(phi_) < 1E-12 ) return aBoost; 00052 00053 TMatrixD tmpboost(4,4); 00054 00055 // Lorentz boost to frame where the collision is head-on 00056 // phi is the half crossing angle in the plane ZS 00057 // alpha is the angle to the S axis from the X axis in the XY plane 00058 00059 double calpha = std::cos(alpha_); 00060 double salpha = std::sin(alpha_); 00061 double cphi = std::cos(phi_); 00062 double sphi = std::sin(phi_); 00063 double tphi = sphi/cphi; 00064 tmpboost(0,0) = 1./cphi; 00065 tmpboost(0,1) = - calpha*sphi; 00066 tmpboost(0,2) = - tphi*sphi; 00067 tmpboost(0,3) = - salpha*sphi; 00068 tmpboost(1,0) = - calpha*tphi; 00069 tmpboost(1,1) = 1.; 00070 tmpboost(1,2) = calpha*tphi; 00071 tmpboost(1,3) = 0.; 00072 tmpboost(2,0) = 0.; 00073 tmpboost(2,1) = -calpha*sphi; 00074 tmpboost(2,2) = cphi; 00075 tmpboost(2,3) = - salpha*sphi; 00076 tmpboost(3,0) = - salpha*tphi; 00077 tmpboost(3,1) = 0.; 00078 tmpboost(3,2) = salpha*tphi; 00079 tmpboost(3,3) = 1.; 00080 00081 tmpboost.Invert(); 00082 aBoost = new TMatrixD(tmpboost); 00083 00084 return aBoost; 00085 00086 }
double BetaFuncPrimaryVertexGenerator::alpha_ [private] |
Definition at line 40 of file BetaFuncPrimaryVertexGenerator.h.
Referenced by inverseLorentzBoost().
double BetaFuncPrimaryVertexGenerator::fbetastar [private] |
double BetaFuncPrimaryVertexGenerator::femittance [private] |
double BetaFuncPrimaryVertexGenerator::fSigmaZ [private] |
double BetaFuncPrimaryVertexGenerator::fX0 [private] |
Definition at line 38 of file BetaFuncPrimaryVertexGenerator.h.
Referenced by BetaFuncPrimaryVertexGenerator(), and generate().
double BetaFuncPrimaryVertexGenerator::fY0 [private] |
Definition at line 38 of file BetaFuncPrimaryVertexGenerator.h.
Referenced by BetaFuncPrimaryVertexGenerator(), and generate().
double BetaFuncPrimaryVertexGenerator::fZ0 [private] |
Definition at line 38 of file BetaFuncPrimaryVertexGenerator.h.
Referenced by BetaFuncPrimaryVertexGenerator(), and generate().
double BetaFuncPrimaryVertexGenerator::phi_ [private] |
Definition at line 40 of file BetaFuncPrimaryVertexGenerator.h.
Referenced by inverseLorentzBoost().