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 // $Id: FlatEvtVtxGenerator.cc,v 1.5 2009/05/25 12:46:04 fabiocos Exp $
3 
6 
8 
9 #include "CLHEP/Random/RandFlat.h"
10 #include "CLHEP/Units/GlobalSystemOfUnits.h"
11 #include "CLHEP/Units/GlobalPhysicalConstants.h"
12 //#include "CLHEP/Vector/ThreeVector.h"
13 #include "HepMC/SimpleVector.h"
14 
17 {
18 
19  fRandom = new CLHEP::RandFlat(getEngine()) ;
20 
21  fMinX = p.getParameter<double>("MinX")*cm;
22  fMinY = p.getParameter<double>("MinY")*cm;
23  fMinZ = p.getParameter<double>("MinZ")*cm;
24  fMaxX = p.getParameter<double>("MaxX")*cm;
25  fMaxY = p.getParameter<double>("MaxY")*cm;
26  fMaxZ = p.getParameter<double>("MaxZ")*cm;
27  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light;
28 
29  if (fMinX > fMaxX) {
30  throw cms::Exception("Configuration")
31  << "Error in FlatEvtVtxGenerator: "
32  << "MinX is greater than MaxX";
33  }
34  if (fMinY > fMaxY) {
35  throw cms::Exception("Configuration")
36  << "Error in FlatEvtVtxGenerator: "
37  << "MinY is greater than MaxY";
38  }
39  if (fMinZ > fMaxZ) {
40  throw cms::Exception("Configuration")
41  << "Error in FlatEvtVtxGenerator: "
42  << "MinZ is greater than MaxZ";
43  }
44 }
45 
47 {
48  delete fRandom;
49 }
50 
51 
52 //Hep3Vector * FlatEvtVtxGenerator::newVertex() {
53 HepMC::FourVector* FlatEvtVtxGenerator::newVertex() {
54  double aX,aY,aZ;
55  aX = fRandom->fire(fMinX,fMaxX) ;
56  aY = fRandom->fire(fMinY,fMaxY) ;
57  aZ = fRandom->fire(fMinZ,fMaxZ) ;
58 
59  //if (fVertex == 0) fVertex = new CLHEP::Hep3Vector;
60  //fVertex->set(aX,aY,aZ);
61  if ( fVertex == 0 ) fVertex = new HepMC::FourVector() ;
62  fVertex->set(aX,aY,aZ,fTimeOffset);
63 
64  return fVertex;
65 }
66 
68 {
69  fMinX = min;
70 }
71 
73 {
74  fMinY = min;
75 }
76 
78 {
79  fMinZ = min;
80 }
81 
83 {
84  fMaxX = max;
85 }
86 
88 {
89  fMaxY = max;
90 }
91 
93 {
94  fMaxZ = max;
95 }
T getParameter(std::string const &) const
CLHEP::HepRandomEngine & getEngine()
HepMC::FourVector * fVertex
#define min(a, b)
Definition: mlp_lapack.h:161
void minZ(double m=0.0)
set min in Z in cm
virtual HepMC::FourVector * newVertex()
return a new event vertex
const T & max(const T &a, const T &b)
CLHEP::RandFlat * fRandom
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
void maxY(double m=0)
set max in Y in cm