CMS 3D CMS Logo

HLLHCEvtVtxGenerator.h
Go to the documentation of this file.
1 #ifndef IOMC_HLLHCEvtVtxGenerator_H
2 #define IOMC_HLLHCEvtVtxGenerator_H
3 
19 
20 #include <string>
21 
22 namespace CLHEP {
23  class RandFlat;
24 }
25 
26 namespace edm {
28 }
29 
31 public:
33 
36 
39 
40  ~HLLHCEvtVtxGenerator() override = default;
41 
42  void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
43 
44  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
45 
47  HepMC::FourVector newVertex(CLHEP::HepRandomEngine*) const override;
48 
49  TMatrixD const* GetInvLorentzBoost() const override { return nullptr; };
50 
51 private:
52  // Configurable parameters
53  double fMeanX, fMeanY, fMeanZ, fTimeOffset_c_light; //spatial and time offset for mean collision
54  double fEProton; // proton beam energy
55  double fCrossingAngle; // crossing angle
56  double fCrabFrequency; // crab cavity frequency
57  bool fRF800; // 800 MHz RF?
58  double fBetaCrossingPlane; // beta crossing plane (m)
59  double fBetaSeparationPlane; // beta separation plane (m)
60  double fHorizontalEmittance; // horizontal emittance
61  double fVerticalEmittance; // vertical emittance
62  double fBunchLength; // bunch length
63  double fCrabbingAngleCrossing; // crabbing angle crossing
64  double fCrabbingAngleSeparation; // crabbing angle separation
65 
66  // Parameters inferred from configurables
67  double gamma; // beam configurations
68  double beta;
69  double betagamma;
70  double oncc; // ratio of crabbing angle to crossing angle
71  double epsx; // normalized crossing emittance
72  double epss; // normalized separation emittance
73  double sigx; // size in x
74  double phiCR; // crossing angle * crab frequency
75 
76  //width for y plane
77  double sigma(double z, double epsilon, double beta, double betagamma) const;
78 
79  //density with crabbing
80  double integrandCC(double x, double z, double t) const;
81 
82  // 4D intensity
83  double intensity(double x, double y, double z, double t) const;
84 
85  // Read from DB
86  bool readDB_;
87  void update(const edm::EventSetup& iEventSetup);
90 };
91 
92 #endif
HLLHCEvtVtxGenerator & operator=(const HLLHCEvtVtxGenerator &rhs)=delete
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const override
return a new event vertex
HLLHCEvtVtxGenerator(const edm::ParameterSet &p)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double sigma(double z, double epsilon, double beta, double betagamma) const
edm::ESWatcher< SimBeamSpotHLLHCObjectsRcd > parameterWatcher_
double intensity(double x, double y, double z, double t) const
edm::ESGetToken< SimBeamSpotHLLHCObjects, SimBeamSpotHLLHCObjectsRcd > beamToken_
void update(const edm::EventSetup &iEventSetup)
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
HLT enums.
double integrandCC(double x, double z, double t) const
TMatrixD const * GetInvLorentzBoost() const override
~HLLHCEvtVtxGenerator() override=default