Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4LowEnergyBremsstrahlung Class Reference

#include <G4LowEnergyBremsstrahlung.hh>

Inheritance diagram for G4LowEnergyBremsstrahlung:
G4eLowEnergyLoss G4RDVeLowEnergyLoss G4VContinuousDiscreteProcess G4VProcess

Public Member Functions

 G4LowEnergyBremsstrahlung (const G4String &processName="LowEnBrem")
 
 ~G4LowEnergyBremsstrahlung ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &particleType)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step)
 
void SetCutForLowEnSecPhotons (G4double cut)
 
void SetAngularGenerator (G4RDVBremAngularDistribution *distribution)
 
void SetAngularGenerator (const G4String &name)
 
void PrintInfoDefinition ()
 
- Public Member Functions inherited from G4eLowEnergyLoss
 G4eLowEnergyLoss (const G4String &)
 
 ~G4eLowEnergyLoss ()
 
void BuildDEDXTable (const G4ParticleDefinition &aParticleType)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &Step)
 
void ActivateFluorescence (G4bool val)
 
G4bool Fluorescence () const
 
- Public Member Functions inherited from G4RDVeLowEnergyLoss
 G4RDVeLowEnergyLoss (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4RDVeLowEnergyLoss (G4RDVeLowEnergyLoss &)
 
virtual ~G4RDVeLowEnergyLoss ()
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4eLowEnergyLoss
virtual std::vector
< G4DynamicParticle * > * 
DeexciteAtom (const G4MaterialCutsCouple *, G4double, G4double)
 
- Protected Member Functions inherited from G4RDVeLowEnergyLoss
G4double GetLossWithFluct (const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple, G4double MeanLoss, G4double step)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4eLowEnergyLoss
static void SetNbOfProcesses (G4int nb)
 
static void PlusNbOfProcesses ()
 
static void MinusNbOfProcesses ()
 
static G4int GetNbOfProcesses ()
 
static void SetLowerBoundEloss (G4double val)
 
static void SetUpperBoundEloss (G4double val)
 
static void SetNbinEloss (G4int nb)
 
static G4double GetLowerBoundEloss ()
 
static G4double GetUpperBoundEloss ()
 
static G4int GetNbinEloss ()
 
- Static Public Member Functions inherited from G4RDVeLowEnergyLoss
static void SetRndmStep (G4bool value)
 
static void SetEnlossFluc (G4bool value)
 
static void SetStepFunction (G4double c1, G4double c2)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Static Protected Member Functions inherited from G4RDVeLowEnergyLoss
static G4PhysicsTableBuildRangeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildLabTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildProperTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffATable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffBTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffCTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildInverseRangeTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
- Protected Attributes inherited from G4eLowEnergyLoss
G4PhysicsTabletheLossTable
 
G4double MinKineticEnergy
 
G4double Charge
 
G4double lastCharge
 
- Protected Attributes inherited from G4RDVeLowEnergyLoss
const G4MateriallastMaterial
 
G4int imat
 
G4double f1Fluct
 
G4double f2Fluct
 
G4double e1Fluct
 
G4double e2Fluct
 
G4double rateFluct
 
G4double ipotFluct
 
G4double e1LogFluct
 
G4double e2LogFluct
 
G4double ipotLogFluct
 
const G4int nmaxCont1
 
const G4int nmaxCont2
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 
- Static Protected Attributes inherited from G4eLowEnergyLoss
static G4PhysicsTabletheDEDXElectronTable = 0
 
static G4PhysicsTabletheDEDXPositronTable = 0
 
static G4PhysicsTabletheRangeElectronTable = 0
 
static G4PhysicsTabletheRangePositronTable = 0
 
static G4PhysicsTabletheInverseRangeElectronTable = 0
 
static G4PhysicsTabletheInverseRangePositronTable = 0
 
static G4PhysicsTabletheLabTimeElectronTable = 0
 
static G4PhysicsTabletheLabTimePositronTable = 0
 
static G4PhysicsTabletheProperTimeElectronTable = 0
 
static G4PhysicsTabletheProperTimePositronTable = 0
 
static G4int NbOfProcesses = 2
 
static G4int CounterOfElectronProcess = 0
 
static G4int CounterOfPositronProcess = 0
 
static G4PhysicsTable ** RecorderOfElectronProcess
 
static G4PhysicsTable ** RecorderOfPositronProcess
 
- Static Protected Attributes inherited from G4RDVeLowEnergyLoss
static G4double ParticleMass
 
static G4double taulow
 
static G4double tauhigh
 
static G4double ltaulow
 
static G4double ltauhigh
 
static G4double dRoverRange = 20*perCent
 
static G4double finalRange = 200*micrometer
 
static G4double c1lim = dRoverRange
 
static G4double c2lim = 2.*(1.-dRoverRange)*finalRange
 
static G4double c3lim = -(1.-dRoverRange)*finalRange*finalRange
 
static G4bool rndmStepFlag = false
 
static G4bool EnlossFlucFlag = true
 

Detailed Description

Definition at line 70 of file G4LowEnergyBremsstrahlung.hh.

Constructor & Destructor Documentation

G4LowEnergyBremsstrahlung::G4LowEnergyBremsstrahlung ( const G4String processName = "LowEnBrem")

Definition at line 85 of file G4LowEnergyBremsstrahlung.cc.

References G4VProcess::verboseLevel.

86  : G4eLowEnergyLoss(nam),
87  crossSectionHandler(0),
88  theMeanFreePath(0),
89  energySpectrum(0)
90 {
91  cutForPhotons = 0.;
92  verboseLevel = 0;
93  generatorName = "tsai";
94  angularDistribution = new G4RDModifiedTsai("TsaiGenerator"); // default generator
95 // angularDistribution->PrintGeneratorInformation();
96  TsaiAngularDistribution = new G4RDModifiedTsai("TsaiGenerator");
97 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4eLowEnergyLoss(const G4String &)
G4LowEnergyBremsstrahlung::~G4LowEnergyBremsstrahlung ( )

Definition at line 116 of file G4LowEnergyBremsstrahlung.cc.

117 {
118  if(crossSectionHandler) delete crossSectionHandler;
119  if(energySpectrum) delete energySpectrum;
120  if(theMeanFreePath) delete theMeanFreePath;
121  delete angularDistribution;
122  delete TsaiAngularDistribution;
123  energyBins.clear();
124 }

Member Function Documentation

void G4LowEnergyBremsstrahlung::BuildPhysicsTable ( const G4ParticleDefinition particleType)
virtual

Reimplemented from G4VProcess.

Definition at line 127 of file G4LowEnergyBremsstrahlung.cc.

References G4eLowEnergyLoss::BuildDEDXTable(), G4eLowEnergyLoss::CounterOfElectronProcess, G4eLowEnergyLoss::CounterOfPositronProcess, G4Electron::Electron(), G4cout, G4endl, G4eLowEnergyLoss::GetLowerBoundEloss(), G4eLowEnergyLoss::GetNbinEloss(), G4VProcess::GetProcessName(), G4eLowEnergyLoss::GetUpperBoundEloss(), G4RDVCrossSectionHandler::Initialise(), G4RDVCrossSectionHandler::LoadShellData(), G4RDVEnergySpectrum::PrintData(), G4RDVCrossSectionHandler::PrintData(), PrintInfoDefinition(), G4eLowEnergyLoss::RecorderOfElectronProcess, G4eLowEnergyLoss::RecorderOfPositronProcess, G4VProcess::verboseLevel, and test::x.

128 {
129  if(verboseLevel > 0) {
130  G4cout << "G4LowEnergyBremsstrahlung::BuildPhysicsTable start"
131  << G4endl;
132  }
133 
134  cutForSecondaryPhotons.clear();
135 
136  // Create and fill BremsstrahlungParameters once
137  if( energySpectrum != 0 ) delete energySpectrum;
138  energyBins.clear();
139  for(size_t i=0; i<15; i++) {
140  G4double x = 0.1*((G4double)i);
141  if(i == 0) x = 0.01;
142  if(i == 10) x = 0.95;
143  if(i == 11) x = 0.97;
144  if(i == 12) x = 0.99;
145  if(i == 13) x = 0.995;
146  if(i == 14) x = 1.0;
147  energyBins.push_back(x);
148  }
149  const G4String dataName("/brem/br-sp.dat");
150  energySpectrum = new G4RDeBremsstrahlungSpectrum(energyBins,dataName);
151 
152  if(verboseLevel > 0) {
153  G4cout << "G4LowEnergyBremsstrahlungSpectrum is initialized"
154  << G4endl;
155  }
156 
157  // Create and fill G4RDCrossSectionHandler once
158 
159  if( crossSectionHandler != 0 ) delete crossSectionHandler;
160  G4RDVDataSetAlgorithm* interpolation = new G4RDLogLogInterpolation();
161  G4double lowKineticEnergy = GetLowerBoundEloss();
162  G4double highKineticEnergy = GetUpperBoundEloss();
163  G4int totBin = GetNbinEloss();
164  crossSectionHandler = new G4RDBremsstrahlungCrossSectionHandler(energySpectrum, interpolation);
165  crossSectionHandler->Initialise(0,lowKineticEnergy, highKineticEnergy, totBin);
166  crossSectionHandler->LoadShellData("brem/br-cs-");
167 
168  if (verboseLevel > 0) {
170  << " is created; Cross section data: "
171  << G4endl;
172  crossSectionHandler->PrintData();
173  G4cout << "Parameters: "
174  << G4endl;
175  energySpectrum->PrintData();
176  }
177 
178  // Build loss table for Bremsstrahlung
179 
180  BuildLossTable(aParticleType);
181 
182  if(verboseLevel > 0) {
183  G4cout << "The loss table is built"
184  << G4endl;
185  }
186 
187  if (&aParticleType==G4Electron::Electron()) {
188 
189  RecorderOfElectronProcess[CounterOfElectronProcess] = (*this).theLossTable;
192 
193  } else {
194 
195  RecorderOfPositronProcess[CounterOfPositronProcess] = (*this).theLossTable;
197  }
198 
199  // Build mean free path data using cut values
200 
201  if( theMeanFreePath != 0 ) delete theMeanFreePath;
202  theMeanFreePath = crossSectionHandler->
203  BuildMeanFreePathForMaterials(&cutForSecondaryPhotons);
204 
205  if(verboseLevel > 0) {
206  G4cout << "The MeanFreePath table is built"
207  << G4endl;
208  }
209 
210  // Build common DEDX table for all ionisation processes
211 
212  BuildDEDXTable(aParticleType);
213 
214  if(verboseLevel > 0) {
215  G4cout << "G4LowEnergyBremsstrahlung::BuildPhysicsTable end"
216  << G4endl;
217  }
218 
219 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static G4int CounterOfElectronProcess
void Initialise(G4RDVDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
void LoadShellData(const G4String &dataFile)
static G4double GetLowerBoundEloss()
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void BuildDEDXTable(const G4ParticleDefinition &aParticleType)
virtual void PrintData() const =0
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4PhysicsTable ** RecorderOfPositronProcess
static G4double GetUpperBoundEloss()
static G4PhysicsTable ** RecorderOfElectronProcess
static G4int GetNbinEloss()
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static G4int CounterOfPositronProcess
double G4double
Definition: G4Types.hh:76
G4double G4LowEnergyBremsstrahlung::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Implements G4eLowEnergyLoss.

Definition at line 387 of file G4LowEnergyBremsstrahlung.cc.

References G4RDVEMDataSet::FindValue(), G4RDVEMDataSet::GetComponent(), G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), and NotForced.

390 {
391  *cond = NotForced;
392  G4int index = (track.GetMaterialCutsCouple())->GetIndex();
393  const G4RDVEMDataSet* data = theMeanFreePath->GetComponent(index);
394  G4double meanFreePath = data->FindValue(track.GetKineticEnergy());
395  return meanFreePath;
396 }
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual const G4RDVEMDataSet * GetComponent(G4int componentId) const =0
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
double G4double
Definition: G4Types.hh:76
const XML_Char const XML_Char * data
G4bool G4LowEnergyBremsstrahlung::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4eLowEnergyLoss.

Definition at line 381 of file G4LowEnergyBremsstrahlung.cc.

References G4Electron::Electron().

382 {
383  return ( (&particle == G4Electron::Electron()) );
384 }
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4VParticleChange * G4LowEnergyBremsstrahlung::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4eLowEnergyLoss.

Definition at line 302 of file G4LowEnergyBremsstrahlung.cc.

References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, python.hepunit::electron_mass_c2, G4UniformRand, G4Gamma::Gamma(), G4MaterialCutsCouple::GetIndex(), G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), G4Track::GetMomentumDirection(), G4ParticleChange::Initialize(), python.hepunit::keV, G4RDVBremAngularDistribution::PolarAngle(), G4VContinuousDiscreteProcess::PostStepDoIt(), G4ParticleChange::ProposeEnergy(), G4ParticleChange::ProposeMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), G4RDVEnergySpectrum::SampleEnergy(), G4RDVCrossSectionHandler::SelectRandomAtom(), G4VParticleChange::SetNumberOfSecondaries(), python.hepunit::twopi, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

304 {
306 
307  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
308  G4double kineticEnergy = track.GetKineticEnergy();
309  G4int index = couple->GetIndex();
310  G4double tCut = cutForSecondaryPhotons[index];
311 
312  // Control limits
313  if(tCut >= kineticEnergy)
314  return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
315 
316  G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy);
317 
318  G4double tGamma = energySpectrum->SampleEnergy(Z, tCut, kineticEnergy, kineticEnergy);
319  G4double totalEnergy = kineticEnergy + electron_mass_c2;
320  G4double finalEnergy = kineticEnergy - tGamma; // electron/positron final energy
321  G4double theta = 0;
322 
323  if((kineticEnergy < 1*MeV && kineticEnergy > 1*keV && generatorName == "2bn")){
324  theta = angularDistribution->PolarAngle(kineticEnergy,finalEnergy,Z);
325  }else{
326  theta = TsaiAngularDistribution->PolarAngle(kineticEnergy,finalEnergy,Z);
327  }
328 
329  G4double phi = twopi * G4UniformRand();
330  G4double dirZ = std::cos(theta);
331  G4double sinTheta = std::sqrt(1. - dirZ*dirZ);
332  G4double dirX = sinTheta*std::cos(phi);
333  G4double dirY = sinTheta*std::sin(phi);
334 
335  G4ThreeVector gammaDirection (dirX, dirY, dirZ);
336  G4ThreeVector electronDirection = track.GetMomentumDirection();
337 
338  //
339  // Update the incident particle
340  //
341  gammaDirection.rotateUz(electronDirection);
342 
343  // Kinematic problem
344  if (finalEnergy < 0.) {
345  tGamma += finalEnergy;
346  finalEnergy = 0.0;
347  }
348 
349  G4double momentum = std::sqrt((totalEnergy + electron_mass_c2)*kineticEnergy);
350 
351  G4double finalX = momentum*electronDirection.x() - tGamma*gammaDirection.x();
352  G4double finalY = momentum*electronDirection.y() - tGamma*gammaDirection.y();
353  G4double finalZ = momentum*electronDirection.z() - tGamma*gammaDirection.z();
354 
356  G4double norm = 1./std::sqrt(finalX*finalX + finalY*finalY + finalZ*finalZ);
357  aParticleChange.ProposeMomentumDirection(finalX*norm, finalY*norm, finalZ*norm);
358  aParticleChange.ProposeEnergy( finalEnergy );
359 
360  // create G4DynamicParticle object for the gamma
362  gammaDirection, tGamma);
364 
365  return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
366 }
virtual G4double PolarAngle(const G4double initial_energy, const G4double final_energy, const G4int Z)=0
double x() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
int G4int
Definition: G4Types.hh:78
double z() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:87
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
float electron_mass_c2
Definition: hepunit.py:274
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
virtual void Initialize(const G4Track &)
const G4ThreeVector & GetMomentumDirection() const
void SetNumberOfSecondaries(G4int totSecondaries)
double y() const
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void AddSecondary(G4Track *aSecondary)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
virtual G4double SampleEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
void G4LowEnergyBremsstrahlung::PrintInfoDefinition ( )

