CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BetafuncEvtVtxGenerator.h
Go to the documentation of this file.
1 #ifndef IOMC_BetafuncEvtVtxGenerator_H
2 #define IOMC_BetafuncEvtVtxGenerator_H
3 
4 // $Id: BetafuncEvtVtxGenerator.h,v 1.9 2013/02/27 18:41:06 wmtan Exp $
5 /*
6 ________________________________________________________________________
7 
8  BetafuncEvtVtxGenerator
9 
10  Smear vertex according to the Beta function on the transverse plane
11  and a Gaussian on the z axis. It allows the beam to have a crossing
12  angle (dx/dz and dy/dz).
13 
14  Based on GaussEvtVtxGenerator.h
15  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
16 
17  FERMILAB
18  2006
19 ________________________________________________________________________
20 */
21 
25 
26 namespace CLHEP {
27  class RandGaussQ;
28 }
29 
31 {
32 public:
34  virtual ~BetafuncEvtVtxGenerator();
35 
36  virtual void beginRun(const edm::Run & , const edm::EventSetup&) override;
37  virtual void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
38 
40  //virtual CLHEP::Hep3Vector * newVertex();
41  virtual HepMC::FourVector* newVertex() ;
42 
43  virtual TMatrixD* GetInvLorentzBoost();
44 
45 
47  void sigmaZ(double s=1.0);
48 
50  void X0(double m=0) { fX0=m; }
52  void Y0(double m=0) { fY0=m; }
54  void Z0(double m=0) { fZ0=m; }
55 
57  void Phi(double m=0) { phi_=m; }
59  void Alpha(double m=0) { alpha_=m; }
60 
62  void betastar(double m=0) { fbetastar=m; }
64  void emittance(double m=0) { femittance=m; }
65 
67  double BetaFunction(double z, double z0);
68 
69 private:
74 
75 private:
76 
77  bool readDB_;
78 
79  double alpha_, phi_;
80  //TMatrixD boost_;
81 
82  double fX0, fY0, fZ0;
83  double fSigmaZ;
84  //double fdxdz, fdydz;
86  // double falpha;
87  double fTimeOffset;
88 
89  CLHEP::RandGaussQ* fRandom ;
90 
91  void update(const edm::EventSetup& iEventSetup);
93 };
94 
95 #endif
void sigmaZ(double s=1.0)
set resolution in Z in cm
void Phi(double m=0)
set half crossing angle
void X0(double m=0)
set mean in X in cm
virtual HepMC::FourVector * newVertex()
return a new event vertex
BetafuncEvtVtxGenerator(const edm::ParameterSet &p)
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()
void betastar(double m=0)
set beta_star
void Y0(double m=0)
set mean in Y in cm
void update(const edm::EventSetup &iEventSetup)
edm::ESWatcher< SimBeamSpotObjectsRcd > parameterWatcher_
void Z0(double m=0)
set mean in Z in cm
void Alpha(double m=0)
angle between crossing plane and horizontal plane
BetafuncEvtVtxGenerator & operator=(const BetafuncEvtVtxGenerator &rhs)
virtual void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
void emittance(double m=0)
emittance (no the normalized)
Definition: Run.h:36