CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

virtual TMatrixD * GetInvLorentzBoost ()
 
 HLLHCEvtVtxGenerator (const edm::ParameterSet &p)
 
virtual HepMC::FourVector * newVertex (CLHEP::HepRandomEngine *)
 return a new event vertex More...
 
virtual ~HLLHCEvtVtxGenerator ()
 
- Public Member Functions inherited from BaseEvtVtxGenerator
 BaseEvtVtxGenerator (const edm::ParameterSet &)
 
virtual HepMC::FourVector * lastVertex ()
 
virtual void produce (edm::Event &, const edm::EventSetup &) override
 
virtual ~BaseEvtVtxGenerator ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector
< edm::ProductResolverIndex >
const & 
indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, std::unordered_multimap< std::string, edm::ProductResolverIndex > const &iIndicies, std::string const &moduleLabel)
 
virtual ~ProducerBase () noexcept(false)
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector
< ProductResolverIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

 HLLHCEvtVtxGenerator (const HLLHCEvtVtxGenerator &p)
 
double integrandCC (double x, double z, double t) const
 
double intensity (double x, double y, double z, double t) const
 
HLLHCEvtVtxGeneratoroperator= (const HLLHCEvtVtxGenerator &rhs)
 
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::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
- Protected Attributes inherited from BaseEvtVtxGenerator
TMatrixD * boost_
 
double fTimeOffset
 
HepMC::FourVector * fVertex
 

Detailed Description

Definition at line 26 of file HLLHCEvtVtxGenerator.h.

Constructor & Destructor Documentation

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

Definition at line 53 of file HLLHCEvtVtxGenerator.cc.

References phi().

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)))),
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),
77  phiCR(oncc*phi)
78 
79 {
80 
81 }
T getParameter(std::string const &) const
BaseEvtVtxGenerator(const edm::ParameterSet &)
T sqrt(T t)
Definition: SSEVec.h:18
HLLHCEvtVtxGenerator::~HLLHCEvtVtxGenerator ( )
virtual

Definition at line 83 of file HLLHCEvtVtxGenerator.cc.

84 {
85 }
HLLHCEvtVtxGenerator::HLLHCEvtVtxGenerator ( const HLLHCEvtVtxGenerator p)
private

Copy constructor

Member Function Documentation

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

Definition at line 30 of file HLLHCEvtVtxGenerator.cc.

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

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 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual TMatrixD* HLLHCEvtVtxGenerator::GetInvLorentzBoost ( )
inlinevirtual

Implements BaseEvtVtxGenerator.

Definition at line 39 of file HLLHCEvtVtxGenerator.h.

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

Definition at line 166 of file HLLHCEvtVtxGenerator.cc.

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

Referenced by intensity().

