CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 
26 
27 namespace CLHEP {
28  class HepRandomEngine;
29 }
30 
32 public:
38  ~BetafuncEvtVtxGenerator() override;
39 
40  void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
41 
43  //virtual CLHEP::Hep3Vector * newVertex();
44  HepMC::FourVector newVertex(CLHEP::HepRandomEngine*) const override;
45 
46  TMatrixD const* GetInvLorentzBoost() const override;
47 
49  void sigmaZ(double s = 1.0);
50 
52  void X0(double m = 0) { fX0 = m; }
54  void Y0(double m = 0) { fY0 = m; }
56  void Z0(double m = 0) { fZ0 = m; }
57 
59  void betastar(double m = 0) { fbetastar = m; }
61  void emittance(double m = 0) { femittance = m; }
62 
64  double BetaFunction(double z, double z0) const;
65 
66 private:
67  void setBoost(double alpha, double phi);
68 
69 private:
70  bool readDB_;
71 
72  double fX0, fY0, fZ0;
73  double fSigmaZ;
74  //double fdxdz, fdydz;
76  // double falpha;
77  double fTimeOffset;
78 
79  TMatrixD boost_;
80 
81  void update(const edm::EventSetup& iEventSetup);
84 };
85 
86 #endif
float alpha
Definition: AMPTWrapper.h:105
void sigmaZ(double s=1.0)
set resolution in Z in cm
void X0(double m=0)
set mean in X in cm
BetafuncEvtVtxGenerator & operator=(const BetafuncEvtVtxGenerator &rhs)=delete
TMatrixD const * GetInvLorentzBoost() const override
BetafuncEvtVtxGenerator(const edm::ParameterSet &p)
void setBoost(double alpha, double phi)
void betastar(double m=0)
set beta_star
void Y0(double m=0)
set mean in Y in cm
edm::ESGetToken< SimBeamSpotObjects, SimBeamSpotObjectsRcd > beamToken_
void update(const edm::EventSetup &iEventSetup)
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const override
return a new event vertex
edm::ESWatcher< SimBeamSpotObjectsRcd > parameterWatcher_
void Z0(double m=0)
set mean in Z in cm
double BetaFunction(double z, double z0) const
beta function
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
void emittance(double m=0)
emittance (no the normalized)