#include <IOMC/EventVertexGenerators/interface/BeamProfileVtxGenerator.h>
Public Member Functions | |
void | beamPos (double m=0) |
set mean in Z in cm | |
BeamProfileVtxGenerator (const edm::ParameterSet &p) | |
void | eta (double m=0) |
set eta | |
virtual TMatrixD * | GetInvLorentzBoost () |
void | meanX (double m=0) |
set mean in X in cm | |
void | meanY (double m=0) |
set mean in Y in cm | |
virtual HepMC::FourVector * | newVertex () |
return a new event vertex | |
void | phi (double m=0) |
set phi in radian | |
void | setType (bool m=true) |
set type | |
void | sigmaX (double s=1.0) |
set resolution in X in cm | |
void | sigmaY (double s=1.0) |
set resolution in Y in cm | |
virtual | ~BeamProfileVtxGenerator () |
Private Member Functions | |
BeamProfileVtxGenerator (const BeamProfileVtxGenerator &p) | |
Copy constructor. | |
BeamProfileVtxGenerator & | operator= (const BeamProfileVtxGenerator &rhs) |
Copy assignment operator. | |
Private Attributes | |
std::vector< double > | fdistn |
double | fEta |
bool | ffile |
double | fMeanX |
double | fMeanY |
double | fMeanZ |
double | fPhi |
CLHEP::HepRandom * | fRandom |
double | fSigmaX |
double | fSigmaY |
double | fTheta |
double | fTimeOffset |
bool | fType |
int | nBinx |
int | nBiny |
Definition at line 21 of file BeamProfileVtxGenerator.h.
BeamProfileVtxGenerator::BeamProfileVtxGenerator | ( | const edm::ParameterSet & | p | ) |
Definition at line 20 of file BeamProfileVtxGenerator.cc.
References beamPos(), eta(), fdistn, fEta, ffile, file, fMeanX, fMeanY, fMeanZ, fPhi, fSigmaX, fSigmaY, fTheta, fTimeOffset, fType, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), i, in, prof2calltree::last, meanX(), meanY(), nBinx, nBiny, phi(), setType(), sigmaX(), sigmaY(), and sum().
00020 : 00021 BaseEvtVtxGenerator(p), fRandom(0) { 00022 00023 meanX(p.getUntrackedParameter<double>("BeamMeanX",0.0)*cm); 00024 meanY(p.getUntrackedParameter<double>("BeamMeanY",0.0)*cm); 00025 beamPos(p.getUntrackedParameter<double>("BeamPosition",0.0)*cm); 00026 sigmaX(p.getUntrackedParameter<double>("BeamSigmaX",0.0)*cm); 00027 sigmaY(p.getUntrackedParameter<double>("BeamSigmaY",0.0)*cm); 00028 double fMinEta = p.getUntrackedParameter<double>("MinEta",-5.5); 00029 double fMaxEta = p.getUntrackedParameter<double>("MaxEta",5.5); 00030 double fMinPhi = p.getUntrackedParameter<double>("MinPhi",-3.14159265358979323846); 00031 double fMaxPhi = p.getUntrackedParameter<double>("MaxPhi", 3.14159265358979323846); 00032 eta(0.5*(fMinEta+fMaxEta)); 00033 phi(0.5*(fMinPhi+fMaxPhi)); 00034 nBinx = p.getUntrackedParameter<int>("BinX",50); 00035 nBiny = p.getUntrackedParameter<int>("BinY",50); 00036 ffile = p.getUntrackedParameter<bool>("UseFile",false); 00037 fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light; 00038 00039 if (ffile) { 00040 std::string file = p.getUntrackedParameter<std::string>("File","beam.profile"); 00041 ifstream is(file.c_str(), std::ios::in); 00042 if (is) { 00043 double elem,sum=0; 00044 while (!is.eof()) { 00045 is >> elem; 00046 fdistn.push_back(elem); 00047 sum += elem; 00048 } 00049 if (((int)(fdistn.size())) >= nBinx*nBiny) { 00050 double last = 0; 00051 for (unsigned int i=0; i<fdistn.size(); i++) { 00052 fdistn[i] /= sum; 00053 fdistn[i] += last; 00054 last = fdistn[i]; 00055 } 00056 setType(false); 00057 } else { 00058 ffile = false; 00059 } 00060 } else { 00061 ffile = false; 00062 } 00063 } 00064 if (!ffile) { 00065 setType(p.getUntrackedParameter<bool>("GaussianProfile",true)); 00066 } 00067 00068 edm::LogInfo("BeamProfileVtxGenerator") << "BeamProfileVtxGenerator: with " 00069 << "beam along eta = " << fEta 00070 << " (Theta = " << fTheta/deg 00071 << ") phi = " << fPhi/deg 00072 << " centred at (" << fMeanX << ", " 00073 << fMeanY << ", " << fMeanZ << ") " 00074 << "and spread (" << fSigmaX << ", " 00075 << fSigmaY << ") of type Gaussian = " 00076 << fType << " use file " << ffile; 00077 if (ffile) 00078 edm::LogInfo("BeamProfileVtxGenerator") << "With " << nBinx << " bins " 00079 << " along X and " << nBiny 00080 << " bins along Y"; 00081 }
BeamProfileVtxGenerator::~BeamProfileVtxGenerator | ( | ) | [virtual] |
Definition at line 83 of file BeamProfileVtxGenerator.cc.
References fRandom.
00083 { 00084 delete fRandom; 00085 }
BeamProfileVtxGenerator::BeamProfileVtxGenerator | ( | const BeamProfileVtxGenerator & | p | ) | [private] |
Copy constructor.
void BeamProfileVtxGenerator::beamPos | ( | double | m = 0 |
) | [inline] |
set mean in Z in cm
Definition at line 46 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator().
void BeamProfileVtxGenerator::eta | ( | double | m = 0 |
) |
set eta
Definition at line 151 of file BeamProfileVtxGenerator.cc.
References funct::exp(), fEta, and fTheta.
Referenced by BeamProfileVtxGenerator().
virtual TMatrixD* BeamProfileVtxGenerator::GetInvLorentzBoost | ( | ) | [inline, virtual] |
void BeamProfileVtxGenerator::meanX | ( | double | m = 0 |
) | [inline] |
set mean in X in cm
Definition at line 42 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator().
void BeamProfileVtxGenerator::meanY | ( | double | m = 0 |
) | [inline] |
set mean in Y in cm
Definition at line 44 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator().
HepMC::FourVector * BeamProfileVtxGenerator::newVertex | ( | ) | [virtual] |
return a new event vertex
Implements BaseEvtVtxGenerator.
Definition at line 89 of file BeamProfileVtxGenerator.cc.
References funct::cos(), fdistn, ffile, fMeanX, fMeanY, fMeanZ, fPhi, fRandom, fSigmaX, fSigmaY, fTheta, fTimeOffset, fType, BaseEvtVtxGenerator::fVertex, i, nBinx, nBiny, r1, and funct::sin().
00089 { 00090 double aX, aY; 00091 if (ffile) { 00092 double r1 = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire(); 00093 int ixy = 0, ix, iy; 00094 for (unsigned int i=0; i<fdistn.size(); i++) { 00095 if (r1 > fdistn[i]) ixy = i+1; 00096 } 00097 if (ixy >= nBinx*nBiny) { 00098 ix = nBinx-1; iy = nBiny-1; 00099 } else { 00100 ix = ixy%nBinx; iy = (ixy-ix)/nBinx; 00101 } 00102 aX = 0.5*(2*ix-nBinx+2*(dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire())*fSigmaX + fMeanX ; 00103 aY = 0.5*(2*iy-nBiny+2*(dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire())*fSigmaY + fMeanY ; 00104 } else { 00105 if (fType) { 00106 aX = fSigmaX*(dynamic_cast<CLHEP::RandGaussQ*>(fRandom))->fire() +fMeanX; 00107 aY = fSigmaY*(dynamic_cast<CLHEP::RandGaussQ*>(fRandom))->fire() +fMeanY; 00108 } else { 00109 aX = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire(-0.5*fSigmaX,0.5*fSigmaX) + fMeanX ; 00110 aY = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire(-0.5*fSigmaY,0.5*fSigmaY) + fMeanY; 00111 } 00112 } 00113 double xp = -aX*cos(fTheta)*cos(fPhi) +aY*sin(fPhi) +fMeanZ*sin(fTheta)*cos(fPhi); 00114 double yp = -aX*cos(fTheta)*sin(fPhi) -aY*cos(fPhi) +fMeanZ*sin(fTheta)*sin(fPhi); 00115 double zp = aX*sin(fTheta) +fMeanZ*cos(fTheta); 00116 00117 //if (fVertex == 0) fVertex = new CLHEP::Hep3Vector; 00118 //fVertex->set(xp, yp, zp); 00119 if (fVertex == 0 ) fVertex = new HepMC::FourVector() ; 00120 fVertex->set(xp, yp, zp, fTimeOffset ); 00121 00122 // LogDebug("BeamProfileVtxGenerator") << "BeamProfileVtxGenerator: Vertex created " 00123 // << "at " << *fVertex; 00124 return fVertex; 00125 }
BeamProfileVtxGenerator& BeamProfileVtxGenerator::operator= | ( | const BeamProfileVtxGenerator & | rhs | ) | [private] |
Copy assignment operator.
void BeamProfileVtxGenerator::phi | ( | double | m = 0 |
) | [inline] |
set phi in radian
Definition at line 51 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator().
set type
Definition at line 156 of file BeamProfileVtxGenerator.cc.
References fRandom, fType, and BaseEvtVtxGenerator::getEngine().
Referenced by BeamProfileVtxGenerator().
00156 { 00157 00158 fType = s; 00159 delete fRandom; 00160 00161 if (fType == true) 00162 fRandom = new CLHEP::RandGaussQ(getEngine()); 00163 else 00164 fRandom = new CLHEP::RandFlat(getEngine()); 00165 }
void BeamProfileVtxGenerator::sigmaX | ( | double | s = 1.0 |
) |
set resolution in X in cm
Definition at line 127 of file BeamProfileVtxGenerator.cc.
References fSigmaX.
Referenced by BeamProfileVtxGenerator().
00127 { 00128 00129 if (s>=0) { 00130 fSigmaX = s; 00131 } else { 00132 edm::LogWarning("BeamProfileVtxGenerator") << "Warning BeamProfileVtxGenerator:" 00133 << " Illegal resolution in X " << s 00134 << "- set to default value 0 cm"; 00135 fSigmaX = 0; 00136 } 00137 }
void BeamProfileVtxGenerator::sigmaY | ( | double | s = 1.0 |
) |
set resolution in Y in cm
Definition at line 139 of file BeamProfileVtxGenerator.cc.
References fSigmaY.
Referenced by BeamProfileVtxGenerator().
00139 { 00140 00141 if (s>=0) { 00142 fSigmaY = s; 00143 } else { 00144 edm::LogWarning("BeamProfileVtxGenerator") << "Warning BeamProfileVtxGenerator:" 00145 << " Illegal resolution in Y " << s 00146 << "- set to default value 0 cm"; 00147 fSigmaY = 0; 00148 } 00149 }
std::vector<double> BeamProfileVtxGenerator::fdistn [private] |
Definition at line 66 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), and newVertex().
double BeamProfileVtxGenerator::fEta [private] |
Definition at line 63 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), and eta().
bool BeamProfileVtxGenerator::ffile [private] |
Definition at line 64 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), and newVertex().
double BeamProfileVtxGenerator::fMeanX [private] |
Definition at line 62 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), meanX(), and newVertex().
double BeamProfileVtxGenerator::fMeanY [private] |
Definition at line 62 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), meanY(), and newVertex().
double BeamProfileVtxGenerator::fMeanZ [private] |
Definition at line 62 of file BeamProfileVtxGenerator.h.
Referenced by beamPos(), BeamProfileVtxGenerator(), and newVertex().
double BeamProfileVtxGenerator::fPhi [private] |
Definition at line 63 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), newVertex(), and phi().
CLHEP::HepRandom* BeamProfileVtxGenerator::fRandom [private] |
Definition at line 67 of file BeamProfileVtxGenerator.h.
Referenced by newVertex(), setType(), and ~BeamProfileVtxGenerator().
double BeamProfileVtxGenerator::fSigmaX [private] |
Definition at line 61 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), newVertex(), and sigmaX().
double BeamProfileVtxGenerator::fSigmaY [private] |
Definition at line 61 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), newVertex(), and sigmaY().
double BeamProfileVtxGenerator::fTheta [private] |
Definition at line 63 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), eta(), and newVertex().
double BeamProfileVtxGenerator::fTimeOffset [private] |
Reimplemented from BaseEvtVtxGenerator.
Definition at line 68 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), and newVertex().
bool BeamProfileVtxGenerator::fType [private] |
Definition at line 64 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), newVertex(), and setType().
int BeamProfileVtxGenerator::nBinx [private] |
Definition at line 65 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), and newVertex().
int BeamProfileVtxGenerator::nBiny [private] |
Definition at line 65 of file BeamProfileVtxGenerator.h.
Referenced by BeamProfileVtxGenerator(), and newVertex().