CMS 3D CMS Logo

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 
27 
28 namespace CLHEP {
29  class HepRandomEngine;
30 }
31 
33 public:
39  ~BetafuncEvtVtxGenerator() override = default;
40 
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
42 
43  void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
44 
46  //virtual CLHEP::Hep3Vector * newVertex();
47  HepMC::FourVector newVertex(CLHEP::HepRandomEngine*) const override;
48 
49  TMatrixD const* GetInvLorentzBoost() const override;
50 
52  void sigmaZ(double s = 1.0);
53 
55  void X0(double m = 0) { fX0 = m; }
57  void Y0(double m = 0) { fY0 = m; }
59  void Z0(double m = 0) { fZ0 = 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) const;
68 
69 private:
70  void setBoost(double alpha, double phi);
71 
72 private:
73  bool readDB_;
74 
75  double fX0, fY0, fZ0;
76  double fSigmaZ;
77  //double fdxdz, fdydz;
79  // double falpha;
80  double fTimeOffset;
81 
82  TMatrixD boost_;
83 
84  void update(const edm::EventSetup& iEventSetup);
87 };
88 
89 #endif
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
double BetaFunction(double z, double z0) const
beta function
BetafuncEvtVtxGenerator(const edm::ParameterSet &p)
void setBoost(double alpha, double phi)
TMatrixD const * GetInvLorentzBoost() const override
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)
~BetafuncEvtVtxGenerator() override=default
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
void emittance(double m=0)
emittance (no the normalized)