CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HLLHCEvtVtxGenerator Class Reference

#include <HLLHCEvtVtxGenerator.h>

Inheritance diagram for HLLHCEvtVtxGenerator:
BaseEvtVtxGenerator edm::stream::EDProducer<>

Public Member Functions

TMatrixD const * GetInvLorentzBoost () const override
 
 HLLHCEvtVtxGenerator (const edm::ParameterSet &p)
 
 HLLHCEvtVtxGenerator (const HLLHCEvtVtxGenerator &p)=delete
 
HepMC::FourVector newVertex (CLHEP::HepRandomEngine *) const override
 return a new event vertex More...
 
HLLHCEvtVtxGeneratoroperator= (const HLLHCEvtVtxGenerator &rhs)=delete
 
 ~HLLHCEvtVtxGenerator () override
 
- Public Member Functions inherited from BaseEvtVtxGenerator
 BaseEvtVtxGenerator (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~BaseEvtVtxGenerator () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

double integrandCC (double x, double z, double t) const
 
double intensity (double x, double y, double z, double t) const
 
double sigma (double z, double epsilon, double beta, double betagamma) const
 

Private Attributes

const double alphax
 
const double alphay
 
const double beta
 
const double betagamma
 
const double bets
 
const double betx
 
const double epss
 
const double epssn
 
const double epsx
 
const double epsxn
 
const double fMeanX
 
const double fMeanY
 
const double fMeanZ
 
const double fTimeOffset
 
const double gamma
 
const double momeV
 
const double oncc
 
const double phi
 
const double phiCR
 
const bool RF800
 
const double sigs
 
const double sigx
 
const double wcc
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T...>
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T...>
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 26 of file HLLHCEvtVtxGenerator.h.

Constructor & Destructor Documentation

HLLHCEvtVtxGenerator::HLLHCEvtVtxGenerator ( const edm::ParameterSet p)

Definition at line 51 of file HLLHCEvtVtxGenerator.cc.

53  fMeanX(p.getParameter<double>("MeanXIncm") * cm),
54  fMeanY(p.getParameter<double>("MeanYIncm") * cm),
55  fMeanZ(p.getParameter<double>("MeanZIncm") * cm),
56  fTimeOffset(p.getParameter<double>("TimeOffsetInns") * ns * c_light),
57  momeV(p.getParameter<double>("EprotonInGeV") * 1e9),
58  gamma(momeV / pmass + 1.0),
59  beta(std::sqrt((1.0 - 1.0 / gamma) * ((1.0 + 1.0 / gamma)))),
60  betagamma(beta * gamma),
61  phi(p.getParameter<double>("CrossingAngleInurad") * 1e-6),
62  wcc(p.getParameter<double>("CrabFrequencyInMHz") * 1e6),
63  RF800(p.getParameter<bool>("RF800")),
64  betx(p.getParameter<double>("BetaCrossingPlaneInm")),
65  bets(p.getParameter<double>("BetaSeparationPlaneInm")),
66  epsxn(p.getParameter<double>("HorizontalEmittance")),
67  epssn(p.getParameter<double>("VerticalEmittance")),
68  sigs(p.getParameter<double>("BunchLengthInm")),
69  alphax(p.getParameter<double>("CrabbingAngleCrossingInurad") * 1e-6),
70  alphay(p.getParameter<double>("CrabbingAngleSeparationInurad") * 1e-6),
71  oncc(alphax / phi),
72  epsx(epsxn / (betagamma)),
73  epss(epsx),
74  sigx(std::sqrt(epsx * betx)),
75  phiCR(oncc * phi)
76 
77 {}
BaseEvtVtxGenerator(const edm::ParameterSet &)
T sqrt(T t)
Definition: SSEVec.h:19
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HLLHCEvtVtxGenerator::HLLHCEvtVtxGenerator ( const HLLHCEvtVtxGenerator p)
delete

Copy constructor

HLLHCEvtVtxGenerator::~HLLHCEvtVtxGenerator ( )
override

Definition at line 79 of file HLLHCEvtVtxGenerator.cc.

79 {}

Member Function Documentation

void HLLHCEvtVtxGenerator::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 29 of file HLLHCEvtVtxGenerator.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and submitPVResolutionJobs::desc.

29  {
31  desc.add<double>("MeanXIncm", 0.0);
32  desc.add<double>("MeanYIncm", 0.0);
33  desc.add<double>("MeanZIncm", 0.0);
34  desc.add<double>("TimeOffsetInns", 0.0);
35  desc.add<double>("EprotonInGeV", 7000.0);
36  desc.add<double>("CrossingAngleInurad", 510.0);
37  desc.add<double>("CrabbingAngleCrossingInurad", 380.0);
38  desc.add<double>("CrabbingAngleSeparationInurad", 0.0);
39  desc.add<double>("CrabFrequencyInMHz", 400.0);
40  desc.add<bool>("RF800", false);
41  desc.add<double>("BetaCrossingPlaneInm", 0.20);
42  desc.add<double>("BetaSeparationPlaneInm", 0.20);
43  desc.add<double>("HorizontalEmittance", 2.5e-06);
44  desc.add<double>("VerticalEmittance", 2.05e-06);
45  desc.add<double>("BunchLengthInm", 0.09);
46  desc.add<edm::InputTag>("src");
47  desc.add<bool>("readDB");
48  descriptions.add("HLLHCEvtVtxGenerator", desc);
49 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
TMatrixD const* HLLHCEvtVtxGenerator::GetInvLorentzBoost ( ) const
inlineoverridevirtual

This method - and the comment - is a left-over from COBRA-OSCAR time : return the last generated event vertex. If no vertex has been generated yet, a NULL pointer is returned.

Implements BaseEvtVtxGenerator.

Definition at line 43 of file HLLHCEvtVtxGenerator.h.

43 { return nullptr; };
double HLLHCEvtVtxGenerator::integrandCC ( double  x,
double  z,
double  t 
) const
private

Definition at line 146 of file HLLHCEvtVtxGenerator.cc.

References bets, betx, funct::cos(), funct::exp(), isotrackApplyRegressor::k, relval_steps::k2, local_c_light, phi, phiCR, mps_fire::result, RF800, alignCSCRings::s, sigs, sigx, funct::sin(), mathSSE::sqrt(), wcc, and z.

Referenced by intensity().

146  {
147  constexpr double local_c_light = 2.99792458e8;
148 
149  const double k = wcc / local_c_light * two_pi;
150  const double k2 = k * k;
151  const double cos = std::cos(phi / 2.0);
152  const double sin = std::sin(phi / 2.0);
153  const double cos2 = cos * cos;
154  const double sin2 = sin * sin;
155 
156  const double sigx2 = sigx * sigx;
157  const double sigmax2 = sigx2 * (1 + z * z / (betx * betx));
158 
159  const double sigs2 = sigs * sigs;
160 
161  constexpr double factorRMSgauss4 =
162  1. / sqrt2 / gamma34 * gamma14; // # Factor to take rms sigma as input of the supergaussian
163  constexpr double NormFactorGauss4 = sqrt2to5 * gamma54 * gamma54;
164 
165  const double sinCR = std::sin(phiCR / 2.0);
166  const double sinCR2 = sinCR * sinCR;
167 
168  double result = -1.0;
169 
170  if (!RF800) {
171  const double norm = 2.0 / (two_pi * sigs2);
172  const double cosks = std::cos(k * z);
173  const double sinkct = std::sin(k * ct);
174  result = norm *
175  std::exp(-ct * ct / sigs2 - z * z * cos2 / sigs2 -
176  1.0 / (4 * k2 * sigmax2) *
177  (
178  //-4*cosks*cosks * sinkct*sinkct * sinCR2 // comes from integral over x
179  -8 * z * k * std::sin(k * z) * std::cos(k * ct) * sin * sinCR + 2 * sinCR2 -
180  std::cos(2 * k * (z - ct)) * sinCR2 - std::cos(2 * k * (z + ct)) * sinCR2 +
181  4 * k2 * z * z * sin2) -
182  x * x * (cos2 / sigmax2 + sin2 / sigs2) // contribution from x integrand
183  + x * ct * sin / sigs2 // contribution from x integrand
184  + 2 * x * cos * cosks * sinkct * sinCR / k / sigmax2 // contribution from x integrand
185  //+(2*ct/k)*np.cos(k*s)*np.sin(k*ct) *(sin*sinCR)/(sigs2*cos) # small term
186  //+ct**2*(sin2/sigs4)/(cos2/sigmax2) # small term
187  ) /
188  (1.0 + (z * z) / (betx * betx)) / std::sqrt(1.0 + (z * z) / (bets * bets));
189 
190  } else {
191  const double norm = 2.0 / (NormFactorGauss4 * sigs2 * factorRMSgauss4);
192  const double sigs4 = sigs2 * sigs2 * factorRMSgauss4 * factorRMSgauss4;
193  const double cosks = std::cos(k * z);
194  const double sinct = std::sin(k * ct);
195  result =
196  norm *
197  std::exp(-ct * ct * ct * ct / sigs4 - z * z * z * z * cos2 * cos2 / sigs4 - 6 * ct * ct * z * z * cos2 / sigs4 -
198  sin2 / (4 * k2 * sigmax2) *
199  (2 + 4 * k2 * z * z - std::cos(2 * k * (z - ct)) - std::cos(2 * k * (z + ct)) -
200  8 * k * s * std::cos(k * ct) * std::sin(k * z) - 4 * cosks * cosks * sinct * sinct)) /
201  std::sqrt(1 + z * z / (betx * betx)) / std::sqrt(1 + z * z / (bets * bets));
202  }
203 
204  return result;
205 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
tuple result
Definition: mps_fire.py:311
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static constexpr float local_c_light
Definition: pTFrom2Stubs.cc:9
double HLLHCEvtVtxGenerator::intensity ( double  x,
double  y,
double  z,
double  t 
) const
private

Definition at line 126 of file HLLHCEvtVtxGenerator.cc.

References alphay, betagamma, bets, c, funct::cos(), PVValHelper::dy, epssn, funct::exp(), integrandCC(), sigma(), funct::sin(), submitPVValidationJobs::t, and wcc.

Referenced by newVertex().

126  {
127  //---c in m/s --- remember t is already in meters
128  constexpr double c = 2.99792458e+8; // m/s
129 
130  const double sigmay = sigma(z, epssn, bets, betagamma);
131 
132  const double alphay_mod = alphay * std::cos(wcc * (z - t) / c);
133 
134  const double cay = std::cos(alphay_mod);
135  const double say = std::sin(alphay_mod);
136 
137  const double dy = -(z - t) * say - y * cay;
138 
139  const double xzt_density = integrandCC(x, z, t);
140 
141  const double norm = two_pi * sigmay;
142 
143  return (std::exp(-dy * dy / (sigmay * sigmay)) * xzt_density / norm);
144 }
const edm::EventSetup & c
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double integrandCC(double x, double z, double t) const
double sigma(double z, double epsilon, double beta, double betagamma) const
HepMC::FourVector HLLHCEvtVtxGenerator::newVertex ( CLHEP::HepRandomEngine *  engine) const
overridevirtual

return a new event vertex

Implements BaseEvtVtxGenerator.

Definition at line 81 of file HLLHCEvtVtxGenerator.cc.

References betagamma, bets, betx, submitPVResolutionJobs::count, epssn, epsxn, fMeanX, fMeanY, fMeanZ, fTimeOffset, mps_fire::i, intensity(), sigma(), sigs, submitPVValidationJobs::t, x, y, and z.

81  {
82  double imax = intensity(0., 0., 0., 0.);
83 
84  double x(0.), y(0.), z(0.), t(0.), i(0.);
85 
86  int count = 0;
87 
88  auto shoot = [&]() { return CLHEP::RandFlat::shoot(engine); };
89 
90  do {
91  z = (shoot() - 0.5) * 6.0 * sigs;
92  t = (shoot() - 0.5) * 6.0 * sigs;
93  x = (shoot() - 0.5) * 12.0 * sigma(0.0, epsxn, betx, betagamma);
94  y = (shoot() - 0.5) * 12.0 * sigma(0.0, epssn, bets, betagamma);
95 
96  i = intensity(x, y, z, t);
97 
98  if (i > imax)
99  edm::LogError("Too large intensity") << "i>imax : " << i << " > " << imax << endl;
100  ++count;
101  } while ((i < imax * shoot()) && count < 10000);
102 
103  if (count > 9999)
104  edm::LogError("Too many tries ") << " count : " << count << endl;
105 
106  //---convert to mm
107  x *= 1000.0;
108  y *= 1000.0;
109  z *= 1000.0;
110  t *= 1000.0;
111 
112  x += fMeanX;
113  y += fMeanY;
114  z += fMeanZ;
115  t += fTimeOffset;
116 
117  return HepMC::FourVector(x, y, z, t);
118 }
Log< level::Error, false > LogError
double intensity(double x, double y, double z, double t) const
double sigma(double z, double epsilon, double beta, double betagamma) const
HLLHCEvtVtxGenerator& HLLHCEvtVtxGenerator::operator= ( const HLLHCEvtVtxGenerator rhs)
delete

Copy assignment operator

double HLLHCEvtVtxGenerator::sigma ( double  z,
double  epsilon,
double  beta,
double  betagamma 
) const
private

Definition at line 120 of file HLLHCEvtVtxGenerator.cc.

References mathSSE::sqrt().

Referenced by intensity(), and newVertex().

120  {
121  double sigma = std::sqrt(epsilon * (beta + z * z / beta) / betagamma);
122 
123  return sigma;
124 }
T sqrt(T t)
Definition: SSEVec.h:19
double sigma(double z, double epsilon, double beta, double betagamma) const

Member Data Documentation

const double HLLHCEvtVtxGenerator::alphax
private

Definition at line 80 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::alphay
private

Definition at line 83 of file HLLHCEvtVtxGenerator.h.

Referenced by intensity().

const double HLLHCEvtVtxGenerator::beta
private

Definition at line 52 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::betagamma
private

Definition at line 53 of file HLLHCEvtVtxGenerator.h.

Referenced by intensity(), and newVertex().

const double HLLHCEvtVtxGenerator::bets
private

Definition at line 68 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC(), intensity(), and newVertex().

const double HLLHCEvtVtxGenerator::betx
private

Definition at line 65 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC(), and newVertex().

const double HLLHCEvtVtxGenerator::epss
private

Definition at line 92 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::epssn
private

Definition at line 74 of file HLLHCEvtVtxGenerator.h.

Referenced by intensity(), and newVertex().

const double HLLHCEvtVtxGenerator::epsx
private

Definition at line 89 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::epsxn
private

Definition at line 71 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fMeanX
private

Definition at line 43 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fMeanY
private

Definition at line 43 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fMeanZ
private

Definition at line 43 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fTimeOffset
private

Definition at line 43 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::gamma
private

Definition at line 51 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::momeV
private

Definition at line 50 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::oncc
private

Definition at line 86 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::phi
private
const double HLLHCEvtVtxGenerator::phiCR
private

Definition at line 98 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC().

const bool HLLHCEvtVtxGenerator::RF800
private

Definition at line 62 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC().

const double HLLHCEvtVtxGenerator::sigs
private

Definition at line 77 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC(), and newVertex().

const double HLLHCEvtVtxGenerator::sigx
private

Definition at line 95 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC().

const double HLLHCEvtVtxGenerator::wcc
private

Definition at line 59 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC(), and intensity().