CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FlatEvtVtxGenerator.cc
Go to the documentation of this file.
1 
2 
5 
7 
8 #include "CLHEP/Random/RandFlat.h"
9 #include "CLHEP/Units/GlobalSystemOfUnits.h"
10 #include "CLHEP/Units/GlobalPhysicalConstants.h"
11 //#include "CLHEP/Vector/ThreeVector.h"
12 #include "HepMC/SimpleVector.h"
13 
16 {
17  fMinX = p.getParameter<double>("MinX")*cm;
18  fMinY = p.getParameter<double>("MinY")*cm;
19  fMinZ = p.getParameter<double>("MinZ")*cm;
20  fMaxX = p.getParameter<double>("MaxX")*cm;
21  fMaxY = p.getParameter<double>("MaxY")*cm;
22  fMaxZ = p.getParameter<double>("MaxZ")*cm;
23  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light;
24 
25  if (fMinX > fMaxX) {
26  throw cms::Exception("Configuration")
27  << "Error in FlatEvtVtxGenerator: "
28  << "MinX is greater than MaxX";
29  }
30  if (fMinY > fMaxY) {
31  throw cms::Exception("Configuration")
32  << "Error in FlatEvtVtxGenerator: "
33  << "MinY is greater than MaxY";
34  }
35  if (fMinZ > fMaxZ) {
36  throw cms::Exception("Configuration")
37  << "Error in FlatEvtVtxGenerator: "
38  << "MinZ is greater than MaxZ";
39  }
40 }
41 
43 {
44 }
45 
46 //Hep3Vector * FlatEvtVtxGenerator::newVertex() {
47 HepMC::FourVector* FlatEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) {
48  double aX,aY,aZ;
49  aX = CLHEP::RandFlat::shoot(engine, fMinX, fMaxX);
50  aY = CLHEP::RandFlat::shoot(engine, fMinY, fMaxY);
51  aZ = CLHEP::RandFlat::shoot(engine, fMinZ, fMaxZ);
52 
53  //if (fVertex == 0) fVertex = new CLHEP::Hep3Vector;
54  //fVertex->set(aX,aY,aZ);
55  if ( fVertex == 0 ) fVertex = new HepMC::FourVector() ;
56  fVertex->set(aX,aY,aZ,fTimeOffset);
57 
58  return fVertex;
59 }
60 
62 {
63  fMinX = min;
64 }
65 
67 {
68  fMinY = min;
69 }
70 
72 {
73  fMinZ = min;
74 }
75 
77 {
78  fMaxX = max;
79 }
80 
82 {
83  fMaxY = max;
84 }
85 
87 {
88  fMaxZ = max;
89 }
T getParameter(std::string const &) const
HepMC::FourVector * fVertex
void minZ(double m=0.0)
set min in Z in cm
const T & max(const T &a, const T &b)
void maxZ(double m=0)
set max in Z in cm
FlatEvtVtxGenerator(const edm::ParameterSet &p)
void minX(double m=0.0)
set min in X in cm
void minY(double m=0.0)
set min in Y in cm
void maxX(double m=0)
set max in X in cm
virtual HepMC::FourVector * newVertex(CLHEP::HepRandomEngine *)
return a new event vertex
void maxY(double m=0)
set max in Y in cm