Definition at line 369 of file G4LowEnergyBremsstrahlung.cc.

References G4cout, G4endl, and G4VProcess::GetProcessName().

Referenced by BuildPhysicsTable().

370 {
371  G4String comments = "Total cross sections from EEDL database.";
372  comments += "\n Gamma energy sampled from a parameterised formula.";
373  comments += "\n Implementation of the continuous dE/dx part.";
374  comments += "\n At present it can be used for electrons ";
375  comments += "in the energy range [250eV,100GeV].";
376  comments += "\n The process must work with G4LowEnergyIonisation.";
377 
378  G4cout << G4endl << GetProcessName() << ": " << comments << G4endl;
379 }
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
void G4LowEnergyBremsstrahlung::SetAngularGenerator ( G4RDVBremAngularDistribution distribution)

Definition at line 403 of file G4LowEnergyBremsstrahlung.cc.

References G4RDVBremAngularDistribution::PrintGeneratorInformation().

404 {
405  angularDistribution = distribution;
406  angularDistribution->PrintGeneratorInformation();
407 }
virtual void PrintGeneratorInformation() const =0
void G4LowEnergyBremsstrahlung::SetAngularGenerator ( const G4String name)

Definition at line 409 of file G4LowEnergyBremsstrahlung.cc.

References FatalException, G4Exception(), and G4RDVBremAngularDistribution::PrintGeneratorInformation().

410 {
411  if (name == "tsai")
412  {
413  delete angularDistribution;
414  angularDistribution = new G4RDModifiedTsai("TsaiGenerator");
415  generatorName = name;
416  }
417  else if (name == "2bn")
418  {
419  delete angularDistribution;
420  angularDistribution = new G4RDGenerator2BN("2BNGenerator");
421  generatorName = name;
422  }
423  else if (name == "2bs")
424  {
425  delete angularDistribution;
426  angularDistribution = new G4RDGenerator2BS("2BSGenerator");
427  generatorName = name;
428  }
429  else
430  {
431  G4Exception("G4LowEnergyBremsstrahlung::SetAngularGenerator()",
432  "InvalidSetup", FatalException, "Generator does not exist!");
433  }
434 
435  angularDistribution->PrintGeneratorInformation();
436 }
const XML_Char * name
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void PrintGeneratorInformation() const =0
void G4LowEnergyBremsstrahlung::SetCutForLowEnSecPhotons ( G4double  cut)

Definition at line 398 of file G4LowEnergyBremsstrahlung.cc.

399 {
400  cutForPhotons = cut;
401 }

The documentation for this class was generated from the following files: