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 /*
5 ________________________________________________________________________
6 
7  BetafuncEvtVtxGenerator
8 
9  Smear vertex according to the Beta function on the transverse plane
10  and a Gaussian on the z axis. It allows the beam to have a crossing
11  angle (dx/dz and dy/dz).
12 
13  Based on GaussEvtVtxGenerator.h
14  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
15 
16  FERMILAB
17  2006
18 ________________________________________________________________________
19 */
20 
24 
25 namespace CLHEP {
26  class HepRandomEngine;
27 }
28 
30 {
31 public:
33  virtual ~BetafuncEvtVtxGenerator();
34 
35  virtual void beginRun(const edm::Run & , const edm::EventSetup&) override;
36  virtual void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
37 
39  //virtual CLHEP::Hep3Vector * newVertex();
40  virtual HepMC::FourVector* newVertex(CLHEP::HepRandomEngine*) ;
41 
42  virtual TMatrixD* GetInvLorentzBoost();
43 
44 
46  void sigmaZ(double s=1.0);
47 
49  void X0(double m=0) { fX0=m; }
51  void Y0(double m=0) { fY0=m; }
53  void Z0(double m=0) { fZ0=m; }
54 
56  void Phi(double m=0) { phi_=m; }
58  void Alpha(double m=0) { alpha_=m; }
59 
61  void betastar(double m=0) { fbetastar=m; }
63  void emittance(double m=0) { femittance=m; }
64 
66  double BetaFunction(double z, double z0);
67 
68 private:
73 
74 private:
75 
76  bool readDB_;
77 
78  double alpha_, phi_;
79  //TMatrixD boost_;
80 
81  double fX0, fY0, fZ0;
82  double fSigmaZ;
83  //double fdxdz, fdydz;
85  // double falpha;
86  double fTimeOffset;
87 
88  void update(const edm::EventSetup& iEventSetup);
90 };
91 
92 #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
BetafuncEvtVtxGenerator(const edm::ParameterSet &p)
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)
virtual HepMC::FourVector * newVertex(CLHEP::HepRandomEngine *)
return a new event vertex
Definition: Run.h:43