168  {
169  constexpr double local_c_light = 2.99792458e8;
170 
171  const double k = wcc/local_c_light*two_pi;
172  const double k2 = k*k;
173  const double cos = std::cos(phi/2.0);
174  const double sin = std::sin(phi/2.0);
175  const double cos2 = cos*cos;
176  const double sin2 = sin*sin;
177 
178  const double sigx2 = sigx*sigx;
179  const double sigmax2=sigx2*(1+z*z/(betx*betx));
180 
181  const double sigs2 = sigs*sigs;
182 
183  constexpr double factorRMSgauss4 = 1./sqrt2/gamma34 * gamma14; // # Factor to take rms sigma as input of the supergaussian
184  constexpr double NormFactorGauss4 = sqrt2to5 * gamma54 * gamma54;
185 
186  const double sinCR = std::sin(phiCR/2.0);
187  const double sinCR2 = sinCR*sinCR;
188 
189  double result = -1.0;
190 
191  if( !RF800 ) {
192  const double norm =2.0/(two_pi*sigs2);
193  const double cosks = std::cos(k*z);
194  const double sinkct = std::sin(k*ct);
195  result = norm*std::exp(-ct*ct/sigs2
196  -z*z*cos2/sigs2
197  -1.0/(4*k2*sigmax2)*(
198  //-4*cosks*cosks * sinkct*sinkct * sinCR2 // comes from integral over x
199  -8*z*k*std::sin(k*z)*std::cos(k*ct) * sin * sinCR
200  +2 * sinCR2
201  -std::cos(2*k*(z-ct)) * sinCR2
202  -std::cos(2*k*(z+ct)) * sinCR2
203  +4*k2*z*z *sin2
204  )
205  - x*x*(cos2/sigmax2 + sin2/sigs2) // contribution from x integrand
206  + x*ct*sin/sigs2 // contribution from x integrand
207  + 2*x*cos*cosks*sinkct*sinCR/k/sigmax2 // contribution from x integrand
208  //+(2*ct/k)*np.cos(k*s)*np.sin(k*ct) *(sin*sinCR)/(sigs2*cos) # small term
209  //+ct**2*(sin2/sigs4)/(cos2/sigmax2) # small term
210  )/(1.0+(z*z)/(betx*betx))/std::sqrt(1.0+(z*z)/(bets*bets));
211 
212  } else {
213 
214  const double norm = 2.0/(NormFactorGauss4*sigs2*factorRMSgauss4);
215  const double sigs4=sigs2*sigs2*factorRMSgauss4*factorRMSgauss4;
216  const double cosks = std::cos(k*z);
217  const double sinct = std::sin(k*ct);
218  result = norm*std::exp(
219  -ct*ct*ct*ct/sigs4
220  -z*z*z*z*cos2*cos2/sigs4
221  -6*ct*ct*z*z*cos2/sigs4
222  -sin2/(4*k2*sigmax2)*(
223  2
224  +4*k2*z*z
225  -std::cos(2*k*(z-ct))
226  -std::cos(2*k*(z+ct))
227  -8*k*s*std::cos(k*ct)*std::sin(k*z)
228  -4 * cosks*cosks * sinct*sinct)
229  )/std::sqrt(1+z*z/(betx*betx))/std::sqrt(1+z*z/(bets*bets));
230  }
231 
232  return result;
233 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define constexpr
tuple result
Definition: mps_fire.py:84
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double HLLHCEvtVtxGenerator::intensity ( double  x,
double  y,
double  z,
double  t 
) const
private

Definition at line 142 of file HLLHCEvtVtxGenerator.cc.

References alphay, betagamma, bets, EnergyCorrector::c, constexpr, funct::cos(), epssn, create_public_lumi_plots::exp, integrandCC(), sigma(), funct::sin(), t, and wcc.

Referenced by newVertex().

145  {
146  //---c in m/s --- remember t is already in meters
147  constexpr double c= 2.99792458e+8; // m/s
148 
149  const double sigmay=sigma(z,epssn,bets,betagamma);
150 
151  const double alphay_mod=alphay*std::cos(wcc*(z-t)/c);
152 
153  const double cay=std::cos(alphay_mod);
154  const double say=std::sin(alphay_mod);
155 
156  const double dy=-(z-t)*say-y*cay;
157 
158  const double xzt_density= integrandCC(x,z,t);
159 
160  const double norm = two_pi*sigmay;
161 
162  return ( std::exp(-dy*dy/(sigmay*sigmay))*
163  xzt_density/norm );
164 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define constexpr
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)
virtual

return a new event vertex

Implements BaseEvtVtxGenerator.

Definition at line 88 of file HLLHCEvtVtxGenerator.cc.

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

88  {
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  if(fVertex == 0)
115  fVertex = new HepMC::FourVector();
116 
117  //---convert to mm
118  x*=1000.0;
119  y*=1000.0;
120  z*=1000.0;
121  t*=1000.0;
122 
123  x+=fMeanX;
124  y+=fMeanY;
125  z+=fMeanZ;
126  t+=fTimeOffset;
127 
128  fVertex->set( x, y, z, t );
129 
130  return fVertex;
131 }
int i
Definition: DBlmapReader.cc:9
HepMC::FourVector * fVertex
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)
private

Copy assignment operator

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

Definition at line 135 of file HLLHCEvtVtxGenerator.cc.

References mathSSE::sqrt().

Referenced by intensity(), and newVertex().

136 {
138 
139  return sigma;
140 }
T sqrt(T t)
Definition: SSEVec.h:18
double sigma(double z, double epsilon, double beta, double betagamma) const

Member Data Documentation

const double HLLHCEvtVtxGenerator::alphax
private

Definition at line 82 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::alphay
private

Definition at line 85 of file HLLHCEvtVtxGenerator.h.

Referenced by intensity().

const double HLLHCEvtVtxGenerator::beta
private

Definition at line 54 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::betagamma
private

Definition at line 55 of file HLLHCEvtVtxGenerator.h.

Referenced by intensity(), and newVertex().

const double HLLHCEvtVtxGenerator::bets
private

Definition at line 70 of file HLLHCEvtVtxGenerator.h.

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

const double HLLHCEvtVtxGenerator::betx
private

Definition at line 67 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC(), and newVertex().

const double HLLHCEvtVtxGenerator::epss
private

Definition at line 94 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::epssn
private

Definition at line 76 of file HLLHCEvtVtxGenerator.h.

Referenced by intensity(), and newVertex().

const double HLLHCEvtVtxGenerator::epsx
private

Definition at line 91 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::epsxn
private

Definition at line 73 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fMeanX
private

Definition at line 49 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fMeanY
private

Definition at line 49 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fMeanZ
private

Definition at line 49 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::fTimeOffset
private

Definition at line 49 of file HLLHCEvtVtxGenerator.h.

Referenced by newVertex().

const double HLLHCEvtVtxGenerator::gamma
private

Definition at line 53 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::momeV
private

Definition at line 52 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::oncc
private

Definition at line 88 of file HLLHCEvtVtxGenerator.h.

const double HLLHCEvtVtxGenerator::phi
private

Definition at line 58 of file HLLHCEvtVtxGenerator.h.

Referenced by Particle.Particle::__str__(), and integrandCC().

const double HLLHCEvtVtxGenerator::phiCR
private

Definition at line 100 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC().

const bool HLLHCEvtVtxGenerator::RF800
private

Definition at line 64 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC().

const double HLLHCEvtVtxGenerator::sigs
private

Definition at line 79 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC(), and newVertex().

const double HLLHCEvtVtxGenerator::sigx
private

Definition at line 97 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC().

const double HLLHCEvtVtxGenerator::wcc
private

Definition at line 61 of file HLLHCEvtVtxGenerator.h.

Referenced by integrandCC(), and intensity().