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