CMS 3D CMS Logo

HLLHCEvtVtxGenerator.cc
Go to the documentation of this file.
1 // $Id: HLLHCEvtVtxGenerator_Fix.cc, v 1.0 2015/03/15 10:32:11 Exp $
2 
5 
10 
11 
12 #include "CLHEP/Random/RandFlat.h"
13 #include "CLHEP/Units/GlobalSystemOfUnits.h"
14 #include "CLHEP/Units/GlobalPhysicalConstants.h"
15 #include "HepMC/SimpleVector.h"
16 
17 using namespace std;
18 
19 namespace {
20 
21  constexpr double pmass = 0.9382720813e9; // eV
22  constexpr double gamma34 = 1.22541670246517764513; // Gamma(3/4)
23  constexpr double gamma14 = 3.62560990822190831193; // Gamma(1/4)
24  constexpr double gamma54 = 0.90640247705547798267; // Gamma(5/4)
25  constexpr double sqrt2 = 1.41421356237;
26  constexpr double sqrt2to5 = 5.65685424949;
27  constexpr double two_pi=2.0*M_PI;
28 }
29 
31 {
33  desc.add<double>("MeanXIncm",0.0);
34  desc.add<double>("MeanYIncm",0.0);
35  desc.add<double>("MeanZIncm",0.0);
36  desc.add<double>("TimeOffsetInns",0.0);
37  desc.add<double>("EprotonInGeV", 7000.0);
38  desc.add<double>("CrossingAngleInurad", 510.0);
39  desc.add<double>("CrabbingAngleCrossingInurad", 380.0);
40  desc.add<double>("CrabbingAngleSeparationInurad", 0.0);
41  desc.add<double>("CrabFrequencyInMHz", 400.0);
42  desc.add<bool>("RF800", false);
43  desc.add<double>("BetaCrossingPlaneInm", 0.20);
44  desc.add<double>("BetaSeparationPlaneInm", 0.20);
45  desc.add<double>("HorizontalEmittance", 2.5e-06);
46  desc.add<double>("VerticalEmittance", 2.05e-06);
47  desc.add<double>("BunchLengthInm", 0.09);
48  desc.add<edm::InputTag>("src");
49  desc.add<bool>("readDB");
50  descriptions.add("HLLHCEvtVtxGenerator",desc);
51 }
52 
55  fMeanX(p.getParameter<double>("MeanXIncm")*cm),
56  fMeanY(p.getParameter<double>("MeanYIncm")*cm),
57  fMeanZ(p.getParameter<double>("MeanZIncm")*cm),
58  fTimeOffset(p.getParameter<double>("TimeOffsetInns")*ns*c_light),
59  momeV(p.getParameter<double>("EprotonInGeV")*1e9),
60  gamma(momeV/pmass + 1.0),
61  beta(std::sqrt((1.0 - 1.0/gamma)*((1.0 + 1.0/gamma)))),
62  betagamma(beta*gamma),
63  phi(p.getParameter<double>("CrossingAngleInurad")*1e-6),
64  wcc(p.getParameter<double>("CrabFrequencyInMHz")*1e6),
65  RF800(p.getParameter<bool>("RF800")),
66  betx(p.getParameter<double>("BetaCrossingPlaneInm")),
67  bets(p.getParameter<double>("BetaSeparationPlaneInm")),
68  epsxn(p.getParameter<double>("HorizontalEmittance")),
69  epssn(p.getParameter<double>("VerticalEmittance")),
70  sigs(p.getParameter<double>("BunchLengthInm")),
71  alphax(p.getParameter<double>("CrabbingAngleCrossingInurad")*1e-6),
72  alphay(p.getParameter<double>("CrabbingAngleSeparationInurad")*1e-6),
73  oncc(alphax/phi),
74  epsx(epsxn/(betagamma)),
75  epss(epsx),
76  sigx(std::sqrt(epsx*betx)),
77  phiCR(oncc*phi)
78 
79 {
80 
81 }
82 
84 {
85 }
86 
87 
88 HepMC::FourVector HLLHCEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) const {
89 
90  double imax=intensity(0.,0.,0.,0.);
91 
92  double x(0.),y(0.),z(0.),t(0.),i(0.);
93 
94  int count=0;
95 
96  auto shoot = [&](){ return CLHEP::RandFlat::shoot(engine); };
97 
98  do {
99  z=(shoot()-0.5)*6.0*sigs;
100  t=(shoot()-0.5)*6.0*sigs;
101  x=(shoot()-0.5)*12.0*sigma(0.0,epsxn,betx,betagamma);
102  y=(shoot()-0.5)*12.0*sigma(0.0,epssn,bets,betagamma);
103 
104  i=intensity(x,y,z,t);
105 
106  if (i>imax) edm::LogError("Too large intensity")
107  << "i>imax : "<<i<<" > "<<imax<<endl;
108  ++count;
109  } while ((i<imax*shoot())&&count<10000);
110 
111  if (count>9999) edm::LogError("Too many tries ")
112  << " count : "<<count<<endl;
113 
114  //---convert to mm
115  x*=1000.0;
116  y*=1000.0;
117  z*=1000.0;
118  t*=1000.0;
119 
120  x+=fMeanX;
121  y+=fMeanY;
122  z+=fMeanZ;
123  t+=fTimeOffset;
124 
125  return HepMC::FourVector( x, y, z, t );
126 }
127 
128 
129 
130 double HLLHCEvtVtxGenerator::sigma(double z, double epsilon, double beta, double betagamma) const
131 {
132  double sigma=std::sqrt(epsilon*(beta+z*z/beta)/betagamma);
133 
134  return sigma;
135 }
136 
138  double y,
139  double z,
140  double t) const {
141  //---c in m/s --- remember t is already in meters
142  constexpr double c= 2.99792458e+8; // m/s
143 
144  const double sigmay=sigma(z,epssn,bets,betagamma);
145 
146  const double alphay_mod=alphay*std::cos(wcc*(z-t)/c);
147 
148  const double cay=std::cos(alphay_mod);
149  const double say=std::sin(alphay_mod);
150 
151  const double dy=-(z-t)*say-y*cay;
152 
153  const double xzt_density= integrandCC(x,z,t);
154 
155  const double norm = two_pi*sigmay;
156 
157  return ( std::exp(-dy*dy/(sigmay*sigmay))*
158  xzt_density/norm );
159 }
160 
162  double z,
163  double ct) const {
164  constexpr double local_c_light = 2.99792458e8;
165 
166  const double k = wcc/local_c_light*two_pi;
167  const double k2 = k*k;
168  const double cos = std::cos(phi/2.0);
169  const double sin = std::sin(phi/2.0);
170  const double cos2 = cos*cos;
171  const double sin2 = sin*sin;
172 
173  const double sigx2 = sigx*sigx;
174  const double sigmax2=sigx2*(1+z*z/(betx*betx));
175 
176  const double sigs2 = sigs*sigs;
177 
178  constexpr double factorRMSgauss4 = 1./sqrt2/gamma34 * gamma14; // # Factor to take rms sigma as input of the supergaussian
179  constexpr double NormFactorGauss4 = sqrt2to5 * gamma54 * gamma54;
180 
181  const double sinCR = std::sin(phiCR/2.0);
182  const double sinCR2 = sinCR*sinCR;
183 
184  double result = -1.0;
185 
186  if( !RF800 ) {
187  const double norm =2.0/(two_pi*sigs2);
188  const double cosks = std::cos(k*z);
189  const double sinkct = std::sin(k*ct);
190  result = norm*std::exp(-ct*ct/sigs2
191  -z*z*cos2/sigs2
192  -1.0/(4*k2*sigmax2)*(
193  //-4*cosks*cosks * sinkct*sinkct * sinCR2 // comes from integral over x
194  -8*z*k*std::sin(k*z)*std::cos(k*ct) * sin * sinCR
195  +2 * sinCR2
196  -std::cos(2*k*(z-ct)) * sinCR2
197  -std::cos(2*k*(z+ct)) * sinCR2
198  +4*k2*z*z *sin2
199  )
200  - x*x*(cos2/sigmax2 + sin2/sigs2) // contribution from x integrand
201  + x*ct*sin/sigs2 // contribution from x integrand
202  + 2*x*cos*cosks*sinkct*sinCR/k/sigmax2 // contribution from x integrand
203  //+(2*ct/k)*np.cos(k*s)*np.sin(k*ct) *(sin*sinCR)/(sigs2*cos) # small term
204  //+ct**2*(sin2/sigs4)/(cos2/sigmax2) # small term
205  )/(1.0+(z*z)/(betx*betx))/std::sqrt(1.0+(z*z)/(bets*bets));
206 
207  } else {
208 
209  const double norm = 2.0/(NormFactorGauss4*sigs2*factorRMSgauss4);
210  const double sigs4=sigs2*sigs2*factorRMSgauss4*factorRMSgauss4;
211  const double cosks = std::cos(k*z);
212  const double sinct = std::sin(k*ct);
213  result = norm*std::exp(
214  -ct*ct*ct*ct/sigs4
215  -z*z*z*z*cos2*cos2/sigs4
216  -6*ct*ct*z*z*cos2/sigs4
217  -sin2/(4*k2*sigmax2)*(
218  2
219  +4*k2*z*z
220  -std::cos(2*k*(z-ct))
221  -std::cos(2*k*(z+ct))
222  -8*k*s*std::cos(k*ct)*std::sin(k*z)
223  -4 * cosks*cosks * sinct*sinct)
224  )/std::sqrt(1+z*z/(betx*betx))/std::sqrt(1+z*z/(bets*bets));
225  }
226 
227  return result;
228 }
229 
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const override
return a new event vertex
HLLHCEvtVtxGenerator(const edm::ParameterSet &p)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define constexpr
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double integrandCC(double x, double z, double t) const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define M_PI
int k[5][pyjets_maxn]
double intensity(double x, double y, double z, double t) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double sigma(double z, double epsilon, double beta, double betagamma) const