00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00031
00032
00043
00044
00045
00046 #include <CLHEP/Random/RandGauss.h>
00047
00048 #include "G4ios.hh"
00049 #include "G4Event.hh"
00050 #include "G4Gamma.hh"
00051 #include "G4ParticleGun.hh"
00052 #include "G4Run.hh"
00053 #include "G4ThreeVector.hh"
00054
00055 #include "generators/MJGeneratorTUNLFELMessenger.hh"
00056 #include "io/MJLogger.hh"
00057
00058
00059
00060 #include "generators/MJGeneratorTUNLFEL.hh"
00061
00062
00063
00064 MJGeneratorTUNLFEL::MJGeneratorTUNLFEL()
00065 {
00066 fGeneratorName = "TUNLFEL";
00067 fG4Messenger = new MJGeneratorTUNLFELMessenger(this);
00068 fParticleGun = new G4ParticleGun(1);
00069 SetDefaults();
00070 }
00071
00072
00073
00074 MJGeneratorTUNLFEL::MJGeneratorTUNLFEL(const MJGeneratorTUNLFEL & other)
00075 {;}
00076
00077
00078
00079 MJGeneratorTUNLFEL::~MJGeneratorTUNLFEL()
00080 {
00081 delete fG4Messenger;
00082 delete fParticleGun;
00083 }
00084
00085
00086
00087 void MJGeneratorTUNLFEL::BeginOfRunAction(G4Run const *run)
00088 {
00089 if(fMajorSigma < fMinorSigma){
00090 G4double temp = fMajorSigma;
00091 fMajorSigma = fMinorSigma;
00092 fMinorSigma = temp;
00093 MJLog(warning) << "Minor sigma > Major Sigma, values swapped." << endlog;
00094 }
00095 if(MJLogger::GetSeverity() <= MJLogger::routine)
00096 Dump();
00097 }
00098
00099
00100
00101 void MJGeneratorTUNLFEL::Dump()
00102 {
00103 G4cout << " TUNL FEL Parameters :" << G4endl;
00104 G4cout << " -----------------------------------------" << G4endl;
00105 G4cout << " Rho (degrees): " << fRho/deg << G4endl;
00106 G4cout << " Major sigma (cm): " << fMajorSigma/cm << G4endl;
00107 G4cout << " Minor sigma (cm): " << fMinorSigma/cm << G4endl;
00108 G4cout << " Mean Energy (MeV): " << fMeanEnergy/MeV << G4endl;
00109 G4cout << " Energy sigma (keV): " << fEnergySigma/keV << G4endl;
00110 G4cout << " Origin (cm) : " << fOrigin/cm << G4endl;
00111 G4cout << "-------------------------------------------" << G4endl << G4endl;
00112 }
00113
00114
00115
00116 void MJGeneratorTUNLFEL::EndOfRunAction(G4Run const *run)
00117 {;}
00118
00119
00120
00121 void MJGeneratorTUNLFEL::GeneratePrimaryVertex(G4Event *event)
00122 {
00123 fParticleGun->SetParticleDefinition(G4Gamma::GammaDefinition());
00124 fParticleGun->SetParticleMomentumDirection(fDirection);
00125 fCurrentPosition.set(RandGauss::shoot(0.0, fMajorSigma),
00126 RandGauss::shoot(0.0, fMinorSigma), 0.0);
00127 fCurrentPosition.rotateZ(fRho);
00128 fCurrentPosition += fOrigin;
00129 fParticleGun->SetParticlePosition(fCurrentPosition);
00130 fCurrentEnergy = RandGauss::shoot(fMeanEnergy, fEnergySigma);
00131 fParticleGun->SetParticleEnergy(fCurrentEnergy);
00132 fParticleGun->GeneratePrimaryVertex(event);
00133 }
00134
00135
00136
00137 void MJGeneratorTUNLFEL::SetDefaults()
00138 {
00139 fDirection.set(0.0, 0.0, -1.0);
00140 fMajorSigma = fMinorSigma = 0.6 * cm;
00141 fOrigin.set(0.0, 0.0, 10.0 * m);
00142 fRho = 0.0;
00143 fMeanEnergy = 2.2 * MeV;
00144 fEnergySigma = 22.0 * keV;
00145 }
00146
00147
00148