CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  BaseEvtVtxGenerator(p), fRandom(0) {
21 
22  meanX(p.getParameter<double>("BeamMeanX")*cm);
23  meanY(p.getParameter<double>("BeamMeanY")*cm);
24  beamPos(p.getParameter<double>("BeamPosition")*cm);
25  sigmaX(p.getParameter<double>("BeamSigmaX")*cm);
26  sigmaY(p.getParameter<double>("BeamSigmaY")*cm);
27  double fMinEta = p.getParameter<double>("MinEta");
28  double fMaxEta = p.getParameter<double>("MaxEta");
29  double fMinPhi = p.getParameter<double>("MinPhi");
30  double fMaxPhi = p.getParameter<double>("MaxPhi");
31  eta(0.5*(fMinEta+fMaxEta));
32  phi(0.5*(fMinPhi+fMaxPhi));
33  psi(p.getParameter<double>("Psi"));
34  nBinx = p.getParameter<int>("BinX");
35  nBiny = p.getParameter<int>("BinY");
36  ffile = p.getParameter<bool>("UseFile");
37  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light;
38 
39  if (ffile) {
41  std::ifstream is(file.c_str(), std::ios::in);
42  if (is) {
43  double elem,sum=0;
44  while (!is.eof()) {
45  is >> elem;
46  fdistn.push_back(elem);
47  sum += elem;
48  }
49  if (((int)(fdistn.size())) >= nBinx*nBiny) {
50  double last = 0;
51  for (unsigned int i=0; i<fdistn.size(); i++) {
52  fdistn[i] /= sum;
53  fdistn[i] += last;
54  last = fdistn[i];
55  }
56  setType(false);
57  } else {
58  ffile = false;
59  }
60  } else {
61  ffile = false;
62  }
63  }
64  if (!ffile) {
65  setType(p.getParameter<bool>("GaussianProfile"));
66  }
67 
68  edm::LogInfo("VertexGenerator") << "BeamProfileVtxGenerator: with "
69  << "beam along eta = " << fEta
70  << " (Theta = " << fTheta/deg
71  << ") phi = " << fPhi/deg
72  << ") psi = " << fPsi/deg
73  << " centred at (" << fMeanX << ", "
74  << fMeanY << ", " << fMeanZ << ") "
75  << "and spread (" << fSigmaX << ", "
76  << fSigmaY << ") of type Gaussian = "
77  << fType << " use file " << ffile;
78  if (ffile)
79  edm::LogInfo("VertexGenerator") << "With " << nBinx << " bins "
80  << " along X and " << nBiny
81  << " bins along Y";
82 }
83 
85  delete fRandom;
86 }
87 
88 
89 //Hep3Vector * BeamProfileVtxGenerator::newVertex() {
90 HepMC::FourVector* BeamProfileVtxGenerator::newVertex() {
91  double aX, aY;
92  if (ffile) {
93  double r1 = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire();
94  int ixy = 0, ix, iy;
95  for (unsigned int i=0; i<fdistn.size(); i++) {
96  if (r1 > fdistn[i]) ixy = i+1;
97  }
98  if (ixy >= nBinx*nBiny) {
99  ix = nBinx-1; iy = nBiny-1;
100  } else {
101  ix = ixy%nBinx; iy = (ixy-ix)/nBinx;
102  }
103  aX = 0.5*(2*ix-nBinx+2*(dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire())*fSigmaX + fMeanX ;
104  aY = 0.5*(2*iy-nBiny+2*(dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire())*fSigmaY + fMeanY ;
105  } else {
106  if (fType) {
107  aX = fSigmaX*(dynamic_cast<CLHEP::RandGaussQ*>(fRandom))->fire() +fMeanX;
108  aY = fSigmaY*(dynamic_cast<CLHEP::RandGaussQ*>(fRandom))->fire() +fMeanY;
109  } else {
110  aX = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire(-0.5*fSigmaX,0.5*fSigmaX) + fMeanX ;
111  aY = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire(-0.5*fSigmaY,0.5*fSigmaY) + fMeanY;
112  }
113  }
114 
115  double xp, yp, zp ;
116  if( 2.*M_PI < fabs(fPsi) )
117  {
118  xp = -aX*cos(fTheta)*cos(fPhi) +aY*sin(fPhi) +fMeanZ*sin(fTheta)*cos(fPhi);
119  yp = -aX*cos(fTheta)*sin(fPhi) -aY*cos(fPhi) +fMeanZ*sin(fTheta)*sin(fPhi);
120  zp = aX*sin(fTheta) +fMeanZ*cos(fTheta);
121  }
122  else // for endcap testbeam, use 3rd angle psi
123  {
124  const HepGeom::Vector3D<double> av ( aX, aY, fMeanZ ) ;
125 
126 /*
127  static const double kRadToDeg ( 180./M_PI ) ;
128  std::cout<<"theta = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
129  <<fTheta*kRadToDeg<<", phi="<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
130  << fPhi*kRadToDeg <<", PSI="<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
131  <<fPsi*kRadToDeg<<std::endl ;
132 */
133  const HepGeom::RotateZ3D R1 ( fPhi - M_PI ) ;
134  const HepGeom::Point3D<double> xUnit ( 0,1,0 ) ;
135  const HepGeom::Point3D<double> zUnit ( 0,0,1 ) ;
136  const HepGeom::Transform3D RXRZ ( HepGeom::Rotate3D( -fTheta, R1*xUnit )*R1 ) ;
137  const HepGeom::Transform3D TRF ( HepGeom::Rotate3D(fPsi,RXRZ*zUnit)*RXRZ ) ;
138 /*
139  std::cout<<"\n\n$$$$$$$$$$$Transform="
140  <<" thetaZ = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
141  <<TRF.getRotation().thetaZ()*kRadToDeg
142  <<", phiZ = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
143  <<TRF.getRotation().phiZ()*kRadToDeg
144  <<", thetaY = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
145  <<TRF.getRotation().thetaY()*kRadToDeg
146  <<", phiY = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
147  <<TRF.getRotation().phiY()*kRadToDeg
148  <<", thetaX = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
149  <<TRF.getRotation().thetaX()*kRadToDeg
150  <<", phiX = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
151  <<TRF.getRotation().phiX()*kRadToDeg
152  <<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
153  <<TRF.getTranslation()<<std::endl ;
154 */
155  const HepGeom::Vector3D<double> pv ( TRF*av ) ;
156 
157  xp = pv.x() ;
158  yp = pv.y() ;
159  zp = pv.z() ;
160  }
161  //if (fVertex == 0) fVertex = new CLHEP::Hep3Vector;
162  //fVertex->set(xp, yp, zp);
163  if (fVertex == 0 ) fVertex = new HepMC::FourVector() ;
164  fVertex->set(xp, yp, zp, fTimeOffset );
165 
166  LogDebug("VertexGenerator") << "BeamProfileVtxGenerator: Vertex created "
167  << "at (" << xp << ", " << yp << ", "
168  << zp << ", " << fTimeOffset << ")";
169  return fVertex;
170 }
171 
173 
174  if (s>=0) {
175  fSigmaX = s;
176  } else {
177  edm::LogWarning("VertexGenerator") << "Warning BeamProfileVtxGenerator:"
178  << " Illegal resolution in X " << s
179  << "- set to default value 0 cm";
180  fSigmaX = 0;
181  }
182 }
183 
185 
186  if (s>=0) {
187  fSigmaY = s;
188  } else {
189  edm::LogWarning("VertexGenerator") << "Warning BeamProfileVtxGenerator:"
190  << " Illegal resolution in Y " << s
191  << "- set to default value 0 cm";
192  fSigmaY = 0;
193  }
194 }
195 
197  fEta = s;
198  fTheta = 2.0*atan(exp(-fEta));
199 }
200 
202 
203  fType = s;
204  delete fRandom;
205 
206  if (fType == true)
207  fRandom = new CLHEP::RandGaussQ(getEngine());
208  else
209  fRandom = new CLHEP::RandFlat(getEngine());
210 }
#define LogDebug(id)
void meanY(double m=0)
set mean in Y in cm
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual HepMC::FourVector * newVertex()
return a new event vertex
CLHEP::HepRandomEngine & getEngine()
HepMC::FourVector * fVertex
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
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
Definition: BFit3D.cc:3
BeamProfileVtxGenerator(const edm::ParameterSet &p)
void sigmaY(double s=1.0)
set resolution in Y in cm