CMS 3D CMS Logo

BeamProfileVtxGenerator.cc
Go to the documentation of this file.
1 
2 
5 
8 
9 #include "CLHEP/Geometry/Transform3D.h"
10 #include "CLHEP/Random/RandFlat.h"
11 #include "CLHEP/Random/RandGaussQ.h"
12 #include "CLHEP/Units/GlobalSystemOfUnits.h"
13 #include "CLHEP/Units/GlobalPhysicalConstants.h"
14 #include "HepMC/SimpleVector.h"
15 
16 #include <fstream>
17 #include <string>
18 
20  meanX(p.getParameter<double>("BeamMeanX") * cm);
21  meanY(p.getParameter<double>("BeamMeanY") * cm);
22  beamPos(p.getParameter<double>("BeamPosition") * cm);
23  sigmaX(p.getParameter<double>("BeamSigmaX") * cm);
24  sigmaY(p.getParameter<double>("BeamSigmaY") * cm);
25  double fMinEta = p.getParameter<double>("MinEta");
26  double fMaxEta = p.getParameter<double>("MaxEta");
27  double fMinPhi = p.getParameter<double>("MinPhi");
28  double fMaxPhi = p.getParameter<double>("MaxPhi");
29  eta(0.5 * (fMinEta + fMaxEta));
30  phi(0.5 * (fMinPhi + fMaxPhi));
31  psi(p.getParameter<double>("Psi"));
32  nBinx = p.getParameter<int>("BinX");
33  nBiny = p.getParameter<int>("BinY");
34  ffile = p.getParameter<bool>("UseFile");
35  fTimeOffset = p.getParameter<double>("TimeOffset") * ns * c_light;
36 
37  if (ffile) {
38  std::string file = p.getParameter<std::string>("File");
39  std::ifstream is(file.c_str(), std::ios::in);
40  if (is) {
41  double elem, sum = 0;
42  while (!is.eof()) {
43  is >> elem;
44  fdistn.push_back(elem);
45  sum += elem;
46  }
47  if (((int)(fdistn.size())) >= nBinx * nBiny) {
48  double last = 0;
49  for (unsigned int i = 0; i < fdistn.size(); i++) {
50  fdistn[i] /= sum;
51  fdistn[i] += last;
52  last = fdistn[i];
53  }
54  setType(false);
55  } else {
56  ffile = false;
57  }
58  } else {
59  ffile = false;
60  }
61  }
62  if (!ffile) {
63  setType(p.getParameter<bool>("GaussianProfile"));
64  }
65 
66  edm::LogInfo("VertexGenerator") << "BeamProfileVtxGenerator: with "
67  << "beam along eta = " << fEta << " (Theta = " << fTheta / deg
68  << ") phi = " << fPhi / deg << ") psi = " << fPsi / deg << " centred at (" << fMeanX
69  << ", " << fMeanY << ", " << fMeanZ << ") "
70  << "and spread (" << fSigmaX << ", " << fSigmaY << ") of type Gaussian = " << fType
71  << " use file " << ffile;
72  if (ffile)
73  edm::LogInfo("VertexGenerator") << "With " << nBinx << " bins "
74  << " along X and " << nBiny << " bins along Y";
75 }
76 
78 
79 //Hep3Vector * BeamProfileVtxGenerator::newVertex() {
80 HepMC::FourVector BeamProfileVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) const {
81  double aX, aY;
82  if (ffile) {
83  double r1 = engine->flat();
84  int ixy = 0, ix, iy;
85  for (unsigned int i = 0; i < fdistn.size(); i++) {
86  if (r1 > fdistn[i])
87  ixy = i + 1;
88  }
89  if (ixy >= nBinx * nBiny) {
90  ix = nBinx - 1;
91  iy = nBiny - 1;
92  } else {
93  ix = ixy % nBinx;
94  iy = (ixy - ix) / nBinx;
95  }
96  aX = 0.5 * (2 * ix - nBinx + 2 * engine->flat()) * fSigmaX + fMeanX;
97  aY = 0.5 * (2 * iy - nBiny + 2 * engine->flat()) * fSigmaY + fMeanY;
98  } else {
99  if (fType) {
100  aX = fSigmaX * CLHEP::RandGaussQ::shoot(engine) + fMeanX;
101  aY = fSigmaY * CLHEP::RandGaussQ::shoot(engine) + fMeanY;
102  } else {
103  aX = CLHEP::RandFlat::shoot(engine, -0.5 * fSigmaX, 0.5 * fSigmaX) + fMeanX;
104  aY = CLHEP::RandFlat::shoot(engine, -0.5 * fSigmaY, 0.5 * fSigmaY) + fMeanY;
105  }
106  }
107 
108  double xp, yp, zp;
109  if (2. * M_PI < fabs(fPsi)) {
110  xp = -aX * cos(fTheta) * cos(fPhi) + aY * sin(fPhi) + fMeanZ * sin(fTheta) * cos(fPhi);
111  yp = -aX * cos(fTheta) * sin(fPhi) - aY * cos(fPhi) + fMeanZ * sin(fTheta) * sin(fPhi);
112  zp = aX * sin(fTheta) + fMeanZ * cos(fTheta);
113  } else // for endcap testbeam, use 3rd angle psi
114  {
115  const HepGeom::Vector3D<double> av(aX, aY, fMeanZ);
116 
117  /*
118  static const double kRadToDeg ( 180./M_PI ) ;
119  std::cout<<"theta = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
120  <<fTheta*kRadToDeg<<", phi="<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
121  << fPhi*kRadToDeg <<", PSI="<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
122  <<fPsi*kRadToDeg<<std::endl ;
123 */
124  const HepGeom::RotateZ3D R1(fPhi - M_PI);
125  const HepGeom::Point3D<double> xUnit(0, 1, 0);
126  const HepGeom::Point3D<double> zUnit(0, 0, 1);
127  const HepGeom::Transform3D RXRZ(HepGeom::Rotate3D(-fTheta, R1 * xUnit) * R1);
128  const HepGeom::Transform3D TRF(HepGeom::Rotate3D(fPsi, RXRZ * zUnit) * RXRZ);
129  /*
130  std::cout<<"\n\n$$$$$$$$$$$Transform="
131  <<" thetaZ = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
132  <<TRF.getRotation().thetaZ()*kRadToDeg
133  <<", phiZ = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
134  <<TRF.getRotation().phiZ()*kRadToDeg
135  <<", thetaY = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
136  <<TRF.getRotation().thetaY()*kRadToDeg
137  <<", phiY = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
138  <<TRF.getRotation().phiY()*kRadToDeg
139  <<", thetaX = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
140  <<TRF.getRotation().thetaX()*kRadToDeg
141  <<", phiX = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
142  <<TRF.getRotation().phiX()*kRadToDeg
143  <<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
144  <<TRF.getTranslation()<<std::endl ;
145 */
146  const HepGeom::Vector3D<double> pv(TRF * av);
147 
148  xp = pv.x();
149  yp = pv.y();
150  zp = pv.z();
151  }
152 
153  LogDebug("VertexGenerator") << "BeamProfileVtxGenerator: Vertex created "
154  << "at (" << xp << ", " << yp << ", " << zp << ", " << fTimeOffset << ")";
155  return HepMC::FourVector(xp, yp, zp, fTimeOffset);
156  ;
157 }
158 
160  if (s >= 0) {
161  fSigmaX = s;
162  } else {
163  edm::LogWarning("VertexGenerator") << "Warning BeamProfileVtxGenerator:"
164  << " Illegal resolution in X " << s << "- set to default value 0 cm";
165  fSigmaX = 0;
166  }
167 }
168 
170  if (s >= 0) {
171  fSigmaY = s;
172  } else {
173  edm::LogWarning("VertexGenerator") << "Warning BeamProfileVtxGenerator:"
174  << " Illegal resolution in Y " << s << "- set to default value 0 cm";
175  fSigmaY = 0;
176  }
177 }
178 
180  fEta = s;
181  fTheta = 2.0 * atan(exp(-fEta));
182 }
183 
void meanY(double m=0)
set mean in Y in cm
void setType(bool m=true)
set type
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void meanX(double m=0)
set mean in X in cm
void psi(double m=999)
set psi in radian
std::vector< double > fdistn
void sigmaX(double s=1.0)
set resolution in X in cm
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
def pv(vc)
Definition: MetAnalyzer.py:7
void beamPos(double m=0)
set mean in Z in cm
void phi(double m=0)
set phi in radian
void eta(double m=0)
set eta
#define M_PI
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const override
return a new event vertex
Log< level::Info, false > LogInfo
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
BeamProfileVtxGenerator(const edm::ParameterSet &p)
void sigmaY(double s=1.0)
set resolution in Y in cm
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)
Log< level::Warning, false > LogWarning
#define LogDebug(id)